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_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/DYN/wet_dry.F90 @ 11048

Last change on this file since 11048 was 10499, checked in by deazer, 5 years ago

Fix ticket #2154

  • Property svn:keywords set to Id
File size: 19.7 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      REWIND( numnam_ref )              ! Namelist namwad in reference namelist : Parameters for Wetting/Drying
82      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
83905   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namwad in reference namelist', .TRUE.) 
84      REWIND( numnam_cfg )              ! Namelist namwad in configuration namelist : Parameters for Wetting/Drying
85      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
86906   IF( ios >  0 )   CALL ctl_nam ( ios , 'namwad in configuration namelist', .TRUE. )
87      IF(lwm) WRITE ( numond, namwad )
88      !
89      IF( rn_wd_sbcfra>=1 )   CALL ctl_stop( 'STOP', 'rn_wd_sbcfra >=1 : must be < 1' )
90      IF(lwp) THEN                  ! control print
91         WRITE(numout,*)
92         WRITE(numout,*) 'wad_init : Wetting and drying initialization through namelist read'
93         WRITE(numout,*) '~~~~~~~~'
94         WRITE(numout,*) '   Namelist namwad'
95         WRITE(numout,*) '      Logical for Iter Lim wd option   ln_wd_il     = ', ln_wd_il
96         WRITE(numout,*) '      Logical for Dir. Lim wd option   ln_wd_dl     = ', ln_wd_dl
97         WRITE(numout,*) '      Depth at which wet/drying starts rn_wdmin0    = ', rn_wdmin0
98         WRITE(numout,*) '      Minimum wet depth on dried cells rn_wdmin1    = ', rn_wdmin1
99         WRITE(numout,*) '      Tolerance of min wet depth       rn_wdmin2    = ', rn_wdmin2
100         WRITE(numout,*) '      land elevation threshold         rn_wdld      = ', rn_wdld
101         WRITE(numout,*) '      Max iteration for W/D limiter    nn_wdit      = ', nn_wdit
102         WRITE(numout,*) '      T => baroclinic u,v=0 at dry pts: ln_wd_dl_bc = ', ln_wd_dl_bc     
103         WRITE(numout,*) '      use a ramp for rwd limiter:  ln_wd_dl_rwd_rmp = ', ln_wd_dl_rmp
104         WRITE(numout,*) '      cut off depth sbc for wd   rn_wd_sbcdep       = ', rn_wd_sbcdep
105         WRITE(numout,*) '      fraction to start sbc wgt rn_wd_sbcfra        = ', rn_wd_sbcfra
106      ENDIF
107      IF( .NOT. ln_read_cfg ) THEN
108         IF(lwp) WRITE(numout,*) '      No configuration file so seting ssh_ref to zero  '
109         ssh_ref=0._wp
110      ENDIF
111
112      r_rn_wdmin1 = 1 / rn_wdmin1
113      ll_wd = .FALSE.
114      IF( ln_wd_il .OR. ln_wd_dl ) THEN
115         ll_wd = .TRUE.
116         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
117         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
118         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
119      ENDIF
120      !
121   END SUBROUTINE wad_init
122
123
124   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
125      !!----------------------------------------------------------------------
126      !!                  ***  ROUTINE wad_lmt  ***
127      !!                   
128      !! ** Purpose :   generate flux limiters for wetting/drying
129      !!
130      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
131      !!
132      !! ** Action  : - calculate flux limiter and W/D flag
133      !!----------------------------------------------------------------------
134      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1        !!gm DOCTOR names: should start with p !
135      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sshemp
136      REAL(wp)                , INTENT(in   ) ::   z2dt
137      !
138      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
139      INTEGER  ::   jflag               ! local scalar
140      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
141      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
142      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
143      REAL(wp) ::   ztmp                ! local scalars
144      REAL(wp),  DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv   ! W/D flux limiters
145      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxp  ,  zflxn    ! local 2D workspace
146      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu  ,  zflxv    ! local 2D workspace
147      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu1 , zflxv1    ! local 2D workspace
148      !!----------------------------------------------------------------------
149      IF( ln_timing )   CALL timing_start('wad_lmt')      !
150      !
151      DO jk = 1, jpkm1
152         un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:) 
153         vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:) 
154      END DO
155      jflag  = 0
156      zdepwd = 50._wp      ! maximum depth on which that W/D could possibly happen
157      !
158      zflxp(:,:)   = 0._wp
159      zflxn(:,:)   = 0._wp
160      zflxu(:,:)   = 0._wp
161      zflxv(:,:)   = 0._wp
162      !
163      zwdlmtu(:,:) = 1._wp
164      zwdlmtv(:,:) = 1._wp
165      !
166      DO jk = 1, jpkm1     ! Horizontal Flux in u and v direction
167         zflxu(:,:) = zflxu(:,:) + e3u_n(:,:,jk) * un(:,:,jk) * umask(:,:,jk)
168         zflxv(:,:) = zflxv(:,:) + e3v_n(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk)
169      END DO
170      zflxu(:,:) = zflxu(:,:) * e2u(:,:)
171      zflxv(:,:) = zflxv(:,:) * e1v(:,:)
172      !
173      wdmask(:,:) = 1._wp
174      DO jj = 2, jpj
175         DO ji = 2, jpi 
176            !
177            IF( tmask(ji,jj,1)        < 0.5_wp )   CYCLE    ! we don't care about land cells
178            IF( ht_0(ji,jj) - ssh_ref > zdepwd )   CYCLE    ! and cells which are unlikely to dry
179            !
180            zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
181               &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
182            zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
183               &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
184            !
185            zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1
186            IF( zdep2 <= 0._wp ) THEN     ! add more safty, but not necessary
187               sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
188               IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
189               IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
190               IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
191               IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
192               wdmask(ji,jj) = 0._wp
193            END IF
194         END DO
195      END DO
196      !
197      !           ! HPG limiter from jholt
198      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
199      !jth assume don't need a lbc_lnk here
200      DO jj = 1, jpjm1
201         DO ji = 1, jpim1
202            wdrampu(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji+1,jj) )
203            wdrampv(ji,jj) = MIN( wdramp(ji,jj) , wdramp(ji,jj+1) )
204         END DO
205      END DO
206      !           ! end HPG limiter
207      !
208      !
209      DO jk1 = 1, nn_wdit + 1      !==  start limiter iterations  ==!
210         !
211         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
212         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
213         jflag = 0     ! flag indicating if any further iterations are needed
214         !
215         DO jj = 2, jpj
216            DO ji = 2, jpi 
217               IF( tmask(ji, jj, 1) < 0.5_wp )   CYCLE
218               IF( ht_0(ji,jj)      > zdepwd )   CYCLE
219               !
220               ztmp = e1e2t(ji,jj)
221               !
222               zzflxp = MAX( zflxu1(ji,jj) , 0._wp ) - MIN( zflxu1(ji-1,jj  ) , 0._wp)   &
223                  &   + MAX( zflxv1(ji,jj) , 0._wp ) - MIN( zflxv1(ji,  jj-1) , 0._wp) 
224               zzflxn = MIN( zflxu1(ji,jj) , 0._wp ) - MAX( zflxu1(ji-1,jj  ) , 0._wp)   &
225                  &   + MIN( zflxv1(ji,jj) , 0._wp ) - MAX( zflxv1(ji,  jj-1) , 0._wp) 
226               !
227               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
228               zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)
229               !
230               IF( zdep1 > zdep2 ) THEN
231                  wdmask(ji, jj) = 0._wp
232                  zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
233                  !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
234                  ! flag if the limiter has been used but stop flagging if the only
235                  ! changes have zeroed the coefficient since further iterations will
236                  ! not change anything
237                  IF( zcoef > 0._wp ) THEN   ;   jflag = 1 
238                  ELSE                       ;   zcoef = 0._wp
239                  ENDIF
240                  IF( jk1 > nn_wdit )   zcoef = 0._wp
241                  IF( zflxu1(ji  ,jj  ) > 0._wp )   zwdlmtu(ji  ,jj  ) = zcoef
242                  IF( zflxu1(ji-1,jj  ) < 0._wp )   zwdlmtu(ji-1,jj  ) = zcoef
243                  IF( zflxv1(ji  ,jj  ) > 0._wp )   zwdlmtv(ji  ,jj  ) = zcoef
244                  IF( zflxv1(ji  ,jj-1) < 0._wp )   zwdlmtv(ji  ,jj-1) = zcoef
245               ENDIF
246            END DO
247         END DO
248         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. )
249         !
250         CALL mpp_max('wet_dry', jflag)   !max over the global domain
251         !
252         IF( jflag == 0 )   EXIT
253         !
254      END DO  ! jk1 loop
255      !
256      DO jk = 1, jpkm1
257         un(:,:,jk) = un(:,:,jk) * zwdlmtu(:,:) 
258         vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:,:) 
259      END DO
260      un_b(:,:) = un_b(:,:) * zwdlmtu(:, :)
261      vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :)
262      !
263!!gm TO BE SUPPRESSED ?  these lbc_lnk are useless since zwdlmtu and zwdlmtv are defined everywhere !
264      CALL lbc_lnk_multi( 'wet_dry', un  , 'U', -1., vn  , 'V', -1. )
265      CALL lbc_lnk_multi( 'wet_dry', un_b, 'U', -1., vn_b, 'V', -1. )
266!!gm
267      !
268      IF(jflag == 1 .AND. lwp)   WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
269      !
270      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
271      !
272      IF( ln_timing )   CALL timing_stop('wad_lmt')      !
273      !
274   END SUBROUTINE wad_lmt
275
276
277   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
278      !!----------------------------------------------------------------------
279      !!                  ***  ROUTINE wad_lmt  ***
280      !!                   
281      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
282      !!
283      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
284      !!
285      !! ** Action  : - calculate flux limiter and W/D flag
286      !!----------------------------------------------------------------------
287      REAL(wp)                , INTENT(in   ) ::   rdtbt    ! ocean time-step index
288      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
289      !
290      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
291      INTEGER  ::   jflag               ! local integer
292      REAL(wp) ::   z2dt
293      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
294      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
295      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
296      REAL(wp) ::   ztmp                ! local scalars
297      REAL(wp), DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
298      REAL(wp), DIMENSION(jpi,jpj) ::   zflxp,  zflxn            ! local 2D workspace
299      REAL(wp), DIMENSION(jpi,jpj) ::   zflxu1, zflxv1           ! local 2D workspace
300      !!----------------------------------------------------------------------
301      IF( ln_timing )   CALL timing_start('wad_lmt_bt')      !
302      !
303      jflag  = 0
304      zdepwd = 50._wp   ! maximum depth that ocean cells can have W/D processes
305      !
306      z2dt = rdtbt   
307      !
308      zflxp(:,:)   = 0._wp
309      zflxn(:,:)   = 0._wp
310      zwdlmtu(:,:) = 1._wp
311      zwdlmtv(:,:) = 1._wp
312      !
313      DO jj = 2, jpj      ! Horizontal Flux in u and v direction
314         DO ji = 2, jpi 
315            !
316            IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
317            IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
318            !
319            zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
320               &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
321            zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
322               &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
323            !
324            zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
325            IF( zdep2 <= 0._wp ) THEN  !add more safety, but not necessary
326              sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
327              IF( zflxu(ji  ,jj  ) > 0._wp)   zwdlmtu(ji  ,jj  ) = 0._wp
328              IF( zflxu(ji-1,jj  ) < 0._wp)   zwdlmtu(ji-1,jj  ) = 0._wp
329              IF( zflxv(ji  ,jj  ) > 0._wp)   zwdlmtv(ji  ,jj  ) = 0._wp
330              IF( zflxv(ji  ,jj-1) < 0._wp)   zwdlmtv(ji  ,jj-1) = 0._wp 
331            ENDIF
332         END DO
333      END DO
334      !
335      DO jk1 = 1, nn_wdit + 1      !! start limiter iterations
336         !
337         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
338         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
339         jflag = 0     ! flag indicating if any further iterations are needed
340         !
341         DO jj = 2, jpj
342            DO ji = 2, jpi 
343               !
344               IF( tmask(ji, jj, 1 ) < 0.5_wp )   CYCLE
345               IF( ht_0(ji,jj)       > zdepwd )   CYCLE
346               !
347               ztmp = e1e2t(ji,jj)
348               !
349               zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp)   &
350                  &   + max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
351               zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp)   &
352                  &   + min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
353         
354               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
355               zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
356         
357               IF(zdep1 > zdep2) THEN
358                 zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
359                 !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
360                 ! flag if the limiter has been used but stop flagging if the only
361                 ! changes have zeroed the coefficient since further iterations will
362                 ! not change anything
363                 IF( zcoef > 0._wp ) THEN
364                    jflag = 1 
365                 ELSE
366                    zcoef = 0._wp
367                 ENDIF
368                 IF(jk1 > nn_wdit) zcoef = 0._wp
369                 IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
370                 IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
371                 IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
372                 IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
373               END IF
374            END DO ! ji loop
375         END DO  ! jj loop
376         !
377         CALL lbc_lnk_multi( 'wet_dry', zwdlmtu, 'U', 1., zwdlmtv, 'V', 1. )
378         !
379         CALL mpp_max('wet_dry', jflag)   !max over the global domain
380         !
381         IF(jflag == 0)   EXIT
382         !
383      END DO  ! jk1 loop
384      !
385      zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
386      zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
387      !
388!!gm THIS lbc_lnk is useless since it is already done at the end of the jk1-loop
389      CALL lbc_lnk_multi( 'wet_dry', zflxu, 'U', -1., zflxv, 'V', -1. )
390!!gm end
391      !
392      IF( jflag == 1 .AND. lwp )   WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
393      !
394      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
395      !
396      IF( ln_timing )   CALL timing_stop('wad_lmt_bt')      !
397      !
398   END SUBROUTINE wad_lmt_bt
399
400   !!==============================================================================
401END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.