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/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/DOM/dommsk.F90 @ 7412

Last change on this file since 7412 was 7412, checked in by lovato, 8 years ago

Merge dev_NOC_CMCC_merge_2016 into branch

  • Property svn:keywords set to Id
File size: 19.0 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 mbathy 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   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   dom_msk        : compute land/ocean mask
23   !!----------------------------------------------------------------------
24   USE oce             ! ocean dynamics and tracers
25   USE dom_oce         ! ocean space and time domain
26   USE bdy_oce     
27   USE in_out_manager  ! I/O manager
28   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
29   USE lib_mpp         !
30   USE iom
31   USE wrk_nemo        ! Memory allocation
32   USE timing          ! Timing
33
34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC   dom_msk    ! routine called by inidom.F90
38
39   !                            !!* Namelist namlbc : lateral boundary condition *
40   REAL(wp)        :: rn_shlat   ! type of lateral boundary condition on velocity
41   LOGICAL, PUBLIC :: ln_vorlat  !  consistency of vorticity boundary condition
42   !                                            with analytical eqs.
43
44   !! * Substitutions
45#  include "vectopt_loop_substitute.h90"
46   !!----------------------------------------------------------------------
47   !! NEMO/OPA 3.2 , LODYC-IPSL  (2009)
48   !! $Id$
49   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
50   !!----------------------------------------------------------------------
51CONTAINS
52
53   SUBROUTINE dom_msk
54      !!---------------------------------------------------------------------
55      !!                 ***  ROUTINE dom_msk  ***
56      !!
57      !! ** Purpose :   Compute land/ocean mask arrays at tracer points, hori-
58      !!      zontal velocity points (u & v), vorticity points (f) points.
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.
62      !!      mbathy equals 0 over continental T-point
63      !!      and the number of ocean level over the ocean.
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)
75      !!                  and mbathy(ji+1,jj) and mbathy(ji+1,jj+1) >= jk.
76      !!      tmask_i : interior ocean mask at t-point, i.e. excluding duplicated
77      !!                rows/lines due to cyclic or North Fold boundaries as well
78      !!                as MPP halos.
79      !!
80      !!        The lateral friction is set through the value of fmask along
81      !!      the coast and topography. This value is defined by rn_shlat, a
82      !!      namelist parameter:
83      !!         rn_shlat = 0, free slip  (no shear along the coast)
84      !!         rn_shlat = 2, no slip  (specified zero velocity at the coast)
85      !!         0 < rn_shlat < 2, partial slip   | non-linear velocity profile
86      !!         2 < rn_shlat, strong slip        | in the lateral boundary layer
87      !!
88      !!      N.B. If nperio not equal to 0, the land/ocean mask arrays
89      !!      are defined with the proper value at lateral domain boundaries.
90      !!
91      !!      In case of open boundaries (ln_bdy=T):
92      !!        - tmask is set to 1 on the points to be computed bay the open
93      !!          boundaries routines.
94      !!
95      !! ** Action :   tmask    : land/ocean mask at t-point (=0. or 1.)
96      !!               umask    : land/ocean mask at u-point (=0. or 1.)
97      !!               vmask    : land/ocean mask at v-point (=0. or 1.)
98      !!               fmask    : land/ocean mask at f-point (=0. or 1.)
99      !!                          =rn_shlat along lateral boundaries
100      !!               tmask_i  : interior ocean mask
101      !!----------------------------------------------------------------------
102      INTEGER  ::   ji, jj, jk               ! dummy loop indices
103      INTEGER  ::   iif, iil, ii0, ii1, ii   ! local integers
104      INTEGER  ::   ijf, ijl, ij0, ij1       !   -       -
105      INTEGER  ::   ios, inum
106      INTEGER  ::   isrow                    ! index for ORCA1 starting row
107      INTEGER , POINTER, DIMENSION(:,:) ::  imsk
108      REAL(wp), POINTER, DIMENSION(:,:) ::  zwf
109      !!
110      NAMELIST/namlbc/ rn_shlat, ln_vorlat
111      NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file,         &
112         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
113         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
114         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
115         &             cn_ice_lim, nn_ice_lim_dta,                           &
116         &             rn_ice_tem, rn_ice_sal, rn_ice_age,                 &
117         &             ln_vol, nn_volctl, nn_rimwidth, nb_jpk_bdy
118      !!---------------------------------------------------------------------
119      !
120      IF( nn_timing == 1 )  CALL timing_start('dom_msk')
121      !
122      CALL wrk_alloc( jpi, jpj, imsk )
123      CALL wrk_alloc( jpi, jpj, zwf  )
124      !
125      REWIND( numnam_ref )              ! Namelist namlbc in reference namelist : Lateral momentum boundary condition
126      READ  ( numnam_ref, namlbc, IOSTAT = ios, ERR = 901 )
127901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namlbc in reference namelist', lwp )
128
129      REWIND( numnam_cfg )              ! Namelist namlbc in configuration namelist : Lateral momentum boundary condition
130      READ  ( numnam_cfg, namlbc, IOSTAT = ios, ERR = 902 )
131902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namlbc in configuration namelist', lwp )
132      IF(lwm) WRITE ( numond, namlbc )
133     
134      IF(lwp) THEN                  ! control print
135         WRITE(numout,*)
136         WRITE(numout,*) 'dommsk : ocean mask '
137         WRITE(numout,*) '~~~~~~'
138         WRITE(numout,*) '   Namelist namlbc'
139         WRITE(numout,*) '      lateral momentum boundary cond.    rn_shlat  = ',rn_shlat
140         WRITE(numout,*) '      consistency with analytical form   ln_vorlat = ',ln_vorlat 
141      ENDIF
142
143      IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  free-slip '
144      ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  no-slip '
145      ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  partial-slip '
146      ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ocean lateral  strong-slip '
147      ELSE
148         WRITE(ctmp1,*) ' rn_shlat is negative = ', rn_shlat
149         CALL ctl_stop( ctmp1 )
150      ENDIF
151
152      ! 1. Ocean/land mask at t-point (computed from mbathy)
153      ! -----------------------------
154      ! N.B. tmask has already the right boundary conditions since mbathy is ok
155      !
156      tmask(:,:,:) = 0._wp
157      DO jk = 1, jpk
158         DO jj = 1, jpj
159            DO ji = 1, jpi
160               IF( REAL( mbathy(ji,jj) - jk, wp ) + 0.1_wp >= 0._wp )   tmask(ji,jj,jk) = 1._wp
161            END DO 
162         END DO 
163      END DO 
164     
165      REWIND( numnam_ref )              ! Namelist nambdy in reference namelist :Unstructured open boundaries 
166      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 903)
167903   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in reference namelist', lwp )
168
169      REWIND( numnam_cfg )              ! Namelist nambdy in configuration namelist :Unstructured open boundaries
170      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 904 )
171904   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in configuration namelist', lwp )
172      IF(lwm) WRITE ( numond, nambdy )
173
174     IF( ln_bdy .AND. ln_mask_file ) THEN ! correct for bdy mask
175         CALL iom_open( cn_mask_file, inum )
176         CALL iom_get ( inum, jpdom_data, 'bdy_msk', bdytmask(:,:) )
177         CALL iom_close( inum )
178
179         ! Mask corrections
180         ! ----------------
181         DO jk = 1, jpkm1
182            DO jj = 1, jpj
183               DO ji = 1, jpi
184                  tmask(ji,jj,jk) = tmask(ji,jj,jk) * bdytmask(ji,jj)
185               END DO
186            END DO
187         END DO
188      ENDIF
189
190      ! (ISF) define barotropic mask and mask the ice shelf point
191      ssmask(:,:)=tmask(:,:,1) ! at this stage ice shelf is not masked
192     
193      DO jk = 1, jpk
194         DO jj = 1, jpj
195            DO ji = 1, jpi
196               IF( REAL( misfdep(ji,jj) - jk, wp ) - 0.1_wp >= 0._wp )   THEN
197                  tmask(ji,jj,jk) = 0._wp
198               END IF
199            END DO 
200         END DO 
201      END DO 
202
203      ! Interior domain mask (used for global sum)
204      ! --------------------
205      tmask_i(:,:) = ssmask(:,:)            ! (ISH) tmask_i = 1 even on the ice shelf
206
207      tmask_h(:,:) = 1._wp                 ! 0 on the halo and 1 elsewhere
208      iif = jpreci                         ! ???
209      iil = nlci - jpreci + 1
210      ijf = jprecj                         ! ???
211      ijl = nlcj - jprecj + 1
212
213      tmask_h( 1 :iif,   :   ) = 0._wp      ! first columns
214      tmask_h(iil:jpi,   :   ) = 0._wp      ! last  columns (including mpp extra columns)
215      tmask_h(   :   , 1 :ijf) = 0._wp      ! first rows
216      tmask_h(   :   ,ijl:jpj) = 0._wp      ! last  rows (including mpp extra rows)
217
218      ! north fold mask
219      ! ---------------
220      tpol(1:jpiglo) = 1._wp 
221      fpol(1:jpiglo) = 1._wp
222      IF( jperio == 3 .OR. jperio == 4 ) THEN      ! T-point pivot
223         tpol(jpiglo/2+1:jpiglo) = 0._wp
224         fpol(     1    :jpiglo) = 0._wp
225         IF( mjg(nlej) == jpjglo ) THEN                  ! only half of the nlcj-1 row
226            DO ji = iif+1, iil-1
227               tmask_h(ji,nlej-1) = tmask_h(ji,nlej-1) * tpol(mig(ji))
228            END DO
229         ENDIF
230      ENDIF
231     
232      tmask_i(:,:) = tmask_i(:,:) * tmask_h(:,:)
233
234      IF( jperio == 5 .OR. jperio == 6 ) THEN      ! F-point pivot
235         tpol(     1    :jpiglo) = 0._wp
236         fpol(jpiglo/2+1:jpiglo) = 0._wp
237      ENDIF
238
239      ! 2. Ocean/land mask at u-,  v-, and z-points (computed from tmask)
240      ! -------------------------------------------
241      DO jk = 1, jpk
242         DO jj = 1, jpjm1
243            DO ji = 1, fs_jpim1   ! vector loop
244               umask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)
245               vmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji  ,jj+1,jk)
246            END DO
247            DO ji = 1, jpim1      ! NO vector opt.
248               fmask(ji,jj,jk) = tmask(ji,jj  ,jk) * tmask(ji+1,jj  ,jk)   &
249                  &            * tmask(ji,jj+1,jk) * tmask(ji+1,jj+1,jk)
250            END DO
251         END DO
252      END DO
253      ! (ISF) MIN(1,SUM(umask)) is here to check if you have effectively at least 1 wet cell at u point
254      DO jj = 1, jpjm1
255         DO ji = 1, fs_jpim1   ! vector loop
256            ssumask(ji,jj)  = ssmask(ji,jj) * ssmask(ji+1,jj  )  * MIN(1._wp,SUM(umask(ji,jj,:)))
257            ssvmask(ji,jj)  = ssmask(ji,jj) * ssmask(ji  ,jj+1)  * MIN(1._wp,SUM(vmask(ji,jj,:)))
258         END DO
259         DO ji = 1, jpim1      ! NO vector opt.
260            ssfmask(ji,jj) =  ssmask(ji,jj  ) * ssmask(ji+1,jj  )   &
261               &            * ssmask(ji,jj+1) * ssmask(ji+1,jj+1) * MIN(1._wp,SUM(fmask(ji,jj,:)))
262         END DO
263      END DO
264      CALL lbc_lnk( umask  , 'U', 1._wp )      ! Lateral boundary conditions
265      CALL lbc_lnk( vmask  , 'V', 1._wp )
266      CALL lbc_lnk( fmask  , 'F', 1._wp )
267      CALL lbc_lnk( ssumask, 'U', 1._wp )      ! Lateral boundary conditions
268      CALL lbc_lnk( ssvmask, 'V', 1._wp )
269      CALL lbc_lnk( ssfmask, 'F', 1._wp )
270
271      ! 3. Ocean/land mask at wu-, wv- and w points
272      !----------------------------------------------
273      wmask (:,:,1) = tmask(:,:,1)     ! surface
274      wumask(:,:,1) = umask(:,:,1)
275      wvmask(:,:,1) = vmask(:,:,1)
276      DO jk = 2, jpk                   ! interior values
277         wmask (:,:,jk) = tmask(:,:,jk) * tmask(:,:,jk-1)
278         wumask(:,:,jk) = umask(:,:,jk) * umask(:,:,jk-1)   
279         wvmask(:,:,jk) = vmask(:,:,jk) * vmask(:,:,jk-1)
280      END DO
281
282      ! Lateral boundary conditions on velocity (modify fmask)
283      ! ---------------------------------------     
284      DO jk = 1, jpk
285         zwf(:,:) = fmask(:,:,jk)         
286         DO jj = 2, jpjm1
287            DO ji = fs_2, fs_jpim1   ! vector opt.
288               IF( fmask(ji,jj,jk) == 0._wp ) THEN
289                  fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1),   &
290                     &                                           zwf(ji-1,jj), zwf(ji,jj-1)  )  )
291               ENDIF
292            END DO
293         END DO
294         DO jj = 2, jpjm1
295            IF( fmask(1,jj,jk) == 0._wp ) THEN
296               fmask(1  ,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
297            ENDIF
298            IF( fmask(jpi,jj,jk) == 0._wp ) THEN
299               fmask(jpi,jj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
300            ENDIF
301         END DO         
302         DO ji = 2, jpim1
303            IF( fmask(ji,1,jk) == 0._wp ) THEN
304               fmask(ji, 1 ,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
305            ENDIF
306            IF( fmask(ji,jpj,jk) == 0._wp ) THEN
307               fmask(ji,jpj,jk) = rn_shlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
308            ENDIF
309         END DO
310      END DO
311      !
312      IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN   ! ORCA_R2 configuration
313         !                                                 ! Increased lateral friction near of some straits
314         !                                ! Gibraltar strait  : partial slip (fmask=0.5)
315         ij0 = 101   ;   ij1 = 101
316         ii0 = 139   ;   ii1 = 140   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  0.5_wp
317         ij0 = 102   ;   ij1 = 102
318         ii0 = 139   ;   ii1 = 140   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  0.5_wp
319         !
320         !                                ! Bab el Mandeb : partial slip (fmask=1)
321         ij0 =  87   ;   ij1 =  88
322         ii0 = 160   ;   ii1 = 160   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  1._wp
323         ij0 =  88   ;   ij1 =  88
324         ii0 = 159   ;   ii1 = 159   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) =  1._wp
325         !
326         !                                ! Danish straits  : strong slip (fmask > 2)
327! We keep this as an example but it is instable in this case
328!         ij0 = 115   ;   ij1 = 115
329!         ii0 = 145   ;   ii1 = 146   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) = 4._wp
330!         ij0 = 116   ;   ij1 = 116
331!         ii0 = 145   ;   ii1 = 146   ;   fmask( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , 1:jpk ) = 4._wp
332         !
333      ENDIF
334      !
335      IF( cp_cfg == "orca" .AND. jp_cfg == 1 ) THEN   ! ORCA R1 configuration
336         !                                                 ! Increased lateral friction near of some straits
337         ! This dirty section will be suppressed by simplification process:
338         ! all this will come back in input files
339         ! Currently these hard-wired indices relate to configuration with
340         ! extend grid (jpjglo=332)
341         !
342         isrow = 332 - jpjglo
343         !
344         IF(lwp) WRITE(numout,*)
345         IF(lwp) WRITE(numout,*) '   orca_r1: increase friction near the following straits : '
346         IF(lwp) WRITE(numout,*) '      Gibraltar '
347         ii0 = 282           ;   ii1 = 283        ! Gibraltar Strait
348         ij0 = 241 - isrow   ;   ij1 = 241 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
349
350         IF(lwp) WRITE(numout,*) '      Bhosporus '
351         ii0 = 314           ;   ii1 = 315        ! Bhosporus Strait
352         ij0 = 248 - isrow   ;   ij1 = 248 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
353
354         IF(lwp) WRITE(numout,*) '      Makassar (Top) '
355         ii0 =  48           ;   ii1 =  48        ! Makassar Strait (Top)
356         ij0 = 189 - isrow   ;   ij1 = 190 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
357
358         IF(lwp) WRITE(numout,*) '      Lombok '
359         ii0 =  44           ;   ii1 =  44        ! Lombok Strait
360         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
361
362         IF(lwp) WRITE(numout,*) '      Ombai '
363         ii0 =  53           ;   ii1 =  53        ! Ombai Strait
364         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
365
366         IF(lwp) WRITE(numout,*) '      Timor Passage '
367         ii0 =  56           ;   ii1 =  56        ! Timor Passage
368         ij0 = 164 - isrow   ;   ij1 = 165 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 2._wp 
369
370         IF(lwp) WRITE(numout,*) '      West Halmahera '
371         ii0 =  58           ;   ii1 =  58        ! West Halmahera Strait
372         ij0 = 181 - isrow   ;   ij1 = 182 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
373
374         IF(lwp) WRITE(numout,*) '      East Halmahera '
375         ii0 =  55           ;   ii1 =  55        ! East Halmahera Strait
376         ij0 = 181 - isrow   ;   ij1 = 182 - isrow   ;   fmask( mi0(ii0):mi1(ii1),mj0(ij0):mj1(ij1),1:jpk ) = 3._wp 
377         !
378      ENDIF
379      !
380      CALL lbc_lnk( fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask
381      !
382      ! CAUTION : The fmask may be further modified in dyn_vor_init ( dynvor.F90 )
383      !
384      CALL wrk_dealloc( jpi, jpj, imsk )
385      CALL wrk_dealloc( jpi, jpj, zwf  )
386      !
387      IF( nn_timing == 1 )  CALL timing_stop('dom_msk')
388      !
389   END SUBROUTINE dom_msk
390   
391   !!======================================================================
392END MODULE dommsk
Note: See TracBrowser for help on using the repository browser.