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

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 6004

Last change on this file since 6004 was 6004, checked in by gm, 8 years ago

#1613: vvl by default, step III: Merge with the trunk (free surface simplification) (see wiki)

  • Property svn:keywords set to Id
File size: 12.7 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
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 "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      INTEGER, INTENT(in) ::   kt   ! time step
73      !
74      INTEGER  ::   jk            ! dummy loop indice
75      REAL(wp) ::   z2dt, zcoef   ! local scalars
76      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhdiv   ! 2D workspace
77      !!----------------------------------------------------------------------
78      !
79      IF( nn_timing == 1 )   CALL timing_start('ssh_nxt')
80      !
81      CALL wrk_alloc( jpi,jpj,   zhdiv ) 
82      !
83      IF( kt == nit000 ) THEN
84         IF(lwp) WRITE(numout,*)
85         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
86         IF(lwp) WRITE(numout,*) '~~~~~~~ '
87      ENDIF
88      !
89      CALL div_hor( kt )                              ! Horizontal divergence
90      !
91      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
92      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
93
94      !                                           !------------------------------!
95      !                                           !   After Sea Surface Height   !
96      !                                           !------------------------------!
97      zhdiv(:,:) = 0._wp
98      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
99        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
100      END DO
101      !                                                ! Sea surface elevation time stepping
102      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
103      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
104      !
105      zcoef = 0.5_wp * r1_rau0
106      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
107
108      IF ( .NOT.ln_dynspg_ts ) THEN
109         ! These lines are not necessary with time splitting since
110         ! boundary condition on sea level is set during ts loop
111# if defined key_agrif
112         CALL agrif_ssh( kt )
113# endif
114# if defined key_bdy
115         IF( lk_bdy ) THEN
116            CALL lbc_lnk( ssha, 'T', 1. )    ! Not sure that's necessary
117            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
118         ENDIF
119# endif
120      ENDIF
121
122#if defined key_asminc
123      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN     ! Include the IAU weighted SSH increment
124         CALL ssh_asm_inc( kt )
125         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
126      ENDIF
127#endif
128      !                                           !------------------------------!
129      !                                           !           outputs            !
130      !                                           !------------------------------!
131      !
132      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
133      !
134      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
135      !
136      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
137      !
138   END SUBROUTINE ssh_nxt
139
140   
141   SUBROUTINE wzv( kt )
142      !!----------------------------------------------------------------------
143      !!                ***  ROUTINE wzv  ***
144      !!                   
145      !! ** Purpose :   compute the now vertical velocity
146      !!
147      !! ** Method  : - Using the incompressibility hypothesis, the vertical
148      !!      velocity is computed by integrating the horizontal divergence 
149      !!      from the bottom to the surface minus the scale factor evolution.
150      !!        The boundary conditions are w=0 at the bottom (no flux) and.
151      !!
152      !! ** action  :   wn      : now vertical velocity
153      !!
154      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
155      !!----------------------------------------------------------------------
156      INTEGER, INTENT(in) ::   kt   ! time step
157      !
158      INTEGER  ::   ji, jj, jk   ! dummy loop indices
159      REAL(wp) ::   z1_2dt       ! local scalars
160      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
161      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
162      !!----------------------------------------------------------------------
163      !
164      IF( nn_timing == 1 )   CALL timing_start('wzv')
165      !
166      IF( kt == nit000 ) THEN
167         IF(lwp) WRITE(numout,*)
168         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
169         IF(lwp) WRITE(numout,*) '~~~~~ '
170         !
171         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
172      ENDIF
173      !                                           !------------------------------!
174      !                                           !     Now Vertical Velocity    !
175      !                                           !------------------------------!
176      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
177      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
178      !
179      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
180         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
181         !
182         DO jk = 1, jpkm1
183            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
184            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
185            DO jj = 2, jpjm1
186               DO ji = fs_2, fs_jpim1   ! vector opt.
187                  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) )
188               END DO
189            END DO
190         END DO
191         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
192         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
193         !                             ! Same question holds for hdivn. Perhaps just for security
194         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
195            ! computation of w
196            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
197               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
198         END DO
199         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
200         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
201      ELSE   ! z_star and linear free surface cases
202         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
203            ! computation of w
204            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
205               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
206         END DO
207      ENDIF
208
209#if defined key_bdy
210      IF( lk_bdy ) THEN
211         DO jk = 1, jpkm1
212            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
213         END DO
214      ENDIF
215#endif
216      !
217      IF( nn_timing == 1 )  CALL timing_stop('wzv')
218      !
219   END SUBROUTINE wzv
220
221
222   SUBROUTINE ssh_swp( kt )
223      !!----------------------------------------------------------------------
224      !!                    ***  ROUTINE ssh_nxt  ***
225      !!
226      !! ** Purpose :   achieve the sea surface  height time stepping by
227      !!              applying Asselin time filter and swapping the arrays
228      !!              ssha  already computed in ssh_nxt 
229      !!
230      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
231      !!              from the filter, see Leclair and Madec 2010) and swap :
232      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
233      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
234      !!                sshn = ssha
235      !!
236      !! ** action  : - sshb, sshn   : before & now sea surface height
237      !!                               ready for the next time step
238      !!
239      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
240      !!----------------------------------------------------------------------
241      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
242      !
243      REAL(wp) ::   zcoef   ! local scalar
244      !!----------------------------------------------------------------------
245      !
246      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
247      !
248      IF( kt == nit000 ) THEN
249         IF(lwp) WRITE(numout,*)
250         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
251         IF(lwp) WRITE(numout,*) '~~~~~~~ '
252      ENDIF
253      !              !==  Euler time-stepping: no filter, just swap  ==!
254      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
255         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN
256         sshb(:,:) = sshn(:,:)                              ! before <-- now
257         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
258         !
259      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
260         !                                                  ! before <-- now filtered
261         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
262         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
263            zcoef = atfp * rdt * r1_rau0
264            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
265               &                             -    rnf_b(:,:) + rnf   (:,:)   &
266               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
267         ENDIF
268         sshn(:,:) = ssha(:,:)                              ! now <-- after
269      ENDIF
270      !
271      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
272      !
273      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
274      !
275   END SUBROUTINE ssh_swp
276
277   !!======================================================================
278END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.