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_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DYN/sshwzv.F90 @ 11822

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

Branch 2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps. Sette tested updates to branch to align with trunk changes between 10721 and 11740. Sette tests are passing but results differ from branch before these changes (except for GYRE_PISCES and VORTEX) and branch results already differed from trunk because of algorithmic fixes. Will need more checks to confirm correctness.

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