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

source: branches/2017/dev_METO_2017/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 8868

Last change on this file since 8868 was 8868, checked in by timgraham, 6 years ago

Merged dev_r8789_sbc into branch

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