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/UKMO/dev_r5518_sshinc_with_VVL/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/UKMO/dev_r5518_sshinc_with_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 6526

Last change on this file since 6526 was 6526, checked in by drew, 8 years ago

Add an ifdef for key_vvl to isolate this change to VVL.

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