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 @ 8870

Last change on this file since 8870 was 8870, checked in by deazer, 6 years ago

Changed WAD option names to Iterative and Directional
Removed old Diagnostics
Updated Domain CFG to allow domain generation with ref height for wad cases
Cleaned up TEST_CASES/cfg.txt file (need to not include WAD2 etc)
TEST caaes run ok
SETTE runs OK
AMM15 5 level runs OK

File size: 20.0 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
37                                                                     !  (can include negative depths)
38   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdramp, wdrampu, wdrampv ! for hpg limiting
39
40   LOGICAL,  PUBLIC  ::   ln_wd_il    !: Wetting/drying il activation switch (T:on,F:off)
41   LOGICAL,  PUBLIC  ::   ln_wd_dl    !: Wetting/drying dl activation switch (T:on,F:off)
42   REAL(wp), PUBLIC  ::   rn_wdmin0   !: depth at which wetting/drying starts
43   REAL(wp), PUBLIC  ::   rn_wdmin1   !: 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
46                                           !: will be considered
47   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
48   LOGICAL,  PUBLIC  ::   ln_wd_dl_bc !: DL scheme: True implies 3D velocities are set to the barotropic values at points
49                                      !: where the flow is from wet points on less than half the barotropic sub-steps 
50   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)     
51   REAL(wp), PUBLIC  ::   rn_ssh_ref  !: height of z=0 with respect to the geoid;
52   REAL(wp), PUBLIC  ::   rn_ht_0     !: the height at which ht_0 = 0   
53
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, ln_wd_dl, rn_wdmin1, rn_wdmin2, rn_wdld, nn_wdit, ln_wd_dl_bc, &
72                     &  rn_ht_0,ln_wd_dl_rmp
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      !
105      IF(ln_wd_il .OR. ln_wd_dl) THEN
106         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
107         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
108         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
109      ENDIF
110      !
111   END SUBROUTINE wad_init
112
113
114   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
115      !!----------------------------------------------------------------------
116      !!                  ***  ROUTINE wad_lmt  ***
117      !!                   
118      !! ** Purpose :   generate flux limiters for wetting/drying
119      !!
120      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
121      !!
122      !! ** Action  : - calculate flux limiter and W/D flag
123      !!----------------------------------------------------------------------
124      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1
125      REAL(wp), DIMENSION(:,:), INTENT(in)    ::   sshemp
126      REAL(wp), INTENT(in) :: z2dt
127      !
128      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
129      INTEGER  ::   jflag               ! local scalar
130      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
131      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
132      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
133      REAL(wp) ::   ztmp                ! local scalars
134      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
135      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
136      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu,  zflxv            ! local 2D workspace
137      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
138      !!----------------------------------------------------------------------
139      !
140
141      IF( nn_timing == 1 )  CALL timing_start('wad_lmt')
142
143      IF(ln_wd_il) THEN
144
145        CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
146        CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
147        !
148       
149        !IF(lwp) WRITE(numout,*)
150        !IF(lwp) WRITE(numout,*) 'wad_lmt : wetting/drying limiters and velocity limiting'
151       
152        jflag  = 0
153        zdepwd = 50._wp   !maximum depth on which that W/D could possibly happen
154
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        ! Horizontal Flux in u and v direction
165        DO jk = 1, jpkm1 
166           DO jj = 1, jpj
167              DO ji = 1, jpi
168                 zflxu(ji,jj) = zflxu(ji,jj) + e3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
169                 zflxv(ji,jj) = zflxv(ji,jj) + e3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
170              END DO 
171           END DO 
172        END DO
173       
174        zflxu(:,:) = zflxu(:,:) * e2u(:,:)
175        zflxv(:,:) = zflxv(:,:) * e1v(:,:)
176       
177        wdmask(:,:) = 1
178        DO jj = 2, jpj
179           DO ji = 2, jpi 
180
181             IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE   ! we don't care about land cells
182             IF( ht_0(ji,jj)-rn_ssh_ref > zdepwd )      CYCLE   ! and cells which are unlikely to dry
183
184              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
185                           & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
186              zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
187                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
188
189              zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1
190              IF(zdep2 .le. 0._wp) THEN  !add more safty, but not necessary
191                sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
192                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
193                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
194                IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
195                IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
196                wdmask(ji,jj) = 0._wp
197              END IF
198           ENDDO
199        END DO
200
201
202! slwa
203! HPG limiter from jholt
204      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
205!      write(6,*)'wdramp ',wdramp(10,10),wdramp(10,11)
206!jth assume don't need a lbc_lnk here
207       DO jj = 1, jpjm1
208         DO ji = 1, jpim1
209           wdrampu(ji,jj) = min(wdramp(ji,jj),wdramp(ji+1,jj))
210           wdrampv(ji,jj) = min(wdramp(ji,jj),wdramp(ji,jj+1))
211         ENDDO
212       ENDDO
213        !wdrampu(:,:)=1.0_wp
214        !wdrampv(:,:)=1.0_wp
215! end HPG limiter
216
217
218     
219        !! start limiter iterations
220        DO jk1 = 1, nn_wdit + 1
221       
222         
223           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
224           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
225           jflag = 0     ! flag indicating if any further iterations are needed
226         
227           DO jj = 2, jpj
228              DO ji = 2, jpi 
229       
230                 IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE
231                 IF( ht_0(ji,jj) > zdepwd )      CYCLE
232       
233                 ztmp = e1e2t(ji,jj)
234
235                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
236                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
237                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
238                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
239         
240                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
241                 zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)
242         
243                 IF( zdep1 > zdep2 ) THEN
244                   wdmask(ji, jj) = 0
245                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
246                   !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
247                   ! flag if the limiter has been used but stop flagging if the only
248                   ! changes have zeroed the coefficient since further iterations will
249                   ! not change anything
250                   IF( zcoef > 0._wp ) THEN
251                      jflag = 1 
252                   ELSE
253                      zcoef = 0._wp
254                   ENDIF
255                   IF(jk1 > nn_wdit) zcoef = 0._wp
256                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
257                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
258                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
259                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
260                 END IF
261              END DO ! ji loop
262           END DO  ! jj loop
263
264           CALL lbc_lnk( zwdlmtu, 'U', 1. )
265           CALL lbc_lnk( zwdlmtv, 'V', 1. )
266
267           IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
268
269           IF(jflag == 0) EXIT
270         
271        END DO  ! jk1 loop
272       
273        DO jk = 1, jpkm1
274          un(:,:,jk) = un(:,:,jk) * zwdlmtu(:, :) 
275          vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:, :) 
276        END DO
277
278        CALL lbc_lnk( un, 'U', -1. )
279        CALL lbc_lnk( vn, 'V', -1. )
280      !
281        un_b(:,:) = un_b(:,:) * zwdlmtu(:, :)
282        vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :)
283        CALL lbc_lnk( un_b, 'U', -1. )
284        CALL lbc_lnk( vn_b, 'V', -1. )
285       
286        IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
287       
288        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
289        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
290        !
291        !
292        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
293        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
294        !
295      ENDIF
296      !
297      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
298      !
299   END SUBROUTINE wad_lmt
300
301
302   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
303      !!----------------------------------------------------------------------
304      !!                  ***  ROUTINE wad_lmt  ***
305      !!                   
306      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
307      !!
308      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
309      !!
310      !! ** Action  : - calculate flux limiter and W/D flag
311      !!----------------------------------------------------------------------
312      REAL(wp), INTENT(in)    ::   rdtbt    ! ocean time-step index
313      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
314      !
315      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
316      INTEGER  ::   jflag               ! local scalar
317      REAL(wp) ::   z2dt
318      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
319      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
320      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
321      REAL(wp) ::   ztmp                ! local scalars
322      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
323      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
324      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
325      !!----------------------------------------------------------------------
326      !
327      IF( nn_timing == 1 )  CALL timing_start('wad_lmt_bt')
328
329      IF(ln_wd_il) THEN
330
331        CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
332        CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
333        !
334       
335        !IF(lwp) WRITE(numout,*)
336        !IF(lwp) WRITE(numout,*) 'wad_lmt_bt : wetting/drying limiters and velocity limiting'
337       
338        jflag  = 0
339        zdepwd = 50._wp   !maximum depth that ocean cells can have W/D processes
340
341        z2dt = rdtbt   
342       
343        zflxp(:,:)   = 0._wp
344        zflxn(:,:)   = 0._wp
345
346        zwdlmtu(:,:)  = 1._wp
347        zwdlmtv(:,:)  = 1._wp
348       
349        ! Horizontal Flux in u and v direction
350       
351        DO jj = 2, jpj
352           DO ji = 2, jpi 
353
354             IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
355             IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
356
357              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
358                           & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
359              zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
360                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
361
362              zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
363              IF(zdep2 .le. 0._wp) THEN  !add more safety, but not necessary
364                sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
365                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
366                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
367                IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
368                IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
369              END IF
370           ENDDO
371        END DO
372
373     
374        !! start limiter iterations
375        DO jk1 = 1, nn_wdit + 1
376       
377         
378           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
379           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
380           jflag = 0     ! flag indicating if any further iterations are needed
381         
382           DO jj = 2, jpj
383              DO ji = 2, jpi 
384       
385                 IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE
386                 IF( ht_0(ji,jj) > zdepwd )      CYCLE
387       
388                 ztmp = e1e2t(ji,jj)
389
390                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
391                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
392                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
393                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
394         
395                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
396                 zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
397         
398                 IF(zdep1 > zdep2) THEN
399                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
400                   !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
401                   ! flag if the limiter has been used but stop flagging if the only
402                   ! changes have zeroed the coefficient since further iterations will
403                   ! not change anything
404                   IF( zcoef > 0._wp ) THEN
405                      jflag = 1 
406                   ELSE
407                      zcoef = 0._wp
408                   ENDIF
409                   IF(jk1 > nn_wdit) zcoef = 0._wp
410                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
411                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
412                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
413                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
414                 END IF
415              END DO ! ji loop
416           END DO  ! jj loop
417
418           CALL lbc_lnk( zwdlmtu, 'U', 1. )
419           CALL lbc_lnk( zwdlmtv, 'V', 1. )
420
421           IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
422
423           IF(jflag == 0) EXIT
424         
425        END DO  ! jk1 loop
426       
427        zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
428        zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
429
430        CALL lbc_lnk( zflxu, 'U', -1. )
431        CALL lbc_lnk( zflxv, 'V', -1. )
432       
433        IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
434       
435        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
436        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
437        !
438        !
439        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
440        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
441        !
442      END IF
443
444      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
445   END SUBROUTINE wad_lmt_bt
446
447   !!==============================================================================
448END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.