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

source: branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 8741

Last change on this file since 8741 was 8741, checked in by jchanut, 6 years ago

AGRIF + vvl Main changes - #1965

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