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

source: branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 4263

Last change on this file since 4263 was 4263, checked in by acc, 11 years ago

Branch 2013/dev_r3858_NOC_ZTC, #863. Fixes to allow successful compilation of ORCA2_LIM with key_dynspg_ts

  • Property svn:keywords set to Id
File size: 14.2 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
[3865]10   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
[3]11   !!----------------------------------------------------------------------
[1438]12
[3]13   !!----------------------------------------------------------------------
[3865]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)
[1438]23   USE iom             ! I/O library
[3865]24   USE restart         ! only for lrst_oce
[3]25   USE in_out_manager  ! I/O manager
[258]26   USE prtctl          ! Print control
[592]27   USE phycst
28   USE lbclnk          ! ocean lateral boundary condition (or mpp link)
[2715]29   USE lib_mpp         ! MPP library
[1241]30   USE obc_par         ! open boundary cond. parameter
31   USE obc_oce
[2528]32   USE bdy_oce
[4254]33   USE bdy_par         
34   USE bdydyn2d        ! bdy_ssh routine
[2715]35   USE diaar5, ONLY:   lk_diaar5
[1482]36   USE iom
[4178]37   USE sbcrnf, ONLY: h_rnf, nk_rnf, sbc_rnf_div   ! River runoff
[4263]38   USE dynspg_oce,  ONLY: un_b, vn_b, un_adv, vn_adv
39   USE dynspg_ts,   ONLY: ln_bt_fw
[4178]40   USE dynspg_oce, ONLY: lk_dynspg_ts
[2528]41#if defined key_agrif
42   USE agrif_opa_update
[2486]43   USE agrif_opa_interp
[2528]44#endif
45#if defined key_asminc   
46   USE asminc          ! Assimilation increment
47#endif
[3294]48   USE wrk_nemo        ! Memory Allocation
49   USE timing          ! Timing
[592]50
[3]51   IMPLICIT NONE
52   PRIVATE
53
[1438]54   PUBLIC   ssh_nxt    ! called by step.F90
[4246]55   PUBLIC   wzv        ! called by step.F90
[3865]56   PUBLIC   ssh_swp    ! called by step.F90
[3]57
58   !! * Substitutions
59#  include "domzgr_substitute.h90"
[1438]60#  include "vectopt_loop_substitute.h90"
[3]61   !!----------------------------------------------------------------------
[2528]62   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[888]63   !! $Id$
[2715]64   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[592]65   !!----------------------------------------------------------------------
[3]66CONTAINS
67
[3865]68   SUBROUTINE ssh_nxt( kt )
[3]69      !!----------------------------------------------------------------------
[3865]70      !!                ***  ROUTINE ssh_nxt  ***
[1438]71      !!                   
[3865]72      !! ** Purpose :   compute the after ssh (ssha)
[3]73      !!
[3865]74      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
75      !!      is computed by integrating the horizontal divergence and multiply by
76      !!      by the time step.
[3]77      !!
[1438]78      !! ** action  :   ssha    : after sea surface height
[2528]79      !!
80      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[3]81      !!----------------------------------------------------------------------
[2715]82      !
[3865]83      REAL(wp), POINTER, DIMENSION(:,:  ) ::  zhdiv
84      INTEGER, INTENT(in) ::   kt                      ! time step
85      !
86      INTEGER             ::   jk                      ! dummy loop indice
87      REAL(wp)            ::   z2dt, z1_rau0           ! local scalars
[3]88      !!----------------------------------------------------------------------
[3294]89      !
[3865]90      IF( nn_timing == 1 )  CALL timing_start('ssh_nxt')
[3294]91      !
[3865]92      CALL wrk_alloc( jpi, jpj, zhdiv ) 
[3294]93      !
[3]94      IF( kt == nit000 ) THEN
[2528]95         !
[3]96         IF(lwp) WRITE(numout,*)
[3865]97         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
[1438]98         IF(lwp) WRITE(numout,*) '~~~~~~~ '
99         !
[3]100      ENDIF
[2528]101      !
102      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
103      !
[2715]104      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
105      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
[3]106
[1438]107      !                                           !------------------------------!
108      !                                           !   After Sea Surface Height   !
109      !                                           !------------------------------!
[2715]110      zhdiv(:,:) = 0._wp
[1438]111      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
[3865]112        zhdiv(:,:) = zhdiv(:,:) + fse3t_n(:,:,jk) * hdivn(:,:,jk)
[1438]113      END DO
114      !                                                ! Sea surface elevation time stepping
[2528]115      ! In forward Euler time stepping case, the same formulation as in the leap-frog case can be used
116      ! because emp_b field is initialized with the vlaues of emp field. Hence, 0.5 * ( emp + emp_b ) = emp
[3896]117      z1_rau0 = 0.5_wp * r1_rau0
[2715]118      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * tmask(:,:,1)
[1438]119
[2486]120#if defined key_agrif
[2715]121      CALL agrif_ssh( kt )
[2486]122#endif
[1438]123#if defined key_obc
[2528]124      IF( Agrif_Root() ) THEN
[1438]125         ssha(:,:) = ssha(:,:) * obctmsk(:,:)
[2715]126         CALL lbc_lnk( ssha, 'T', 1. )                 ! absolutly compulsory !! (jmm)
[1438]127      ENDIF
128#endif
[2528]129#if defined key_bdy
[4254]130      ! bg jchanut tschanges
131      ! These lines are not necessary with time splitting since
132      ! boundary condition on sea level is set during ts loop
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
[4254]138      ! end jchanut tschanges
[3764]139#if defined key_asminc
140      !                                                ! Include the IAU weighted SSH increment
141      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
142         CALL ssh_asm_inc( kt )
143         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
144      ENDIF
145#endif
[3865]146
147      !                                           !------------------------------!
148      !                                           !           outputs            !
149      !                                           !------------------------------!
150      CALL iom_put( "ssh" , sshn                  )   ! sea surface height
151      CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height
152      !
153      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
154      !
155      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
156      !
157      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
158      !
159   END SUBROUTINE ssh_nxt
160
161   
[4246]162   SUBROUTINE wzv( kt )
[3865]163      !!----------------------------------------------------------------------
[4246]164      !!                ***  ROUTINE wzv  ***
[3865]165      !!                   
166      !! ** Purpose :   compute the now vertical velocity
167      !!
168      !! ** Method  : - Using the incompressibility hypothesis, the vertical
169      !!      velocity is computed by integrating the horizontal divergence 
170      !!      from the bottom to the surface minus the scale factor evolution.
171      !!        The boundary conditions are w=0 at the bottom (no flux) and.
172      !!
173      !! ** action  :   wn      : now vertical velocity
174      !!
175      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
176      !!----------------------------------------------------------------------
177      !
178      INTEGER, INTENT(in) ::   kt           ! time step
[4246]179      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
180      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
[3865]181      !
182      INTEGER             ::   ji, jj, jk   ! dummy loop indices
183      REAL(wp)            ::   z1_2dt       ! local scalars
184      !!----------------------------------------------------------------------
185     
[4246]186      IF( nn_timing == 1 )  CALL timing_start('wzv')
[3865]187      !
188      IF( kt == nit000 ) THEN
189         !
190         IF(lwp) WRITE(numout,*)
[4246]191         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
[4178]192         IF(lwp) WRITE(numout,*) '~~~~~ '
[3865]193         !
194         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
195         !
196      ENDIF
197      !                                           !------------------------------!
198      !                                           !     Now Vertical Velocity    !
199      !                                           !------------------------------!
200      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
201      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
202      !
203      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
[4246]204         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
205         !
[3865]206         DO jk = 1, jpkm1
207            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
208            ! - ML - note: computation allready done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
209            DO jj = 2, jpjm1
210               DO ji = fs_2, fs_jpim1   ! vector opt.
211                  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) )
212               END DO
[592]213            END DO
214         END DO
[3865]215         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
216         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
217         !                             ! Same question holds for hdivn. Perhaps just for security
218         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
219            ! computation of w
220            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)                    &
221               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
222         END DO
223         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
[4246]224         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
[3865]225      ELSE   ! z_star and linear free surface cases
226         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
227            ! computation of w
228            wn(:,:,jk) = wn(:,:,jk+1) - (   fse3t_n(:,:,jk) * hdivn(:,:,jk)                                   &
229               &                          + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
230         END DO
[1438]231      ENDIF
[592]232
[2528]233#if defined key_bdy
234         wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
235#endif
[4178]236      !
[2528]237      !                                           !------------------------------!
238      !                                           !           outputs            !
239      !                                           !------------------------------!
[3865]240      CALL iom_put( "woce", wn )   ! vertical velocity
[2528]241      IF( lk_diaar5 ) THEN                            ! vertical mass transport & its square value
[4178]242         CALL wrk_alloc( jpi, jpj, z2d ) 
[4246]243         CALL wrk_alloc( jpi, jpj, jpk, z3d ) 
[2528]244         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
[3865]245         z2d(:,:) = rau0 * e12t(:,:)
[1756]246         DO jk = 1, jpk
247            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:)
248         END DO
[2528]249         CALL iom_put( "w_masstr" , z3d                     ) 
250         CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
[4178]251         CALL wrk_dealloc( jpi, jpj, z2d  ) 
[4246]252         CALL wrk_dealloc( jpi, jpj, jpk, z3d ) 
[1756]253      ENDIF
[1438]254      !
[4246]255      IF( nn_timing == 1 )  CALL timing_stop('wzv')
[592]256
257
[4246]258   END SUBROUTINE wzv
[3865]259
260   SUBROUTINE ssh_swp( kt )
[1438]261      !!----------------------------------------------------------------------
262      !!                    ***  ROUTINE ssh_nxt  ***
263      !!
264      !! ** Purpose :   achieve the sea surface  height time stepping by
265      !!              applying Asselin time filter and swapping the arrays
[3865]266      !!              ssha  already computed in ssh_nxt 
[1438]267      !!
[2528]268      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
269      !!              from the filter, see Leclair and Madec 2010) and swap :
270      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
271      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
272      !!                sshn = ssha
[1438]273      !!
274      !! ** action  : - sshb, sshn   : before & now sea surface height
275      !!                               ready for the next time step
[2528]276      !!
277      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[1438]278      !!----------------------------------------------------------------------
[2528]279      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[1438]280      !!----------------------------------------------------------------------
[3294]281      !
[3865]282      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
[3294]283      !
[1438]284      IF( kt == nit000 ) THEN
285         IF(lwp) WRITE(numout,*)
[3865]286         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
[1438]287         IF(lwp) WRITE(numout,*) '~~~~~~~ '
288      ENDIF
[592]289
[4250]290      IF( ( neuler == 0 .AND. kt == nit000 ) .OR. ln_bt_fw ) THEN    !** Euler time-stepping at first time-step : no filter
291         sshb(:,:) = sshn(:,:)                           ! before <-- now
292         sshn(:,:) = ssha(:,:)                           ! now    <-- after  (before already = now)
[3865]293      ELSE                                         !** Leap-Frog time-stepping: Asselin filter + swap
294         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )     ! before <-- now filtered
295         IF( lk_vvl ) sshb(:,:) = sshb(:,:) - atfp * rdt / rau0 * ( emp_b(:,:) - emp(:,:) ) * tmask(:,:,1)
296         sshn(:,:) = ssha(:,:)                           ! now <-- after
[1438]297      ENDIF
298      !
[2528]299      ! Update velocity at AGRIF zoom boundaries
[2486]300#if defined key_agrif
[2528]301      IF ( .NOT.Agrif_Root() ) CALL Agrif_Update_Dyn( kt )
[2486]302#endif
[1438]303      !
[2528]304      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
305      !
[3865]306      IF( nn_timing == 1 )  CALL timing_stop('ssh_swp')
[3294]307      !
[3865]308   END SUBROUTINE ssh_swp
[3]309
310   !!======================================================================
[1565]311END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.