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

source: branches/UKMO/dev_r8864_restart_date/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 9235

Last change on this file since 9235 was 9235, checked in by davestorkey, 6 years ago

UKMO/dev_r8864_restart_date : clear SVN keywords.

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