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/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90 @ 9168

Last change on this file since 9168 was 9168, checked in by gm, 6 years ago

dev_merge_2017: OPA_SRC & CONFIG: remove useless warning when reading namelist_cfg

  • Property svn:keywords set to Id
File size: 14.4 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
[6140]9   !!             -   ! 1996-05  (G. Madec)  mask computed from tmask
[2528]10   !!            8.0  ! 1997-02  (G. Madec)  mesh information put in domhgr.F
[7646]11   !!            8.1  ! 1997-07  (G. Madec)  modification of kbat and fmask
[1566]12   !!             -   ! 1998-05  (G. Roullet)  free surface
[2528]13   !!            8.2  ! 2000-03  (G. Madec)  no slip accurate
[1566]14   !!             -   ! 2001-09  (J.-M. Molines)  Open boundaries
15   !!   NEMO     1.0  ! 2002-08  (G. Madec)  F90: Free form and module
16   !!             -   ! 2005-11  (V. Garnier) Surface pressure gradient organization
17   !!            3.2  ! 2009-07  (R. Benshila) Suppression of rigid-lid option
[5836]18   !!            3.6  ! 2015-05  (P. Mathiot) ISF: add wmask,wumask and wvmask
[7646]19   !!            4.0  ! 2016-06  (G. Madec, S. Flavoni)  domain configuration / user defined interface
[1566]20   !!----------------------------------------------------------------------
[3]21
22   !!----------------------------------------------------------------------
[7646]23   !!   dom_msk       : compute land/ocean mask
[3]24   !!----------------------------------------------------------------------
[7646]25   USE oce            ! ocean dynamics and tracers
26   USE dom_oce        ! ocean space and time domain
27   USE usrdef_fmask   ! user defined fmask
[9124]28   USE bdy_oce        ! open boundary
29   !
[7646]30   USE in_out_manager ! I/O manager
31   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
32   USE lib_mpp        ! Massively Parallel Processing library
[3]33
34   IMPLICIT NONE
35   PRIVATE
36
[6140]37   PUBLIC   dom_msk    ! routine called by inidom.F90
[3]38
[1601]39   !                            !!* Namelist namlbc : lateral boundary condition *
[4147]40   REAL(wp)        :: rn_shlat   ! type of lateral boundary condition on velocity
41   LOGICAL, PUBLIC :: ln_vorlat  !  consistency of vorticity boundary condition
[3294]42   !                                            with analytical eqs.
[2715]43
[3]44   !! * Substitutions
45#  include "vectopt_loop_substitute.h90"
[1566]46   !!----------------------------------------------------------------------
47   !! NEMO/OPA 3.2 , LODYC-IPSL  (2009)
[7753]48   !! $Id$
[2528]49   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[1566]50   !!----------------------------------------------------------------------
[3]51CONTAINS
[2715]52
[7646]53   SUBROUTINE dom_msk( k_top, k_bot )
[3]54      !!---------------------------------------------------------------------
55      !!                 ***  ROUTINE dom_msk  ***
56      !!
57      !! ** Purpose :   Compute land/ocean mask arrays at tracer points, hori-
[6125]58      !!      zontal velocity points (u & v), vorticity points (f) points.
[3]59      !!
[7646]60      !! ** Method  :   The ocean/land mask  at t-point is deduced from ko_top
61      !!      and ko_bot, the indices of the fist and last ocean t-levels which
62      !!      are either defined in usrdef_zgr or read in zgr_read.
63      !!                The velocity masks (umask, vmask, wmask, wumask, wvmask)
64      !!      are deduced from a product of the two neighboring tmask.
65      !!                The vorticity mask (fmask) is deduced from tmask taking
66      !!      into account the choice of lateral boundary condition (rn_shlat) :
[1601]67      !!         rn_shlat = 0, free slip  (no shear along the coast)
68      !!         rn_shlat = 2, no slip  (specified zero velocity at the coast)
69      !!         0 < rn_shlat < 2, partial slip   | non-linear velocity profile
70      !!         2 < rn_shlat, strong slip        | in the lateral boundary layer
[3]71      !!
[7646]72      !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated
73      !!                rows/lines due to cyclic or North Fold boundaries as well
74      !!                as MPP halos.
75      !!      tmask_h : halo mask at t-point, i.e. excluding duplicated rows/lines
76      !!                due to cyclic or North Fold boundaries as well as MPP halos.
[3]77      !!
[7646]78      !! ** Action :   tmask, umask, vmask, wmask, wumask, wvmask : land/ocean mask
79      !!                         at t-, u-, v- w, wu-, and wv-points (=0. or 1.)
80      !!               fmask   : land/ocean mask at f-point (=0., or =1., or
81      !!                         =rn_shlat along lateral boundaries)
82      !!               tmask_i : interior ocean mask
83      !!               tmask_h : halo mask
84      !!               ssmask , ssumask, ssvmask, ssfmask : 2D ocean mask
[3]85      !!----------------------------------------------------------------------
[7646]86      INTEGER, DIMENSION(:,:), INTENT(in) ::   k_top, k_bot   ! first and last ocean level
87      !
88      INTEGER  ::   ji, jj, jk     ! dummy loop indices
89      INTEGER  ::   iif, iil       ! local integers
90      INTEGER  ::   ijf, ijl       !   -       -
91      INTEGER  ::   iktop, ikbot   !   -       -
92      INTEGER  ::   ios, inum
[9019]93      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zwf   ! 2D workspace
[1601]94      !!
[3294]95      NAMELIST/namlbc/ rn_shlat, ln_vorlat
[7646]96      NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file,         &
97         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
98         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
99         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
100         &             cn_ice_lim, nn_ice_lim_dta,                           &
101         &             rn_ice_tem, rn_ice_sal, rn_ice_age,                 &
102         &             ln_vol, nn_volctl, nn_rimwidth, nb_jpk_bdy
[3]103      !!---------------------------------------------------------------------
[3294]104      !
[4147]105      REWIND( numnam_ref )              ! Namelist namlbc in reference namelist : Lateral momentum boundary condition
106      READ  ( numnam_ref, namlbc, IOSTAT = ios, ERR = 901 )
[9168]107901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namlbc in reference namelist', lwp )
[4147]108      REWIND( numnam_cfg )              ! Namelist namlbc in configuration namelist : Lateral momentum boundary condition
109      READ  ( numnam_cfg, namlbc, IOSTAT = ios, ERR = 902 )
[9168]110902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namlbc in configuration namelist', lwp )
[4624]111      IF(lwm) WRITE ( numond, namlbc )
[1566]112     
113      IF(lwp) THEN                  ! control print
[3]114         WRITE(numout,*)
115         WRITE(numout,*) 'dommsk : ocean mask '
116         WRITE(numout,*) '~~~~~~'
[1566]117         WRITE(numout,*) '   Namelist namlbc'
[3294]118         WRITE(numout,*) '      lateral momentum boundary cond.    rn_shlat  = ',rn_shlat
119         WRITE(numout,*) '      consistency with analytical form   ln_vorlat = ',ln_vorlat 
[3]120      ENDIF
121
[2528]122      IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  free-slip '
[1601]123      ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  no-slip '
124      ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  partial-slip '
125      ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  strong-slip '
126      ELSE
127         WRITE(ctmp1,*) ' rn_shlat is negative = ', rn_shlat
128         CALL ctl_stop( ctmp1 )
[3]129      ENDIF
130
[7646]131
132      !  Ocean/land mask at t-point  (computed from ko_top and ko_bot)
133      ! ----------------------------
[1566]134      !
[7753]135      tmask(:,:,:) = 0._wp
[7646]136      DO jj = 1, jpj
137         DO ji = 1, jpi
138            iktop = k_top(ji,jj)
139            ikbot = k_bot(ji,jj)
140            IF( iktop /= 0 ) THEN       ! water in the column
141               tmask(ji,jj,iktop:ikbot  ) = 1._wp
142            ENDIF
[3]143         END DO 
144      END DO 
[7646]145!SF  add here lbc_lnk: bug not still understood : cause now domain configuration is read !
146!!gm I don't understand why... 
147      CALL lbc_lnk( tmask  , 'T', 1._wp )      ! Lateral boundary conditions
[3]148
[7646]149     ! Mask corrections for bdy (read in mppini2)
150      REWIND( numnam_ref )              ! Namelist nambdy in reference namelist :Unstructured open boundaries
151      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 903)
[9168]152903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy in reference namelist', lwp )
[7646]153      REWIND( numnam_cfg )              ! Namelist nambdy in configuration namelist :Unstructured open boundaries
154      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 904 )
[9168]155904   IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy in configuration namelist', lwp )
[7646]156      ! ------------------------
157      IF ( ln_bdy .AND. ln_mask_file ) THEN
158         DO jk = 1, jpkm1
159            DO jj = 1, jpj
160               DO ji = 1, jpi
161                  tmask(ji,jj,jk) = tmask(ji,jj,jk) * bdytmask(ji,jj)
162               END DO
[291]163            END DO
[7646]164         END DO
[3]165      ENDIF
[7646]166         
167      !  Ocean/land mask at u-, v-, and f-points   (computed from tmask)
168      ! ----------------------------------------
169      ! NB: at this point, fmask is designed for free slip lateral boundary condition
[3]170      DO jk = 1, jpk
171         DO jj = 1, jpjm1
172            DO ji = 1, fs_jpim1   ! vector loop
173               umask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)
174               vmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji  ,jj+1,jk)
[1271]175            END DO
[1694]176            DO ji = 1, jpim1      ! NO vector opt.
[3]177               fmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)   &
[62]178                  &            * tmask(ji,jj+1,jk) * tmask(ji+1,jj+1,jk)
[3]179            END DO
180         END DO
181      END DO
[9090]182      CALL lbc_lnk_multi( umask, 'U', 1., vmask, 'V', 1., fmask, 'F', 1. )      ! Lateral boundary conditions
[7646]183 
184      ! Ocean/land mask at wu-, wv- and w points    (computed from tmask)
185      !-----------------------------------------
[7753]186      wmask (:,:,1) = tmask(:,:,1)     ! surface
187      wumask(:,:,1) = umask(:,:,1)
188      wvmask(:,:,1) = vmask(:,:,1)
[5836]189      DO jk = 2, jpk                   ! interior values
[7753]190         wmask (:,:,jk) = tmask(:,:,jk) * tmask(:,:,jk-1)
191         wumask(:,:,jk) = umask(:,:,jk) * umask(:,:,jk-1)   
192         wvmask(:,:,jk) = vmask(:,:,jk) * vmask(:,:,jk-1)
[5120]193      END DO
[3]194
[7646]195
196      ! Ocean/land column mask at t-, u-, and v-points   (i.e. at least 1 wet cell in the vertical)
197      ! ----------------------------------------------
198      ssmask (:,:) = MAXVAL( tmask(:,:,:), DIM=3 )
199      ssumask(:,:) = MAXVAL( umask(:,:,:), DIM=3 )
200      ssvmask(:,:) = MAXVAL( vmask(:,:,:), DIM=3 )
201
202
203      ! Interior domain mask  (used for global sum)
204      ! --------------------
205      !
[9019]206      iif = nn_hls   ;   iil = nlci - nn_hls + 1
207      ijf = nn_hls   ;   ijl = nlcj - nn_hls + 1
[7646]208      !
209      !                          ! halo mask : 0 on the halo and 1 elsewhere
[7753]210      tmask_h(:,:) = 1._wp                 
[7646]211      tmask_h( 1 :iif,   :   ) = 0._wp      ! first columns
212      tmask_h(iil:jpi,   :   ) = 0._wp      ! last  columns (including mpp extra columns)
213      tmask_h(   :   , 1 :ijf) = 0._wp      ! first rows
214      tmask_h(   :   ,ijl:jpj) = 0._wp      ! last  rows (including mpp extra rows)
215      !
216      !                          ! north fold mask
217      tpol(1:jpiglo) = 1._wp 
218      fpol(1:jpiglo) = 1._wp
219      IF( jperio == 3 .OR. jperio == 4 ) THEN      ! T-point pivot
220         tpol(jpiglo/2+1:jpiglo) = 0._wp
221         fpol(     1    :jpiglo) = 0._wp
222         IF( mjg(nlej) == jpjglo ) THEN                  ! only half of the nlcj-1 row for tmask_h
223            DO ji = iif+1, iil-1
224               tmask_h(ji,nlej-1) = tmask_h(ji,nlej-1) * tpol(mig(ji))
225            END DO
226         ENDIF
227      ENDIF
228      !
229      IF( jperio == 5 .OR. jperio == 6 ) THEN      ! F-point pivot
230         tpol(     1    :jpiglo) = 0._wp
231         fpol(jpiglo/2+1:jpiglo) = 0._wp
232      ENDIF
233      !
234      !                          ! interior mask : 2D ocean mask x halo mask
[7753]235      tmask_i(:,:) = ssmask(:,:) * tmask_h(:,:)
[7646]236
237
[3]238      ! Lateral boundary conditions on velocity (modify fmask)
[7646]239      ! --------------------------------------- 
240      IF( rn_shlat /= 0 ) THEN      ! Not free-slip lateral boundary condition
241         !
[9019]242         ALLOCATE( zwf(jpi,jpj) )
[7646]243         !
244         DO jk = 1, jpk
[7753]245            zwf(:,:) = fmask(:,:,jk)         
[7646]246            DO jj = 2, jpjm1
247               DO ji = fs_2, fs_jpim1   ! vector opt.
248                  IF( fmask(ji,jj,jk) == 0._wp ) THEN
249                     fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1),   &
250                        &                                           zwf(ji-1,jj), zwf(ji,jj-1)  )  )
251                  ENDIF
252               END DO
253            END DO
254            DO jj = 2, jpjm1
255               IF( fmask(1,jj,jk) == 0._wp ) THEN
256                  fmask(1  ,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
[3]257               ENDIF
[7646]258               IF( fmask(jpi,jj,jk) == 0._wp ) THEN
259                  fmask(jpi,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
260               ENDIF
261            END DO         
262            DO ji = 2, jpim1
263               IF( fmask(ji,1,jk) == 0._wp ) THEN
264                  fmask(ji, 1 ,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
265               ENDIF
266               IF( fmask(ji,jpj,jk) == 0._wp ) THEN
267                  fmask(ji,jpj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
268               ENDIF
[3]269            END DO
[9023]270#if defined key_agrif 
[9090]271            IF( .NOT. AGRIF_Root() ) THEN
272               IF ((nbondi ==  1).OR.(nbondi == 2)) fmask(nlci-1 , :     ,jk) = 0.e0      ! east
273               IF ((nbondi == -1).OR.(nbondi == 2)) fmask(1      , :     ,jk) = 0.e0      ! west
274               IF ((nbondj ==  1).OR.(nbondj == 2)) fmask(:      ,nlcj-1 ,jk) = 0.e0      ! north
275               IF ((nbondj == -1).OR.(nbondj == 2)) fmask(:      ,1      ,jk) = 0.e0      ! south
276            ENDIF 
[9023]277#endif
[3]278         END DO
[5836]279         !
[9019]280         DEALLOCATE( zwf )
[5836]281         !
[7646]282         CALL lbc_lnk( fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask
[3]283         !
[7646]284         ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 ) depending on ln_vorlat
[5506]285         !
[2528]286      ENDIF
[7646]287     
288      ! User defined alteration of fmask (use to reduce ocean transport in specified straits)
289      ! --------------------------------
[2528]290      !
[7646]291      CALL usr_def_fmask( cn_cfg, nn_cfg, fmask )
[6140]292      !
[3]293   END SUBROUTINE dom_msk
294   
295   !!======================================================================
296END MODULE dommsk
Note: See TracBrowser for help on using the repository browser.