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.
Changeset 13463 for NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/DYN/sshwzv.F90 – NEMO

Ignore:
Timestamp:
2020-09-14T17:40:34+02:00 (4 years ago)
Author:
andmirek
Message:

Ticket #2195:update to trunk 13461

Location:
NEMO/branches/2019/dev_r11351_fldread_with_XIOS
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r11351_fldread_with_XIOS

    • Property svn:externals
      •  

        old new  
        33^/utils/build/mk@HEAD         mk 
        44^/utils/tools@HEAD            tools 
        5 ^/vendors/AGRIF/dev@HEAD      ext/AGRIF 
         5^/vendors/AGRIF/dev_r12970_AGRIF_CMEMS      ext/AGRIF 
        66^/vendors/FCM@HEAD            ext/FCM 
        77^/vendors/IOIPSL@HEAD         ext/IOIPSL 
         8 
         9# SETTE 
         10^/utils/CI/sette@13382        sette 
  • NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/DYN/sshwzv.F90

    r11293 r13463  
    99   !!             -   !  2010-09  (D.Storkey and E.O'Dea) bug fixes for BDY module 
    1010   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work 
     11   !!            4.0  !  2018-12  (A. Coward) add mixed implicit/explicit advection 
     12   !!            4.1  !  2019-08  (A. Coward, D. Storkey) Rename ssh_nxt -> ssh_atf. Now only does time filtering. 
    1113   !!---------------------------------------------------------------------- 
    1214 
    1315   !!---------------------------------------------------------------------- 
    1416   !!   ssh_nxt       : after ssh 
    15    !!   ssh_swp       : filter ans swap the ssh arrays 
     17   !!   ssh_atf       : time filter the ssh arrays 
    1618   !!   wzv           : compute now vertical velocity 
    1719   !!---------------------------------------------------------------------- 
    1820   USE oce            ! ocean dynamics and tracers variables 
     21   USE isf_oce        ! ice shelf 
    1922   USE dom_oce        ! ocean space and time domain variables  
    2023   USE sbc_oce        ! surface boundary condition: ocean 
     
    2528   USE bdydyn2d       ! bdy_ssh routine 
    2629#if defined key_agrif 
     30   USE agrif_oce 
    2731   USE agrif_oce_interp 
    2832#endif 
     
    4347   PUBLIC   wzv        ! called by step.F90 
    4448   PUBLIC   wAimp      ! called by step.F90 
    45    PUBLIC   ssh_swp    ! called by step.F90 
     49   PUBLIC   ssh_atf    ! called by step.F90 
    4650 
    4751   !! * Substitutions 
    48 #  include "vectopt_loop_substitute.h90" 
     52#  include "do_loop_substitute.h90" 
     53#  include "domzgr_substitute.h90" 
     54 
    4955   !!---------------------------------------------------------------------- 
    5056   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    5460CONTAINS 
    5561 
    56    SUBROUTINE ssh_nxt( kt ) 
     62   SUBROUTINE ssh_nxt( kt, Kbb, Kmm, pssh, Kaa ) 
    5763      !!---------------------------------------------------------------------- 
    5864      !!                ***  ROUTINE ssh_nxt  *** 
    5965      !!                    
    60       !! ** Purpose :   compute the after ssh (ssha) 
     66      !! ** Purpose :   compute the after ssh (ssh(Kaa)) 
    6167      !! 
    6268      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment 
     
    6470      !!      by the time step. 
    6571      !! 
    66       !! ** action  :   ssha, after sea surface height 
     72      !! ** action  :   ssh(:,:,Kaa), after sea surface height 
    6773      !! 
    6874      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling. 
    6975      !!---------------------------------------------------------------------- 
    70       INTEGER, INTENT(in) ::   kt   ! time step 
     76      INTEGER                         , INTENT(in   ) ::   kt             ! time step 
     77      INTEGER                         , INTENT(in   ) ::   Kbb, Kmm, Kaa  ! time level index 
     78      REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) ::   pssh           ! sea-surface height 
    7179      !  
    72       INTEGER  ::   jk            ! dummy loop indice 
    73       REAL(wp) ::   z2dt, zcoef   ! local scalars 
     80      INTEGER  ::   jk      ! dummy loop index 
     81      REAL(wp) ::   zcoef   ! local scalar 
    7482      REAL(wp), DIMENSION(jpi,jpj) ::   zhdiv   ! 2D workspace 
    7583      !!---------------------------------------------------------------------- 
     
    8391      ENDIF 
    8492      ! 
    85       z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog) 
    86       IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt 
    87       zcoef = 0.5_wp * r1_rau0 
     93      zcoef = 0.5_wp * r1_rho0 
    8894 
    8995      !                                           !------------------------------! 
     
    9197      !                                           !------------------------------! 
    9298      IF(ln_wd_il) THEN 
    93          CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt) 
    94       ENDIF 
    95  
    96       CALL div_hor( kt )                               ! Horizontal divergence 
     99         CALL wad_lmt(pssh(:,:,Kbb), zcoef * (emp_b(:,:) + emp(:,:)), rDt, Kmm, uu, vv ) 
     100      ENDIF 
     101 
     102      CALL div_hor( kt, Kbb, Kmm )                     ! Horizontal divergence 
    97103      ! 
    98104      zhdiv(:,:) = 0._wp 
    99105      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports 
    100         zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk) 
     106        zhdiv(:,:) = zhdiv(:,:) + e3t(:,:,jk,Kmm) * hdiv(:,:,jk) 
    101107      END DO 
    102108      !                                                ! Sea surface elevation time stepping 
     
    104110      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations. 
    105111      !  
    106       ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:) 
     112      pssh(:,:,Kaa) = (  pssh(:,:,Kbb) - rDt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:) 
    107113      ! 
    108114#if defined key_agrif 
     115      Kbb_a = Kbb   ;   Kmm_a = Kmm   ;   Krhs_a = Kaa 
    109116      CALL agrif_ssh( kt ) 
    110117#endif 
     
    112119      IF ( .NOT.ln_dynspg_ts ) THEN 
    113120         IF( ln_bdy ) THEN 
    114             CALL lbc_lnk( 'sshwzv', ssha, 'T', 1. )    ! Not sure that's necessary 
    115             CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries 
     121            CALL lbc_lnk( 'sshwzv', pssh(:,:,Kaa), 'T', 1.0_wp )    ! Not sure that's necessary 
     122            CALL bdy_ssh( pssh(:,:,Kaa) )             ! Duplicate sea level across open boundaries 
    116123         ENDIF 
    117124      ENDIF 
     
    120127      !                                           !------------------------------! 
    121128      ! 
    122       IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask ) 
     129      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pssh(:,:,Kaa), clinfo1=' pssh(:,:,Kaa)  - : ', mask1=tmask ) 
    123130      ! 
    124131      IF( ln_timing )   CALL timing_stop('ssh_nxt') 
     
    127134 
    128135    
    129    SUBROUTINE wzv( kt ) 
     136   SUBROUTINE wzv( kt, Kbb, Kmm, Kaa, pww ) 
    130137      !!---------------------------------------------------------------------- 
    131138      !!                ***  ROUTINE wzv  *** 
     
    138145      !!        The boundary conditions are w=0 at the bottom (no flux) and. 
    139146      !! 
    140       !! ** action  :   wn      : now vertical velocity 
     147      !! ** action  :   pww      : now vertical velocity 
    141148      !! 
    142149      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling. 
    143150      !!---------------------------------------------------------------------- 
    144       INTEGER, INTENT(in) ::   kt   ! time step 
     151      INTEGER                         , INTENT(in)    ::   kt             ! time step 
     152      INTEGER                         , INTENT(in)    ::   Kbb, Kmm, Kaa  ! time level indices 
     153      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   pww            ! vertical velocity at Kmm 
    145154      ! 
    146155      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    147       REAL(wp) ::   z1_2dt       ! local scalars 
    148156      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zhdiv 
    149157      !!---------------------------------------------------------------------- 
     
    156164         IF(lwp) WRITE(numout,*) '~~~~~ ' 
    157165         ! 
    158          wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all) 
     166         pww(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all) 
    159167      ENDIF 
    160168      !                                           !------------------------------! 
    161169      !                                           !     Now Vertical Velocity    ! 
    162170      !                                           !------------------------------! 
    163       z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog) 
    164       IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt 
    165       ! 
    166       IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases 
     171      ! 
     172      !                                               !===============================! 
     173      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      !==  z_tilde and layer cases  ==! 
     174         !                                            !===============================! 
    167175         ALLOCATE( zhdiv(jpi,jpj,jpk) )  
    168176         ! 
     
    170178            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t) 
    171179            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way) 
    172             DO jj = 2, jpjm1 
    173                DO ji = fs_2, fs_jpim1   ! vector opt. 
    174                   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) ) 
    175                END DO 
    176             END DO 
    177          END DO 
    178          CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions" 
     180            DO_2D( 0, 0, 0, 0 ) 
     181               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) ) 
     182            END_2D 
     183         END DO 
     184         CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.0_wp)  ! - ML - Perhaps not necessary: not used for horizontal "connexions" 
    179185         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells? 
    180          !                             ! Same question holds for hdivn. Perhaps just for security 
     186         !                             ! Same question holds for hdiv. Perhaps just for security 
    181187         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence 
    182188            ! computation of w 
    183             wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    & 
    184                &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk) 
    185          END DO 
    186          !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0 
     189            pww(:,:,jk) = pww(:,:,jk+1) - (   e3t(:,:,jk,Kmm) * hdiv(:,:,jk)   & 
     190               &                            +                  zhdiv(:,:,jk)   & 
     191               &                            + r1_Dt * (  e3t(:,:,jk,Kaa)       & 
     192               &                                       - e3t(:,:,jk,Kbb) )   ) * tmask(:,:,jk) 
     193         END DO 
     194         !          IF( ln_vvl_layer ) pww(:,:,:) = 0.e0 
    187195         DEALLOCATE( zhdiv )  
    188       ELSE   ! z_star and linear free surface cases 
     196         !                                            !=================================! 
     197      ELSEIF( ln_linssh )   THEN                      !==  linear free surface cases  ==! 
     198         !                                            !=================================! 
     199         DO jk = jpkm1, 1, -1                               ! integrate from the bottom the hor. divergence 
     200            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)  ) * tmask(:,:,jk) 
     201         END DO 
     202         !                                            !==========================================! 
     203      ELSE                                            !==  Quasi-Eulerian vertical coordinate  ==!   ('key_qco') 
     204         !                                            !==========================================! 
    189205         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence 
    190             ! computation of w 
    191             wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 & 
    192                &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk) 
     206            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)                 & 
     207               &                            + r1_Dt * (  e3t(:,:,jk,Kaa)        & 
     208               &                                       - e3t(:,:,jk,Kbb)  )   ) * tmask(:,:,jk) 
    193209         END DO 
    194210      ENDIF 
     
    196212      IF( ln_bdy ) THEN 
    197213         DO jk = 1, jpkm1 
    198             wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:) 
    199          END DO 
    200       ENDIF 
    201       ! 
    202 #if defined key_agrif  
    203       IF( .NOT. AGRIF_Root() ) THEN  
    204          IF ((nbondi ==  1).OR.(nbondi == 2)) wn(nlci-1 , :     ,:) = 0.e0      ! east  
    205          IF ((nbondi == -1).OR.(nbondi == 2)) wn(2      , :     ,:) = 0.e0      ! west  
    206          IF ((nbondj ==  1).OR.(nbondj == 2)) wn(:      ,nlcj-1 ,:) = 0.e0      ! north  
    207          IF ((nbondj == -1).OR.(nbondj == 2)) wn(:      ,2      ,:) = 0.e0      ! south  
     214            pww(:,:,jk) = pww(:,:,jk) * bdytmask(:,:) 
     215         END DO 
     216      ENDIF 
     217      ! 
     218#if defined key_agrif 
     219      IF( .NOT. AGRIF_Root() ) THEN 
     220         ! 
     221         ! Mask vertical velocity at first/last columns/row  
     222         ! inside computational domain (cosmetic)  
     223         DO jk = 1, jpkm1 
     224            IF( lk_west ) THEN                             ! --- West --- ! 
     225               DO ji = mi0(2+nn_hls), mi1(2+nn_hls) 
     226                  DO jj = 1, jpj 
     227                     pww(ji,jj,jk) = 0._wp  
     228                  END DO 
     229               END DO 
     230            ENDIF 
     231            IF( lk_east ) THEN                             ! --- East --- ! 
     232               DO ji = mi0(jpiglo-1-nn_hls), mi1(jpiglo-1-nn_hls) 
     233                  DO jj = 1, jpj 
     234                     pww(ji,jj,jk) = 0._wp 
     235                  END DO 
     236               END DO 
     237            ENDIF 
     238            IF( lk_south ) THEN                            ! --- South --- ! 
     239               DO jj = mj0(2+nn_hls), mj1(2+nn_hls) 
     240                  DO ji = 1, jpi 
     241                     pww(ji,jj,jk) = 0._wp 
     242                  END DO 
     243               END DO 
     244            ENDIF 
     245            IF( lk_north ) THEN                            ! --- North --- ! 
     246               DO jj = mj0(jpjglo-1-nn_hls), mj1(jpjglo-1-nn_hls) 
     247                  DO ji = 1, jpi 
     248                     pww(ji,jj,jk) = 0._wp 
     249                  END DO 
     250               END DO 
     251            ENDIF 
     252            ! 
     253         END DO 
     254         ! 
    208255      ENDIF  
    209 #endif  
     256#endif 
    210257      ! 
    211258      IF( ln_timing )   CALL timing_stop('wzv') 
     
    214261 
    215262 
    216    SUBROUTINE ssh_swp( kt ) 
    217       !!---------------------------------------------------------------------- 
    218       !!                    ***  ROUTINE ssh_nxt  *** 
    219       !! 
    220       !! ** Purpose :   achieve the sea surface  height time stepping by  
    221       !!              applying Asselin time filter and swapping the arrays 
    222       !!              ssha  already computed in ssh_nxt   
     263   SUBROUTINE ssh_atf( kt, Kbb, Kmm, Kaa, pssh, pssh_f ) 
     264      !!---------------------------------------------------------------------- 
     265      !!                    ***  ROUTINE ssh_atf  *** 
     266      !! 
     267      !! ** Purpose :   Apply Asselin time filter to now SSH. 
    223268      !! 
    224269      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing 
    225270      !!              from the filter, see Leclair and Madec 2010) and swap : 
    226       !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha ) 
    227       !!                            - atfp * rdt * ( emp_b - emp ) / rau0 
    228       !!                sshn = ssha 
    229       !! 
    230       !! ** action  : - sshb, sshn   : before & now sea surface height 
    231       !!                               ready for the next time step 
     271      !!                pssh(:,:,Kmm) = pssh(:,:,Kaa) + rn_atfp * ( pssh(:,:,Kbb) -2 pssh(:,:,Kmm) + pssh(:,:,Kaa) ) 
     272      !!                            - rn_atfp * rn_Dt * ( emp_b - emp ) / rho0 
     273      !! 
     274      !! ** action  : - pssh(:,:,Kmm) time filtered 
    232275      !! 
    233276      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling. 
    234277      !!---------------------------------------------------------------------- 
    235       INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     278      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index 
     279      INTEGER                         , INTENT(in   ) ::   Kbb, Kmm, Kaa  ! ocean time level indices 
     280      REAL(wp), DIMENSION(jpi,jpj,jpt)          , TARGET, INTENT(inout) ::   pssh           ! SSH field 
     281      REAL(wp), DIMENSION(jpi,jpj    ), OPTIONAL, TARGET, INTENT(  out) ::   pssh_f         ! filtered SSH field 
    236282      ! 
    237283      REAL(wp) ::   zcoef   ! local scalar 
    238       !!---------------------------------------------------------------------- 
    239       ! 
    240       IF( ln_timing )   CALL timing_start('ssh_swp') 
     284      REAL(wp), POINTER, DIMENSION(:,:) ::   zssh   ! pointer for filtered SSH  
     285      !!---------------------------------------------------------------------- 
     286      ! 
     287      IF( ln_timing )   CALL timing_start('ssh_atf') 
    241288      ! 
    242289      IF( kt == nit000 ) THEN 
    243290         IF(lwp) WRITE(numout,*) 
    244          IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height' 
     291         IF(lwp) WRITE(numout,*) 'ssh_atf : Asselin time filter of sea surface height' 
    245292         IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
    246293      ENDIF 
    247294      !              !==  Euler time-stepping: no filter, just swap  ==! 
    248       IF ( neuler == 0 .AND. kt == nit000 ) THEN 
    249          sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now) 
    250          ! 
    251       ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==! 
    252          !                                                  ! before <-- now filtered 
    253          sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) ) 
    254          IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed 
    255             zcoef = atfp * rdt * r1_rau0 
    256             sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   & 
    257                &                             -    rnf_b(:,:) + rnf   (:,:)   & 
    258                &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:) 
     295      IF ( .NOT.( l_1st_euler ) ) THEN   ! Only do time filtering for leapfrog timesteps 
     296         IF( PRESENT( pssh_f ) ) THEN   ;   zssh => pssh_f 
     297         ELSE                           ;   zssh => pssh(:,:,Kmm) 
    259298         ENDIF 
    260          sshn(:,:) = ssha(:,:)                              ! now <-- after 
    261       ENDIF 
    262       ! 
    263       IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask ) 
    264       ! 
    265       IF( ln_timing )   CALL timing_stop('ssh_swp') 
    266       ! 
    267    END SUBROUTINE ssh_swp 
    268  
    269    SUBROUTINE wAimp( kt ) 
     299         !                                                  ! filtered "now" field 
     300         pssh(:,:,Kmm) = pssh(:,:,Kmm) + rn_atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) ) 
     301         IF( .NOT.ln_linssh ) THEN                          ! "now" <-- with forcing removed 
     302            zcoef = rn_atfp * rn_Dt * r1_rho0 
     303            pssh(:,:,Kmm) = pssh(:,:,Kmm) - zcoef * (     emp_b(:,:) - emp   (:,:)   & 
     304               &                             - rnf_b(:,:)        + rnf   (:,:)       & 
     305               &                             + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   & 
     306               &                             + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:) 
     307 
     308            ! ice sheet coupling 
     309            IF ( ln_isf .AND. ln_isfcpl .AND. kt == nit000+1) pssh(:,:,Kbb) = pssh(:,:,Kbb) - rn_atfp * rn_Dt * ( risfcpl_ssh(:,:) - 0.0 ) * ssmask(:,:) 
     310 
     311         ENDIF 
     312      ENDIF 
     313      ! 
     314      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab2d_1=pssh(:,:,Kmm), clinfo1=' pssh(:,:,Kmm)  - : ', mask1=tmask ) 
     315      ! 
     316      IF( ln_timing )   CALL timing_stop('ssh_atf') 
     317      ! 
     318   END SUBROUTINE ssh_atf 
     319 
     320    
     321   SUBROUTINE wAimp( kt, Kmm ) 
    270322      !!---------------------------------------------------------------------- 
    271323      !!                ***  ROUTINE wAimp  *** 
     
    276328      !! ** Method  : -  
    277329      !! 
    278       !! ** action  :   wn      : now vertical velocity (to be handled explicitly) 
     330      !! ** action  :   ww      : now vertical velocity (to be handled explicitly) 
    279331      !!            :   wi      : now vertical velocity (for implicit treatment) 
    280332      !! 
    281       !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling. 
     333      !! Reference  : Shchepetkin, A. F. (2015): An adaptive, Courant-number-dependent 
     334      !!              implicit scheme for vertical advection in oceanic modeling.  
     335      !!              Ocean Modelling, 91, 38-69. 
    282336      !!---------------------------------------------------------------------- 
    283337      INTEGER, INTENT(in) ::   kt   ! time step 
     338      INTEGER, INTENT(in) ::   Kmm  ! time level index 
    284339      ! 
    285340      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    286       REAL(wp)             ::   zCu, zcff, z1_e3t                     ! local scalars 
     341      REAL(wp)             ::   zCu, zcff, z1_e3t, zdt                ! local scalars 
    287342      REAL(wp) , PARAMETER ::   Cu_min = 0.15_wp                      ! local parameters 
    288       REAL(wp) , PARAMETER ::   Cu_max = 0.27                         ! local parameters 
     343      REAL(wp) , PARAMETER ::   Cu_max = 0.30_wp                      ! local parameters 
    289344      REAL(wp) , PARAMETER ::   Cu_cut = 2._wp*Cu_max - Cu_min        ! local parameters 
    290345      REAL(wp) , PARAMETER ::   Fcu    = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters 
     
    300355      ENDIF 
    301356      ! 
    302       ! 
    303       DO jk = 1, jpkm1            ! calculate Courant numbers 
    304          DO jj = 2, jpjm1 
    305             DO ji = 2, fs_jpim1   ! vector opt. 
    306                z1_e3t = 1._wp / e3t_n(ji,jj,jk) 
    307                Cu_adv(ji,jj,jk) = 2._wp * rdt * ( ( MAX( wn(ji,jj,jk) , 0._wp ) - MIN( wn(ji,jj,jk+1) , 0._wp ) )   &  ! 2*rdt and not r2dt (for restartability) 
    308                   &                             + ( MAX( e2u(ji  ,jj)*e3u_n(ji  ,jj,jk)*un(ji  ,jj,jk), 0._wp ) -   & 
    309                   &                                 MIN( e2u(ji-1,jj)*e3u_n(ji-1,jj,jk)*un(ji-1,jj,jk), 0._wp ) )   & 
    310                   &                               * r1_e1e2t(ji,jj)                                                 & 
    311                   &                             + ( MAX( e1v(ji,jj  )*e3v_n(ji,jj  ,jk)*vn(ji,jj  ,jk), 0._wp ) -   & 
    312                   &                                 MIN( e1v(ji,jj-1)*e3v_n(ji,jj-1,jk)*vn(ji,jj-1,jk), 0._wp ) )   & 
    313                   &                               * r1_e1e2t(ji,jj)                                                 & 
    314                   &                             ) * z1_e3t 
    315             END DO 
    316          END DO 
    317       END DO 
    318       CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. ) 
     357      ! Calculate Courant numbers 
     358      zdt = 2._wp * rn_Dt                            ! 2*rn_Dt and not rDt (for restartability) 
     359      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
     360         DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     361            z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) 
     362            Cu_adv(ji,jj,jk) =   zdt *                                                         & 
     363               &  ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )            & 
     364               &  + ( MAX( e2u(ji  ,jj) * e3u(ji  ,jj,jk,Kmm)                                  & 
     365               &                        * uu (ji  ,jj,jk,Kmm) + un_td(ji  ,jj,jk), 0._wp ) -   & 
     366               &      MIN( e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm)                                  & 
     367               &                        * uu (ji-1,jj,jk,Kmm) + un_td(ji-1,jj,jk), 0._wp ) )   & 
     368               &                               * r1_e1e2t(ji,jj)                                                                     & 
     369               &  + ( MAX( e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kmm)                                  & 
     370               &                        * vv (ji,jj  ,jk,Kmm) + vn_td(ji,jj  ,jk), 0._wp ) -   & 
     371               &      MIN( e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm)                                  & 
     372               &                        * vv (ji,jj-1,jk,Kmm) + vn_td(ji,jj-1,jk), 0._wp ) )   & 
     373               &                               * r1_e1e2t(ji,jj)                                                                     & 
     374               &                             ) * z1_e3t 
     375         END_3D 
     376      ELSE 
     377         DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     378            z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) 
     379            Cu_adv(ji,jj,jk) =   zdt *                                                      & 
     380               &  ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )         & 
     381               &                             + ( MAX( e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kmm)*uu(ji  ,jj,jk,Kmm), 0._wp ) -   & 
     382               &                                 MIN( e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kmm)*uu(ji-1,jj,jk,Kmm), 0._wp ) )   & 
     383               &                               * r1_e1e2t(ji,jj)                                                 & 
     384               &                             + ( MAX( e1v(ji,jj  )*e3v(ji,jj  ,jk,Kmm)*vv(ji,jj  ,jk,Kmm), 0._wp ) -   & 
     385               &                                 MIN( e1v(ji,jj-1)*e3v(ji,jj-1,jk,Kmm)*vv(ji,jj-1,jk,Kmm), 0._wp ) )   & 
     386               &                               * r1_e1e2t(ji,jj)                                                 & 
     387               &                             ) * z1_e3t 
     388         END_3D 
     389      ENDIF 
     390      CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1.0_wp ) 
    319391      ! 
    320392      CALL iom_put("Courant",Cu_adv) 
    321393      ! 
    322394      IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN       ! Quick check if any breaches anywhere 
    323          DO jk = jpkm1, 2, -1                           ! or scan Courant criterion and partition 
    324             DO jj = 1, jpj                              ! w where necessary 
    325                DO ji = 1, jpi 
    326                   ! 
    327                   zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) ) 
     395         DO_3DS( 1, 1, 1, 1, jpkm1, 2, -1 ) 
     396            ! 
     397            zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) ) 
    328398! alt: 
    329 !                  IF ( wn(ji,jj,jk) > 0._wp ) THEN  
     399!                  IF ( ww(ji,jj,jk) > 0._wp ) THEN  
    330400!                     zCu =  Cu_adv(ji,jj,jk)  
    331401!                  ELSE 
    332402!                     zCu =  Cu_adv(ji,jj,jk-1) 
    333403!                  ENDIF  
    334                   ! 
    335                   IF( zCu <= Cu_min ) THEN              !<-- Fully explicit 
    336                      zcff = 0._wp 
    337                   ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit 
    338                      zcff = ( zCu - Cu_min )**2 
    339                      zcff = zcff / ( Fcu + zcff ) 
    340                   ELSE                                  !<-- Mostly implicit 
    341                      zcff = ( zCu - Cu_max )/ zCu 
    342                   ENDIF 
    343                   zcff = MIN(1._wp, zcff) 
    344                   ! 
    345                   wi(ji,jj,jk) =           zcff   * wn(ji,jj,jk) 
    346                   wn(ji,jj,jk) = ( 1._wp - zcff ) * wn(ji,jj,jk) 
    347                   ! 
    348                   Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient 
    349                END DO 
    350             END DO 
    351          END DO 
     404            ! 
     405            IF( zCu <= Cu_min ) THEN              !<-- Fully explicit 
     406               zcff = 0._wp 
     407            ELSEIF( zCu < Cu_cut ) THEN           !<-- Mixed explicit 
     408               zcff = ( zCu - Cu_min )**2 
     409               zcff = zcff / ( Fcu + zcff ) 
     410            ELSE                                  !<-- Mostly implicit 
     411               zcff = ( zCu - Cu_max )/ zCu 
     412            ENDIF 
     413            zcff = MIN(1._wp, zcff) 
     414            ! 
     415            wi(ji,jj,jk) =           zcff   * ww(ji,jj,jk) 
     416            ww(ji,jj,jk) = ( 1._wp - zcff ) * ww(ji,jj,jk) 
     417            ! 
     418            Cu_adv(ji,jj,jk) = zcff               ! Reuse array to output coefficient below and in stp_ctl 
     419         END_3D 
    352420         Cu_adv(:,:,1) = 0._wp  
    353421      ELSE 
    354422         ! Fully explicit everywhere 
    355          Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient 
     423         Cu_adv(:,:,:) = 0._wp                          ! Reuse array to output coefficient below and in stp_ctl 
    356424         wi    (:,:,:) = 0._wp 
    357425      ENDIF 
    358426      CALL iom_put("wimp",wi)  
    359427      CALL iom_put("wi_cff",Cu_adv) 
    360       CALL iom_put("wexp",wn) 
     428      CALL iom_put("wexp",ww) 
    361429      ! 
    362430      IF( ln_timing )   CALL timing_stop('wAimp') 
Note: See TracChangeset for help on using the changeset viewer.