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

source: branches/UKMO/r6232_tracer_advection/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 9295

Last change on this file since 9295 was 9295, checked in by jcastill, 6 years ago

Remove svn keywords

File size: 12.8 KB
RevLine 
[1565]1MODULE sshwzv   
[3]2   !!==============================================================================
[1438]3   !!                       ***  MODULE  sshwzv  ***
4   !! Ocean dynamics : sea surface height and vertical velocity
[3]5   !!==============================================================================
[1438]6   !! History :  3.1  !  2009-02  (G. Madec, M. Leclair)  Original code
[2528]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
[4292]10   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
[3]11   !!----------------------------------------------------------------------
[1438]12
[3]13   !!----------------------------------------------------------------------
[4292]14   !!   ssh_nxt        : after ssh
15   !!   ssh_swp        : filter ans swap the ssh arrays
16   !!   wzv            : compute now vertical velocity
[1438]17   !!----------------------------------------------------------------------
[3]18   USE oce             ! ocean dynamics and tracers variables
19   USE dom_oce         ! ocean space and time domain variables
[888]20   USE sbc_oce         ! surface boundary condition: ocean
21   USE domvvl          ! Variable volume
[1565]22   USE divcur          ! hor. divergence and curl      (div & cur routines)
[4292]23   USE restart         ! only for lrst_oce
[3]24   USE in_out_manager  ! I/O manager
[258]25   USE prtctl          ! Print control
[592]26   USE phycst
27   USE lbclnk          ! ocean lateral boundary condition (or mpp link)
[2715]28   USE lib_mpp         ! MPP library
[2528]29   USE bdy_oce
[4292]30   USE bdy_par         
31   USE bdydyn2d        ! bdy_ssh routine
[2528]32#if defined key_agrif
[2486]33   USE agrif_opa_interp
[2528]34#endif
35#if defined key_asminc   
36   USE asminc          ! Assimilation increment
37#endif
[3294]38   USE wrk_nemo        ! Memory Allocation
39   USE timing          ! Timing
[592]40
[3]41   IMPLICIT NONE
42   PRIVATE
43
[1438]44   PUBLIC   ssh_nxt    ! called by step.F90
[4292]45   PUBLIC   wzv        ! called by step.F90
46   PUBLIC   ssh_swp    ! called by step.F90
[3]47
48   !! * Substitutions
49#  include "domzgr_substitute.h90"
[1438]50#  include "vectopt_loop_substitute.h90"
[3]51   !!----------------------------------------------------------------------
[2528]52   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[888]53   !! $Id$
[2715]54   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[592]55   !!----------------------------------------------------------------------
[3]56CONTAINS
57
[4292]58   SUBROUTINE ssh_nxt( kt )
[3]59      !!----------------------------------------------------------------------
[4292]60      !!                ***  ROUTINE ssh_nxt  ***
[1438]61      !!                   
[4292]62      !! ** Purpose :   compute the after ssh (ssha)
[3]63      !!
[4292]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.
[3]67      !!
[1438]68      !! ** action  :   ssha    : after sea surface height
[2528]69      !!
70      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[3]71      !!----------------------------------------------------------------------
[2715]72      !
[4292]73      REAL(wp), POINTER, DIMENSION(:,:  ) ::  zhdiv
74      INTEGER, INTENT(in) ::   kt                      ! time step
75      !
76      INTEGER             ::   jk                      ! dummy loop indice
77      REAL(wp)            ::   z2dt, z1_rau0           ! local scalars
[3]78      !!----------------------------------------------------------------------
[3294]79      !
[4292]80      IF( nn_timing == 1 )  CALL timing_start('ssh_nxt')
[3294]81      !
[4292]82      CALL wrk_alloc( jpi, jpj, zhdiv ) 
[3294]83      !
[3]84      IF( kt == nit000 ) THEN
[2528]85         !
[3]86         IF(lwp) WRITE(numout,*)
[4292]87         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
[1438]88         IF(lwp) WRITE(numout,*) '~~~~~~~ '
89         !
[3]90      ENDIF
[2528]91      !
92      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
93      !
[2715]94      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
95      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
[3]96
[1438]97      !                                           !------------------------------!
98      !                                           !   After Sea Surface Height   !
99      !                                           !------------------------------!
[2715]100      zhdiv(:,:) = 0._wp
[1438]101      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
[4292]102        zhdiv(:,:) = zhdiv(:,:) + fse3t_n(:,:,jk) * hdivn(:,:,jk)
[1438]103      END DO
104      !                                                ! Sea surface elevation time stepping
[4338]105      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
106      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
107      !
[4292]108      z1_rau0 = 0.5_wp * r1_rau0
[4990]109      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
[1438]110
[4486]111#if ! defined key_dynspg_ts
112      ! These lines are not necessary with time splitting since
113      ! boundary condition on sea level is set during ts loop
[2486]114#if defined key_agrif
[2715]115      CALL agrif_ssh( kt )
[2486]116#endif
[2528]117#if defined key_bdy
[4292]118      IF (lk_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
[2528]122#endif
[4486]123#endif
124
[3764]125#if defined key_asminc
126      !                                                ! Include the IAU weighted SSH increment
127      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
128         CALL ssh_asm_inc( kt )
129         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
130      ENDIF
131#endif
[4292]132
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      !
162      INTEGER, INTENT(in) ::   kt           ! time step
163      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
164      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
165      !
166      INTEGER             ::   ji, jj, jk   ! dummy loop indices
167      REAL(wp)            ::   z1_2dt       ! local scalars
168      !!----------------------------------------------------------------------
169     
170      IF( nn_timing == 1 )  CALL timing_start('wzv')
171      !
172      IF( kt == nit000 ) THEN
173         !
174         IF(lwp) WRITE(numout,*)
175         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
176         IF(lwp) WRITE(numout,*) '~~~~~ '
177         !
178         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
179         !
180      ENDIF
181      !                                           !------------------------------!
182      !                                           !     Now Vertical Velocity    !
183      !                                           !------------------------------!
184      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
185      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
186      !
187      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
188         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
189         !
190         DO jk = 1, jpkm1
191            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
[4338]192            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
[4292]193            DO jj = 2, jpjm1
194               DO ji = fs_2, fs_jpim1   ! vector opt.
195                  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) )
196               END DO
[592]197            END DO
198         END DO
[4292]199         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
200         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
201         !                             ! Same question holds for hdivn. Perhaps just for security
202         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
203            ! computation of w
204            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)                    &
205               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
206         END DO
207         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
208         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
209      ELSE   ! z_star and linear free surface cases
210         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
211            ! computation of w
212            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk)                                   &
213               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
214         END DO
[1438]215      ENDIF
[592]216
[2528]217#if defined key_bdy
[4327]218      IF (lk_bdy) THEN
219         DO jk = 1, jpkm1
220            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
221         END DO
222      ENDIF
[2528]223#endif
[4292]224      !
225      IF( nn_timing == 1 )  CALL timing_stop('wzv')
[592]226
227
[4292]228   END SUBROUTINE wzv
229
230   SUBROUTINE ssh_swp( kt )
[1438]231      !!----------------------------------------------------------------------
232      !!                    ***  ROUTINE ssh_nxt  ***
233      !!
234      !! ** Purpose :   achieve the sea surface  height time stepping by
235      !!              applying Asselin time filter and swapping the arrays
[4292]236      !!              ssha  already computed in ssh_nxt 
[1438]237      !!
[2528]238      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
239      !!              from the filter, see Leclair and Madec 2010) and swap :
240      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
241      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
242      !!                sshn = ssha
[1438]243      !!
244      !! ** action  : - sshb, sshn   : before & now sea surface height
245      !!                               ready for the next time step
[2528]246      !!
247      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[1438]248      !!----------------------------------------------------------------------
[2528]249      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[1438]250      !!----------------------------------------------------------------------
[3294]251      !
[4292]252      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
[3294]253      !
[1438]254      IF( kt == nit000 ) THEN
255         IF(lwp) WRITE(numout,*)
[4292]256         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
[1438]257         IF(lwp) WRITE(numout,*) '~~~~~~~ '
258      ENDIF
[592]259
[4292]260# if defined key_dynspg_ts
261      IF( ( neuler == 0 .AND. kt == nit000 ) .OR. ln_bt_fw ) THEN    !** Euler time-stepping: no filter
262# else
263      IF ( neuler == 0 .AND. kt == nit000 ) THEN   !** Euler time-stepping at first time-step : no filter
264#endif
265         sshb(:,:) = sshn(:,:)                           ! before <-- now
266         sshn(:,:) = ssha(:,:)                           ! now    <-- after  (before already = now)
267      ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
268         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )     ! before <-- now filtered
[5628]269         IF( lk_vvl ) sshb(:,:) = sshb(:,:) - atfp * rdt / rau0 * ( emp_b(:,:)    - emp(:,:)    &
270                                &                                 - rnf_b(:,:)    + rnf(:,:)    &
271                                &                                 + fwfisf_b(:,:) - fwfisf(:,:) ) * ssmask(:,:)
[4292]272         sshn(:,:) = ssha(:,:)                           ! now <-- after
[1438]273      ENDIF
274      !
[2528]275      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
276      !
[4292]277      IF( nn_timing == 1 )  CALL timing_stop('ssh_swp')
[3294]278      !
[4292]279   END SUBROUTINE ssh_swp
[3]280
281   !!======================================================================
[1565]282END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.