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

source: branches/UKMO/r6232_collate_bgc_diagnostics/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 11134

Last change on this file since 11134 was 11134, checked in by jcastill, 5 years ago

Full set of changes as in the original branch

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