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 NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/DYN – NEMO

source: NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/DYN/sshwzv.F90 @ 9863

Last change on this file since 9863 was 9863, checked in by gm, 6 years ago

#1911 (ENHANCE-04): simplified implementation of the Euler stepping at nit000

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