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

source: NEMO/trunk/src/OCE/DYN/sshwzv.F90 @ 10907

Last change on this file since 10907 was 10907, checked in by clem, 5 years ago

solve ticket #2275

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