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

source: branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 5842

Last change on this file since 5842 was 5842, checked in by hliu, 8 years ago

Wetting and Drying update based on r:5803

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