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_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2017/dev_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 8010

Last change on this file since 8010 was 8010, checked in by jchanut, 7 years ago

AGRIF vvl add on

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