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 NEMO/branches/2019/fix_sn_cfctl_ticket2328/src/OCE/DYN – NEMO

source: NEMO/branches/2019/fix_sn_cfctl_ticket2328/src/OCE/DYN/sshwzv.F90 @ 11872

Last change on this file since 11872 was 11872, checked in by acc, 4 years ago

Branch 2019/fix_sn_cfctl_ticket2328. See #2328. Replacement of ln_ctl and activation of full functionality with
sn_cfctl structure. These changes rename structure components l_mppout and l_mpptop as l_prtctl and l_prttrc
and introduce l_glochk to activate former ln_ctl code in stpctl.F90 to perform global location of min and max
checks. Also added is l_allon which can be used to activate all output (much like the former ln_ctl). If l_allon
is .false. then l_config decides whether or not the suboptions are used.

   sn_cfctl%l_glochk = .FALSE.    ! Range sanity checks are local (F) or global (T). Set T for debugging only
   sn_cfctl%l_allon  = .FALSE.    ! IF T activate all options. If F deactivate all unless l_config is T
   sn_cfctl%l_config = .TRUE.     ! IF .true. then control which reports are written with the remaining options

Note, these changes pass SETTE tests but all references to ln_ctl need to be removed from the sette scripts.

  • Property svn:keywords set to Id
