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

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 9040

Last change on this file since 9040 was 9023, checked in by timgraham, 7 years ago

Merged METO_MERCATOR branch and resolved all conflicts in OPA_SRC

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