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/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DYN/wet_dry.F90 @ 9090

Last change on this file since 9090 was 9090, checked in by flavoni, 6 years ago

change lbc_lnk in lbc_lnk_multi

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