File size: 18.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
[11414]11   !!            4.0  !  2018-12  (A. Coward) add mixed implicit/explicit advection
[3]12   !!----------------------------------------------------------------------
[1438]13
[3]14   !!----------------------------------------------------------------------
[6140]15   !!   ssh_nxt       : after ssh
16   !!   ssh_swp       : filter ans swap the ssh arrays
17   !!   wzv           : compute now vertical velocity
[1438]18   !!----------------------------------------------------------------------
[6140]19   USE oce            ! ocean dynamics and tracers variables
20   USE dom_oce        ! ocean space and time domain variables
21   USE sbc_oce        ! surface boundary condition: ocean
22   USE domvvl         ! Variable volume
23   USE divhor         ! horizontal divergence
24   USE phycst         ! physical constants
[9019]25   USE bdy_oce , ONLY : ln_bdy, bdytmask   ! Open BounDarY
[6140]26   USE bdydyn2d       ! bdy_ssh routine
[2528]27#if defined key_agrif
[9570]28   USE agrif_oce_interp
[2528]29#endif
[6140]30   !
[10364]31   USE iom 
[6140]32   USE in_out_manager ! I/O manager
33   USE restart        ! only for lrst_oce
34   USE prtctl         ! Print control
35   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
36   USE lib_mpp        ! MPP library
37   USE timing         ! Timing
[9023]38   USE wet_dry        ! Wetting/Drying flux limiting
[592]39
[3]40   IMPLICIT NONE
41   PRIVATE
42
[1438]43   PUBLIC   ssh_nxt    ! called by step.F90
[4292]44   PUBLIC   wzv        ! called by step.F90
[10364]45   PUBLIC   wAimp      ! called by step.F90
[4292]46   PUBLIC   ssh_swp    ! called by step.F90
[3]47
48   !! * Substitutions
[1438]49#  include "vectopt_loop_substitute.h90"
[3]50   !!----------------------------------------------------------------------
[9598]51   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[888]52   !! $Id$
[10068]53   !! Software governed by the CeCILL license (see ./LICENSE)
[592]54   !!----------------------------------------------------------------------
[3]55CONTAINS
56
[4292]57   SUBROUTINE ssh_nxt( kt )
[3]58      !!----------------------------------------------------------------------
[4292]59      !!                ***  ROUTINE ssh_nxt  ***
[1438]60      !!                   
[4292]61      !! ** Purpose :   compute the after ssh (ssha)
[3]62      !!
[4292]63      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
64      !!      is computed by integrating the horizontal divergence and multiply by
65      !!      by the time step.
[3]66      !!
[5836]67      !! ** action  :   ssha, after sea surface height
[2528]68      !!
69      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[3]70      !!----------------------------------------------------------------------
[5836]71      INTEGER, INTENT(in) ::   kt   ! time step
[4292]72      !
[7753]73      INTEGER  ::   jk            ! dummy loop indice
[5836]74      REAL(wp) ::   z2dt, zcoef   ! local scalars
[9019]75      REAL(wp), DIMENSION(jpi,jpj) ::   zhdiv   ! 2D workspace
[3]76      !!----------------------------------------------------------------------
[3294]77      !
[9019]78      IF( ln_timing )   CALL timing_start('ssh_nxt')
[3294]79      !
[3]80      IF( kt == nit000 ) THEN
81         IF(lwp) WRITE(numout,*)
[4292]82         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
[1438]83         IF(lwp) WRITE(numout,*) '~~~~~~~ '
[3]84      ENDIF
[2528]85      !
[7646]86      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
[2715]87      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
[7646]88      zcoef = 0.5_wp * r1_rau0
[3]89
[1438]90      !                                           !------------------------------!
91      !                                           !   After Sea Surface Height   !
92      !                                           !------------------------------!
[9023]93      IF(ln_wd_il) THEN
[7753]94         CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
95      ENDIF
[7646]96
[7753]97      CALL div_hor( kt )                               ! Horizontal divergence
[7646]98      !
[7753]99      zhdiv(:,:) = 0._wp
[1438]100      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
[7753]101        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
[1438]102      END DO
103      !                                                ! Sea surface elevation time stepping
[4338]104      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
105      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
106      !
[7753]107      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
[9023]108      !
109#if defined key_agrif
110      CALL agrif_ssh( kt )
111#endif
112      !
[5930]113      IF ( .NOT.ln_dynspg_ts ) THEN
[7646]114         IF( ln_bdy ) THEN
[10425]115            CALL lbc_lnk( 'sshwzv', ssha, 'T', 1. )    ! Not sure that's necessary
[5930]116            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
117         ENDIF
[4292]118      ENDIF
119      !                                           !------------------------------!
120      !                                           !           outputs            !
121      !                                           !------------------------------!
122      !
[11872]123      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask )
[4292]124      !
[9019]125      IF( ln_timing )   CALL timing_stop('ssh_nxt')
[4292]126      !
127   END SUBROUTINE ssh_nxt
128
129   
130   SUBROUTINE wzv( kt )
131      !!----------------------------------------------------------------------
132      !!                ***  ROUTINE wzv  ***
133      !!                   
134      !! ** Purpose :   compute the now vertical velocity
135      !!
136      !! ** Method  : - Using the incompressibility hypothesis, the vertical
137      !!      velocity is computed by integrating the horizontal divergence 
138      !!      from the bottom to the surface minus the scale factor evolution.
139      !!        The boundary conditions are w=0 at the bottom (no flux) and.
140      !!
141      !! ** action  :   wn      : now vertical velocity
142      !!
143      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
144      !!----------------------------------------------------------------------
[5836]145      INTEGER, INTENT(in) ::   kt   ! time step
[4292]146      !
[5836]147      INTEGER  ::   ji, jj, jk   ! dummy loop indices
148      REAL(wp) ::   z1_2dt       ! local scalars
[9019]149      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdiv
[4292]150      !!----------------------------------------------------------------------
151      !
[9019]152      IF( ln_timing )   CALL timing_start('wzv')
[5836]153      !
[4292]154      IF( kt == nit000 ) THEN
155         IF(lwp) WRITE(numout,*)
156         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
157         IF(lwp) WRITE(numout,*) '~~~~~ '
158         !
[7753]159         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
[4292]160      ENDIF
161      !                                           !------------------------------!
162      !                                           !     Now Vertical Velocity    !
163      !                                           !------------------------------!
164      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
165      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
166      !
167      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
[9019]168         ALLOCATE( zhdiv(jpi,jpj,jpk) ) 
[4292]169         !
170         DO jk = 1, jpkm1
171            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
[4338]172            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
[4292]173            DO jj = 2, jpjm1
174               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]175                  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]176               END DO
[592]177            END DO
178         END DO
[10425]179         CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
[4292]180         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
181         !                             ! Same question holds for hdivn. Perhaps just for security
182         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
183            ! computation of w
[7753]184            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
185               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
[4292]186         END DO
187         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
[9019]188         DEALLOCATE( zhdiv ) 
[4292]189      ELSE   ! z_star and linear free surface cases
190         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
[7753]191            ! computation of w
192            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
193               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
[4292]194         END DO
[1438]195      ENDIF
[592]196
[7646]197      IF( ln_bdy ) THEN
[4327]198         DO jk = 1, jpkm1
[7753]199            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
[4327]200         END DO
201      ENDIF
[4292]202      !
[9023]203#if defined key_agrif 
204      IF( .NOT. AGRIF_Root() ) THEN
205         IF ((nbondi ==  1).OR.(nbondi == 2)) wn(nlci-1 , :     ,:) = 0.e0      ! east
206         IF ((nbondi == -1).OR.(nbondi == 2)) wn(2      , :     ,:) = 0.e0      ! west
207         IF ((nbondj ==  1).OR.(nbondj == 2)) wn(:      ,nlcj-1 ,:) = 0.e0      ! north
208         IF ((nbondj == -1).OR.(nbondj == 2)) wn(:      ,2      ,:) = 0.e0      ! south
209      ENDIF 
210#endif 
[5836]211      !
[9124]212      IF( ln_timing )   CALL timing_stop('wzv')
[9023]213      !
[5836]214   END SUBROUTINE wzv
[592]215
216
[4292]217   SUBROUTINE ssh_swp( kt )
[1438]218      !!----------------------------------------------------------------------
219      !!                    ***  ROUTINE ssh_nxt  ***
220      !!
221      !! ** Purpose :   achieve the sea surface  height time stepping by
222      !!              applying Asselin time filter and swapping the arrays
[4292]223      !!              ssha  already computed in ssh_nxt 
[1438]224      !!
[2528]225      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
226      !!              from the filter, see Leclair and Madec 2010) and swap :
227      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
228      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
229      !!                sshn = ssha
[1438]230      !!
231      !! ** action  : - sshb, sshn   : before & now sea surface height
232      !!                               ready for the next time step
[2528]233      !!
234      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[1438]235      !!----------------------------------------------------------------------
[2528]236      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[6140]237      !
238      REAL(wp) ::   zcoef   ! local scalar
[1438]239      !!----------------------------------------------------------------------
[3294]240      !
[9124]241      IF( ln_timing )   CALL timing_start('ssh_swp')
[3294]242      !
[1438]243      IF( kt == nit000 ) THEN
244         IF(lwp) WRITE(numout,*)
[4292]245         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
[1438]246         IF(lwp) WRITE(numout,*) '~~~~~~~ '
247      ENDIF
[6140]248      !              !==  Euler time-stepping: no filter, just swap  ==!
[9239]249      IF ( neuler == 0 .AND. kt == nit000 ) THEN
[7753]250         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
[5836]251         !
[6140]252      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
253         !                                                  ! before <-- now filtered
[7753]254         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
[6140]255         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
256            zcoef = atfp * rdt * r1_rau0
[7753]257            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
258               &                             -    rnf_b(:,:) + rnf   (:,:)   &
259               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
[6140]260         ENDIF
[7753]261         sshn(:,:) = ssha(:,:)                              ! now <-- after
[1438]262      ENDIF
263      !
[11872]264      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask )
[2528]265      !
[9019]266      IF( ln_timing )   CALL timing_stop('ssh_swp')
[3294]267      !
[4292]268   END SUBROUTINE ssh_swp
[3]269
[10364]270   SUBROUTINE wAimp( kt )
271      !!----------------------------------------------------------------------
272      !!                ***  ROUTINE wAimp  ***
273      !!                   
274      !! ** Purpose :   compute the Courant number and partition vertical velocity
275      !!                if a proportion needs to be treated implicitly
276      !!
277      !! ** Method  : -
278      !!
279      !! ** action  :   wn      : now vertical velocity (to be handled explicitly)
280      !!            :   wi      : now vertical velocity (for implicit treatment)
281      !!
[11414]282      !! Reference  : Shchepetkin, A. F. (2015): An adaptive, Courant-number-dependent
283      !!              implicit scheme for vertical advection in oceanic modeling.
284      !!              Ocean Modelling, 91, 38-69.
[10364]285      !!----------------------------------------------------------------------
286      INTEGER, INTENT(in) ::   kt   ! time step
287      !
288      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[11293]289      REAL(wp)             ::   zCu, zcff, z1_e3t                     ! local scalars
[10364]290      REAL(wp) , PARAMETER ::   Cu_min = 0.15_wp                      ! local parameters
[11407]291      REAL(wp) , PARAMETER ::   Cu_max = 0.30_wp                      ! local parameters
[10364]292      REAL(wp) , PARAMETER ::   Cu_cut = 2._wp*Cu_max - Cu_min        ! local parameters
293      REAL(wp) , PARAMETER ::   Fcu    = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters
294      !!----------------------------------------------------------------------
295      !
296      IF( ln_timing )   CALL timing_start('wAimp')
297      !
298      IF( kt == nit000 ) THEN
299         IF(lwp) WRITE(numout,*)
300         IF(lwp) WRITE(numout,*) 'wAimp : Courant number-based partitioning of now vertical velocity '
301         IF(lwp) WRITE(numout,*) '~~~~~ '
[11293]302         wi(:,:,:) = 0._wp
[10364]303      ENDIF
304      !
[11414]305      ! Calculate Courant numbers
306      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
307         DO jk = 1, jpkm1
308            DO jj = 2, jpjm1
309               DO ji = 2, fs_jpim1   ! vector opt.
310                  z1_e3t = 1._wp / e3t_n(ji,jj,jk)
311                  ! 2*rdt and not r2dt (for restartability)
312                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )                       & 
313                     &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk) + un_td(ji  ,jj,jk), 0._wp ) -   &
314                     &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk) + un_td(ji-1,jj,jk), 0._wp ) )   &
315                     &                               * r1_e1e2t(ji,jj)                                                                     &
316                     &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk) + vn_td(ji,jj  ,jk), 0._wp ) -   &
317                     &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk) + vn_td(ji,jj-1,jk), 0._wp ) )   &
318                     &                               * r1_e1e2t(ji,jj)                                                                     &
319                     &                             ) * z1_e3t
320               END DO
[10364]321            END DO
322         END DO
[11414]323      ELSE
324         DO jk = 1, jpkm1
325            DO jj = 2, jpjm1
326               DO ji = 2, fs_jpim1   ! vector opt.
327                  z1_e3t = 1._wp / e3t_n(ji,jj,jk)
328                  ! 2*rdt and not r2dt (for restartability)
329                  Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )   & 
330                     &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk), 0._wp ) -   &
331                     &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk), 0._wp ) )   &
332                     &                               * r1_e1e2t(ji,jj)                                                 &
333                     &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk), 0._wp ) -   &
334                     &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk), 0._wp ) )   &
335                     &                               * r1_e1e2t(ji,jj)                                                 &
336                     &                             ) * z1_e3t
337               END DO
338            END DO
339         END DO
340      ENDIF
[10907]341      CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. )
[10364]342      !
343      CALL iom_put("Courant",Cu_adv)
344      !
345      IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN       ! Quick check if any breaches anywhere
[11293]346         DO jk = jpkm1, 2, -1                           ! or scan Courant criterion and partition
[10907]347            DO jj = 1, jpj                              ! w where necessary
348               DO ji = 1, jpi
[10364]349                  !
[11293]350                  zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) )
351! alt:
352!                  IF ( wn(ji,jj,jk) > 0._wp ) THEN
353!                     zCu =  Cu_adv(ji,jj,jk)
354!                  ELSE
355!                     zCu =  Cu_adv(ji,jj,jk-1)
356!                  ENDIF
[10364]357                  !
[10907]358                  IF( zCu <= Cu_min ) THEN              !<-- Fully explicit
[10364]359                     zcff = 0._wp
360                  ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit
361                     zcff = ( zCu - Cu_min )**2
362                     zcff = zcff / ( Fcu + zcff )
363                  ELSE                                  !<-- Mostly implicit
364                     zcff = ( zCu - Cu_max )/ zCu
365                  ENDIF
366                  zcff = MIN(1._wp, zcff)
367                  !
368                  wi(ji,jj,jk) =           zcff   * wn(ji,jj,jk)
369                  wn(ji,jj,jk) = ( 1._wp - zcff ) * wn(ji,jj,jk)
370                  !
[11407]371                  Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient below and in stp_ctl
[10364]372               END DO
373            END DO
374         END DO
[11293]375         Cu_adv(:,:,1) = 0._wp 
[10364]376      ELSE
377         ! Fully explicit everywhere
[11407]378         Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient below and in stp_ctl
[10907]379         wi    (:,:,:) = 0._wp
[10364]380      ENDIF
381      CALL iom_put("wimp",wi) 
382      CALL iom_put("wi_cff",Cu_adv)
383      CALL iom_put("wexp",wn)
384      !
385      IF( ln_timing )   CALL timing_stop('wAimp')
386      !
387   END SUBROUTINE wAimp
[3]388   !!======================================================================
[1565]389END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.