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 NEMO/trunk/src/OCE/DOM – NEMO

source: NEMO/trunk/src/OCE/DOM/dommsk.F90 @ 15717

Last change on this file since 15717 was 15556, checked in by jchanut, 3 years ago

#2638: Add the possibility to read bottom levels at U/V/F points in the mesh file. Store fe3mask (i.e. fmask as it is prior updating it for lateral boundary conditions). All is this is only needed to ensure a correct update of parent grid variables with AgRIF. This also anticipates the possible use of coarsened meshes.

  • Property svn:keywords set to Id
File size: 11.7 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
[13237]20   !!            4.x  ! 2020-02  (G. Madec, S. Techene) introduce ssh to h0 ratio
[1566]21   !!----------------------------------------------------------------------
[3]22
23   !!----------------------------------------------------------------------
[7646]24   !!   dom_msk       : compute land/ocean mask
[3]25   !!----------------------------------------------------------------------
[7646]26   USE oce            ! ocean dynamics and tracers
27   USE dom_oce        ! ocean space and time domain
[13286]28   USE domutl         !
[7646]29   USE usrdef_fmask   ! user defined fmask
[9124]30   USE bdy_oce        ! open boundary
31   !
[7646]32   USE in_out_manager ! I/O manager
[9600]33   USE iom            ! IOM library
[7646]34   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
35   USE lib_mpp        ! Massively Parallel Processing library
[3]36
37   IMPLICIT NONE
38   PRIVATE
39
[6140]40   PUBLIC   dom_msk    ! routine called by inidom.F90
[3]41
[1601]42   !                            !!* Namelist namlbc : lateral boundary condition *
[4147]43   REAL(wp)        :: rn_shlat   ! type of lateral boundary condition on velocity
44   LOGICAL, PUBLIC :: ln_vorlat  !  consistency of vorticity boundary condition
[3294]45   !                                            with analytical eqs.
[2715]46
[3]47   !! * Substitutions
[12377]48#  include "do_loop_substitute.h90"
[1566]49   !!----------------------------------------------------------------------
[9598]50   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[7753]51   !! $Id$
[10068]52   !! Software governed by the CeCILL license (see ./LICENSE)
[1566]53   !!----------------------------------------------------------------------
[3]54CONTAINS
[2715]55
[7646]56   SUBROUTINE dom_msk( k_top, k_bot )
[3]57      !!---------------------------------------------------------------------
58      !!                 ***  ROUTINE dom_msk  ***
59      !!
60      !! ** Purpose :   Compute land/ocean mask arrays at tracer points, hori-
[6125]61      !!      zontal velocity points (u & v), vorticity points (f) points.
[3]62      !!
[7646]63      !! ** Method  :   The ocean/land mask  at t-point is deduced from ko_top
64      !!      and ko_bot, the indices of the fist and last ocean t-levels which
65      !!      are either defined in usrdef_zgr or read in zgr_read.
66      !!                The velocity masks (umask, vmask, wmask, wumask, wvmask)
67      !!      are deduced from a product of the two neighboring tmask.
68      !!                The vorticity mask (fmask) is deduced from tmask taking
69      !!      into account the choice of lateral boundary condition (rn_shlat) :
[1601]70      !!         rn_shlat = 0, free slip  (no shear along the coast)
71      !!         rn_shlat = 2, no slip  (specified zero velocity at the coast)
72      !!         0 < rn_shlat < 2, partial slip   | non-linear velocity profile
73      !!         2 < rn_shlat, strong slip        | in the lateral boundary layer
[3]74      !!
[7646]75      !! ** Action :   tmask, umask, vmask, wmask, wumask, wvmask : land/ocean mask
76      !!                         at t-, u-, v- w, wu-, and wv-points (=0. or 1.)
77      !!               fmask   : land/ocean mask at f-point (=0., or =1., or
78      !!                         =rn_shlat along lateral boundaries)
[15145]79      !!               ssmask , ssumask, ssvmask, ssfmask : 2D ocean mask, i.e. at least 1 wet cell in the vertical
80      !!               tmask_h : halo mask at t-point, i.e. excluding all duplicated rows/lines
81      !!                         due to cyclic or North Fold boundaries as well as MPP halos.
82      !!               tmask_i : ssmask * tmask_h
[3]83      !!----------------------------------------------------------------------
[7646]84      INTEGER, DIMENSION(:,:), INTENT(in) ::   k_top, k_bot   ! first and last ocean level
85      !
86      INTEGER  ::   ji, jj, jk     ! dummy loop indices
87      INTEGER  ::   iktop, ikbot   !   -       -
88      INTEGER  ::   ios, inum
[1601]89      !!
[3294]90      NAMELIST/namlbc/ rn_shlat, ln_vorlat
[7646]91      NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file,         &
92         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
93         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
94         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
[9657]95         &             cn_ice, nn_ice_dta,                                     &
[11536]96         &             ln_vol, nn_volctl, nn_rimwidth
[3]97      !!---------------------------------------------------------------------
[3294]98      !
[4147]99      READ  ( numnam_ref, namlbc, IOSTAT = ios, ERR = 901 )
[11536]100901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namlbc in reference namelist' )
[4147]101      READ  ( numnam_cfg, namlbc, IOSTAT = ios, ERR = 902 )
[11536]102902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namlbc in configuration namelist' )
[4624]103      IF(lwm) WRITE ( numond, namlbc )
[1566]104     
105      IF(lwp) THEN                  ! control print
[3]106         WRITE(numout,*)
107         WRITE(numout,*) 'dommsk : ocean mask '
108         WRITE(numout,*) '~~~~~~'
[1566]109         WRITE(numout,*) '   Namelist namlbc'
[3294]110         WRITE(numout,*) '      lateral momentum boundary cond.    rn_shlat  = ',rn_shlat
111         WRITE(numout,*) '      consistency with analytical form   ln_vorlat = ',ln_vorlat 
[3]112      ENDIF
[9169]113      !
114      IF(lwp) WRITE(numout,*)
115      IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  free-slip'
116      ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  no-slip'
117      ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  partial-slip'
118      ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  strong-slip'
[1601]119      ELSE
[9527]120         CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' )
[3]121      ENDIF
122
[7646]123      !  Ocean/land mask at t-point  (computed from ko_top and ko_bot)
124      ! ----------------------------
[1566]125      !
[7753]126      tmask(:,:,:) = 0._wp
[13305]127      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[12377]128         iktop = k_top(ji,jj)
129         ikbot = k_bot(ji,jj)
130         IF( iktop /= 0 ) THEN       ! water in the column
[13286]131            tmask(ji,jj,iktop:ikbot) = 1._wp
[12377]132         ENDIF
133      END_2D
[11233]134      !
[13286]135      ! Mask corrections for bdy (read in mppini2)
[7646]136      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 903)
[11536]137903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy in reference namelist' )
[7646]138      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 904 )
[11536]139904   IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy in configuration namelist' )
[7646]140      ! ------------------------
141      IF ( ln_bdy .AND. ln_mask_file ) THEN
[9600]142         CALL iom_open( cn_mask_file, inum )
[13286]143         CALL iom_get ( inum, jpdom_global, 'bdy_msk', bdytmask(:,:) )
[9600]144         CALL iom_close( inum )
[13295]145         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[12377]146            tmask(ji,jj,jk) = tmask(ji,jj,jk) * bdytmask(ji,jj)
147         END_3D
[3]148      ENDIF
[7646]149         
150      !  Ocean/land mask at u-, v-, and f-points   (computed from tmask)
151      ! ----------------------------------------
152      ! NB: at this point, fmask is designed for free slip lateral boundary condition
[13295]153      DO_3D( 0, 0, 0, 0, 1, jpk )
[13286]154         umask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)
155         vmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji  ,jj+1,jk)
156         fmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)   &
157            &            * tmask(ji,jj+1,jk) * tmask(ji+1,jj+1,jk)
158      END_3D
[15556]159      !
160      ! In case of a coarsened grid, account her for possibly aditionnal 
161      ! masked points; these have been read in the mesh file and stored in mbku, mbkv, mbkf
162      DO_2D( 0, 0, 0, 0 )
163         IF (mbku(ji,jj)<=1 ) umask(ji,jj,:) = 0._wp
164         IF (mbkv(ji,jj)<=1 ) vmask(ji,jj,:) = 0._wp
165         IF (mbkf(ji,jj)<=1 ) fmask(ji,jj,:) = 0._wp
166         IF ( MAXVAL(umask(ji,jj,:))/=0._wp )  umask(ji,jj,mbku(ji,jj)+1:jpk) = 0._wp
167         IF ( MAXVAL(vmask(ji,jj,:))/=0._wp )  vmask(ji,jj,mbkv(ji,jj)+1:jpk) = 0._wp
168         IF ( MAXVAL(fmask(ji,jj,:))/=0._wp )  fmask(ji,jj,mbkf(ji,jj)+1:jpk) = 0._wp
169      END_2D
[14433]170      CALL lbc_lnk( 'dommsk', umask, 'U', 1.0_wp, vmask, 'V', 1.0_wp, fmask, 'F', 1.0_wp )      ! Lateral boundary conditions
[7646]171 
172      ! Ocean/land mask at wu-, wv- and w points    (computed from tmask)
173      !-----------------------------------------
[7753]174      wmask (:,:,1) = tmask(:,:,1)     ! surface
175      wumask(:,:,1) = umask(:,:,1)
176      wvmask(:,:,1) = vmask(:,:,1)
[5836]177      DO jk = 2, jpk                   ! interior values
[7753]178         wmask (:,:,jk) = tmask(:,:,jk) * tmask(:,:,jk-1)
179         wumask(:,:,jk) = umask(:,:,jk) * umask(:,:,jk-1)   
180         wvmask(:,:,jk) = vmask(:,:,jk) * vmask(:,:,jk-1)
[5120]181      END DO
[3]182
[7646]183      ! Ocean/land column mask at t-, u-, and v-points   (i.e. at least 1 wet cell in the vertical)
184      ! ----------------------------------------------
185      ssmask (:,:) = MAXVAL( tmask(:,:,:), DIM=3 )
186      ssumask(:,:) = MAXVAL( umask(:,:,:), DIM=3 )
187      ssvmask(:,:) = MAXVAL( vmask(:,:,:), DIM=3 )
[13237]188      ssfmask(:,:) = MAXVAL( fmask(:,:,:), DIM=3 )
[14053]189      IF( lk_SWE ) THEN      ! Shallow Water Eq. case : redefine ssfmask
[15014]190         DO_2D( 0, 0, 0, 0 )
[14053]191            ssfmask(ji,jj) = MAX(  ssmask(ji,jj+1), ssmask(ji+1,jj+1),  & 
192               &                   ssmask(ji,jj  ), ssmask(ji+1,jj  )   )
193         END_2D
194         CALL lbc_lnk( 'dommsk', ssfmask, 'F', 1.0_wp )
195      ENDIF
196      fe3mask(:,:,:) = fmask(:,:,:)
[7646]197
198      ! Interior domain mask  (used for global sum)
199      ! --------------------
200      !
[13286]201      CALL dom_uniq( tmask_h, 'T' )
[7646]202      !
203      !                          ! interior mask : 2D ocean mask x halo mask
[7753]204      tmask_i(:,:) = ssmask(:,:) * tmask_h(:,:)
[7646]205
[3]206      ! Lateral boundary conditions on velocity (modify fmask)
[7646]207      ! --------------------------------------- 
[15014]208      IF( rn_shlat /= 0._wp ) THEN      ! Not free-slip lateral boundary condition
[7646]209         !
[15014]210         DO_3D( 0, 0, 0, 0, 1, jpk )
211            IF( fmask(ji,jj,jk) == 0._wp ) THEN
212               fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), &
213                  &                                           vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) )
214            ENDIF
215         END_3D
[10425]216         CALL lbc_lnk( 'dommsk', fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask
[3]217         !
[7646]218         ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 ) depending on ln_vorlat
[5506]219         !
[2528]220      ENDIF
[7646]221     
222      ! User defined alteration of fmask (use to reduce ocean transport in specified straits)
223      ! --------------------------------
[2528]224      !
[7646]225      CALL usr_def_fmask( cn_cfg, nn_cfg, fmask )
[6140]226      !
[3]227   END SUBROUTINE dom_msk
228   
229   !!======================================================================
230END MODULE dommsk
Note: See TracBrowser for help on using the repository browser.