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/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90 @ 5956

Last change on this file since 5956 was 5956, checked in by mathiot, 8 years ago

ISF : merged trunk (5936) into branch

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