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 @ 13305

Last change on this file since 13305 was 13305, checked in by acc, 4 years ago

Trunk changes required to avoid issues with the outer halo in ORCA2_ICE_PISCES,
REPRO_8_4 tests with nn_hls=2. These changes ensure that tmask and output
from sbc_blk are set correctly in the outer halo. Failure to set valid
values in the outer halo can generate NaNs? which lead to OOB errors in the
XRGB lookup table used for the TRC optics.See #2366 for details. With these
changes all variants of the ORCA2_ICE_PISCES SETTE test will complete. There
are still differences between the 1 and 2 halo width runs but running with:
no land suppression; partial land suppression or full land suppression does
not alter either set of results. Likewise setting ln_nnogather either true
or false does not alter results. Differences in run.stat start after 140
timesteps and differences in tracer.stat start after 60 timesteps between
the different halo width sets. Equivalent tests with ln_icebergs = F show no
differences in run.stat but halo-width dependent differences in tracer.stat
persist (now after 64 timesteps).

  • Property svn:keywords set to Id
File size: 11.9 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
10   !!            8.0  ! 1997-02  (G. Madec)  mesh information put in domhgr.F
11   !!            8.1  ! 1997-07  (G. Madec)  modification of kbat and fmask
12   !!             -   ! 1998-05  (G. Roullet)  free surface
13   !!            8.2  ! 2000-03  (G. Madec)  no slip accurate
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
18   !!            3.6  ! 2015-05  (P. Mathiot) ISF: add wmask,wumask and wvmask
19   !!            4.0  ! 2016-06  (G. Madec, S. Flavoni)  domain configuration / user defined interface
20   !!            4.x  ! 2020-02  (G. Madec, S. Techene) introduce ssh to h0 ratio
21   !!----------------------------------------------------------------------
22
23   !!----------------------------------------------------------------------
24   !!   dom_msk       : compute land/ocean mask
25   !!----------------------------------------------------------------------
26   USE oce            ! ocean dynamics and tracers
27   USE dom_oce        ! ocean space and time domain
28   USE domutl         !
29   USE usrdef_fmask   ! user defined fmask
30   USE bdy_oce        ! open boundary
31   !
32   USE in_out_manager ! I/O manager
33   USE iom            ! IOM library
34   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
35   USE lib_mpp        ! Massively Parallel Processing library
36
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   dom_msk    ! routine called by inidom.F90
41
42   !                            !!* Namelist namlbc : lateral boundary condition *
43   REAL(wp)        :: rn_shlat   ! type of lateral boundary condition on velocity
44   LOGICAL, PUBLIC :: ln_vorlat  !  consistency of vorticity boundary condition
45   !                                            with analytical eqs.
46
47   !! * Substitutions
48#  include "do_loop_substitute.h90"
49   !!----------------------------------------------------------------------
50   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
51   !! $Id$
52   !! Software governed by the CeCILL license (see ./LICENSE)
53   !!----------------------------------------------------------------------
54CONTAINS
55
56   SUBROUTINE dom_msk( k_top, k_bot )
57      !!---------------------------------------------------------------------
58      !!                 ***  ROUTINE dom_msk  ***
59      !!
60      !! ** Purpose :   Compute land/ocean mask arrays at tracer points, hori-
61      !!      zontal velocity points (u & v), vorticity points (f) points.
62      !!
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) :
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
74      !!
75      !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated
76      !!                rows/lines due to cyclic or North Fold boundaries as well
77      !!                as MPP halos.
78      !!      tmask_h : halo mask at t-point, i.e. excluding duplicated rows/lines
79      !!                due to cyclic or North Fold boundaries as well as MPP halos.
80      !!
81      !! ** Action :   tmask, umask, vmask, wmask, wumask, wvmask : land/ocean mask
82      !!                         at t-, u-, v- w, wu-, and wv-points (=0. or 1.)
83      !!               fmask   : land/ocean mask at f-point (=0., or =1., or
84      !!                         =rn_shlat along lateral boundaries)
85      !!               tmask_i : interior ocean mask
86      !!               tmask_h : halo mask
87      !!               ssmask , ssumask, ssvmask, ssfmask : 2D ocean mask
88      !!----------------------------------------------------------------------
89      INTEGER, DIMENSION(:,:), INTENT(in) ::   k_top, k_bot   ! first and last ocean level
90      !
91      INTEGER  ::   ji, jj, jk     ! dummy loop indices
92      INTEGER  ::   iktop, ikbot   !   -       -
93      INTEGER  ::   ios, inum
94      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zwf   ! 2D workspace
95      !!
96      NAMELIST/namlbc/ rn_shlat, ln_vorlat
97      NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file,         &
98         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
99         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
100         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
101         &             cn_ice, nn_ice_dta,                                     &
102         &             ln_vol, nn_volctl, nn_rimwidth
103      !!---------------------------------------------------------------------
104      !
105      READ  ( numnam_ref, namlbc, IOSTAT = ios, ERR = 901 )
106901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namlbc in reference namelist' )
107      READ  ( numnam_cfg, namlbc, IOSTAT = ios, ERR = 902 )
108902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namlbc in configuration namelist' )
109      IF(lwm) WRITE ( numond, namlbc )
110     
111      IF(lwp) THEN                  ! control print
112         WRITE(numout,*)
113         WRITE(numout,*) 'dommsk : ocean mask '
114         WRITE(numout,*) '~~~~~~'
115         WRITE(numout,*) '   Namelist namlbc'
116         WRITE(numout,*) '      lateral momentum boundary cond.    rn_shlat  = ',rn_shlat
117         WRITE(numout,*) '      consistency with analytical form   ln_vorlat = ',ln_vorlat 
118      ENDIF
119      !
120      IF(lwp) WRITE(numout,*)
121      IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  free-slip'
122      ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  no-slip'
123      ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  partial-slip'
124      ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  strong-slip'
125      ELSE
126         CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' )
127      ENDIF
128
129      !  Ocean/land mask at t-point  (computed from ko_top and ko_bot)
130      ! ----------------------------
131      !
132      tmask(:,:,:) = 0._wp
133      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
134         iktop = k_top(ji,jj)
135         ikbot = k_bot(ji,jj)
136         IF( iktop /= 0 ) THEN       ! water in the column
137            tmask(ji,jj,iktop:ikbot) = 1._wp
138         ENDIF
139      END_2D
140      !
141      ! Mask corrections for bdy (read in mppini2)
142      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 903)
143903   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy in reference namelist' )
144      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 904 )
145904   IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy in configuration namelist' )
146      ! ------------------------
147      IF ( ln_bdy .AND. ln_mask_file ) THEN
148         CALL iom_open( cn_mask_file, inum )
149         CALL iom_get ( inum, jpdom_global, 'bdy_msk', bdytmask(:,:) )
150         CALL iom_close( inum )
151         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
152            tmask(ji,jj,jk) = tmask(ji,jj,jk) * bdytmask(ji,jj)
153         END_3D
154      ENDIF
155         
156      !  Ocean/land mask at u-, v-, and f-points   (computed from tmask)
157      ! ----------------------------------------
158      ! NB: at this point, fmask is designed for free slip lateral boundary condition
159      DO_3D( 0, 0, 0, 0, 1, jpk )
160         umask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)
161         vmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji  ,jj+1,jk)
162         fmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)   &
163            &            * tmask(ji,jj+1,jk) * tmask(ji+1,jj+1,jk)
164      END_3D
165      CALL lbc_lnk_multi( 'dommsk', umask, 'U', 1.0_wp, vmask, 'V', 1.0_wp, fmask, 'F', 1.0_wp )      ! Lateral boundary conditions
166 
167      ! Ocean/land mask at wu-, wv- and w points    (computed from tmask)
168      !-----------------------------------------
169      wmask (:,:,1) = tmask(:,:,1)     ! surface
170      wumask(:,:,1) = umask(:,:,1)
171      wvmask(:,:,1) = vmask(:,:,1)
172      DO jk = 2, jpk                   ! interior values
173         wmask (:,:,jk) = tmask(:,:,jk) * tmask(:,:,jk-1)
174         wumask(:,:,jk) = umask(:,:,jk) * umask(:,:,jk-1)   
175         wvmask(:,:,jk) = vmask(:,:,jk) * vmask(:,:,jk-1)
176      END DO
177
178      ! Ocean/land column mask at t-, u-, and v-points   (i.e. at least 1 wet cell in the vertical)
179      ! ----------------------------------------------
180      ssmask (:,:) = MAXVAL( tmask(:,:,:), DIM=3 )
181      ssumask(:,:) = MAXVAL( umask(:,:,:), DIM=3 )
182      ssvmask(:,:) = MAXVAL( vmask(:,:,:), DIM=3 )
183      ssfmask(:,:) = MAXVAL( fmask(:,:,:), DIM=3 )
184
185      ! Interior domain mask  (used for global sum)
186      ! --------------------
187      !
188      CALL dom_uniq( tmask_h, 'T' )
189      !
190      !                          ! interior mask : 2D ocean mask x halo mask
191      tmask_i(:,:) = ssmask(:,:) * tmask_h(:,:)
192
193      ! Lateral boundary conditions on velocity (modify fmask)
194      ! --------------------------------------- 
195      IF( rn_shlat /= 0 ) THEN      ! Not free-slip lateral boundary condition
196         !
197         ALLOCATE( zwf(jpi,jpj) )
198         !
199         DO jk = 1, jpk
200            zwf(:,:) = fmask(:,:,jk)         
201            DO_2D( 0, 0, 0, 0 )
202               IF( fmask(ji,jj,jk) == 0._wp ) THEN
203                  fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1),   &
204                     &                                           zwf(ji-1,jj), zwf(ji,jj-1)  )  )
205               ENDIF
206            END_2D
207            DO jj = 2, jpjm1
208               IF( fmask(1,jj,jk) == 0._wp ) THEN
209                  fmask(1  ,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
210               ENDIF
211               IF( fmask(jpi,jj,jk) == 0._wp ) THEN
212                  fmask(jpi,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
213               ENDIF
214            END DO         
215            DO ji = 2, jpim1
216               IF( fmask(ji,1,jk) == 0._wp ) THEN
217                  fmask(ji, 1 ,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
218               ENDIF
219               IF( fmask(ji,jpj,jk) == 0._wp ) THEN
220                  fmask(ji,jpj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
221               ENDIF
222            END DO
223         END DO
224         !
225         DEALLOCATE( zwf )
226         !
227         CALL lbc_lnk( 'dommsk', fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask
228         !
229         ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 ) depending on ln_vorlat
230         !
231      ENDIF
232     
233      ! User defined alteration of fmask (use to reduce ocean transport in specified straits)
234      ! --------------------------------
235      !
236      CALL usr_def_fmask( cn_cfg, nn_cfg, fmask )
237      !
238   END SUBROUTINE dom_msk
239   
240   !!======================================================================
241END MODULE dommsk
Note: See TracBrowser for help on using the repository browser.