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

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

Modified implicit bfriciton to only be done in WAD cases
This aids bit comparison tests to show WAD code does not effect CONFIGS,
unless you want implicit bfriction.

File size: 19.9 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  ::   rn_wdmin2   !: tolerance of minimum water depth on dried cells
44   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying will be considered
45   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
46   LOGICAL,  PUBLIC  ::   ln_wd_dl_bc !: DL scheme: True implies 3D velocities are set to the barotropic values at points
47                                      !: where the flow is from wet points on less than half the barotropic sub-steps 
48   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)     
49   REAL(wp), PUBLIC  ::   rn_ssh_ref  !: height of z=0 with respect to the geoid;
50   REAL(wp), PUBLIC  ::   rn_ht_0     !: the height at which ht_0 = 0   
51
52
53   PUBLIC   wad_init                  ! initialisation routine called by step.F90
54   PUBLIC   wad_lmt                   ! routine called by sshwzv.F90
55   PUBLIC   wad_lmt_bt                ! routine called by dynspg_ts.F90
56
57   !! * Substitutions
58#  include "vectopt_loop_substitute.h90"
59CONTAINS
60
61   SUBROUTINE wad_init
62      !!----------------------------------------------------------------------
63      !!                     ***  ROUTINE wad_init  ***
64      !!                   
65      !! ** Purpose :   read wetting and drying namelist and print the variables.
66      !!
67      !! ** input   : - namwad namelist
68      !!----------------------------------------------------------------------
69      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, &
70                     &  rn_ht_0,ln_wd_dl_rmp
71
72      INTEGER  ::   ios                 ! Local integer output status for namelist read
73      INTEGER  ::   ierr                ! Local integer status array allocation
74      !!----------------------------------------------------------------------
75
76      REWIND( numnam_ref )              ! Namelist namwad in reference namelist
77                                        ! : Parameters for Wetting/Drying
78      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
79905   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in reference namelist', .TRUE.) 
80      REWIND( numnam_cfg )              ! Namelist namwad in configuration namelist
81                                        ! : Parameters for Wetting/Drying
82      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
83906   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in configuration namelist', .TRUE. )
84      IF(lwm) WRITE ( numond, namwad )
85
86      IF(lwp) THEN                  ! control print
87         WRITE(numout,*)
88         WRITE(numout,*) 'wad_init : Wetting and drying initialization through namelist read'
89         WRITE(numout,*) '~~~~~~~~'
90         WRITE(numout,*) '   Namelist namwad'
91         WRITE(numout,*) '      Logical for Iter Lim wd option ln_wd_il       = ', ln_wd_il
92         WRITE(numout,*) '      Logical for Dir. Lim wd option ln_wd_dl       = ', ln_wd_dl
93         WRITE(numout,*) '      Depth at which wet/drying starts rn_wdmin0    = ', rn_wdmin0
94         WRITE(numout,*) '      Minimum wet depth on dried cells rn_wdmin1    = ', rn_wdmin1
95         WRITE(numout,*) '      Tolerance of min wet depth     rn_wdmin2      = ', rn_wdmin2
96         WRITE(numout,*) '      land elevation threshold         rn_wdld      = ', rn_wdld
97         WRITE(numout,*) '      Max iteration for W/D limiter    nn_wdit      = ', nn_wdit
98         WRITE(numout,*) '      T => baroclinic u,v = 0 at dry pts: ln_wd_dl_bc = ', ln_wd_dl_bc     
99         WRITE(numout,*) '      use a ramp for rwd limiter:      ln_wd_dl_rwd_rmp   = ', ln_wd_dl_rmp
100         WRITE(numout,*) '      the height (z) at which ht_0 = 0:rn_ht_0      = ', rn_ht_0 
101      ENDIF
102      !
103      IF(ln_wd_il .OR. ln_wd_dl) THEN
104         ALLOCATE( wdmask(jpi,jpj),   STAT=ierr )
105         ALLOCATE( wdramp(jpi,jpj), wdrampu(jpi,jpj), wdrampv(jpi,jpj), STAT=ierr ) 
106         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
107      ENDIF
108      !
109   END SUBROUTINE wad_init
110
111
112   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
113      !!----------------------------------------------------------------------
114      !!                  ***  ROUTINE wad_lmt  ***
115      !!                   
116      !! ** Purpose :   generate flux limiters for wetting/drying
117      !!
118      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
119      !!
120      !! ** Action  : - calculate flux limiter and W/D flag
121      !!----------------------------------------------------------------------
122      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1
123      REAL(wp), DIMENSION(:,:), INTENT(in)    ::   sshemp
124      REAL(wp), INTENT(in) :: z2dt
125      !
126      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
127      INTEGER  ::   jflag               ! local scalar
128      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
129      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
130      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
131      REAL(wp) ::   ztmp                ! local scalars
132      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
133      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
134      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu,  zflxv            ! local 2D workspace
135      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
136      !!----------------------------------------------------------------------
137      !
138
139      IF( nn_timing == 1 )  CALL timing_start('wad_lmt')
140
141      IF(ln_wd_il) THEN
142
143        CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
144        CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
145        !
146       
147        !IF(lwp) WRITE(numout,*)
148        !IF(lwp) WRITE(numout,*) 'wad_lmt : wetting/drying limiters and velocity limiting'
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! slwa
201! HPG limiter from jholt
202      wdramp(:,:) = min((ht_0(:,:) + sshb1(:,:) - rn_wdmin1)/(rn_wdmin0 - rn_wdmin1),1.0_wp)
203!      write(6,*)'wdramp ',wdramp(10,10),wdramp(10,11)
204!jth assume don't need a lbc_lnk here
205       DO jj = 1, jpjm1
206         DO ji = 1, jpim1
207           wdrampu(ji,jj) = min(wdramp(ji,jj),wdramp(ji+1,jj))
208           wdrampv(ji,jj) = min(wdramp(ji,jj),wdramp(ji,jj+1))
209         ENDDO
210       ENDDO
211        !wdrampu(:,:)=1.0_wp
212        !wdrampv(:,:)=1.0_wp
213! end HPG limiter
214
215
216     
217        !! start limiter iterations
218        DO jk1 = 1, nn_wdit + 1
219       
220         
221           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
222           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
223           jflag = 0     ! flag indicating if any further iterations are needed
224         
225           DO jj = 2, jpj
226              DO ji = 2, jpi 
227       
228                 IF( tmask(ji, jj, 1) < 0.5_wp ) CYCLE
229                 IF( ht_0(ji,jj) > zdepwd )      CYCLE
230       
231                 ztmp = e1e2t(ji,jj)
232
233                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
234                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
235                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
236                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
237         
238                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
239                 zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)
240         
241                 IF( zdep1 > zdep2 ) THEN
242                   wdmask(ji, jj) = 0
243                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
244                   !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
245                   ! flag if the limiter has been used but stop flagging if the only
246                   ! changes have zeroed the coefficient since further iterations will
247                   ! not change anything
248                   IF( zcoef > 0._wp ) THEN
249                      jflag = 1 
250                   ELSE
251                      zcoef = 0._wp
252                   ENDIF
253                   IF(jk1 > nn_wdit) zcoef = 0._wp
254                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
255                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
256                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
257                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
258                 END IF
259              END DO ! ji loop
260           END DO  ! jj loop
261
262           CALL lbc_lnk( zwdlmtu, 'U', 1. )
263           CALL lbc_lnk( zwdlmtv, 'V', 1. )
264
265           IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
266
267           IF(jflag == 0) EXIT
268         
269        END DO  ! jk1 loop
270       
271        DO jk = 1, jpkm1
272          un(:,:,jk) = un(:,:,jk) * zwdlmtu(:, :) 
273          vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:, :) 
274        END DO
275
276        CALL lbc_lnk( un, 'U', -1. )
277        CALL lbc_lnk( vn, 'V', -1. )
278      !
279        un_b(:,:) = un_b(:,:) * zwdlmtu(:, :)
280        vn_b(:,:) = vn_b(:,:) * zwdlmtv(:, :)
281        CALL lbc_lnk( un_b, 'U', -1. )
282        CALL lbc_lnk( vn_b, 'V', -1. )
283       
284        IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
285       
286        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
287        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
288        !
289        !
290        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
291        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
292        !
293      ENDIF
294      !
295      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
296      !
297   END SUBROUTINE wad_lmt
298
299
300   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
301      !!----------------------------------------------------------------------
302      !!                  ***  ROUTINE wad_lmt  ***
303      !!                   
304      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
305      !!
306      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
307      !!
308      !! ** Action  : - calculate flux limiter and W/D flag
309      !!----------------------------------------------------------------------
310      REAL(wp), INTENT(in)    ::   rdtbt    ! ocean time-step index
311      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
312      !
313      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
314      INTEGER  ::   jflag               ! local scalar
315      REAL(wp) ::   z2dt
316      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
317      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
318      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
319      REAL(wp) ::   ztmp                ! local scalars
320      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
321      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
322      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
323      !!----------------------------------------------------------------------
324      !
325      IF( nn_timing == 1 )  CALL timing_start('wad_lmt_bt')
326
327      IF(ln_wd_il) THEN
328
329        CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
330        CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
331        !
332       
333        !IF(lwp) WRITE(numout,*)
334        !IF(lwp) WRITE(numout,*) 'wad_lmt_bt : wetting/drying limiters and velocity limiting'
335       
336        jflag  = 0
337        zdepwd = 50._wp   !maximum depth that ocean cells can have W/D processes
338
339        z2dt = rdtbt   
340       
341        zflxp(:,:)   = 0._wp
342        zflxn(:,:)   = 0._wp
343
344        zwdlmtu(:,:)  = 1._wp
345        zwdlmtv(:,:)  = 1._wp
346       
347        ! Horizontal Flux in u and v direction
348       
349        DO jj = 2, jpj
350           DO ji = 2, jpi 
351
352             IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE   ! we don't care about land cells
353             IF( ht_0(ji,jj) > zdepwd )      CYCLE   ! and cells which are unlikely to dry
354
355              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
356                           & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
357              zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
358                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
359
360              zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
361              IF(zdep2 .le. 0._wp) THEN  !add more safety, but not necessary
362                sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
363                IF(zflxu(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = 0._wp
364                IF(zflxu(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = 0._wp
365                IF(zflxv(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = 0._wp
366                IF(zflxv(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = 0._wp 
367              END IF
368           ENDDO
369        END DO
370
371     
372        !! start limiter iterations
373        DO jk1 = 1, nn_wdit + 1
374       
375         
376           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
377           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
378           jflag = 0     ! flag indicating if any further iterations are needed
379         
380           DO jj = 2, jpj
381              DO ji = 2, jpi 
382       
383                 IF( tmask(ji, jj, 1 ) < 0.5_wp) CYCLE
384                 IF( ht_0(ji,jj) > zdepwd )      CYCLE
385       
386                 ztmp = e1e2t(ji,jj)
387
388                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
389                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
390                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
391                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
392         
393                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
394                 zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1 - z2dt * zssh_frc(ji,jj)
395         
396                 IF(zdep1 > zdep2) THEN
397                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
398                   !zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zzflxp * z2dt )
399                   ! flag if the limiter has been used but stop flagging if the only
400                   ! changes have zeroed the coefficient since further iterations will
401                   ! not change anything
402                   IF( zcoef > 0._wp ) THEN
403                      jflag = 1 
404                   ELSE
405                      zcoef = 0._wp
406                   ENDIF
407                   IF(jk1 > nn_wdit) zcoef = 0._wp
408                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
409                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
410                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
411                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji,jj-1) = zcoef
412                 END IF
413              END DO ! ji loop
414           END DO  ! jj loop
415
416           CALL lbc_lnk( zwdlmtu, 'U', 1. )
417           CALL lbc_lnk( zwdlmtv, 'V', 1. )
418
419           IF(lk_mpp) CALL mpp_max(jflag)   !max over the global domain
420
421           IF(jflag == 0) EXIT
422         
423        END DO  ! jk1 loop
424       
425        zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
426        zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
427
428        CALL lbc_lnk( zflxu, 'U', -1. )
429        CALL lbc_lnk( zflxv, 'V', -1. )
430       
431        IF(jflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
432       
433        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
434        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
435        !
436        !
437        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
438        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
439        !
440      END IF
441
442      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
443   END SUBROUTINE wad_lmt_bt
444
445   !!==============================================================================
446END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.