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 branches/2018/dev_r8864_nemo_v3_6_ZTILDE/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2018/dev_r8864_nemo_v3_6_ZTILDE/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 9529

Last change on this file since 9529 was 9529, checked in by jchanut, 6 years ago

ztilde stability update

  • Property svn:keywords set to Id
File size: 13.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   !!----------------------------------------------------------------------
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 divcur          ! hor. divergence and curl      (div & cur routines)
23   USE restart         ! only for lrst_oce
24   USE in_out_manager  ! I/O manager
25   USE prtctl          ! Print control
26   USE phycst
27   USE lbclnk          ! ocean lateral boundary condition (or mpp link)
28   USE lib_mpp         ! MPP library
29   USE bdy_oce
30   USE bdy_par         
31   USE bdydyn2d        ! bdy_ssh routine
32#if defined key_agrif
33   USE agrif_opa_interp
34#endif
35#if defined key_asminc   
36   USE asminc          ! Assimilation increment
37#endif
38   USE wrk_nemo        ! Memory Allocation
39   USE timing          ! Timing
40
41   IMPLICIT NONE
42   PRIVATE
43
44   PUBLIC   ssh_nxt    ! called by step.F90
45   PUBLIC   wzv        ! called by step.F90
46   PUBLIC   ssh_swp    ! called by step.F90
47
48   !! * Substitutions
49#  include "domzgr_substitute.h90"
50#  include "vectopt_loop_substitute.h90"
51   !!----------------------------------------------------------------------
52   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
53   !! $Id$
54   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
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      !
73      REAL(wp), POINTER, DIMENSION(:,:  ) ::  zhdiv
74      INTEGER, INTENT(in) ::   kt                      ! time step
75      !
76      INTEGER             ::   jk                      ! dummy loop indice
77      REAL(wp)            ::   z2dt, z1_rau0           ! local scalars
78      !!----------------------------------------------------------------------
79      !
80      IF( nn_timing == 1 )  CALL timing_start('ssh_nxt')
81      !
82      CALL wrk_alloc( jpi, jpj, zhdiv ) 
83      !
84      IF( kt == nit000 ) THEN
85         !
86         IF(lwp) WRITE(numout,*)
87         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
88         IF(lwp) WRITE(numout,*) '~~~~~~~ '
89         !
90      ENDIF
91      !
92      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
93      !
94      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
95      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
96
97      !                                           !------------------------------!
98      !                                           !   After Sea Surface Height   !
99      !                                           !------------------------------!
100      zhdiv(:,:) = 0._wp
101      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
102        zhdiv(:,:) = zhdiv(:,:) + fse3t_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      z1_rau0 = 0.5_wp * r1_rau0
109      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
110
111#if ! defined key_dynspg_ts
112      ! These lines are not necessary with time splitting since
113      ! boundary condition on sea level is set during ts loop
114#if defined key_agrif
115      CALL agrif_ssh( kt )
116#endif
117#if defined key_bdy
118      IF (lk_bdy) THEN
119         CALL lbc_lnk( ssha, 'T', 1. ) ! Not sure that's necessary
120         CALL bdy_ssh( ssha ) ! Duplicate sea level across open boundaries
121      ENDIF
122#endif
123#endif
124
125#if defined key_asminc
126      !                                                ! Include the IAU weighted SSH increment
127      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
128         CALL ssh_asm_inc( kt )
129         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
130      ENDIF
131#endif
132
133      !                                           !------------------------------!
134      !                                           !           outputs            !
135      !                                           !------------------------------!
136      !
137      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
138      !
139      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
140      !
141      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
142      !
143   END SUBROUTINE ssh_nxt
144
145   
146   SUBROUTINE wzv( kt, kcall )
147      !!----------------------------------------------------------------------
148      !!                ***  ROUTINE wzv  ***
149      !!                   
150      !! ** Purpose :   compute the now vertical velocity
151      !!
152      !! ** Method  : - Using the incompressibility hypothesis, the vertical
153      !!      velocity is computed by integrating the horizontal divergence 
154      !!      from the bottom to the surface minus the scale factor evolution.
155      !!        The boundary conditions are w=0 at the bottom (no flux) and.
156      !!
157      !! ** action  :   wn      : now vertical velocity
158      !!
159      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
160      !!----------------------------------------------------------------------
161      !
162      INTEGER, INTENT(in) ::   kt           ! time step
163      INTEGER, INTENT( in ), OPTIONAL     :: kcall   ! optional argument
164      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
165      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
166      !
167      LOGICAL             ::   ll_use_totflx
168      INTEGER             ::   ji, jj, jk   ! dummy loop indices
169      REAL(wp)            ::   z1_2dt       ! local scalars
170      !!----------------------------------------------------------------------
171     
172      IF( nn_timing == 1 )  CALL timing_start('wzv')
173      !
174      ll_use_totflx=.FALSE.
175      IF (( PRESENT(kcall) ).AND.(ln_vvl_ztilde .OR. ln_vvl_layer)) THEN
176         IF ( kcall==2 ) ll_use_totflx=.TRUE.
177      ENDIF
178
179      IF( kt == nit000 ) THEN
180         !
181         IF(lwp) WRITE(numout,*)
182         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
183         IF(lwp) WRITE(numout,*) '~~~~~ '
184         !
185         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
186         !
187      ENDIF
188      !                                           !------------------------------!
189      !                                           !     Now Vertical Velocity    !
190      !                                           !------------------------------!
191      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
192      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
193      !
194      IF (( ln_vvl_ztilde .OR. ln_vvl_layer ).AND.ll_use_totflx) THEN      ! z_tilde and layer cases
195         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
196         !
197         DO jk = 1, jpkm1
198            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
199            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
200            DO jj = 2, jpjm1
201               DO ji = fs_2, fs_jpim1   ! vector opt.
202                  zhdiv(ji,jj,jk) = r1_e12t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) )
203               END DO
204            END DO
205         END DO
206         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
207         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
208         !                             ! Same question holds for hdivn. Perhaps just for security
209         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
210            ! computation of w
211            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)                    &
212               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
213         END DO
214         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
215         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
216      ELSE   ! z_star and linear free surface cases
217         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
218            ! computation of w
219            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk)                                   &
220               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
221         END DO
222      ENDIF
223
224#if defined key_bdy
225      IF (lk_bdy) THEN
226         DO jk = 1, jpkm1
227            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
228         END DO
229      ENDIF
230#endif
231      !
232      IF( nn_timing == 1 )  CALL timing_stop('wzv')
233
234
235   END SUBROUTINE wzv
236
237   SUBROUTINE ssh_swp( kt )
238      !!----------------------------------------------------------------------
239      !!                    ***  ROUTINE ssh_nxt  ***
240      !!
241      !! ** Purpose :   achieve the sea surface  height time stepping by
242      !!              applying Asselin time filter and swapping the arrays
243      !!              ssha  already computed in ssh_nxt 
244      !!
245      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
246      !!              from the filter, see Leclair and Madec 2010) and swap :
247      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
248      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
249      !!                sshn = ssha
250      !!
251      !! ** action  : - sshb, sshn   : before & now sea surface height
252      !!                               ready for the next time step
253      !!
254      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
255      !!----------------------------------------------------------------------
256      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
257      !!----------------------------------------------------------------------
258      !
259      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
260      !
261      IF( kt == nit000 ) THEN
262         IF(lwp) WRITE(numout,*)
263         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
264         IF(lwp) WRITE(numout,*) '~~~~~~~ '
265      ENDIF
266
267# if defined key_dynspg_ts
268      IF( ( neuler == 0 .AND. kt == nit000 ) .OR. ln_bt_fw ) THEN    !** Euler time-stepping: no filter
269# else
270      IF ( neuler == 0 .AND. kt == nit000 ) THEN   !** Euler time-stepping at first time-step : no filter
271#endif
272         sshb(:,:) = sshn(:,:)                           ! before <-- now
273         sshn(:,:) = ssha(:,:)                           ! now    <-- after  (before already = now)
274      ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
275         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )     ! before <-- now filtered
276         IF( lk_vvl ) sshb(:,:) = sshb(:,:) - atfp * rdt / rau0 * ( emp_b(:,:)    - emp(:,:)    &
277                                &                                 - rnf_b(:,:)    + rnf(:,:)    &
278                                &                                 + fwfisf_b(:,:) - fwfisf(:,:) ) * ssmask(:,:)
279         sshn(:,:) = ssha(:,:)                           ! now <-- after
280      ENDIF
281      !
282      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
283      !
284      IF( nn_timing == 1 )  CALL timing_stop('ssh_swp')
285      !
286   END SUBROUTINE ssh_swp
287
288   !!======================================================================
289END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.