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

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DYN/sshwzv.F90 @ 12236

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

Branch 2019/dev_r11943_MERGE_2019. Merge in changes from 2019/fix_sn_cfctl_ticket2328. Fully SETTE tested

  • Property svn:keywords set to Id
File size: 19.8 KB
Line 
1MODULE sshwzv   
2   !!==============================================================================
3   !!                       ***  MODULE  sshwzv  ***
4   !! Ocean dynamics : sea surface height and vertical velocity
5   !!==============================================================================
6   !! History :  3.1  !  2009-02  (G. Madec, M. Leclair)  Original code
7   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  modified LF-RA
8   !!             -   !  2010-05  (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
9   !!             -   !  2010-09  (D.Storkey and E.O'Dea) bug fixes for BDY module
10   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
11   !!            4.0  !  2018-12  (A. Coward) add mixed implicit/explicit advection
12   !!            4.1  !  2019-08  (A. Coward, D. Storkey) Rename ssh_nxt -> ssh_atf. Now only does time filtering.
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   ssh_nxt       : after ssh
17   !!   ssh_atf       : time filter the ssh arrays
18   !!   wzv           : compute now vertical velocity
19   !!----------------------------------------------------------------------
20   USE oce            ! ocean dynamics and tracers variables
21   USE isf_oce        ! ice shelf
22   USE dom_oce        ! ocean space and time domain variables
23   USE sbc_oce        ! surface boundary condition: ocean
24   USE domvvl         ! Variable volume
25   USE divhor         ! horizontal divergence
26   USE phycst         ! physical constants
27   USE bdy_oce , ONLY : ln_bdy, bdytmask   ! Open BounDarY
28   USE bdydyn2d       ! bdy_ssh routine
29#if defined key_agrif
30   USE agrif_oce_interp
31#endif
32   !
33   USE iom 
34   USE in_out_manager ! I/O manager
35   USE restart        ! only for lrst_oce
36   USE prtctl         ! Print control
37   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
38   USE lib_mpp        ! MPP library
39   USE timing         ! Timing
40   USE wet_dry        ! Wetting/Drying flux limiting
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC   ssh_nxt    ! called by step.F90
46   PUBLIC   wzv        ! called by step.F90
47   PUBLIC   wAimp      ! called by step.F90
48   PUBLIC   ssh_atf    ! called by step.F90
49
50   !! * Substitutions
51#  include "vectopt_loop_substitute.h90"
52   !!----------------------------------------------------------------------
53   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
54   !! $Id$
55   !! Software governed by the CeCILL license (see ./LICENSE)
56   !!----------------------------------------------------------------------
57CONTAINS
58
59   SUBROUTINE ssh_nxt( kt, Kbb, Kmm, pssh, Kaa )
60      !!----------------------------------------------------------------------
61      !!                ***  ROUTINE ssh_nxt  ***
62      !!                   
63      !! ** Purpose :   compute the after ssh (ssh(Kaa))
64      !!
65      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
66      !!      is computed by integrating the horizontal divergence and multiply by
67      !!      by the time step.
68      !!
69      !! ** action  :   ssh(:,:,Kaa), after sea surface height
70      !!
71      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
72      !!----------------------------------------------------------------------
73      INTEGER                         , INTENT(in   ) ::   kt             ! time step
74      INTEGER                         , INTENT(in   ) ::   Kbb, Kmm, Kaa  ! time level index
75      REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) ::   pssh           ! sea-surface height
76      !
77      INTEGER  ::   jk            ! dummy loop indice
78      REAL(wp) ::   z2dt, zcoef   ! local scalars
79      REAL(wp), DIMENSION(jpi,jpj) ::   zhdiv   ! 2D workspace
80      !!----------------------------------------------------------------------
81      !
82      IF( ln_timing )   CALL timing_start('ssh_nxt')
83      !
84      IF( kt == nit000 ) THEN
85         IF(lwp) WRITE(numout,*)
86         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
87         IF(lwp) WRITE(numout,*) '~~~~~~~ '
88      ENDIF
89      !
90      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
91      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
92      zcoef = 0.5_wp * r1_rau0
93
94      !                                           !------------------------------!
95      !                                           !   After Sea Surface Height   !
96      !                                           !------------------------------!
97      IF(ln_wd_il) THEN
98         CALL wad_lmt(pssh(:,:,Kbb), zcoef * (emp_b(:,:) + emp(:,:)), z2dt, Kmm, uu, vv )
99      ENDIF
100
101      CALL div_hor( kt, Kbb, Kmm )                     ! Horizontal divergence
102      !
103      zhdiv(:,:) = 0._wp
104      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
105        zhdiv(:,:) = zhdiv(:,:) + e3t(:,:,jk,Kmm) * hdiv(:,:,jk)
106      END DO
107      !                                                ! Sea surface elevation time stepping
108      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
109      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
110      !
111      pssh(:,:,Kaa) = (  pssh(:,:,Kbb) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
112      !
113#if defined key_agrif
114      Kbb_a = Kbb; Kmm_a = Kmm; Krhs_a = Kaa; CALL agrif_ssh( kt )
115#endif
116      !
117      IF ( .NOT.ln_dynspg_ts ) THEN
118         IF( ln_bdy ) THEN
119            CALL lbc_lnk( 'sshwzv', pssh(:,:,Kaa), 'T', 1. )    ! Not sure that's necessary
120            CALL bdy_ssh( pssh(:,:,Kaa) )             ! Duplicate sea level across open boundaries
121         ENDIF
122      ENDIF
123      !                                           !------------------------------!
124      !                                           !           outputs            !
125      !                                           !------------------------------!
126      !
127      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pssh(:,:,Kaa), clinfo1=' pssh(:,:,Kaa)  - : ', mask1=tmask )
128      !
129      IF( ln_timing )   CALL timing_stop('ssh_nxt')
130      !
131   END SUBROUTINE ssh_nxt
132
133   
134   SUBROUTINE wzv( kt, Kbb, Kmm, pww, Kaa )
135      !!----------------------------------------------------------------------
136      !!                ***  ROUTINE wzv  ***
137      !!                   
138      !! ** Purpose :   compute the now vertical velocity
139      !!
140      !! ** Method  : - Using the incompressibility hypothesis, the vertical
141      !!      velocity is computed by integrating the horizontal divergence 
142      !!      from the bottom to the surface minus the scale factor evolution.
143      !!        The boundary conditions are w=0 at the bottom (no flux) and.
144      !!
145      !! ** action  :   pww      : now vertical velocity
146      !!
147      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
148      !!----------------------------------------------------------------------
149      INTEGER                         , INTENT(in)    ::   kt             ! time step
150      INTEGER                         , INTENT(in)    ::   Kbb, Kmm, Kaa  ! time level indices
151      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pww            ! now vertical velocity
152      !
153      INTEGER  ::   ji, jj, jk   ! dummy loop indices
154      REAL(wp) ::   z1_2dt       ! local scalars
155      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdiv
156      !!----------------------------------------------------------------------
157      !
158      IF( ln_timing )   CALL timing_start('wzv')
159      !
160      IF( kt == nit000 ) THEN
161         IF(lwp) WRITE(numout,*)
162         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
163         IF(lwp) WRITE(numout,*) '~~~~~ '
164         !
165         pww(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
166      ENDIF
167      !                                           !------------------------------!
168      !                                           !     Now Vertical Velocity    !
169      !                                           !------------------------------!
170      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
171      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
172      !
173      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
174         ALLOCATE( zhdiv(jpi,jpj,jpk) ) 
175         !
176         DO jk = 1, jpkm1
177            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
178            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
179            DO jj = 2, jpjm1
180               DO ji = fs_2, fs_jpim1   ! vector opt.
181                  zhdiv(ji,jj,jk) = r1_e1e2t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) )
182               END DO
183            END DO
184         END DO
185         CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
186         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
187         !                             ! Same question holds for hdiv. Perhaps just for security
188         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
189            ! computation of w
190            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk) + zhdiv(:,:,jk)    &
191               &                         + z1_2dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )     ) * tmask(:,:,jk)
192         END DO
193         !          IF( ln_vvl_layer ) pww(:,:,:) = 0.e0
194         DEALLOCATE( zhdiv ) 
195      ELSE   ! z_star and linear free surface cases
196         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
197            ! computation of w
198            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)                 &
199               &                         + z1_2dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )  ) * tmask(:,:,jk)
200         END DO
201      ENDIF
202
203      IF( ln_bdy ) THEN
204         DO jk = 1, jpkm1
205            pww(:,:,jk) = pww(:,:,jk) * bdytmask(:,:)
206         END DO
207      ENDIF
208      !
209#if defined key_agrif 
210      IF( .NOT. AGRIF_Root() ) THEN
211         IF ((nbondi ==  1).OR.(nbondi == 2)) pww(nlci-1 , :     ,:) = 0.e0      ! east
212         IF ((nbondi == -1).OR.(nbondi == 2)) pww(2      , :     ,:) = 0.e0      ! west
213         IF ((nbondj ==  1).OR.(nbondj == 2)) pww(:      ,nlcj-1 ,:) = 0.e0      ! north
214         IF ((nbondj == -1).OR.(nbondj == 2)) pww(:      ,2      ,:) = 0.e0      ! south
215      ENDIF 
216#endif 
217      !
218      IF( ln_timing )   CALL timing_stop('wzv')
219      !
220   END SUBROUTINE wzv
221
222
223   SUBROUTINE ssh_atf( kt, Kbb, Kmm, Kaa, pssh )
224      !!----------------------------------------------------------------------
225      !!                    ***  ROUTINE ssh_atf  ***
226      !!
227      !! ** Purpose :   Apply Asselin time filter to now SSH.
228      !!
229      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
230      !!              from the filter, see Leclair and Madec 2010) and swap :
231      !!                pssh(:,:,Kmm) = pssh(:,:,Kaa) + atfp * ( pssh(:,:,Kbb) -2 pssh(:,:,Kmm) + pssh(:,:,Kaa) )
232      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
233      !!
234      !! ** action  : - pssh(:,:,Kmm) time filtered
235      !!
236      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
237      !!----------------------------------------------------------------------
238      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index
239      INTEGER                         , INTENT(in   ) ::   Kbb, Kmm, Kaa  ! ocean time level indices
240      REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) ::   pssh           ! SSH field
241      !
242      REAL(wp) ::   zcoef   ! local scalar
243      !!----------------------------------------------------------------------
244      !
245      IF( ln_timing )   CALL timing_start('ssh_atf')
246      !
247      IF( kt == nit000 ) THEN
248         IF(lwp) WRITE(numout,*)
249         IF(lwp) WRITE(numout,*) 'ssh_atf : Asselin time filter of sea surface height'
250         IF(lwp) WRITE(numout,*) '~~~~~~~ '
251      ENDIF
252      !              !==  Euler time-stepping: no filter, just swap  ==!
253      IF ( .NOT.( neuler == 0 .AND. kt == nit000 ) ) THEN   ! Only do time filtering for leapfrog timesteps
254         !                                                  ! filtered "now" field
255         pssh(:,:,Kmm) = pssh(:,:,Kmm) + atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) )
256         IF( .NOT.ln_linssh ) THEN                          ! "now" <-- with forcing removed
257            zcoef = atfp * rdt * r1_rau0
258            pssh(:,:,Kmm) = pssh(:,:,Kmm) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
259               &                             - rnf_b(:,:)        + rnf   (:,:)       &
260               &                             + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   &
261               &                             + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:)
262
263            ! ice sheet coupling
264            IF ( ln_isf .AND. ln_isfcpl .AND. kt == nit000+1) pssh(:,:,Kbb) = pssh(:,:,Kbb) - atfp * rdt * ( risfcpl_ssh(:,:) - 0.0 ) * ssmask(:,:)
265
266         ENDIF
267      ENDIF
268      !
269      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pssh(:,:,Kmm), clinfo1=' pssh(:,:,Kmm)  - : ', mask1=tmask )
270      !
271      IF( ln_timing )   CALL timing_stop('ssh_atf')
272      !
273   END SUBROUTINE ssh_atf
274
275   SUBROUTINE wAimp( kt, Kmm )
276      !!----------------------------------------------------------------------
277      !!                ***  ROUTINE wAimp  ***
278      !!                   
279      !! ** Purpose :   compute the Courant number and partition vertical velocity
280      !!                if a proportion needs to be treated implicitly
281      !!
282      !! ** Method  : -
283      !!
284      !! ** action  :   ww      : now vertical velocity (to be handled explicitly)
285      !!            :   wi      : now vertical velocity (for implicit treatment)
286      !!
287      !! Reference  : Shchepetkin, A. F. (2015): An adaptive, Courant-number-dependent
288      !!              implicit scheme for vertical advection in oceanic modeling.
289      !!              Ocean Modelling, 91, 38-69.
290      !!----------------------------------------------------------------------
291      INTEGER, INTENT(in) ::   kt   ! time step
292      INTEGER, INTENT(in) ::   Kmm  ! time level index
293      !
294      INTEGER  ::   ji, jj, jk   ! dummy loop indices
295      REAL(wp)             ::   zCu, zcff, z1_e3t                     ! local scalars
296      REAL(wp) , PARAMETER ::   Cu_min = 0.15_wp                      ! local parameters
297      REAL(wp) , PARAMETER ::   Cu_max = 0.30_wp                      ! local parameters
298      REAL(wp) , PARAMETER ::   Cu_cut = 2._wp*Cu_max - Cu_min        ! local parameters
299      REAL(wp) , PARAMETER ::   Fcu    = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters
300      !!----------------------------------------------------------------------
301      !
302      IF( ln_timing )   CALL timing_start('wAimp')
303      !
304      IF( kt == nit000 ) THEN
305         IF(lwp) WRITE(numout,*)
306         IF(lwp) WRITE(numout,*) 'wAimp : Courant number-based partitioning of now vertical velocity '
307         IF(lwp) WRITE(numout,*) '~~~~~ '
308         wi(:,:,:) = 0._wp
309      ENDIF
310      !
311      ! Calculate Courant numbers
312      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
313         DO jk = 1, jpkm1
314            DO jj = 2, jpjm1
315               DO ji = 2, fs_jpim1   ! vector opt.
316                  z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm)
317                  ! 2*rdt and not r2dt (for restartability)
318                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )                       & 
319                     &                             + ( MAX( e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kmm)*uu(ji  ,jj,jk,Kmm) + un_td(ji  ,jj,jk), 0._wp ) -   &
320                     &                                 MIN( e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kmm)*uu(ji-1,jj,jk,Kmm) + un_td(ji-1,jj,jk), 0._wp ) )   &
321                     &                               * r1_e1e2t(ji,jj)                                                                     &
322                     &                             + ( MAX( e1v(ji,jj  )*e3v(ji,jj  ,jk,Kmm)*vv(ji,jj  ,jk,Kmm) + vn_td(ji,jj  ,jk), 0._wp ) -   &
323                     &                                 MIN( e1v(ji,jj-1)*e3v(ji,jj-1,jk,Kmm)*vv(ji,jj-1,jk,Kmm) + vn_td(ji,jj-1,jk), 0._wp ) )   &
324                     &                               * r1_e1e2t(ji,jj)                                                                     &
325                     &                             ) * z1_e3t
326               END DO
327            END DO
328         END DO
329      ELSE
330         DO jk = 1, jpkm1
331            DO jj = 2, jpjm1
332               DO ji = 2, fs_jpim1   ! vector opt.
333                  z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm)
334                  ! 2*rdt and not r2dt (for restartability)
335                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )   & 
336                     &                             + ( MAX( e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kmm)*uu(ji  ,jj,jk,Kmm), 0._wp ) -   &
337                     &                                 MIN( e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kmm)*uu(ji-1,jj,jk,Kmm), 0._wp ) )   &
338                     &                               * r1_e1e2t(ji,jj)                                                 &
339                     &                             + ( MAX( e1v(ji,jj  )*e3v(ji,jj  ,jk,Kmm)*vv(ji,jj  ,jk,Kmm), 0._wp ) -   &
340                     &                                 MIN( e1v(ji,jj-1)*e3v(ji,jj-1,jk,Kmm)*vv(ji,jj-1,jk,Kmm), 0._wp ) )   &
341                     &                               * r1_e1e2t(ji,jj)                                                 &
342                     &                             ) * z1_e3t
343               END DO
344            END DO
345         END DO
346      ENDIF
347      CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. )
348      !
349      CALL iom_put("Courant",Cu_adv)
350      !
351      IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN       ! Quick check if any breaches anywhere
352         DO jk = jpkm1, 2, -1                           ! or scan Courant criterion and partition
353            DO jj = 1, jpj                              ! w where necessary
354               DO ji = 1, jpi
355                  !
356                  zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) )
357! alt:
358!                  IF ( ww(ji,jj,jk) > 0._wp ) THEN
359!                     zCu =  Cu_adv(ji,jj,jk)
360!                  ELSE
361!                     zCu =  Cu_adv(ji,jj,jk-1)
362!                  ENDIF
363                  !
364                  IF( zCu <= Cu_min ) THEN              !<-- Fully explicit
365                     zcff = 0._wp
366                  ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit
367                     zcff = ( zCu - Cu_min )**2
368                     zcff = zcff / ( Fcu + zcff )
369                  ELSE                                  !<-- Mostly implicit
370                     zcff = ( zCu - Cu_max )/ zCu
371                  ENDIF
372                  zcff = MIN(1._wp, zcff)
373                  !
374                  wi(ji,jj,jk) =           zcff   * ww(ji,jj,jk)
375                  ww(ji,jj,jk) = ( 1._wp - zcff ) * ww(ji,jj,jk)
376                  !
377                  Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient below and in stp_ctl
378               END DO
379            END DO
380         END DO
381         Cu_adv(:,:,1) = 0._wp 
382      ELSE
383         ! Fully explicit everywhere
384         Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient below and in stp_ctl
385         wi    (:,:,:) = 0._wp
386      ENDIF
387      CALL iom_put("wimp",wi) 
388      CALL iom_put("wi_cff",Cu_adv)
389      CALL iom_put("wexp",ww)
390      !
391      IF( ln_timing )   CALL timing_stop('wAimp')
392      !
393   END SUBROUTINE wAimp
394   !!======================================================================
395END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.