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/2014/dev_r4650_UKMO14.12_STAND_ALONE_OBSOPER/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2014/dev_r4650_UKMO14.12_STAND_ALONE_OBSOPER/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 5600

Last change on this file since 5600 was 5600, checked in by andrewryan, 9 years ago

merged in latest version of trunk alongside changes to SAO_SRC to be compatible with latest OBS

  • Property svn:keywords set to Id
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
33   USE agrif_opa_update
[2486]34   USE agrif_opa_interp
[2528]35#endif
36#if defined key_asminc   
37   USE asminc          ! Assimilation increment
38#endif
[3294]39   USE wrk_nemo        ! Memory Allocation
40   USE timing          ! Timing
[592]41
[3]42   IMPLICIT NONE
43   PRIVATE
44
[1438]45   PUBLIC   ssh_nxt    ! called by step.F90
[4292]46   PUBLIC   wzv        ! called by step.F90
47   PUBLIC   ssh_swp    ! called by step.F90
[3]48
49   !! * Substitutions
50#  include "domzgr_substitute.h90"
[1438]51#  include "vectopt_loop_substitute.h90"
[3]52   !!----------------------------------------------------------------------
[2528]53   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[888]54   !! $Id$
[2715]55   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[592]56   !!----------------------------------------------------------------------
[3]57CONTAINS
58
[4292]59   SUBROUTINE ssh_nxt( kt )
[3]60      !!----------------------------------------------------------------------
[4292]61      !!                ***  ROUTINE ssh_nxt  ***
[1438]62      !!                   
[4292]63      !! ** Purpose :   compute the after ssh (ssha)
[3]64      !!
[4292]65      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
66      !!      is computed by integrating the horizontal divergence and multiply by
67      !!      by the time step.
[3]68      !!
[1438]69      !! ** action  :   ssha    : after sea surface height
[2528]70      !!
71      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[3]72      !!----------------------------------------------------------------------
[2715]73      !
[4292]74      REAL(wp), POINTER, DIMENSION(:,:  ) ::  zhdiv
75      INTEGER, INTENT(in) ::   kt                      ! time step
76      !
77      INTEGER             ::   jk                      ! dummy loop indice
78      REAL(wp)            ::   z2dt, z1_rau0           ! local scalars
[3]79      !!----------------------------------------------------------------------
[3294]80      !
[4292]81      IF( nn_timing == 1 )  CALL timing_start('ssh_nxt')
[3294]82      !
[4292]83      CALL wrk_alloc( jpi, jpj, zhdiv ) 
[3294]84      !
[3]85      IF( kt == nit000 ) THEN
[2528]86         !
[3]87         IF(lwp) WRITE(numout,*)
[4292]88         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
[1438]89         IF(lwp) WRITE(numout,*) '~~~~~~~ '
90         !
[3]91      ENDIF
[2528]92      !
93      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
94      !
[2715]95      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
96      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
[3]97
[1438]98      !                                           !------------------------------!
99      !                                           !   After Sea Surface Height   !
100      !                                           !------------------------------!
[2715]101      zhdiv(:,:) = 0._wp
[1438]102      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
[4292]103        zhdiv(:,:) = zhdiv(:,:) + fse3t_n(:,:,jk) * hdivn(:,:,jk)
[1438]104      END DO
105      !                                                ! Sea surface elevation time stepping
[4338]106      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
107      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
108      !
[4292]109      z1_rau0 = 0.5_wp * r1_rau0
[5034]110      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
[1438]111
[4486]112#if ! defined key_dynspg_ts
113      ! These lines are not necessary with time splitting since
114      ! boundary condition on sea level is set during ts loop
[2486]115#if defined key_agrif
[2715]116      CALL agrif_ssh( kt )
[2486]117#endif
[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
[4292]225      !
226      IF( nn_timing == 1 )  CALL timing_stop('wzv')
[592]227
228
[4292]229   END SUBROUTINE wzv
230
231   SUBROUTINE ssh_swp( kt )
[1438]232      !!----------------------------------------------------------------------
233      !!                    ***  ROUTINE ssh_nxt  ***
234      !!
235      !! ** Purpose :   achieve the sea surface  height time stepping by
236      !!              applying Asselin time filter and swapping the arrays
[4292]237      !!              ssha  already computed in ssh_nxt 
[1438]238      !!
[2528]239      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
240      !!              from the filter, see Leclair and Madec 2010) and swap :
241      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
242      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
243      !!                sshn = ssha
[1438]244      !!
245      !! ** action  : - sshb, sshn   : before & now sea surface height
246      !!                               ready for the next time step
[2528]247      !!
248      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[1438]249      !!----------------------------------------------------------------------
[2528]250      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[1438]251      !!----------------------------------------------------------------------
[3294]252      !
[4292]253      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
[3294]254      !
[1438]255      IF( kt == nit000 ) THEN
256         IF(lwp) WRITE(numout,*)
[4292]257         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
[1438]258         IF(lwp) WRITE(numout,*) '~~~~~~~ '
259      ENDIF
[592]260
[4292]261# if defined key_dynspg_ts
262      IF( ( neuler == 0 .AND. kt == nit000 ) .OR. ln_bt_fw ) THEN    !** Euler time-stepping: no filter
263# else
264      IF ( neuler == 0 .AND. kt == nit000 ) THEN   !** Euler time-stepping at first time-step : no filter
265#endif
266         sshb(:,:) = sshn(:,:)                           ! before <-- now
267         sshn(:,:) = ssha(:,:)                           ! now    <-- after  (before already = now)
268      ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
269         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )     ! before <-- now filtered
[5600]270         IF( lk_vvl ) sshb(:,:) = sshb(:,:) - atfp * rdt / rau0 * ( emp_b(:,:) - emp(:,:) - rnf_b(:,:) + rnf(:,:) ) * ssmask(:,:)
[4292]271         sshn(:,:) = ssha(:,:)                           ! now <-- after
[1438]272      ENDIF
273      !
[2528]274      ! Update velocity at AGRIF zoom boundaries
[2486]275#if defined key_agrif
[2528]276      IF ( .NOT.Agrif_Root() ) CALL Agrif_Update_Dyn( kt )
[2486]277#endif
[1438]278      !
[2528]279      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
280      !
[4292]281      IF( nn_timing == 1 )  CALL timing_stop('ssh_swp')
[3294]282      !
[4292]283   END SUBROUTINE ssh_swp
[3]284
285   !!======================================================================
[1565]286END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.