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/ENHANCE-02_ISF_nemo_TEST_MERGE/src/OCE/DYN – NEMO

source: NEMO/branches/2019/ENHANCE-02_ISF_nemo_TEST_MERGE/src/OCE/DYN/sshwzv.F90 @ 11970

Last change on this file since 11970 was 11970, checked in by davestorkey, 4 years ago

2019/ENHANCE-02_ISF_nemo_TEST_MERGE : copy changes from Pierre's branch.

  • Property svn:keywords set to Id
File size: 19.1 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   !!----------------------------------------------------------------------
13
14   !!----------------------------------------------------------------------
15   !!   ssh_nxt       : after ssh
16   !!   ssh_swp       : filter ans swap the ssh arrays
17   !!   wzv           : compute now vertical velocity
18   !!----------------------------------------------------------------------
19   USE oce            ! ocean dynamics and tracers variables
20   USE isf            ! ice shelf
21   USE isfutils
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_swp    ! 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 )
60      !!----------------------------------------------------------------------
61      !!                ***  ROUTINE ssh_nxt  ***
62      !!                   
63      !! ** Purpose :   compute the after ssh (ssha)
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  :   ssha, after sea surface height
70      !!
71      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
72      !!----------------------------------------------------------------------
73      INTEGER, INTENT(in) ::   kt   ! time step
74      !
75      INTEGER  ::   jk            ! dummy loop indice
76      REAL(wp) ::   z2dt, zcoef   ! local scalars
77      REAL(wp), DIMENSION(jpi,jpj) ::   zhdiv   ! 2D workspace
78      !!----------------------------------------------------------------------
79      !
80      IF( ln_timing )   CALL timing_start('ssh_nxt')
81      !
82      IF( kt == nit000 ) THEN
83         IF(lwp) WRITE(numout,*)
84         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
85         IF(lwp) WRITE(numout,*) '~~~~~~~ '
86      ENDIF
87      !
88      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
89      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
90      zcoef = 0.5_wp * r1_rau0
91
92      !                                           !------------------------------!
93      !                                           !   After Sea Surface Height   !
94      !                                           !------------------------------!
95      IF(ln_wd_il) THEN
96         CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
97      ENDIF
98
99      CALL div_hor( kt )                               ! Horizontal divergence
100      !
101      zhdiv(:,:) = 0._wp
102      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
103        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
104      END DO
105      !                                                ! Sea surface elevation time stepping
106      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
107      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
108      !
109      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
110      !
111#if defined key_agrif
112      CALL agrif_ssh( kt )
113#endif
114      !
115      IF ( .NOT.ln_dynspg_ts ) THEN
116         IF( ln_bdy ) THEN
117            CALL lbc_lnk( 'sshwzv', ssha, 'T', 1. )    ! Not sure that's necessary
118            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
119         ENDIF
120      ENDIF
121      !                                           !------------------------------!
122      !                                           !           outputs            !
123      !                                           !------------------------------!
124      !
125      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask )
126      !
127      IF( ln_timing )   CALL timing_stop('ssh_nxt')
128      !
129   END SUBROUTINE ssh_nxt
130
131   
132   SUBROUTINE wzv( kt )
133      !!----------------------------------------------------------------------
134      !!                ***  ROUTINE wzv  ***
135      !!                   
136      !! ** Purpose :   compute the now vertical velocity
137      !!
138      !! ** Method  : - Using the incompressibility hypothesis, the vertical
139      !!      velocity is computed by integrating the horizontal divergence 
140      !!      from the bottom to the surface minus the scale factor evolution.
141      !!        The boundary conditions are w=0 at the bottom (no flux) and.
142      !!
143      !! ** action  :   wn      : now vertical velocity
144      !!
145      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
146      !!----------------------------------------------------------------------
147      INTEGER, INTENT(in) ::   kt   ! time step
148      !
149      INTEGER  ::   ji, jj, jk   ! dummy loop indices
150      REAL(wp) ::   z1_2dt       ! local scalars
151      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdiv
152      !!----------------------------------------------------------------------
153      !
154      IF( ln_timing )   CALL timing_start('wzv')
155      !
156      IF( kt == nit000 ) THEN
157         IF(lwp) WRITE(numout,*)
158         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
159         IF(lwp) WRITE(numout,*) '~~~~~ '
160         !
161         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
162      ENDIF
163      !                                           !------------------------------!
164      !                                           !     Now Vertical Velocity    !
165      !                                           !------------------------------!
166      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
167      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
168      !
169      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
170         ALLOCATE( zhdiv(jpi,jpj,jpk) ) 
171         !
172         DO jk = 1, jpkm1
173            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
174            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
175            DO jj = 2, jpjm1
176               DO ji = fs_2, fs_jpim1   ! vector opt.
177                  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) )
178               END DO
179            END DO
180         END DO
181         CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
182         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
183         !                             ! Same question holds for hdivn. Perhaps just for security
184         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
185            ! computation of w
186            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
187               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
188         END DO
189         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
190         DEALLOCATE( zhdiv ) 
191      ELSE   ! z_star and linear free surface cases
192         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
193            ! computation of w
194            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
195               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
196         END DO
197      ENDIF
198
199      IF( ln_bdy ) THEN
200         DO jk = 1, jpkm1
201            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
202         END DO
203      ENDIF
204      !
205#if defined key_agrif 
206      IF( .NOT. AGRIF_Root() ) THEN
207         IF ((nbondi ==  1).OR.(nbondi == 2)) wn(nlci-1 , :     ,:) = 0.e0      ! east
208         IF ((nbondi == -1).OR.(nbondi == 2)) wn(2      , :     ,:) = 0.e0      ! west
209         IF ((nbondj ==  1).OR.(nbondj == 2)) wn(:      ,nlcj-1 ,:) = 0.e0      ! north
210         IF ((nbondj == -1).OR.(nbondj == 2)) wn(:      ,2      ,:) = 0.e0      ! south
211      ENDIF 
212#endif 
213      !
214      IF( ln_timing )   CALL timing_stop('wzv')
215      !
216   END SUBROUTINE wzv
217
218
219   SUBROUTINE ssh_swp( kt )
220      !!----------------------------------------------------------------------
221      !!                    ***  ROUTINE ssh_nxt  ***
222      !!
223      !! ** Purpose :   achieve the sea surface  height time stepping by
224      !!              applying Asselin time filter and swapping the arrays
225      !!              ssha  already computed in ssh_nxt 
226      !!
227      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
228      !!              from the filter, see Leclair and Madec 2010) and swap :
229      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
230      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
231      !!                sshn = ssha
232      !!
233      !! ** action  : - sshb, sshn   : before & now sea surface height
234      !!                               ready for the next time step
235      !!
236      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
237      !!----------------------------------------------------------------------
238      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
239      !
240      REAL(wp) ::   zcoef   ! local scalar
241      !!----------------------------------------------------------------------
242      !
243      IF( ln_timing )   CALL timing_start('ssh_swp')
244      !
245      IF( kt == nit000 ) THEN
246         IF(lwp) WRITE(numout,*)
247         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
248         IF(lwp) WRITE(numout,*) '~~~~~~~ '
249      ENDIF
250      !              !==  Euler time-stepping: no filter, just swap  ==!
251      IF ( neuler == 0 .AND. kt == nit000 ) THEN
252         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
253         !
254      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
255         !                                                  ! before <-- now filtered
256         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
257         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
258            zcoef = atfp * rdt * r1_rau0
259            sshb(:,:) = sshb(:,:) - zcoef * (  emp_b(:,:)        - emp   (:,:)       &
260               &                             - rnf_b(:,:)        + rnf   (:,:)       &
261               &                             + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   &
262               &                             + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:)
263
264            ! ice sheet coupling
265            IF ( ln_isf .AND. ln_isfcpl .AND. kt == nit000+1) sshb(:,:) = sshb(:,:) - atfp * rdt * ( risfcpl_ssh(:,:) - 0.0 ) * ssmask(:,:)
266
267         ENDIF
268
269         sshn(:,:) = ssha(:,:)                              ! now <-- after
270      ENDIF
271      !
272      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask )
273      !
274      IF( ln_timing )   CALL timing_stop('ssh_swp')
275      !
276   END SUBROUTINE ssh_swp
277
278   SUBROUTINE wAimp( kt )
279      !!----------------------------------------------------------------------
280      !!                ***  ROUTINE wAimp  ***
281      !!                   
282      !! ** Purpose :   compute the Courant number and partition vertical velocity
283      !!                if a proportion needs to be treated implicitly
284      !!
285      !! ** Method  : -
286      !!
287      !! ** action  :   wn      : now vertical velocity (to be handled explicitly)
288      !!            :   wi      : now vertical velocity (for implicit treatment)
289      !!
290      !! Reference  : Shchepetkin, A. F. (2015): An adaptive, Courant-number-dependent
291      !!              implicit scheme for vertical advection in oceanic modeling.
292      !!              Ocean Modelling, 91, 38-69.
293      !!----------------------------------------------------------------------
294      INTEGER, INTENT(in) ::   kt   ! time step
295      !
296      INTEGER  ::   ji, jj, jk   ! dummy loop indices
297      REAL(wp)             ::   zCu, zcff, z1_e3t                     ! local scalars
298      REAL(wp) , PARAMETER ::   Cu_min = 0.15_wp                      ! local parameters
299      REAL(wp) , PARAMETER ::   Cu_max = 0.30_wp                      ! local parameters
300      REAL(wp) , PARAMETER ::   Cu_cut = 2._wp*Cu_max - Cu_min        ! local parameters
301      REAL(wp) , PARAMETER ::   Fcu    = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters
302      !!----------------------------------------------------------------------
303      !
304      IF( ln_timing )   CALL timing_start('wAimp')
305      !
306      IF( kt == nit000 ) THEN
307         IF(lwp) WRITE(numout,*)
308         IF(lwp) WRITE(numout,*) 'wAimp : Courant number-based partitioning of now vertical velocity '
309         IF(lwp) WRITE(numout,*) '~~~~~ '
310         wi(:,:,:) = 0._wp
311      ENDIF
312      !
313      ! Calculate Courant numbers
314      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
315         DO jk = 1, jpkm1
316            DO jj = 2, jpjm1
317               DO ji = 2, fs_jpim1   ! vector opt.
318                  z1_e3t = 1._wp / e3t_n(ji,jj,jk)
319                  ! 2*rdt and not r2dt (for restartability)
320                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )                       & 
321                     &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk) + un_td(ji  ,jj,jk), 0._wp ) -   &
322                     &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk) + un_td(ji-1,jj,jk), 0._wp ) )   &
323                     &                               * r1_e1e2t(ji,jj)                                                                     &
324                     &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk) + vn_td(ji,jj  ,jk), 0._wp ) -   &
325                     &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk) + vn_td(ji,jj-1,jk), 0._wp ) )   &
326                     &                               * r1_e1e2t(ji,jj)                                                                     &
327                     &                             ) * z1_e3t
328               END DO
329            END DO
330         END DO
331      ELSE
332         DO jk = 1, jpkm1
333            DO jj = 2, jpjm1
334               DO ji = 2, fs_jpim1   ! vector opt.
335                  z1_e3t = 1._wp / e3t_n(ji,jj,jk)
336                  ! 2*rdt and not r2dt (for restartability)
337                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )   & 
338                     &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk), 0._wp ) -   &
339                     &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk), 0._wp ) )   &
340                     &                               * r1_e1e2t(ji,jj)                                                 &
341                     &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk), 0._wp ) -   &
342                     &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk), 0._wp ) )   &
343                     &                               * r1_e1e2t(ji,jj)                                                 &
344                     &                             ) * z1_e3t
345               END DO
346            END DO
347         END DO
348      ENDIF
349      CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. )
350      !
351      CALL iom_put("Courant",Cu_adv)
352      !
353      IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN       ! Quick check if any breaches anywhere
354         DO jk = jpkm1, 2, -1                           ! or scan Courant criterion and partition
355            DO jj = 1, jpj                              ! w where necessary
356               DO ji = 1, jpi
357                  !
358                  zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) )
359! alt:
360!                  IF ( wn(ji,jj,jk) > 0._wp ) THEN
361!                     zCu =  Cu_adv(ji,jj,jk)
362!                  ELSE
363!                     zCu =  Cu_adv(ji,jj,jk-1)
364!                  ENDIF
365                  !
366                  IF( zCu <= Cu_min ) THEN              !<-- Fully explicit
367                     zcff = 0._wp
368                  ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit
369                     zcff = ( zCu - Cu_min )**2
370                     zcff = zcff / ( Fcu + zcff )
371                  ELSE                                  !<-- Mostly implicit
372                     zcff = ( zCu - Cu_max )/ zCu
373                  ENDIF
374                  zcff = MIN(1._wp, zcff)
375                  !
376                  wi(ji,jj,jk) =           zcff   * wn(ji,jj,jk)
377                  wn(ji,jj,jk) = ( 1._wp - zcff ) * wn(ji,jj,jk)
378                  !
379                  Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient below and in stp_ctl
380               END DO
381            END DO
382         END DO
383         Cu_adv(:,:,1) = 0._wp 
384      ELSE
385         ! Fully explicit everywhere
386         Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient below and in stp_ctl
387         wi    (:,:,:) = 0._wp
388      ENDIF
389      CALL iom_put("wimp",wi) 
390      CALL iom_put("wi_cff",Cu_adv)
391      CALL iom_put("wexp",wn)
392      !
393      IF( ln_timing )   CALL timing_stop('wAimp')
394      !
395   END SUBROUTINE wAimp
396   !!======================================================================
397END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.