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_r12072_MERGE_OPTION2_2019/src/OCE/DYN – NEMO

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/DYN/sshwzv.F90 @ 12166

Last change on this file since 12166 was 12166, checked in by cetlod, 4 years ago

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