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

source: NEMO/trunk/src/OCE/DYN/wet_dry.F90 @ 15102

Last change on this file since 15102 was 14834, checked in by hadcv, 3 years ago

#2600: Merge in dev_r14273_HPC-02_Daley_Tiling

  • 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   !! * Substitutions
34#  include "do_loop_substitute.h90"
35#  include "domzgr_substitute.h90"
36   !!----------------------------------------------------------------------
37   !! critical depths,filters, limiters,and masks for  Wetting and Drying
38   !! ---------------------------------------------------------------------
39
40   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdmask   !: u- and v- limiter
41   !                                                           !  (can include negative depths)
42   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdramp, wdrampu, wdrampv !: for hpg limiting
43
44   LOGICAL,  PUBLIC  ::   ln_wd_il    !: Wetting/drying il activation switch (T:on,F:off)
45   LOGICAL,  PUBLIC  ::   ln_wd_dl    !: Wetting/drying dl activation switch (T:on,F:off)
46   REAL(wp), PUBLIC  ::   rn_wdmin0   !: depth at which wetting/drying starts
47   REAL(wp), PUBLIC  ::   rn_wdmin1   !: minimum water depth on dried cells
48   REAL(wp), PUBLIC  ::   r_rn_wdmin1 !: 1/minimum water depth on dried cells
49   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells
50   REAL(wp), PUBLIC  ::   rn_wd_sbcdep   !: Depth at which to taper sbc fluxes
51   REAL(wp), PUBLIC  ::   rn_wd_sbcfra   !: Fraction of SBC at taper depth
52   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying will be considered
53   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
54   LOGICAL,  PUBLIC  ::   ln_wd_dl_bc !: DL scheme: True implies 3D velocities are set to the barotropic values at points
55                                      !: where the flow is from wet points on less than half the barotropic sub-steps 
56   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)     
57   REAL(wp), PUBLIC  ::   ssh_ref     !: height of z=0 with respect to the geoid;
58
59   LOGICAL,  PUBLIC  ::   ll_wd = .FALSE. !: Wetting/drying activation switch (ln_wd_il or ln_wd_dl) <- default def if wad_init not called
60
61   PUBLIC   wad_init                  ! initialisation routine called by step.F90
62   PUBLIC   wad_lmt                   ! routine called by sshwzv.F90
63   PUBLIC   wad_lmt_bt                ! routine called by dynspg_ts.F90
64
65   !! * Substitutions
66   !!----------------------------------------------------------------------
67CONTAINS
68
69   SUBROUTINE wad_init
70      !!----------------------------------------------------------------------
71      !!                     ***  ROUTINE wad_init  ***
72      !!                   
73      !! ** Purpose :   read wetting and drying namelist and print the variables.
74      !!
75      !! ** input   : - namwad namelist
76      !!----------------------------------------------------------------------
77      INTEGER  ::   ios, ierr   ! Local integer
78      !!
79      NAMELIST/namwad/ ln_wd_il, ln_wd_dl   , rn_wdmin0, rn_wdmin1, rn_wdmin2, rn_wdld,   &
80         &             nn_wdit , ln_wd_dl_bc, ln_wd_dl_rmp, rn_wd_sbcdep,rn_wd_sbcfra
81      !!----------------------------------------------------------------------
82      !
83      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
84905   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namwad in reference namelist' ) 
85      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
86906   IF( ios >  0 )   CALL ctl_nam ( ios , 'namwad in configuration namelist' )
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      IF( ln_wd_il .OR. ln_wd_dl ) THEN
114         ll_wd = .TRUE.
115         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
116         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
117         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
118      ENDIF
119
120      IF( ln_tile .AND. ln_wd_il ) CALL ctl_warn('Tiling has not been tested with ln_wd_il = T')
121      !
122   END SUBROUTINE wad_init
123
124
125   SUBROUTINE wad_lmt( psshb1, psshemp, z2dt, Kmm, puu, pvv )
126      !!----------------------------------------------------------------------
127      !!                  ***  ROUTINE wad_lmt  ***
128      !!                   
129      !! ** Purpose :   generate flux limiters for wetting/drying
130      !!
131      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
132      !!
133      !! ** Action  : - calculate flux limiter and W/D flag
134      !!----------------------------------------------------------------------
135      REAL(wp), DIMENSION(:,:)            , INTENT(inout) ::   psshb1
136      REAL(wp), DIMENSION(:,:)            , INTENT(in   ) ::   psshemp
137      REAL(wp)                            , INTENT(in   ) ::   z2dt
138      INTEGER                             , INTENT(in   ) ::   Kmm       ! time level index
139      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv  ! velocity arrays
140      !
141      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
142      INTEGER  ::   jflag               ! local scalar
143      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
144      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
145      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
146      REAL(wp) ::   ztmp                ! local scalars
147      REAL(wp),  DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv   ! W/D flux limiters
148      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxp  ,  zflxn    ! local 2D workspace
149      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu  ,  zflxv    ! local 2D workspace
150      REAL(wp),  DIMENSION(jpi,jpj) ::   zflxu1 , zflxv1    ! local 2D workspace
151      !!----------------------------------------------------------------------
152      IF( ln_timing )   CALL timing_start('wad_lmt')      !
153      !
154      DO jk = 1, jpkm1
155         puu(:,:,jk,Kmm) = puu(:,:,jk,Kmm) * zwdlmtu(:,:) 
156         pvv(:,:,jk,Kmm) = pvv(:,:,jk,Kmm) * zwdlmtv(:,:) 
157      END DO
158      jflag  = 0
159      zdepwd = 50._wp      ! maximum depth on which that W/D could possibly happen
160      !
161      zflxp(:,:)   = 0._wp
162      zflxn(:,:)   = 0._wp
163      zflxu(:,:)   = 0._wp
164      zflxv(:,:)   = 0._wp
165      !
166      zwdlmtu(:,:) = 1._wp
167      zwdlmtv(:,:) = 1._wp
168      !
169      DO jk = 1, jpkm1     ! Horizontal Flux in u and v direction
170         zflxu(:,:) = zflxu(:,:) + e3u(:,:,jk,Kmm) * puu(:,:,jk,Kmm) * umask(:,:,jk)
171         zflxv(:,:) = zflxv(:,:) + e3v(:,:,jk,Kmm) * pvv(:,:,jk,Kmm) * vmask(:,:,jk)
172      END DO
173      zflxu(:,:) = zflxu(:,:) * e2u(:,:)
174      zflxv(:,:) = zflxv(:,:) * e1v(:,:)
175      !
176      wdmask(:,:) = 1._wp
177      DO_2D( 0, 1, 0, 1 )
178         !
179         IF( tmask(ji,jj,1)        < 0.5_wp )   CYCLE    ! we don't care about land cells
180         IF( ht_0(ji,jj) - ssh_ref > zdepwd )   CYCLE    ! and cells which are unlikely to dry
181         !
182         zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
183            &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
184         zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
185            &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
186         !
187         zdep2 = ht_0(ji,jj) + psshb1(ji,jj) - rn_wdmin1
188         IF( zdep2 <= 0._wp ) THEN     ! add more safty, but not necessary
189            psshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
190            IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
191            IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
192            IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
193            IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
194            wdmask(ji,jj) = 0._wp
195         END IF
196      END_2D
197      !
198      !           ! HPG limiter from jholt
199      wdramp(:,:) = min((ht_0(:,:) + psshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
200      !jth assume don't need a lbc_lnk here
201      DO_2D( 1, 0, 1, 0 )
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_2D
205      !           ! end HPG limiter
206      !
207      !
208      DO jk1 = 1, nn_wdit + 1      !==  start limiter iterations  ==!
209         !
210         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
211         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
212         jflag = 0     ! flag indicating if any further iterations are needed
213         !
214         DO_2D( 0, 1, 0, 1 )
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) + psshb1(ji,jj) - rn_wdmin1 - z2dt * psshemp(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_2D
245         CALL lbc_lnk( 'wet_dry', zwdlmtu, 'U', 1.0_wp, zwdlmtv, 'V', 1.0_wp )
246         !
247         CALL mpp_max('wet_dry', jflag)   !max over the global domain
248         !
249         IF( jflag == 0 )   EXIT
250         !
251      END DO  ! jk1 loop
252      !
253      DO jk = 1, jpkm1
254         puu(:,:,jk,Kmm) = puu(:,:,jk,Kmm) * zwdlmtu(:,:) 
255         pvv(:,:,jk,Kmm) = pvv(:,:,jk,Kmm) * zwdlmtv(:,:) 
256      END DO
257      uu_b(:,:,Kmm) = uu_b(:,:,Kmm) * zwdlmtu(:, :)
258      vv_b(:,:,Kmm) = vv_b(:,:,Kmm) * zwdlmtv(:, :)
259      !
260!!gm TO BE SUPPRESSED ?  these lbc_lnk are useless since zwdlmtu and zwdlmtv are defined everywhere !
261      CALL lbc_lnk( 'wet_dry', puu(:,:,:,Kmm)  , 'U', -1.0_wp, pvv(:,:,:,Kmm)  , 'V', -1.0_wp )
262      CALL lbc_lnk( 'wet_dry', uu_b(:,:,Kmm), 'U', -1.0_wp, vv_b(:,:,Kmm), 'V', -1.0_wp )
263!!gm
264      !
265      IF(jflag == 1 .AND. lwp)   WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
266      !
267      !IF( ln_rnf      )   CALL sbc_rnf_div( hdiv )          ! runoffs (update hdiv field)
268      !
269      IF( ln_timing )   CALL timing_stop('wad_lmt')      !
270      !
271   END SUBROUTINE wad_lmt
272
273
274   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rDt_e )
275      !!----------------------------------------------------------------------
276      !!                  ***  ROUTINE wad_lmt  ***
277      !!                   
278      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
279      !!
280      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
281      !!
282      !! ** Action  : - calculate flux limiter and W/D flag
283      !!----------------------------------------------------------------------
284      REAL(wp)                , INTENT(in   ) ::   rDt_e    ! ocean time-step index
285      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
286      !
287      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
288      INTEGER  ::   jflag               ! local integer
289      REAL(wp) ::   z2dt
290      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
291      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
292      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
293      REAL(wp) ::   ztmp                ! local scalars
294      REAL(wp), DIMENSION(jpi,jpj) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
295      REAL(wp), DIMENSION(jpi,jpj) ::   zflxp,  zflxn            ! local 2D workspace
296      REAL(wp), DIMENSION(jpi,jpj) ::   zflxu1, zflxv1           ! local 2D workspace
297      !!----------------------------------------------------------------------
298      IF( ln_timing )   CALL timing_start('wad_lmt_bt')      !
299      !
300      jflag  = 0
301      zdepwd = 50._wp   ! maximum depth that ocean cells can have W/D processes
302      !
303      z2dt = rDt_e   
304      !
305      zflxp(:,:)   = 0._wp
306      zflxn(:,:)   = 0._wp
307      zwdlmtu(:,:) = 1._wp
308      zwdlmtv(:,:) = 1._wp
309      !
310      DO_2D( 0, 1, 0, 1 )      ! Horizontal Flux in u and v direction
311         !
312         IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
313         IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
314         !
315         zflxp(ji,jj) = MAX( zflxu(ji,jj) , 0._wp ) - MIN( zflxu(ji-1,jj  ) , 0._wp )   &
316            &         + MAX( zflxv(ji,jj) , 0._wp ) - MIN( zflxv(ji,  jj-1) , 0._wp ) 
317         zflxn(ji,jj) = MIN( zflxu(ji,jj) , 0._wp ) - MAX( zflxu(ji-1,jj  ) , 0._wp )   &
318            &         + MIN( zflxv(ji,jj) , 0._wp ) - MAX( zflxv(ji,  jj-1) , 0._wp ) 
319         !
320         zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
321         IF( zdep2 <= 0._wp ) THEN  !add more safety, but not necessary
322           sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
323           IF( zflxu(ji  ,jj  ) > 0._wp)   zwdlmtu(ji  ,jj  ) = 0._wp
324           IF( zflxu(ji-1,jj  ) < 0._wp)   zwdlmtu(ji-1,jj  ) = 0._wp
325           IF( zflxv(ji  ,jj  ) > 0._wp)   zwdlmtv(ji  ,jj  ) = 0._wp
326           IF( zflxv(ji  ,jj-1) < 0._wp)   zwdlmtv(ji  ,jj-1) = 0._wp 
327         ENDIF
328      END_2D
329      !
330      DO jk1 = 1, nn_wdit + 1      !! start limiter iterations
331         !
332         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
333         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
334         jflag = 0     ! flag indicating if any further iterations are needed
335         !
336         DO_2D( 0, 1, 0, 1 )
337            !
338            IF( tmask(ji, jj, 1 ) < 0.5_wp )   CYCLE
339            IF( ht_0(ji,jj)       > zdepwd )   CYCLE
340            !
341            ztmp = e1e2t(ji,jj)
342            !
343            zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp)   &
344               &   + max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
345            zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp)   &
346               &   + min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
347       
348            zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
349            zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
350       
351            IF(zdep1 > zdep2) THEN
352              zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
353              !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
354              ! flag if the limiter has been used but stop flagging if the only
355              ! changes have zeroed the coefficient since further iterations will
356              ! not change anything
357              IF( zcoef > 0._wp ) THEN
358                 jflag = 1 
359              ELSE
360                 zcoef = 0._wp
361              ENDIF
362              IF(jk1 > nn_wdit) zcoef = 0._wp
363              IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
364              IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
365              IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
366              IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
367            END IF
368         END_2D
369         !
370         CALL lbc_lnk( 'wet_dry', zwdlmtu, 'U', 1.0_wp, zwdlmtv, 'V', 1.0_wp )
371         !
372         CALL mpp_max('wet_dry', jflag)   !max over the global domain
373         !
374         IF(jflag == 0)   EXIT
375         !
376      END DO  ! jk1 loop
377      !
378      zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
379      zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
380      !
381!!gm THIS lbc_lnk is useless since it is already done at the end of the jk1-loop
382      CALL lbc_lnk( 'wet_dry', zflxu, 'U', -1.0_wp, zflxv, 'V', -1.0_wp )
383!!gm end
384      !
385      IF( jflag == 1 .AND. lwp )   WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
386      !
387      !IF( ln_rnf      )   CALL sbc_rnf_div( hdiv )          ! runoffs (update hdiv field)
388      !
389      IF( ln_timing )   CALL timing_stop('wad_lmt_bt')      !
390      !
391   END SUBROUTINE wad_lmt_bt
392
393   !!==============================================================================
394END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.