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 @ 7016

Last change on this file since 7016 was 6986, checked in by acc, 8 years ago

Branch dev_r6393_NOC_WAD. Introduced some WAD test cases, tidied and corrected WAD code

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