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

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate_IAUtest/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 12303

Last change on this file since 12303 was 12303, checked in by kingr, 4 years ago

Precalculate ht_n and use when applying SSH increment.

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