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 branches/UKMO/ROMS_WAD_7832/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/UKMO/ROMS_WAD_7832/NEMOGCM/NEMO/OPA_SRC/DYN/wet_dry.F90 @ 8977

Last change on this file since 8977 was 8977, checked in by deazer, 7 years ago

Addresses several issues in the review except for rn_ssh_ref in the TEST cases

Builds and runs ok , little extra bracket dealt with in stpctl also

File size: 19.5 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_lmt    : Compute the horizontal flux limiter and the limited velocity
17   !!                when wetting and drying happens
18   !!----------------------------------------------------------------------
19   USE oce             ! ocean dynamics and tracers
20   USE dom_oce         ! ocean space and time domain
21   USE sbc_oce, ONLY : ln_rnf   ! surface boundary condition: ocean
22   USE sbcrnf          ! river runoff
23   USE in_out_manager  ! I/O manager
24   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
25   USE lib_mpp         ! MPP library
26   USE wrk_nemo        ! Memory Allocation
27   USE timing          ! Timing
28
29   IMPLICIT NONE
30   PRIVATE
31
32   !!----------------------------------------------------------------------
33   !! critical depths,filters, limiters,and masks for  Wetting and Drying
34   !! ---------------------------------------------------------------------
35
36   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdmask         !: u- and v- limiter  (can include negative depths)
37   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdramp, wdrampu, wdrampv !: for hpg limiting
38
39   LOGICAL,  PUBLIC  ::   ln_wd_il    !: Wetting/drying il activation switch (T:on,F:off)
40   LOGICAL,  PUBLIC  ::   ln_wd_dl    !: Wetting/drying dl activation switch (T:on,F:off)
41   REAL(wp), PUBLIC  ::   rn_wdmin0   !: depth at which wetting/drying starts
42   REAL(wp), PUBLIC  ::   rn_wdmin1   !: minimum water depth on dried cells
43   REAL(wp), PUBLIC  ::   r_rn_wdmin1 !: 1/minimum water depth on dried cells
44   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells
45   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying will be considered
46   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
47   LOGICAL,  PUBLIC  ::   ln_wd_dl_bc !: DL scheme: True implies 3D velocities are set to the barotropic values at points
48                                      !: where the flow is from wet points on less than half the barotropic sub-steps 
49   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)     
50   REAL(wp), PUBLIC  ::   rn_ssh_ref  !: height of z=0 with respect to the geoid;
51   REAL(wp), PUBLIC  ::   rn_ht_0     !: the height at which ht_0 = 0   
52
53   LOGICAL,  PUBLIC  ::   ll_wd       !: Wetting/drying activation switch if either ln_wd_il or ln_wd_dl
54
55   PUBLIC   wad_init                  ! initialisation routine called by step.F90
56   PUBLIC   wad_lmt                   ! routine called by sshwzv.F90
57   PUBLIC   wad_lmt_bt                ! routine called by dynspg_ts.F90
58
59   !! * Substitutions
60#  include "vectopt_loop_substitute.h90"
61CONTAINS
62
63   SUBROUTINE wad_init
64      !!----------------------------------------------------------------------
65      !!                     ***  ROUTINE wad_init  ***
66      !!                   
67      !! ** Purpose :   read wetting and drying namelist and print the variables.
68      !!
69      !! ** input   : - namwad namelist
70      !!----------------------------------------------------------------------
71      NAMELIST/namwad/ ln_wd_il, ln_wd_dl, rn_wdmin0, rn_wdmin1, rn_wdmin2, rn_wdld, &
72                     & nn_wdit, ln_wd_dl_bc, ln_wd_dl_rmp, rn_ht_0
73
74      INTEGER  ::   ios                 ! Local integer output status for namelist read
75      INTEGER  ::   ierr                ! Local integer status array allocation
76      !!----------------------------------------------------------------------
77
78      REWIND( numnam_ref )              ! Namelist namwad in reference namelist
79                                        ! : Parameters for Wetting/Drying
80      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
81905   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in reference namelist', .TRUE.) 
82      REWIND( numnam_cfg )              ! Namelist namwad in configuration namelist
83                                        ! : Parameters for Wetting/Drying
84      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
85906   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in configuration namelist', .TRUE. )
86      IF(lwm) WRITE ( numond, namwad )
87
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,*) '      the height (z) at which ht_0=0:  rn_ht_0      = ', rn_ht_0 
103      ENDIF
104      r_rn_wdmin1=1/rn_wdmin1
105      ll_wd = .FALSE.
106      IF(ln_wd_il .OR. ln_wd_dl) THEN
107         ll_wd = .TRUE.
108         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
109         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
110         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
111      ENDIF
112      !
113   END SUBROUTINE wad_init
114
115
116   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
117      !!----------------------------------------------------------------------
118      !!                  ***  ROUTINE wad_lmt  ***
119      !!                   
120      !! ** Purpose :   generate flux limiters for wetting/drying
121      !!
122      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
123      !!
124      !! ** Action  : - calculate flux limiter and W/D flag
125      !!----------------------------------------------------------------------
126      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1
127      REAL(wp), DIMENSION(:,:), INTENT(in)    ::   sshemp
128      REAL(wp), INTENT(in) :: z2dt
129      !
130      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
131      INTEGER  ::   jflag               ! local scalar
132      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
133      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
134      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
135      REAL(wp) ::   ztmp                ! local scalars
136      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
137      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
138      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu,  zflxv            ! local 2D workspace
139      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
140      !!----------------------------------------------------------------------
141      !
142
143      IF( nn_timing == 1 )  CALL timing_start('wad_lmt')
144
145
146      CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
147      CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
148      !
149       
150      jflag  = 0
151      zdepwd = 50._wp   !maximum depth on which that W/D could possibly happen
152
153       
154      zflxp(:,:)   = 0._wp
155      zflxn(:,:)   = 0._wp
156      zflxu(:,:)   = 0._wp
157      zflxv(:,:)   = 0._wp
158
159      zwdlmtu(:,:)  = 1._wp
160      zwdlmtv(:,:)  = 1._wp
161       
162      ! Horizontal Flux in u and v direction
163      DO jk = 1, jpkm1 
164         DO jj = 1, jpj
165            DO ji = 1, jpi
166               zflxu(ji,jj) = zflxu(ji,jj) + e3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
167               zflxv(ji,jj) = zflxv(ji,jj) + e3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
168            END DO 
169         END DO 
170      END DO
171       
172      zflxu(:,:) = zflxu(:,:) * e2u(:,:)
173      zflxv(:,:) = zflxv(:,:) * e1v(:,:)
174       
175      wdmask(:,:) = 1
176      DO jj = 2, jpj
177         DO ji = 2, jpi 
178
179            IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE              ! we don't care about land cells
180            IF( ht_0(ji,jj) - rn_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) + sshb1(ji,jj) - rn_wdmin1
188            IF(zdep2 .le. 0._wp) THEN  !add more safty, but not necessary
189               sshb1(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         ENDDO
197      END DO
198
199
200! HPG limiter from jholt
201      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
202!jth assume don't need a lbc_lnk here
203      DO jj = 1, jpjm1
204         DO ji = 1, jpim1
205            wdrampu(ji,jj) = min(wdramp(ji,jj),wdramp(ji+1,jj))
206            wdrampv(ji,jj) = min(wdramp(ji,jj),wdramp(ji,jj+1))
207         END DO
208      END DO
209! end HPG limiter
210
211
212     
213        !! start limiter iterations
214      DO jk1 = 1, nn_wdit + 1
215       
216         
217         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
218         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
219         jflag = 0     ! flag indicating if any further iterations are needed
220         
221         DO jj = 2, jpj
222            DO ji = 2, jpi 
223       
224               IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE
225               IF( ht_0(ji,jj) > zdepwd )      CYCLE
226       
227               ztmp = e1e2t(ji,jj)
228
229               zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
230                      & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
231               zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
232                      & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
233         
234               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
235               zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)
236         
237               IF( zdep1 > zdep2 ) THEN
238                  wdmask(ji, jj) = 0
239                  zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
240                  !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
241                  ! flag if the limiter has been used but stop flagging if the only
242                  ! changes have zeroed the coefficient since further iterations will
243                  ! not change anything
244                  IF( zcoef > 0._wp ) THEN
245                     jflag = 1 
246                  ELSE
247                     zcoef = 0._wp
248                  ENDIF
249                  IF(jk1 > nn_wdit) zcoef = 0._wp
250                  IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
251                  IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
252                  IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
253                  IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
254               END IF
255            END DO ! ji loop
256         END DO  ! jj loop
257
258         CALL lbc_lnk( zwdlmtu, 'U', 1. )
259         CALL lbc_lnk( zwdlmtv, 'V', 1. )
260
261         IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
262
263         IF(jflag == 0) EXIT
264         
265      END DO  ! jk1 loop
266       
267      DO jk = 1, jpkm1
268         un(:,:,jk) = un(:,:,jk) * zwdlmtu(:, :) 
269         vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:, :) 
270      END DO
271
272      CALL lbc_lnk( un, 'U', -1. )
273      CALL lbc_lnk( vn, 'V', -1. )
274        !
275      un_b(:,:) = un_b(:,:) * zwdlmtu(:, :)
276      vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :)
277      CALL lbc_lnk( un_b, 'U', -1. )
278      CALL lbc_lnk( vn_b, 'V', -1. )
279       
280      IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
281       
282      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
283      !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
284      !
285      !
286      CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
287      CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
288      !
289      !
290      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
291      !
292   END SUBROUTINE wad_lmt
293
294
295   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
296      !!----------------------------------------------------------------------
297      !!                  ***  ROUTINE wad_lmt  ***
298      !!                   
299      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
300      !!
301      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
302      !!
303      !! ** Action  : - calculate flux limiter and W/D flag
304      !!----------------------------------------------------------------------
305      REAL(wp), INTENT(in)    ::   rdtbt    ! ocean time-step index
306      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
307      !
308      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
309      INTEGER  ::   jflag               ! local scalar
310      REAL(wp) ::   z2dt
311      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
312      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
313      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
314      REAL(wp) ::   ztmp                ! local scalars
315      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
316      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
317      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
318      !!----------------------------------------------------------------------
319      !
320      IF( nn_timing == 1 )  CALL timing_start('wad_lmt_bt')
321
322
323      CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
324      CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
325      !
326       
327      !IF(lwp) WRITE(numout,*)
328      !IF(lwp) WRITE(numout,*) 'wad_lmt_bt : wetting/drying limiters and velocity limiting'
329       
330      jflag  = 0
331      zdepwd = 50._wp   !maximum depth that ocean cells can have W/D processes
332
333      z2dt = rdtbt   
334       
335      zflxp(:,:)   = 0._wp
336      zflxn(:,:)   = 0._wp
337
338      zwdlmtu(:,:)  = 1._wp
339      zwdlmtv(:,:)  = 1._wp
340       
341        ! Horizontal Flux in u and v direction
342       
343      DO jj = 2, jpj
344         DO ji = 2, jpi 
345
346           IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
347           IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
348
349            zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
350                         & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
351            zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
352                         & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
353
354            zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
355            IF(zdep2 .le. 0._wp) THEN  !add more safety, but not necessary
356              sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
357              IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
358              IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
359              IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
360              IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
361            END IF
362         ENDDO
363      END DO
364
365     
366      !! start limiter iterations
367      DO jk1 = 1, nn_wdit + 1
368       
369         
370         zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
371         zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
372         jflag = 0     ! flag indicating if any further iterations are needed
373         
374         DO jj = 2, jpj
375            DO ji = 2, jpi 
376       
377               IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE
378               IF( ht_0(ji,jj) > zdepwd )      CYCLE
379       
380               ztmp = e1e2t(ji,jj)
381
382               zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
383                      & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
384               zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
385                      & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
386         
387               zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
388               zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
389         
390               IF(zdep1 > zdep2) THEN
391                 zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
392                 !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
393                 ! flag if the limiter has been used but stop flagging if the only
394                 ! changes have zeroed the coefficient since further iterations will
395                 ! not change anything
396                 IF( zcoef > 0._wp ) THEN
397                    jflag = 1 
398                 ELSE
399                    zcoef = 0._wp
400                 ENDIF
401                 IF(jk1 > nn_wdit) zcoef = 0._wp
402                 IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
403                 IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
404                 IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
405                 IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
406               END IF
407            END DO ! ji loop
408         END DO  ! jj loop
409
410         CALL lbc_lnk( zwdlmtu, 'U', 1. )
411         CALL lbc_lnk( zwdlmtv, 'V', 1. )
412
413         IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
414
415         IF(jflag == 0) EXIT
416         
417      END DO  ! jk1 loop
418       
419      zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
420      zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
421
422      CALL lbc_lnk( zflxu, 'U', -1. )
423      CALL lbc_lnk( zflxv, 'V', -1. )
424       
425      IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
426       
427      !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
428      !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
429      !
430      !
431      CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
432      CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
433      !
434
435      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
436   END SUBROUTINE wad_lmt_bt
437
438   !!==============================================================================
439END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.