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

source: branches/2016/dev_INGV_UKMO_2016/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 7351

Last change on this file since 7351 was 7351, checked in by emanuelaclementi, 7 years ago

ticket #1805 step 3: /2016/dev_INGV_UKMO_2016 aligned to the trunk at revision 7161

  • Property svn:keywords set to Id
File size: 12.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        !
25   USE bdy_par        !
26   USE bdydyn2d       ! bdy_ssh routine
27#if defined key_agrif
28   USE agrif_opa_interp
29#endif
30#if defined key_asminc   
31   USE   asminc       ! Assimilation increment
32#endif
33   !
34   USE in_out_manager ! I/O manager
35   USE restart        ! only for lrst_oce
36   USE prtctl         ! Print control
37   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
38   USE lib_mpp        ! MPP library
39   USE wrk_nemo       ! Memory Allocation
40   USE timing         ! Timing
41   USE wet_dry         ! Wetting/Drying flux limting
42
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC   ssh_nxt    ! called by step.F90
47   PUBLIC   wzv        ! called by step.F90
48   PUBLIC   ssh_swp    ! called by step.F90
49
50   !! * Substitutions
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      INTEGER, INTENT(in) ::   kt   ! time step
74      !
75      INTEGER  ::   jk            ! dummy loop indice
76      REAL(wp) ::   z2dt, zcoef   ! local scalars
77      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhdiv   ! 2D workspace
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         IF(lwp) WRITE(numout,*)
86         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
87         IF(lwp) WRITE(numout,*) '~~~~~~~ '
88      ENDIF
89      !
90      CALL div_hor( kt )                              ! Horizontal divergence
91      !
92      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
93      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
94
95      !                                           !------------------------------!
96      !                                           !   After Sea Surface Height   !
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      zcoef = 0.5_wp * r1_rau0
107
108      IF(ln_wd) CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
109
110      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
111
112      IF ( .NOT.ln_dynspg_ts ) THEN
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      !
136      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
137      !
138      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
139      !
140      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
141      !
142   END SUBROUTINE ssh_nxt
143
144   
145   SUBROUTINE wzv( kt )
146      !!----------------------------------------------------------------------
147      !!                ***  ROUTINE wzv  ***
148      !!                   
149      !! ** Purpose :   compute the now vertical velocity
150      !!
151      !! ** Method  : - Using the incompressibility hypothesis, the vertical
152      !!      velocity is computed by integrating the horizontal divergence 
153      !!      from the bottom to the surface minus the scale factor evolution.
154      !!        The boundary conditions are w=0 at the bottom (no flux) and.
155      !!
156      !! ** action  :   wn      : now vertical velocity
157      !!
158      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
159      !!----------------------------------------------------------------------
160      INTEGER, INTENT(in) ::   kt   ! time step
161      !
162      INTEGER  ::   ji, jj, jk   ! dummy loop indices
163      REAL(wp) ::   z1_2dt       ! local scalars
164      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
165      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
166      !!----------------------------------------------------------------------
167      !
168      IF( nn_timing == 1 )   CALL timing_start('wzv')
169      !
170      IF( kt == nit000 ) THEN
171         IF(lwp) WRITE(numout,*)
172         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
173         IF(lwp) WRITE(numout,*) '~~~~~ '
174         !
175         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
176      ENDIF
177      !                                           !------------------------------!
178      !                                           !     Now Vertical Velocity    !
179      !                                           !------------------------------!
180      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
181      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
182      !
183      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
184         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
185         !
186         DO jk = 1, jpkm1
187            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
188            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
189            DO jj = 2, jpjm1
190               DO ji = fs_2, fs_jpim1   ! vector opt.
191                  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) )
192               END DO
193            END DO
194         END DO
195         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
196         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
197         !                             ! Same question holds for hdivn. Perhaps just for security
198         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
199            ! computation of w
200            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
201               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
202         END DO
203         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
204         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
205      ELSE   ! z_star and linear free surface cases
206         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
207            ! computation of w
208            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
209               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
210         END DO
211      ENDIF
212
213#if defined key_bdy
214      IF( lk_bdy ) THEN
215         DO jk = 1, jpkm1
216            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
217         END DO
218      ENDIF
219#endif
220      !
221      IF( nn_timing == 1 )  CALL timing_stop('wzv')
222      !
223   END SUBROUTINE wzv
224
225
226   SUBROUTINE ssh_swp( kt )
227      !!----------------------------------------------------------------------
228      !!                    ***  ROUTINE ssh_nxt  ***
229      !!
230      !! ** Purpose :   achieve the sea surface  height time stepping by
231      !!              applying Asselin time filter and swapping the arrays
232      !!              ssha  already computed in ssh_nxt 
233      !!
234      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
235      !!              from the filter, see Leclair and Madec 2010) and swap :
236      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
237      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
238      !!                sshn = ssha
239      !!
240      !! ** action  : - sshb, sshn   : before & now sea surface height
241      !!                               ready for the next time step
242      !!
243      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
244      !!----------------------------------------------------------------------
245      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
246      !
247      REAL(wp) ::   zcoef   ! local scalar
248      !!----------------------------------------------------------------------
249      !
250      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
251      !
252      IF( kt == nit000 ) THEN
253         IF(lwp) WRITE(numout,*)
254         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
255         IF(lwp) WRITE(numout,*) '~~~~~~~ '
256      ENDIF
257      !              !==  Euler time-stepping: no filter, just swap  ==!
258      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
259         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN
260         sshb(:,:) = sshn(:,:)                              ! before <-- now
261         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
262         !
263      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
264         !                                                  ! before <-- now filtered
265         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
266         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
267            zcoef = atfp * rdt * r1_rau0
268            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
269               &                             -    rnf_b(:,:) + rnf   (:,:)   &
270               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
271         ENDIF
272         sshn(:,:) = ssha(:,:)                              ! now <-- after
273      ENDIF
274      !
275      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
276      !
277      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
278      !
279   END SUBROUTINE ssh_swp
280
281   !!======================================================================
282END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.