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/2016/dev_CNRS_2016/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2016/dev_CNRS_2016/NEMOGCM/NEMO/OPA_SRC/DYN/wet_dry.F90 @ 7280

Last change on this file since 7280 was 7280, checked in by flavoni, 7 years ago

merge CNRS2016 with aerobulk branch

File size: 16.9 KB
Line 
1MODULE wet_dry
2   !!==============================================================================
3   !!                       ***  MODULE  wet_dry  ***
4   !! Wetting and drying includes initialisation routine and routines to
5   !! compute and apply flux limiters and preserve water depth positivity
6   !! only effects if wetting/drying is on (ln_wd == .true.)
7   !!==============================================================================
8   !! History :  3.6  ! 2014-09  ((H.Liu)  Original code
9   !!                 ! will add the runoff and periodic BC case later
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   wad_lmt    : Compute the horizontal flux limiter and the limited velocity
14   !!                when wetting and drying happens
15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers
17   USE dom_oce         ! ocean space and time domain
18   USE sbc_oce, ONLY : ln_rnf   ! surface boundary condition: ocean
19   USE sbcrnf          ! river runoff
20   USE in_out_manager  ! I/O manager
21   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
22   USE lib_mpp         ! MPP library
23   USE wrk_nemo        ! Memory Allocation
24   USE timing          ! Timing
25
26   IMPLICIT NONE
27   PRIVATE
28
29   !!----------------------------------------------------------------------
30   !! critical depths,filters, limiters,and masks for  Wetting and Drying
31   !! ---------------------------------------------------------------------
32
33   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wduflt, wdvflt !: u- and v- filter
34   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   wdmask         !: u- and v- limiter
35
36   LOGICAL,  PUBLIC  ::   ln_wd       !: Wetting/drying activation switch (T:on,F:off)
37   REAL(wp), PUBLIC  ::   rn_wdmin1   !: minimum water depth on dried cells
38   REAL(wp), PUBLIC  ::   rn_wdmin2   !: tolerrance of minimum water depth on dried cells
39   REAL(wp), PUBLIC  ::   rn_wdld     !: land elevation below which wetting/drying
40                                           !: will be considered
41   INTEGER , PUBLIC  ::   nn_wdit     !: maximum number of iteration for W/D limiter
42
43   PUBLIC   wad_init                  ! initialisation routine called by step.F90
44   PUBLIC   wad_lmt                   ! routine called by sshwzv.F90
45   PUBLIC   wad_lmt_bt                ! routine called by dynspg_ts.F90
46
47   !! * Substitutions
48#  include "vectopt_loop_substitute.h90"
49CONTAINS
50
51   SUBROUTINE wad_init
52      !!----------------------------------------------------------------------
53      !!                     ***  ROUTINE wad_init  ***
54      !!                   
55      !! ** Purpose :   read wetting and drying namelist and print the variables.
56      !!
57      !! ** input   : - namwad namelist
58      !!----------------------------------------------------------------------
59      NAMELIST/namwad/ ln_wd, rn_wdmin1, rn_wdmin2, rn_wdld, nn_wdit
60      INTEGER  ::   ios                 ! Local integer output status for namelist read
61      INTEGER  ::   ierr                ! Local integer status array allocation
62      !!----------------------------------------------------------------------
63
64      REWIND( numnam_ref )              ! Namelist namwad in reference namelist
65                                        ! : Parameters for Wetting/Drying
66      READ  ( numnam_ref, namwad, IOSTAT = ios, ERR = 905)
67905   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in reference namelist', .TRUE.) 
68      REWIND( numnam_cfg )              ! Namelist namwad in configuration namelist
69                                        ! : Parameters for Wetting/Drying
70      READ  ( numnam_cfg, namwad, IOSTAT = ios, ERR = 906)
71906   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namwad in configuration namelist', .TRUE. )
72      IF(lwm) WRITE ( numond, namwad )
73
74      IF(lwp) THEN                  ! control print
75         WRITE(numout,*)
76         WRITE(numout,*) 'wad_init : Wetting and drying initialization through namelist read'
77         WRITE(numout,*) '~~~~~~~~'
78         WRITE(numout,*) '   Namelist namwad'
79         WRITE(numout,*) '      Logical activation                 ln_wd      = ', ln_wd
80         WRITE(numout,*) '      Minimum wet depth on dried cells rn_wdmin1    = ', rn_wdmin1
81         WRITE(numout,*) '      Tolerance of min wet depth     rn_wdmin2      = ', rn_wdmin2
82         WRITE(numout,*) '      land elevation threshold         rn_wdld      = ', rn_wdld
83         WRITE(numout,*) '      Max iteration for W/D limiter    nn_wdit      = ', nn_wdit
84      ENDIF
85      !
86      IF(ln_wd) THEN
87         ALLOCATE( wduflt(jpi,jpj), wdvflt(jpi,jpj), wdmask(jpi,jpj), STAT=ierr )
88         IF( ierr /= 0 ) CALL ctl_stop('STOP', 'wad_init : Array allocation error')
89      ENDIF
90      !
91   END SUBROUTINE wad_init
92
93
94   SUBROUTINE wad_lmt( sshb1, sshemp, z2dt )
95      !!----------------------------------------------------------------------
96      !!                  ***  ROUTINE wad_lmt  ***
97      !!                   
98      !! ** Purpose :   generate flux limiters for wetting/drying
99      !!
100      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
101      !!
102      !! ** Action  : - calculate flux limiter and W/D flag
103      !!----------------------------------------------------------------------
104      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   sshb1
105      REAL(wp), DIMENSION(:,:), INTENT(in)    ::   sshemp
106      REAL(wp), INTENT(in) :: z2dt
107      !
108      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
109      INTEGER  ::   zflag               ! local scalar
110      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
111      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
112      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
113      REAL(wp) ::   ztmp                ! local scalars
114      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
115      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
116      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu,  zflxv            ! local 2D workspace
117      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
118      !!----------------------------------------------------------------------
119      !
120
121      IF( nn_timing == 1 )  CALL timing_start('wad_lmt')
122
123      IF(ln_wd) THEN
124
125        CALL wrk_alloc( jpi,jpj,   zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
126        CALL wrk_alloc( jpi,jpj,   zwdlmtu, zwdlmtv)
127        !
128       
129        !IF(lwp) WRITE(numout,*)
130        !IF(lwp) WRITE(numout,*) 'wad_lmt : wetting/drying limiters and velocity limiting'
131       
132        zflag  = 0
133        zdepwd = 50._wp   !maximum depth on which that W/D could possibly happen
134
135       
136        zflxp(:,:)   = 0._wp
137        zflxn(:,:)   = 0._wp
138        zflxu(:,:)   = 0._wp
139        zflxv(:,:)   = 0._wp
140
141        zwdlmtu(:,:)  = 1._wp
142        zwdlmtv(:,:)  = 1._wp
143       
144        ! Horizontal Flux in u and v direction
145        DO jk = 1, jpkm1 
146           DO jj = 1, jpjm1
147              DO ji = 1, jpim1
148                 zflxu(ji,jj) = zflxu(ji,jj) + e3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
149                 zflxv(ji,jj) = zflxv(ji,jj) + e3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
150              END DO 
151           END DO 
152        END DO
153       
154        zflxu(:,:) = zflxu(:,:) * e2u(:,:)
155        zflxv(:,:) = zflxv(:,:) * e1v(:,:)
156       
157         DO jj = 2, jpjm1
158           DO ji = 2, jpim1 
159
160             IF( tmask(ji,jj,1) == 0._wp  )   CYCLE   ! we don't care about land cells
161             IF( ht_0 (ji,jj)   >  zdepwd )   CYCLE   ! and cells which will unlikely go dried out
162
163              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
164                           & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
165              zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
166                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
167
168              zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1
169              IF(zdep2 < 0._wp) THEN  !add more safty, but not necessary
170                !zdep2 = 0._wp
171                sshb1(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
172              END IF
173           ENDDO
174        END DO
175
176     
177        !! start limiter iterations
178        DO jk1 = 1, nn_wdit + 1
179       
180         
181           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
182           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
183         
184           DO jj = 2, jpjm1
185              DO ji = 2, jpim1 
186       
187                 wdmask(ji,jj) = 0
188                 IF( tmask(ji,jj,1) < 0.5_wp) CYCLE
189                 IF( ht_0(ji,jj) > zdepwd) CYCLE
190       
191                 ztmp = e1e2t(ji,jj)
192
193                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
194                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
195                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
196                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
197         
198                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
199                 zdep2 = ht_0(ji,jj) + sshb1(ji,jj) - rn_wdmin1 - z2dt * sshemp(ji,jj)  ! this one can be moved out of the loop
200         
201                 IF(zdep1 > zdep2) THEN
202                   zflag = 1
203                   wdmask(ji, jj) = 1
204                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
205                   zcoef = max(zcoef, 0._wp)
206                   IF(jk1 > nn_wdit) zcoef = 0._wp
207                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
208                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
209                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
210                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji-1,jj) = zcoef
211                 END IF
212              END DO ! ji loop
213           END DO  ! jj loop
214
215           CALL lbc_lnk( zwdlmtu, 'U', 1. )
216           CALL lbc_lnk( zwdlmtv, 'V', 1. )
217
218           IF(lk_mpp) CALL mpp_max(zflag)   !max over the global domain
219
220           IF(zflag == 0) EXIT
221         
222           zflag = 0     ! flag indicating if any further iteration is needed?
223        END DO  ! jk1 loop
224       
225        DO jk = 1, jpkm1
226          un(:,:,jk) = un(:,:,jk) * zwdlmtu(:, :) 
227          vn(:,:,jk) = vn(:,:,jk) * zwdlmtv(:, :) 
228        END DO
229
230        CALL lbc_lnk( un, 'U', -1. )
231        CALL lbc_lnk( vn, 'V', -1. )
232       
233        IF(zflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt!!!'
234       
235        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
236        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
237        !
238        !
239        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu, zflxv, zflxu1, zflxv1 )
240        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
241         !
242      ENDIF
243      !
244      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
245      !
246   END SUBROUTINE wad_lmt
247
248
249   SUBROUTINE wad_lmt_bt( zflxu, zflxv, sshn_e, zssh_frc, rdtbt )
250      !!----------------------------------------------------------------------
251      !!                  ***  ROUTINE wad_lmt  ***
252      !!                   
253      !! ** Purpose :   limiting flux in the barotropic stepping (dynspg_ts)
254      !!
255      !! ** Method  : - Prevent negative depth occurring (Not ready for Agrif)
256      !!
257      !! ** Action  : - calculate flux limiter and W/D flag
258      !!----------------------------------------------------------------------
259      REAL(wp), INTENT(in)    ::   rdtbt    ! ocean time-step index
260      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   zflxu,  zflxv, sshn_e, zssh_frc 
261      !
262      INTEGER  ::   ji, jj, jk, jk1     ! dummy loop indices
263      INTEGER  ::   zflag         ! local scalar
264      REAL(wp) ::   z2dt
265      REAL(wp) ::   zcoef, zdep1, zdep2 ! local scalars
266      REAL(wp) ::   zzflxp, zzflxn      ! local scalars
267      REAL(wp) ::   zdepwd              ! local scalar, always wet cell depth
268      REAL(wp) ::   ztmp                ! local scalars
269      REAL(wp), POINTER,  DIMENSION(:,:) ::   zwdlmtu, zwdlmtv         !: W/D flux limiters
270      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxp,  zflxn            ! local 2D workspace
271      REAL(wp), POINTER,  DIMENSION(:,:) ::   zflxu1, zflxv1           ! local 2D workspace
272      !!----------------------------------------------------------------------
273      !
274      IF( nn_timing == 1 )  CALL timing_start('wad_lmt_bt')
275
276      IF(ln_wd) THEN
277
278        CALL wrk_alloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
279        CALL wrk_alloc( jpi, jpj, zwdlmtu, zwdlmtv)
280        !
281       
282        !IF(lwp) WRITE(numout,*)
283        !IF(lwp) WRITE(numout,*) 'wad_lmt_bt : wetting/drying limiters and velocity limiting'
284       
285        zflag  = 0
286        zdepwd = 50._wp   !maximum depth that ocean cells can have W/D processes
287
288        z2dt = rdtbt   
289       
290        zflxp(:,:)   = 0._wp
291        zflxn(:,:)   = 0._wp
292        !zflxu(:,:)   = 0._wp
293        !zflxv(:,:)   = 0._wp
294
295        zwdlmtu(:,:)  = 1._wp
296        zwdlmtv(:,:)  = 1._wp
297       
298        ! Horizontal Flux in u and v direction
299       
300        !zflxu(:,:) = zflxu(:,:) * e2u(:,:)
301        !zflxv(:,:) = zflxv(:,:) * e1v(:,:)
302       
303        DO jj = 2, jpjm1
304           DO ji = 2, jpim1 
305
306             IF(tmask(ji,jj,1) < 0.5_wp) CYCLE   ! we don't care about land cells
307             IF(ht_0 (ji,jj)   > zdepwd) CYCLE       ! and cells which will unlikely go dried out
308
309              zflxp(ji,jj) = max(zflxu(ji,jj), 0._wp) - min(zflxu(ji-1,jj),   0._wp) + &
310                           & max(zflxv(ji,jj), 0._wp) - min(zflxv(ji,  jj-1), 0._wp) 
311              zflxn(ji,jj) = min(zflxu(ji,jj), 0._wp) - max(zflxu(ji-1,jj),   0._wp) + &
312                           & min(zflxv(ji,jj), 0._wp) - max(zflxv(ji,  jj-1), 0._wp) 
313
314              zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1
315              IF(zdep2 < 0._wp) THEN  !add more safty, but not necessary
316                !zdep2 = 0._wp
317               sshn_e(ji,jj) = rn_wdmin1 - ht_0(ji,jj)
318              END IF
319           ENDDO
320        END DO
321
322     
323        !! start limiter iterations
324        DO jk1 = 1, nn_wdit + 1
325       
326         
327           zflxu1(:,:) = zflxu(:,:) * zwdlmtu(:,:)
328           zflxv1(:,:) = zflxv(:,:) * zwdlmtv(:,:)
329         
330           DO jj = 2, jpjm1
331              DO ji = 2, jpim1 
332       
333                 wdmask(ji,jj) = 0
334                 IF(tmask(ji,jj,1) < 0.5_wp) CYCLE
335                 IF(ht_0 (ji,jj)   > zdepwd) CYCLE
336       
337                 ztmp = e1e2t(ji,jj)
338
339                 zzflxp = max(zflxu1(ji,jj), 0._wp) - min(zflxu1(ji-1,jj),   0._wp) + &
340                        & max(zflxv1(ji,jj), 0._wp) - min(zflxv1(ji,  jj-1), 0._wp) 
341                 zzflxn = min(zflxu1(ji,jj), 0._wp) - max(zflxu1(ji-1,jj),   0._wp) + &
342                        & min(zflxv1(ji,jj), 0._wp) - max(zflxv1(ji,  jj-1), 0._wp) 
343         
344                 zdep1 = (zzflxp + zzflxn) * z2dt / ztmp
345                 zdep2 = ht_0(ji,jj) + sshn_e(ji,jj) - rn_wdmin1   ! this one can be moved out of the loop
346                 zdep2 = zdep2 - z2dt * zssh_frc(ji,jj)
347         
348                 IF(zdep1 > zdep2) THEN
349                   zflag = 1
350                   !wdmask(ji, jj) = 1
351                   zcoef = ( ( zdep2 - rn_wdmin2 ) * ztmp - zzflxn * z2dt ) / ( zflxp(ji,jj) * z2dt )
352                   zcoef = max(zcoef, 0._wp)
353                   IF(jk1 > nn_wdit) zcoef = 0._wp
354                   IF(zflxu1(ji,  jj) > 0._wp) zwdlmtu(ji  ,jj) = zcoef
355                   IF(zflxu1(ji-1,jj) < 0._wp) zwdlmtu(ji-1,jj) = zcoef
356                   IF(zflxv1(ji,  jj) > 0._wp) zwdlmtv(ji  ,jj) = zcoef
357                   IF(zflxv1(ji,jj-1) < 0._wp) zwdlmtv(ji-1,jj) = zcoef
358                 END IF
359              END DO ! ji loop
360           END DO  ! jj loop
361
362           CALL lbc_lnk( zwdlmtu, 'U', 1. )
363           CALL lbc_lnk( zwdlmtv, 'V', 1. )
364
365           IF(lk_mpp) CALL mpp_max(zflag)   !max over the global domain
366
367           IF(zflag == 0) EXIT
368         
369           zflag = 0     ! flag indicating if any further iteration is needed?
370        END DO  ! jk1 loop
371       
372        zflxu(:,:) = zflxu(:,:) * zwdlmtu(:, :) 
373        zflxv(:,:) = zflxv(:,:) * zwdlmtv(:, :) 
374
375        CALL lbc_lnk( zflxu, 'U', -1. )
376        CALL lbc_lnk( zflxv, 'V', -1. )
377       
378        IF(zflag == 1 .AND. lwp) WRITE(numout,*) 'Need more iterations in wad_lmt_bt!!!'
379       
380        !IF( ln_rnf      )   CALL sbc_rnf_div( hdivn )          ! runoffs (update hdivn field)
381        !IF( nn_cla == 1 )   CALL cla_div    ( kt )             ! Cross Land Advection (update hdivn field)
382        !
383        !
384        CALL wrk_dealloc( jpi, jpj, zflxp, zflxn, zflxu1, zflxv1 )
385        CALL wrk_dealloc( jpi, jpj, zwdlmtu, zwdlmtv)
386         !
387      END IF
388      !
389      IF( nn_timing == 1 )  CALL timing_stop('wad_lmt')
390      !
391   END SUBROUTINE wad_lmt_bt
392
393   !!==============================================================================
394END MODULE wet_dry
Note: See TracBrowser for help on using the repository browser.