New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
dommsk.F90 in branches/NERC/dev_r5840_BDY_MSK/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/NERC/dev_r5840_BDY_MSK/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90 @ 6950

Last change on this file since 6950 was 6950, checked in by jamesharle, 8 years ago

update mppini_2.h90 to account for the bdy_msk (if used) to allow full land supression

  • Property svn:keywords set to Id
File size: 24.3 KB
Line 
1MODULE dommsk
2   !!======================================================================
3   !!                       ***  MODULE dommsk   ***
4   !! Ocean initialization : domain land/sea mask
5   !!======================================================================
6   !! History :  OPA  ! 1987-07  (G. Madec)  Original code
7   !!            6.0  ! 1993-03  (M. Guyon)  symetrical conditions (M. Guyon)
8   !!            7.0  ! 1996-01  (G. Madec)  suppression of common work arrays
9   !!             -   ! 1996-05  (G. Madec)  mask computed from tmask and sup-
10   !!                 !                      pression of the double computation of bmask
11   !!            8.0  ! 1997-02  (G. Madec)  mesh information put in domhgr.F
12   !!            8.1  ! 1997-07  (G. Madec)  modification of mbathy and fmask
13   !!             -   ! 1998-05  (G. Roullet)  free surface
14   !!            8.2  ! 2000-03  (G. Madec)  no slip accurate
15   !!             -   ! 2001-09  (J.-M. Molines)  Open boundaries
16   !!   NEMO     1.0  ! 2002-08  (G. Madec)  F90: Free form and module
17   !!             -   ! 2005-11  (V. Garnier) Surface pressure gradient organization
18   !!            3.2  ! 2009-07  (R. Benshila) Suppression of rigid-lid option
19   !!            3.6  ! 2015-05  (P. Mathiot) ISF: add wmask,wumask and wvmask
20   !!----------------------------------------------------------------------
21
22   !!----------------------------------------------------------------------
23   !!   dom_msk        : compute land/ocean mask
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and tracers
26   USE dom_oce         ! ocean space and time domain
27   USE bdy_oce     
28   USE in_out_manager  ! I/O manager
29   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
30   USE lib_mpp
31   USE iom
32   USE dynspg_oce      ! choice/control of key cpp for surface pressure gradient
33   USE wrk_nemo        ! Memory allocation
34   USE timing          ! Timing
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   dom_msk         ! routine called by inidom.F90
40
41   !                            !!* Namelist namlbc : lateral boundary condition *
42   REAL(wp)        :: rn_shlat   ! type of lateral boundary condition on velocity
43   LOGICAL, PUBLIC :: ln_vorlat  !  consistency of vorticity boundary condition
44   !                                            with analytical eqs.
45
46   !! * Substitutions
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/OPA 3.2 , LODYC-IPSL  (2009)
50   !! $Id$
51   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
52   !!----------------------------------------------------------------------
53CONTAINS
54
55   SUBROUTINE dom_msk
56      !!---------------------------------------------------------------------
57      !!                 ***  ROUTINE dom_msk  ***
58      !!
59      !! ** Purpose :   Compute land/ocean mask arrays at tracer points, hori-
60      !!      zontal velocity points (u & v), vorticity points (f) and baro-
61      !!      tropic stream function  points (b).
62      !!
63      !! ** Method  :   The ocean/land mask is computed from the basin bathy-
64      !!      metry in level (mbathy) which is defined or read in dommba.
65      !!      mbathy equals 0 over continental T-point
66      !!      and the number of ocean level over the ocean.
67      !!
68      !!      At a given position (ji,jj,jk) the ocean/land mask is given by:
69      !!      t-point : 0. IF mbathy( ji ,jj) =< 0
70      !!                1. IF mbathy( ji ,jj) >= jk
71      !!      u-point : 0. IF mbathy( ji ,jj)  or mbathy(ji+1, jj ) =< 0
72      !!                1. IF mbathy( ji ,jj) and mbathy(ji+1, jj ) >= jk.
73      !!      v-point : 0. IF mbathy( ji ,jj)  or mbathy( ji ,jj+1) =< 0
74      !!                1. IF mbathy( ji ,jj) and mbathy( ji ,jj+1) >= jk.
75      !!      f-point : 0. IF mbathy( ji ,jj)  or mbathy( ji ,jj+1)
76      !!                   or mbathy(ji+1,jj)  or mbathy(ji+1,jj+1) =< 0
77      !!                1. IF mbathy( ji ,jj) and mbathy( ji ,jj+1)
78      !!                  and mbathy(ji+1,jj) and mbathy(ji+1,jj+1) >= jk.
79      !!      b-point : the same definition as for f-point of the first ocean
80      !!                level (surface level) but with 0 along coastlines.
81      !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated
82      !!                rows/lines due to cyclic or North Fold boundaries as well
83      !!                as MPP halos.
84      !!
85      !!        The lateral friction is set through the value of fmask along
86      !!      the coast and topography. This value is defined by rn_shlat, a
87      !!      namelist parameter:
88      !!         rn_shlat = 0, free slip  (no shear along the coast)
89      !!         rn_shlat = 2, no slip  (specified zero velocity at the coast)
90      !!         0 < rn_shlat < 2, partial slip   | non-linear velocity profile
91      !!         2 < rn_shlat, strong slip        | in the lateral boundary layer
92      !!
93      !!      N.B. If nperio not equal to 0, the land/ocean mask arrays
94      !!      are defined with the proper value at lateral domain boundaries,
95      !!      but bmask. indeed, bmask defined the domain over which the
96      !!      barotropic stream function is computed. this domain cannot
97      !!      contain identical columns because the matrix associated with
98      !!      the barotropic stream function equation is then no more inverti-
99      !!      ble. therefore bmask is set to 0 along lateral domain boundaries
100      !!      even IF nperio is not zero.
101      !!
102      !!      In case of open boundaries (lk_bdy=T):
103      !!        - tmask is set to 1 on the points to be computed bay the open
104      !!          boundaries routines.
105      !!        - bmask is  set to 0 on the open boundaries.
106      !!
107      !! ** Action :   tmask    : land/ocean mask at t-point (=0. or 1.)
108      !!               umask    : land/ocean mask at u-point (=0. or 1.)
109      !!               vmask    : land/ocean mask at v-point (=0. or 1.)
110      !!               fmask    : land/ocean mask at f-point (=0. or 1.)
111      !!                          =rn_shlat along lateral boundaries
112      !!               bmask    : land/ocean mask at barotropic stream
113      !!                          function point (=0. or 1.) and set to 0 along lateral boundaries
114      !!               tmask_i  : interior ocean mask
115      !!----------------------------------------------------------------------
116      INTEGER  ::   ji, jj, jk               ! dummy loop indices
117      INTEGER  ::   iif, iil, ii0, ii1, ii   ! local integers
118      INTEGER  ::   ijf, ijl, ij0, ij1       !   -       -
119      INTEGER  ::   ios, inum
120      INTEGER  ::   isrow                    ! index for ORCA1 starting row
121      INTEGER , POINTER, DIMENSION(:,:) ::  imsk
122      REAL(wp), POINTER, DIMENSION(:,:) ::  zwf
123      !!
124      NAMELIST/namlbc/ rn_shlat, ln_vorlat
125#if defined key_bdy 
126      NAMELIST/nambdy/ nb_bdy, ln_coords_file, cn_coords_file,                 &
127         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
128         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
129         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
130         &             cn_ice_lim, nn_ice_lim_dta,                           &
131         &             rn_ice_tem, rn_ice_sal, rn_ice_age,                 &
132         &             ln_vol, nn_volctl, nn_rimwidth
133#endif
134      !!---------------------------------------------------------------------
135      !
136      IF( nn_timing == 1 )  CALL timing_start('dom_msk')
137      !
138      CALL wrk_alloc( jpi, jpj, imsk )
139      CALL wrk_alloc( jpi, jpj, zwf  )
140      !
141      REWIND( numnam_ref )              ! Namelist namlbc in reference namelist : Lateral momentum boundary condition
142      READ  ( numnam_ref, namlbc, IOSTAT = ios, ERR = 901 )
143901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namlbc in reference namelist', lwp )
144
145      REWIND( numnam_cfg )              ! Namelist namlbc in configuration namelist : Lateral momentum boundary condition
146      READ  ( numnam_cfg, namlbc, IOSTAT = ios, ERR = 902 )
147902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namlbc in configuration namelist', lwp )
148      IF(lwm) WRITE ( numond, namlbc )
149     
150      IF(lwp) THEN                  ! control print
151         WRITE(numout,*)
152         WRITE(numout,*) 'dommsk : ocean mask '
153         WRITE(numout,*) '~~~~~~'
154         WRITE(numout,*) '   Namelist namlbc'
155         WRITE(numout,*) '      lateral momentum boundary cond.    rn_shlat  = ',rn_shlat
156         WRITE(numout,*) '      consistency with analytical form   ln_vorlat = ',ln_vorlat 
157      ENDIF
158
159      IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  free-slip '
160      ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  no-slip '
161      ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  partial-slip '
162      ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  strong-slip '
163      ELSE
164         WRITE(ctmp1,*) ' rn_shlat is negative = ', rn_shlat
165         CALL ctl_stop( ctmp1 )
166      ENDIF
167
168      ! 1. Ocean/land mask at t-point (computed from mbathy)
169      ! -----------------------------
170      ! N.B. tmask has already the right boundary conditions since mbathy is ok
171      !
172      tmask(:,:,:) = 0._wp
173      DO jk = 1, jpk
174         DO jj = 1, jpj
175            DO ji = 1, jpi
176               IF( REAL( mbathy(ji,jj) - jk, wp ) + 0.1_wp >= 0._wp )   tmask(ji,jj,jk) = 1._wp
177            END DO 
178         END DO 
179      END DO 
180     
181#if defined key_bdy 
182      REWIND( numnam_ref )              ! Namelist nambdy in reference namelist :Unstructured open boundaries 
183      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 901)
184903   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in reference namelist', lwp )
185
186      REWIND( numnam_cfg )              ! Namelist nambdy in configuration namelist :Unstructured open boundaries
187      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 902 )
188904   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in configuration namelist', lwp )
189      IF(lwm) WRITE ( numond, nambdy )
190
191     IF( ln_mask_file ) THEN ! correct for bdy mask
192         CALL iom_open( cn_mask_file, inum )
193         CALL iom_get ( inum, jpdom_data, 'bdy_msk', bdytmask(:,:) )
194         CALL iom_close( inum )
195
196         ! Derive mask on U and V grid from mask on T grid
197         bdyumask(:,:) = 0.e0
198         bdyvmask(:,:) = 0.e0
199         DO jj=1, jpjm1
200            DO ji=1, jpim1
201               bdyumask(ji,jj)=bdytmask(ji,jj)*bdytmask(ji+1, jj )
202               bdyvmask(ji,jj)=bdytmask(ji,jj)*bdytmask(ji  ,jj+1)
203            END DO
204         END DO
205         CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
206
207
208         ! Mask corrections
209         ! ----------------
210         DO jk = 1, jpkm1
211            DO jj = 1, jpj
212               DO ji = 1, jpi
213                  tmask(ji,jj,jk) = tmask(ji,jj,jk) * bdytmask(ji,jj)
214               END DO
215            END DO
216         END DO
217      ENDIF
218#endif
219      ! (ISF) define barotropic mask and mask the ice shelf point
220      ssmask(:,:)=tmask(:,:,1) ! at this stage ice shelf is not masked
221     
222      DO jk = 1, jpk
223         DO jj = 1, jpj
224            DO ji = 1, jpi
225               IF( REAL( misfdep(ji,jj) - jk, wp ) - 0.1_wp >= 0._wp )   THEN
226                  tmask(ji,jj,jk) = 0._wp
227               END IF
228            END DO 
229         END DO 
230      END DO 
231
232      ! Interior domain mask (used for global sum)
233      ! --------------------
234      tmask_i(:,:) = ssmask(:,:)            ! (ISH) tmask_i = 1 even on the ice shelf
235      iif = jpreci                         ! ???
236      iil = nlci - jpreci + 1
237      ijf = jprecj                         ! ???
238      ijl = nlcj - jprecj + 1
239
240      tmask_i( 1 :iif,   :   ) = 0._wp      ! first columns
241      tmask_i(iil:jpi,   :   ) = 0._wp      ! last  columns (including mpp extra columns)
242      tmask_i(   :   , 1 :ijf) = 0._wp      ! first rows
243      tmask_i(   :   ,ijl:jpj) = 0._wp      ! last  rows (including mpp extra rows)
244
245      ! north fold mask
246      ! ---------------
247      tpol(1:jpiglo) = 1._wp 
248      fpol(1:jpiglo) = 1._wp
249      IF( jperio == 3 .OR. jperio == 4 ) THEN      ! T-point pivot
250         tpol(jpiglo/2+1:jpiglo) = 0._wp
251         fpol(     1    :jpiglo) = 0._wp
252         IF( mjg(nlej) == jpjglo ) THEN                  ! only half of the nlcj-1 row
253            DO ji = iif+1, iil-1
254               tmask_i(ji,nlej-1) = tmask_i(ji,nlej-1) * tpol(mig(ji))
255            END DO
256         ENDIF
257      ENDIF
258      IF( jperio == 5 .OR. jperio == 6 ) THEN      ! F-point pivot
259         tpol(     1    :jpiglo) = 0._wp
260         fpol(jpiglo/2+1:jpiglo) = 0._wp
261      ENDIF
262
263      ! 2. Ocean/land mask at u-,  v-, and z-points (computed from tmask)
264      ! -------------------------------------------
265      DO jk = 1, jpk
266         DO jj = 1, jpjm1
267            DO ji = 1, fs_jpim1   ! vector loop
268               umask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)
269               vmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji  ,jj+1,jk)
270            END DO
271            DO ji = 1, jpim1      ! NO vector opt.
272               fmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)   &
273                  &            * tmask(ji,jj+1,jk) * tmask(ji+1,jj+1,jk)
274            END DO
275         END DO
276      END DO
277      ! (ISF) MIN(1,SUM(umask)) is here to check if you have effectively at least 1 wet u point
278      DO jj = 1, jpjm1
279         DO ji = 1, fs_jpim1   ! vector loop
280            umask_i(ji,jj)  = ssmask(ji,jj) * ssmask(ji+1,jj  )  * MIN(1._wp,SUM(umask(ji,jj,:)))
281            vmask_i(ji,jj)  = ssmask(ji,jj) * ssmask(ji  ,jj+1)  * MIN(1._wp,SUM(vmask(ji,jj,:)))
282         END DO
283         DO ji = 1, jpim1      ! NO vector opt.
284            fmask_i(ji,jj) =  ssmask(ji,jj  ) * ssmask(ji+1,jj  )   &
285               &            * ssmask(ji,jj+1) * ssmask(ji+1,jj+1) * MIN(1._wp,SUM(fmask(ji,jj,:)))
286         END DO
287      END DO
288      CALL lbc_lnk( umask, 'U', 1._wp )      ! Lateral boundary conditions
289      CALL lbc_lnk( vmask, 'V', 1._wp )
290      CALL lbc_lnk( fmask, 'F', 1._wp )
291      CALL lbc_lnk( umask_i, 'U', 1._wp )      ! Lateral boundary conditions
292      CALL lbc_lnk( vmask_i, 'V', 1._wp )
293      CALL lbc_lnk( fmask_i, 'F', 1._wp )
294
295      ! 3. Ocean/land mask at wu-, wv- and w points
296      !----------------------------------------------
297      wmask (:,:,1) = tmask(:,:,1)     ! surface
298      wumask(:,:,1) = umask(:,:,1)
299      wvmask(:,:,1) = vmask(:,:,1)
300      DO jk = 2, jpk                   ! interior values
301         wmask (:,:,jk) = tmask(:,:,jk) * tmask(:,:,jk-1)
302         wumask(:,:,jk) = umask(:,:,jk) * umask(:,:,jk-1)   
303         wvmask(:,:,jk) = vmask(:,:,jk) * vmask(:,:,jk-1)
304      END DO
305
306      ! 4. ocean/land mask for the elliptic equation
307      ! --------------------------------------------
308      bmask(:,:) = ssmask(:,:)       ! elliptic equation is written at t-point
309      !
310      !                               ! Boundary conditions
311      !                                    ! cyclic east-west : bmask must be set to 0. on rows 1 and jpi
312      IF( nperio == 1 .OR. nperio == 4 .OR. nperio == 6 ) THEN
313         bmask( 1 ,:) = 0._wp
314         bmask(jpi,:) = 0._wp
315      ENDIF
316      IF( nperio == 2 ) THEN               ! south symmetric :  bmask must be set to 0. on row 1
317         bmask(:, 1 ) = 0._wp
318      ENDIF
319      !                                    ! north fold :
320      IF( nperio == 3 .OR. nperio == 4 ) THEN   ! T-pt pivot : bmask set to 0. on row jpj and on half jpjglo-1 row
321         DO ji = 1, jpi                     
322            ii = ji + nimpp - 1
323            bmask(ji,jpj-1) = bmask(ji,jpj-1) * tpol(ii)
324            bmask(ji,jpj  ) = 0._wp
325         END DO
326      ENDIF
327      IF( nperio == 5 .OR. nperio == 6 ) THEN   ! F-pt pivot and T-pt elliptic eq. : bmask set to 0. on row jpj
328         bmask(:,jpj) = 0._wp
329      ENDIF
330      !
331      IF( lk_mpp ) THEN                    ! mpp specificities
332         !                                      ! bmask is set to zero on the overlap region
333         IF( nbondi /= -1 .AND. nbondi /= 2 )   bmask(  1 :jpreci,:) = 0._wp
334         IF( nbondi /=  1 .AND. nbondi /= 2 )   bmask(nlci:jpi   ,:) = 0._wp
335         IF( nbondj /= -1 .AND. nbondj /= 2 )   bmask(:,  1 :jprecj) = 0._wp
336         IF( nbondj /=  1 .AND. nbondj /= 2 )   bmask(:,nlcj:jpj   ) = 0._wp
337         !
338         IF( npolj == 3 .OR. npolj == 4 ) THEN  ! north fold : bmask must be set to 0. on rows jpj-1 and jpj
339            DO ji = 1, nlci
340               ii = ji + nimpp - 1
341               bmask(ji,nlcj-1) = bmask(ji,nlcj-1) * tpol(ii)
342               bmask(ji,nlcj  ) = 0._wp
343            END DO
344         ENDIF
345         IF( npolj == 5 .OR. npolj == 6 ) THEN  ! F-pt pivot and T-pt elliptic eq. : bmask set to 0. on row jpj
346            DO ji = 1, nlci
347               bmask(ji,nlcj  ) = 0._wp
348            END DO
349         ENDIF
350      ENDIF
351
352      ! Lateral boundary conditions on velocity (modify fmask)
353      ! ---------------------------------------     
354      DO jk = 1, jpk
355         zwf(:,:) = fmask(:,:,jk)         
356         DO jj = 2, jpjm1
357            DO ji = fs_2, fs_jpim1   ! vector opt.
358               IF( fmask(ji,jj,jk) == 0._wp ) THEN
359                  fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1),   &
360                     &                                           zwf(ji-1,jj), zwf(ji,jj-1)  )  )
361               ENDIF
362            END DO
363         END DO
364         DO jj = 2, jpjm1
365            IF( fmask(1,jj,jk) == 0._wp ) THEN
366               fmask(1  ,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
367            ENDIF
368            IF( fmask(jpi,jj,jk) == 0._wp ) THEN
369               fmask(jpi,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
370            ENDIF
371         END DO         
372         DO ji = 2, jpim1
373            IF( fmask(ji,1,jk) == 0._wp ) THEN
374               fmask(ji, 1 ,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
375            ENDIF
376            IF( fmask(ji,jpj,jk) == 0._wp ) THEN
377               fmask(ji,jpj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
378            ENDIF
379         END DO
380      END DO
381      !
382      IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN   ! ORCA_R2 configuration
383         !                                                 ! Increased lateral friction near of some straits
384         !                                ! Gibraltar strait  : partial slip (fmask=0.5)
385         ij0 = 101   ;   ij1 = 101
386         ii0 = 139   ;   ii1 = 140   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  0.5_wp
387         ij0 = 102   ;   ij1 = 102
388         ii0 = 139   ;   ii1 = 140   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  0.5_wp
389         !
390         !                                ! Bab el Mandeb : partial slip (fmask=1)
391         ij0 =  87   ;   ij1 =  88
392         ii0 = 160   ;   ii1 = 160   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  1._wp
393         ij0 =  88   ;   ij1 =  88
394         ii0 = 159   ;   ii1 = 159   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  1._wp
395         !
396         !                                ! Danish straits  : strong slip (fmask > 2)
397! We keep this as an example but it is instable in this case
398!         ij0 = 115   ;   ij1 = 115
399!         ii0 = 145   ;   ii1 = 146   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) = 4._wp
400!         ij0 = 116   ;   ij1 = 116
401!         ii0 = 145   ;   ii1 = 146   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) = 4._wp
402         !
403      ENDIF
404      !
405      IF( cp_cfg == "orca" .AND. jp_cfg == 1 ) THEN   ! ORCA R1 configuration
406         !                                                 ! Increased lateral friction near of some straits
407         ! This dirty section will be suppressed by simplification process:
408         ! all this will come back in input files
409         ! Currently these hard-wired indices relate to configuration with
410         ! extend grid (jpjglo=332)
411         !
412         isrow = 332 - jpjglo
413         !
414         IF(lwp) WRITE(numout,*)
415         IF(lwp) WRITE(numout,*) '   orca_r1: increase friction near the following straits : '
416         IF(lwp) WRITE(numout,*) '      Gibraltar '
417         ii0 = 282           ;   ii1 = 283        ! Gibraltar Strait
418         ij0 = 241 - isrow   ;   ij1 = 241 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
419
420         IF(lwp) WRITE(numout,*) '      Bhosporus '
421         ii0 = 314           ;   ii1 = 315        ! Bhosporus Strait
422         ij0 = 248 - isrow   ;   ij1 = 248 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
423
424         IF(lwp) WRITE(numout,*) '      Makassar (Top) '
425         ii0 =  48           ;   ii1 =  48        ! Makassar Strait (Top)
426         ij0 = 189 - isrow   ;   ij1 = 190 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
427
428         IF(lwp) WRITE(numout,*) '      Lombok '
429         ii0 =  44           ;   ii1 =  44        ! Lombok Strait
430         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
431
432         IF(lwp) WRITE(numout,*) '      Ombai '
433         ii0 =  53           ;   ii1 =  53        ! Ombai Strait
434         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
435
436         IF(lwp) WRITE(numout,*) '      Timor Passage '
437         ii0 =  56           ;   ii1 =  56        ! Timor Passage
438         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
439
440         IF(lwp) WRITE(numout,*) '      West Halmahera '
441         ii0 =  58           ;   ii1 =  58        ! West Halmahera Strait
442         ij0 = 181 - isrow   ;   ij1 = 182 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
443
444         IF(lwp) WRITE(numout,*) '      East Halmahera '
445         ii0 =  55           ;   ii1 =  55        ! East Halmahera Strait
446         ij0 = 181 - isrow   ;   ij1 = 182 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
447         !
448      ENDIF
449      !
450      CALL lbc_lnk( fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask
451
452      ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 )
453           
454      IF( nprint == 1 .AND. lwp ) THEN      ! Control print
455         imsk(:,:) = INT( tmask_i(:,:) )
456         WRITE(numout,*) ' tmask_i : '
457         CALL prihin( imsk(:,:), jpi, jpj, 1, jpi, 1,   &
458               &                           1, jpj, 1, 1, numout)
459         WRITE (numout,*)
460         WRITE (numout,*) ' dommsk: tmask for each level'
461         WRITE (numout,*) ' ----------------------------'
462         DO jk = 1, jpk
463            imsk(:,:) = INT( tmask(:,:,jk) )
464
465            WRITE(numout,*)
466            WRITE(numout,*) ' level = ',jk
467            CALL prihin( imsk(:,:), jpi, jpj, 1, jpi, 1,   &
468               &                              1, jpj, 1, 1, numout)
469         END DO
470         WRITE(numout,*)
471         WRITE(numout,*) ' dom_msk: vmask for each level'
472         WRITE(numout,*) ' -----------------------------'
473         DO jk = 1, jpk
474            imsk(:,:) = INT( vmask(:,:,jk) )
475            WRITE(numout,*)
476            WRITE(numout,*) ' level = ',jk
477            CALL prihin( imsk(:,:), jpi, jpj, 1, jpi, 1,   &
478               &                              1, jpj, 1, 1, numout)
479         END DO
480         WRITE(numout,*)
481         WRITE(numout,*) ' dom_msk: fmask for each level'
482         WRITE(numout,*) ' -----------------------------'
483         DO jk = 1, jpk
484            imsk(:,:) = INT( fmask(:,:,jk) )
485            WRITE(numout,*)
486            WRITE(numout,*) ' level = ',jk
487            CALL prihin( imsk(:,:), jpi, jpj, 1, jpi, 1,   &
488               &                              1, jpj, 1, 1, numout )
489         END DO
490         WRITE(numout,*)
491         WRITE(numout,*) ' dom_msk: bmask '
492         WRITE(numout,*) ' ---------------'
493         WRITE(numout,*)
494         imsk(:,:) = INT( bmask(:,:) )
495         CALL prihin( imsk(:,:), jpi, jpj, 1, jpi, 1,   &
496            &                              1, jpj, 1, 1, numout )
497      ENDIF
498      !
499      CALL wrk_dealloc( jpi, jpj, imsk )
500      CALL wrk_dealloc( jpi, jpj, zwf  )
501      !
502      IF( nn_timing == 1 )  CALL timing_stop('dom_msk')
503      !
504   END SUBROUTINE dom_msk
505   
506   !!======================================================================
507END MODULE dommsk
Note: See TracBrowser for help on using the repository browser.