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

source: branches/2016/dev_r6393_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 6604

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

wet/dry in ssh_nxt

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