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/branches/2019/dev_r11943_MERGE_2019/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DYN/wet_dry.F90 @ 12340

Last change on this file since 12340 was 12340, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. This commit introduces basic do loop macro
substitution to the 2019 option 1, merge branch. These changes have been SETTE
tested. The only addition is the do_loop_substitute.h90 file in the OCE directory but
the macros defined therein are used throughout the code to replace identifiable, 2D-
and 3D- nested loop opening and closing statements with single-line alternatives. Code
indents are also adjusted accordingly.

The following explanation is taken from comments in the new header file:

This header file contains preprocessor definitions and macros used in the do-loop
substitutions introduced between version 4.0 and 4.2. The primary aim of these macros
is to assist in future applications of tiling to improve performance. This is expected
to be achieved by alternative versions of these macros in selected locations. The
initial introduction of these macros simply replaces all identifiable nested 2D- and
3D-loops with single line statements (and adjusts indenting accordingly). Do loops
are identifiable if they comform to either:

DO jk = ....

DO jj = .... DO jj = ...

DO ji = .... DO ji = ...
. OR .
. .

END DO END DO

END DO END DO

END DO

and white-space variants thereof.

Additionally, only loops with recognised jj and ji loops limits are treated; these are:
Lower limits of 1, 2 or fs_2
Upper limits of jpi, jpim1 or fs_jpim1 (for ji) or jpj, jpjm1 or fs_jpjm1 (for jj)

The macro naming convention takes the form: DO_2D_BT_LR where:

B is the Bottom offset from the PE's inner domain;
T is the Top offset from the PE's inner domain;
L is the Left offset from the PE's inner domain;
R is the Right offset from the PE's inner domain

So, given an inner domain of 2,jpim1 and 2,jpjm1, a typical example would replace:

DO jj = 2, jpj

DO ji = 1, jpim1
.
.

END DO

END DO

with:

DO_2D_01_10
.
.
END_2D

similar conventions apply to the 3D loops macros. jk loop limits are retained
through macro arguments and are not restricted. This includes the possibility of
strides for which an extra set of DO_3DS macros are defined.

In the example definition below the inner PE domain is defined by start indices of
(kIs, kJs) and end indices of (kIe, KJe)

#define DO_2D_00_00 DO jj = kJs, kJe ; DO ji = kIs, kIe
#define END_2D END DO ; END DO

TO DO:


Only conventional nested loops have been identified and replaced by this step. There are constructs such as:

DO jk = 2, jpkm1

z2d(:,:) = z2d(:,:) + e3w(:,:,jk,Kmm) * z3d(:,:,jk) * wmask(:,:,jk)

END DO

which may need to be considered.

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