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/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/DYN – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/DYN/sshwzv.F90 @ 8568

Last change on this file since 8568 was 8568, checked in by gm, 7 years ago

#1911 (ENHANCE-09): PART I.2 - _NONE option + remove zts + see associated wiki page

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