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

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 4616

Last change on this file since 4616 was 4616, checked in by gm, 10 years ago

#1260 : see the associated wiki page for explanation

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