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.
wet_dry.F90 in NEMO/branches/2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles/src/OCE/DYN/wet_dry.F90 @ 11671

Last change on this file since 11671 was 11671, checked in by acc, 5 years ago

Branch 2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles. Final, non-substantive changes to complete this branch. These changes remove all REWIND statements on the old namelist fortran units (now character variables for internal files). These changes have been left until last since they are easily repeated via a script and it may be preferable to use the previous revision for merge purposes and reapply these last changes separately. This branch has been fully SETTE tested.

  • Property svn:keywords set to Id
File size: 19.4 KB
Line 
1MODULE wet_dry
2
3   !! includes updates to namelist namwad for diagnostic outputs of ROMS wetting and drying
4
5   !!==============================================================================
6   !!                       ***  MODULE  wet_dry  ***
7   !! Wetting and drying includes initialisation routine and routines to
8   !! compute and apply flux limiters and preserve water depth positivity
9   !! only effects if wetting/drying is on (ln_wd_il == .true. or ln_wd_dl==.true. )
10   !!==============================================================================
11   !! History :  3.6  ! 2014-09  ((H.Liu)  Original code
12   !!                 ! will add the runoff and periodic BC case later
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   wad_init      : initialisation of wetting and drying
17   !!   wad_lmt       : horizontal flux limiter and limited velocity when wetting and drying happens
18   !!   wad_lmt_bt    : same as wad_lmt for the barotropic stepping (dynspg_ts)
19   !!----------------------------------------------------------------------
20   USE oce            ! ocean dynamics and tracers
21   USE dom_oce        ! ocean space and time domain
22   USE sbc_oce  , ONLY: ln_rnf   ! surface boundary condition: ocean
23   USE sbcrnf         ! river runoff
24   !
25   USE in_out_manager ! I/O manager
26   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
27   USE lib_mpp        ! MPP library
28   USE timing         ! timing of the main modules
29
30   IMPLICIT NONE
31   PRIVATE
32
33   !!----------------------------------------------------------------------
34   !! critical depths,filters, limiters,and masks for  Wetting and Drying
35   !! ---------------------------------------------------------------------
36
37   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdmask   !: u- and v- limiter
38   !                                                           !  (can include negative depths)
39   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdramp, wdrampu, wdrampv !: for hpg limiting
40
41   LOGICAL,  PUBLIC  ::   ln_wd_il    !: Wetting/drying il activation switch (T:on,F:off)
42   LOGICAL,  PUBLIC  ::   ln_wd_dl    !: Wetting/drying dl activation switch (T:on,F:off)
43   REAL(wp), PUBLIC  ::   rn_wdmin0   !: depth at which wetting/drying starts
44   REAL(wp), PUBLIC  ::   rn_wdmin1   !: minimum water depth on dried cells
45   REAL(wp), PUBLIC  ::   r_rn_wdmin1 !: 1/minimum water depth on dried cells
46   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells
47   REAL(wp), PUBLIC  ::   rn_wd_sbcdep   !: Depth at which to taper sbc fluxes
48   REAL(wp), PUBLIC  ::   rn_wd_sbcfra   !: Fraction of SBC at taper depth
49   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying will be considered
50   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
51   LOGICAL,  PUBLIC  ::   ln_wd_dl_bc !: DL scheme: True implies 3D velocities are set to the barotropic values at points
52                                      !: where the flow is from wet points on less than half the barotropic sub-steps 
53   LOGICAL,  PUBLIC  ::  ln_wd_dl_rmp !: use a ramp for the dl flux limiter between 2 rn_wdmin1 and rn_wdmin1 (rather than a cut-off at rn_wdmin1)     
54   REAL(wp), PUBLIC  ::   ssh_ref     !: height of z=0 with respect to the geoid;
55
56   LOGICAL,  PUBLIC  ::   ll_wd       !: Wetting/drying activation switch if either ln_wd_il or ln_wd_dl
57
58   PUBLIC   wad_init                  ! initialisation routine called by step.F90
59   PUBLIC   wad_lmt                   ! routine called by sshwzv.F90
60   PUBLIC   wad_lmt_bt                ! routine called by dynspg_ts.F90
61
62   !! * Substitutions
63#  include "vectopt_loop_substitute.h90"
64   !!----------------------------------------------------------------------
65CONTAINS
66
67   SUBROUTINE wad_init
68      !!----------------------------------------------------------------------
69      !!                     ***  ROUTINE wad_init  ***
70      !!                   
71      !! ** Purpose :   read wetting and drying namelist and print the variables.
72      !!
73      !! ** input   : - namwad namelist
74      !!----------------------------------------------------------------------
75      INTEGER  ::   ios, ierr   ! Local integer
76      !!
77      NAMELIST/namwad/ ln_wd_il, ln_wd_dl   , rn_wdmin0, rn_wdmin1, rn_wdmin2, rn_wdld,   &
78         &             nn_wdit , ln_wd_dl_bc, ln_wd_dl_rmp, rn_wd_sbcdep,rn_wd_sbcfra
79      !!----------------------------------------------------------------------
80      !
81      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
82905   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namwad in reference namelist' ) 
83      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
84906   IF( ios >  0 )   CALL ctl_nam ( ios , 'namwad in configuration namelist' )
85      IF(lwm) WRITE ( numond, namwad )
86      !
87      IF( rn_wd_sbcfra>=1 )   CALL ctl_stop( 'STOP', 'rn_wd_sbcfra >=1 : must be < 1' )
88      IF(lwp) THEN                  ! control print
89         WRITE(numout,*)
90         WRITE(numout,*) 'wad_init : Wetting and drying initialization through namelist read'
91         WRITE(numout,*) '~~~~~~~~'
92         WRITE(numout,*) '   Namelist namwad'
93         WRITE(numout,*) '      Logical for Iter Lim wd option   ln_wd_il     = ', ln_wd_il
94         WRITE(numout,*) '      Logical for Dir. Lim wd option   ln_wd_dl     = ', ln_wd_dl
95         WRITE(numout,*) '      Depth at which wet/drying starts rn_wdmin0    = ', rn_wdmin0
96         WRITE(numout,*) '      Minimum wet depth on dried cells rn_wdmin1    = ', rn_wdmin1
97         WRITE(numout,*) '      Tolerance of min wet depth       rn_wdmin2    = ', rn_wdmin2
98         WRITE(numout,*) '      land elevation threshold         rn_wdld      = ', rn_wdld
99         WRITE(numout,*) '      Max iteration for W/D limiter    nn_wdit      = ', nn_wdit
100         WRITE(numout,*) '      T => baroclinic u,v=0 at dry pts: ln_wd_dl_bc = ', ln_wd_dl_bc     
101         WRITE(numout,*) '      use a ramp for rwd limiter:  ln_wd_dl_rwd_rmp = ', ln_wd_dl_rmp
102         WRITE(numout,*) '      cut off depth sbc for wd   rn_wd_sbcdep       = ', rn_wd_sbcdep
103         WRITE(numout,*) '      fraction to start sbc wgt rn_wd_sbcfra        = ', rn_wd_sbcfra
104      ENDIF
105      IF( .NOT. ln_read_cfg ) THEN
106         IF(lwp) WRITE(numout,*) '      No configuration file so seting ssh_ref to zero  '
107         ssh_ref=0._wp
108      ENDIF
109
110      r_rn_wdmin1 = 1 / rn_wdmin1
111      ll_wd = .FALSE.
112      IF( ln_wd_il .OR. ln_wd_dl ) THEN
113         ll_wd = .TRUE.
114         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
115         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
116         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
117      ENDIF
118      !
119   END SUBROUTINE wad_init
120
121
122   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
123      !!----------------------------------------------------------------------
124      !!                  ***  ROUTINE wad_lmt  ***
125      !!                   
126      !! ** Purpose :   generate flux limiters for wetting/drying
127      !!
128      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
129      !!
130      !! ** Action  : - calculate flux limiter and W/D flag
131      !!----------------------------------------------------------------------
132      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1        !!gm DOCTOR names: should start with p !
133      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sshemp
134      REAL(wp)                , INTENT(in   ) ::   z2dt
135      !
136      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
137      INTEGER  ::   jflag               ! local scalar
138      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
139      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
140      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
141      REAL(wp) ::   ztmp                ! local scalars
142      REAL(wp),  DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv   ! W/D flux limiters
143      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxp  ,  zflxn    ! local 2D workspace
144      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu  ,  zflxv    ! local 2D workspace
145      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu1 , zflxv1    ! local 2D workspace
146      !!----------------------------------------------------------------------
147      IF( ln_timing )   CALL timing_start('wad_lmt')      !
148      !
149      DO jk = 1, jpkm1
150         un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:) 
151         vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:) 
152      END DO
153      jflag  = 0
154      zdepwd = 50._wp      ! maximum depth on which that W/D could possibly happen
155      !
156      zflxp(:,:)   = 0._wp
157      zflxn(:,:)   = 0._wp
158      zflxu(:,:)   = 0._wp
159      zflxv(:,:)   = 0._wp
160      !
161      zwdlmtu(:,:) = 1._wp
162      zwdlmtv(:,:) = 1._wp
163      !
164      DO jk = 1, jpkm1     ! Horizontal Flux in u and v direction
165         zflxu(:,:) = zflxu(:,:) + e3u_n(:,:,jk) * un(:,:,jk) * umask(:,:,jk)
166         zflxv(:,:) = zflxv(:,:) + e3v_n(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk)
167      END DO
168      zflxu(:,:) = zflxu(:,:) * e2u(:,:)
169      zflxv(:,:) = zflxv(:,:) * e1v(:,:)
170      !
171      wdmask(:,:) = 1._wp
172      DO jj = 2, jpj
173         DO ji = 2, jpi 
174            !
175            IF( tmask(ji,jj,1)        < 0.5_wp )   CYCLE    ! we don't care about land cells
176            IF( ht_0(ji,jj) - ssh_ref > zdepwd )   CYCLE    ! and cells which are unlikely to dry
177            !
178            zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
179               &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
180            zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
181               &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
182            !
183            zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1
184            IF( zdep2 <= 0._wp ) THEN     ! add more safty, but not necessary
185               sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
186               IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
187               IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
188               IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
189               IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
190               wdmask(ji,jj) = 0._wp
191            END IF
192         END DO
193      END DO
194      !
195      !           ! HPG limiter from jholt
196      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
197      !jth assume don't need a lbc_lnk here
198      DO jj = 1, jpjm1
199         DO ji = 1, jpim1
200            wdrampu(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji+1,jj) )
201            wdrampv(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji,jj+1) )
202         END DO
203      END DO
204      !           ! end HPG limiter
205      !
206      !
207      DO jk1 = 1, nn_wdit + 1      !==  start limiter iterations  ==!
208         !
209         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
210         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
211         jflag = 0     ! flag indicating if any further iterations are needed
212         !
213         DO jj = 2, jpj
214            DO ji = 2, jpi 
215               IF( tmask(ji, jj, 1) < 0.5_wp )   CYCLE
216               IF( ht_0(ji,jj)      > zdepwd )   CYCLE
217               !
218               ztmp = e1e2t(ji,jj)
219               !
220               zzflxp = MAX( zflxu1(ji,jj) , 0._wp ) - MIN( zflxu1(ji-1,jj  ) , 0._wp)   &
221                  &   + MAX( zflxv1(ji,jj) , 0._wp ) - MIN( zflxv1(ji,  jj-1) , 0._wp) 
222               zzflxn = MIN( zflxu1(ji,jj) , 0._wp ) - MAX( zflxu1(ji-1,jj  ) , 0._wp)   &
223                  &   + MIN( zflxv1(ji,jj) , 0._wp ) - MAX( zflxv1(ji,  jj-1) , 0._wp) 
224               !
225               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
226               zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)
227               !
228               IF( zdep1 > zdep2 ) THEN
229                  wdmask(ji, jj) = 0._wp
230                  zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
231                  !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
232                  ! flag if the limiter has been used but stop flagging if the only
233                  ! changes have zeroed the coefficient since further iterations will
234                  ! not change anything
235                  IF( zcoef > 0._wp ) THEN   ;   jflag = 1 
236                  ELSE                       ;   zcoef = 0._wp
237                  ENDIF
238                  IF( jk1 > nn_wdit )   zcoef = 0._wp
239                  IF( zflxu1(ji  ,jj  ) > 0._wp )   zwdlmtu(ji  ,jj  ) = zcoef
240                  IF( zflxu1(ji-1,jj  ) < 0._wp )   zwdlmtu(ji-1,jj  ) = zcoef
241                  IF( zflxv1(ji  ,jj  ) > 0._wp )   zwdlmtv(ji  ,jj  ) = zcoef
242                  IF( zflxv1(ji  ,jj-1) < 0._wp )   zwdlmtv(ji  ,jj-1) = zcoef
243               ENDIF
244            END DO
245         END DO
246         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. )
247         !
248         CALL mpp_max('wet_dry', jflag)   !max over the global domain
249         !
250         IF( jflag == 0 )   EXIT
251         !
252      END DO  ! jk1 loop
253      !
254      DO jk = 1, jpkm1
255         un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:) 
256         vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:) 
257      END DO
258      un_b(:,:) = un_b(:,:) * zwdlmtu(:, :)
259      vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :)
260      !
261!!gm TO BE SUPPRESSED ?  these lbc_lnk are useless since zwdlmtu and zwdlmtv are defined everywhere !
262      CALL lbc_lnk_multi( 'wet_dry', un  , 'U', -1., vn  , 'V', -1. )
263      CALL lbc_lnk_multi( 'wet_dry', un_b, 'U', -1., vn_b, 'V', -1. )
264!!gm
265      !
266      IF(jflag == 1 .AND. lwp)   WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
267      !
268      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
269      !
270      IF( ln_timing )   CALL timing_stop('wad_lmt')      !
271      !
272   END SUBROUTINE wad_lmt
273
274
275   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
276      !!----------------------------------------------------------------------
277      !!                  ***  ROUTINE wad_lmt  ***
278      !!                   
279      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
280      !!
281      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
282      !!
283      !! ** Action  : - calculate flux limiter and W/D flag
284      !!----------------------------------------------------------------------
285      REAL(wp)                , INTENT(in   ) ::   rdtbt    ! ocean time-step index
286      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
287      !
288      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
289      INTEGER  ::   jflag               ! local integer
290      REAL(wp) ::   z2dt
291      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
292      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
293      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
294      REAL(wp) ::   ztmp                ! local scalars
295      REAL(wp), DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
296      REAL(wp), DIMENSION(jpi,jpj) ::   zflxp,  zflxn            ! local 2D workspace
297      REAL(wp), DIMENSION(jpi,jpj) ::   zflxu1, zflxv1           ! local 2D workspace
298      !!----------------------------------------------------------------------
299      IF( ln_timing )   CALL timing_start('wad_lmt_bt')      !
300      !
301      jflag  = 0
302      zdepwd = 50._wp   ! maximum depth that ocean cells can have W/D processes
303      !
304      z2dt = rdtbt   
305      !
306      zflxp(:,:)   = 0._wp
307      zflxn(:,:)   = 0._wp
308      zwdlmtu(:,:) = 1._wp
309      zwdlmtv(:,:) = 1._wp
310      !
311      DO jj = 2, jpj      ! Horizontal Flux in u and v direction
312         DO ji = 2, jpi 
313            !
314            IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
315            IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
316            !
317            zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
318               &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
319            zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
320               &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
321            !
322            zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
323            IF( zdep2 <= 0._wp ) THEN  !add more safety, but not necessary
324              sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
325              IF( zflxu(ji  ,jj  ) > 0._wp)   zwdlmtu(ji  ,jj  ) = 0._wp
326              IF( zflxu(ji-1,jj  ) < 0._wp)   zwdlmtu(ji-1,jj  ) = 0._wp
327              IF( zflxv(ji  ,jj  ) > 0._wp)   zwdlmtv(ji  ,jj  ) = 0._wp
328              IF( zflxv(ji  ,jj-1) < 0._wp)   zwdlmtv(ji  ,jj-1) = 0._wp 
329            ENDIF
330         END DO
331      END DO
332      !
333      DO jk1 = 1, nn_wdit + 1      !! start limiter iterations
334         !
335         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
336         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
337         jflag = 0     ! flag indicating if any further iterations are needed
338         !
339         DO jj = 2, jpj
340            DO ji = 2, jpi 
341               !
342               IF( tmask(ji, jj, 1 ) < 0.5_wp )   CYCLE
343               IF( ht_0(ji,jj)       > zdepwd )   CYCLE
344               !
345               ztmp = e1e2t(ji,jj)
346               !
347               zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp)   &
348                  &   + max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
349               zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp)   &
350                  &   + min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
351         
352               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
353               zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
354         
355               IF(zdep1 > zdep2) THEN
356                 zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
357                 !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
358                 ! flag if the limiter has been used but stop flagging if the only
359                 ! changes have zeroed the coefficient since further iterations will
360                 ! not change anything
361                 IF( zcoef > 0._wp ) THEN
362                    jflag = 1 
363                 ELSE
364                    zcoef = 0._wp
365                 ENDIF
366                 IF(jk1 > nn_wdit) zcoef = 0._wp
367                 IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
368                 IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
369                 IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
370                 IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
371               END IF
372            END DO ! ji loop
373         END DO  ! jj loop
374         !
375         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. )
376         !
377         CALL mpp_max('wet_dry', jflag)   !max over the global domain
378         !
379         IF(jflag == 0)   EXIT
380         !
381      END DO  ! jk1 loop
382      !
383      zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
384      zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
385      !
386!!gm THIS lbc_lnk is useless since it is already done at the end of the jk1-loop
387      CALL lbc_lnk_multi( 'wet_dry', zflxu, 'U', -1., zflxv, 'V', -1. )
388!!gm end
389      !
390      IF( jflag == 1 .AND. lwp )   WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
391      !
392      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
393      !
394      IF( ln_timing )   CALL timing_stop('wad_lmt_bt')      !
395      !
396   END SUBROUTINE wad_lmt_bt
397
398   !!==============================================================================
399END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.