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

source: branches/UKMO/dev_r5107_kara_mld/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 5244

Last change on this file since 5244 was 5244, checked in by davestorkey, 9 years ago

UKMO Kara MLD branch: remove svn keyword property and clear keywords.

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