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.
zdftke.F90 in NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/OCE/ZDF – NEMO

source: NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/OCE/ZDF/zdftke.F90 @ 14601

Last change on this file since 14601 was 14601, checked in by francesca, 3 years ago

[comm_cleanup: ZDF files] - ticket #2607

  • Property svn:keywords set to Id
File size: 52.0 KB
RevLine 
[1531]1MODULE zdftke
[1239]2   !!======================================================================
[1531]3   !!                       ***  MODULE  zdftke  ***
[14072]4   !! Ocean physics:  vertical mixing coefficient computed from the tke
[1239]5   !!                 turbulent closure parameterization
6   !!=====================================================================
[1492]7   !! History :  OPA  !  1991-03  (b. blanke)  Original code
8   !!            7.0  !  1991-11  (G. Madec)   bug fix
9   !!            7.1  !  1992-10  (G. Madec)   new mixing length and eav
10   !!            7.2  !  1993-03  (M. Guyon)   symetrical conditions
11   !!            7.3  !  1994-08  (G. Madec, M. Imbard)  nn_pdl flag
12   !!            7.5  !  1996-01  (G. Madec)   s-coordinates
13   !!            8.0  !  1997-07  (G. Madec)   lbc
14   !!            8.1  !  1999-01  (E. Stretta) new option for the mixing length
15   !!  NEMO      1.0  !  2002-06  (G. Madec) add tke_init routine
16   !!             -   !  2004-10  (C. Ethe )  1D configuration
17   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
18   !!            3.0  !  2008-05  (C. Ethe,  G.Madec) : update TKE physics:
19   !!                 !           - tke penetration (wind steering)
20   !!                 !           - suface condition for tke & mixing length
21   !!                 !           - Langmuir cells
22   !!             -   !  2008-05  (J.-M. Molines, G. Madec)  2D form of avtb
23   !!             -   !  2008-06  (G. Madec)  style + DOCTOR name for namelist parameters
[14072]24   !!             -   !  2008-12  (G. Reffray) stable discretization of the production term
25   !!            3.2  !  2009-06  (G. Madec, S. Masson) TKE restart compatible with key_cpl
[1492]26   !!                 !                                + cleaning of the parameters + bugs correction
[2528]27   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
[5120]28   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
[14072]29   !!            4.0  !  2017-04  (G. Madec)  remove CPP ddm key & avm at t-point only
[13461]30   !!             -   !  2017-05  (G. Madec)  add top/bottom friction as boundary condition
[14007]31   !!            4.2  !  2020-12  (G. Madec, E. Clementi) add wave coupling
32   !                  !           following Couvelard et al., 2019
[1239]33   !!----------------------------------------------------------------------
[9019]34
[1239]35   !!----------------------------------------------------------------------
[3625]36   !!   zdf_tke       : update momentum and tracer Kz from a tke scheme
37   !!   tke_tke       : tke time stepping: update tke at now time step (en)
38   !!   tke_avn       : compute mixing length scale and deduce avm and avt
39   !!   zdf_tke_init  : initialization, namelist read, and parameters control
40   !!   tke_rst       : read/write tke restart in ocean restart file
[1239]41   !!----------------------------------------------------------------------
[2528]42   USE oce            ! ocean: dynamics and active tracers variables
43   USE phycst         ! physical constants
44   USE dom_oce        ! domain: ocean
45   USE domvvl         ! domain: variable volume layer
[1492]46   USE sbc_oce        ! surface boundary condition: ocean
[9019]47   USE zdfdrg         ! vertical physics: top/bottom drag coef.
[2528]48   USE zdfmxl         ! vertical physics: mixed layer
[13007]49#if defined key_si3
50   USE ice, ONLY: hm_i, h_i
51#endif
52#if defined key_cice
53   USE sbc_ice, ONLY: h_i
54#endif
[9019]55   !
[1492]56   USE in_out_manager ! I/O manager
57   USE iom            ! I/O manager library
[2715]58   USE lib_mpp        ! MPP library
[9019]59   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
60   USE prtctl         ! Print control
[14072]61   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
[14007]62   USE sbcwave        ! Surface boundary waves
[1239]63
64   IMPLICIT NONE
65   PRIVATE
66
[2528]67   PUBLIC   zdf_tke        ! routine called in step module
68   PUBLIC   zdf_tke_init   ! routine called in opa module
69   PUBLIC   tke_rst        ! routine called in step module
[1239]70
[4147]71   !                      !!** Namelist  namzdf_tke  **
72   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
[14007]73   LOGICAL  ::   ln_mxhsw  ! mixing length scale surface value as a fonction of wave height
[13472]74   INTEGER  ::   nn_mxlice ! type of scaling under sea-ice (=0/1/2/3)
75   REAL(wp) ::   rn_mxlice ! ice thickness value when scaling under sea-ice
[4147]76   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
77   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
78   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
79   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
[14072]80   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
[4147]81   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
82   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
83   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
84   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
85   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
[9019]86   INTEGER  ::      nn_htau   ! type of tke profile of penetration (=0/1)
[14007]87   INTEGER  ::   nn_bc_surf! surface condition (0/1=Dir/Neum) ! Only applicable for wave coupling
88   INTEGER  ::   nn_bc_bot ! surface condition (0/1=Dir/Neum) ! Only applicable for wave coupling
[9019]89   REAL(wp) ::      rn_efr    ! fraction of TKE surface value which penetrates in the ocean
[4147]90   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
[9019]91   REAL(wp) ::      rn_lc     ! coef to compute vertical velocity of Langmuir cells
[14072]92   INTEGER  ::   nn_eice   ! attenutaion of langmuir & surface wave breaking under ice (=0/1/2/3)
[1239]93
[4147]94   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values)
95   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m]
[2528]96   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3)
97   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3)
[1239]98
[9019]99   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau    ! depth of tke penetration (nn_htau)
100   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl   ! now mixing lenght of dissipation
101   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr   ! now mixing lenght of dissipation
[1492]102
[1239]103   !! * Substitutions
[12377]104#  include "do_loop_substitute.h90"
[13237]105#  include "domzgr_substitute.h90"
[1239]106   !!----------------------------------------------------------------------
[9598]107   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[2528]108   !! $Id$
[10068]109   !! Software governed by the CeCILL license (see ./LICENSE)
[1239]110   !!----------------------------------------------------------------------
111CONTAINS
112
[2715]113   INTEGER FUNCTION zdf_tke_alloc()
114      !!----------------------------------------------------------------------
115      !!                ***  FUNCTION zdf_tke_alloc  ***
116      !!----------------------------------------------------------------------
[9019]117      ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) ,   STAT= zdf_tke_alloc )
118      !
[10425]119      CALL mpp_sum ( 'zdftke', zdf_tke_alloc )
120      IF( zdf_tke_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_alloc: failed to allocate arrays' )
[2715]121      !
122   END FUNCTION zdf_tke_alloc
123
124
[12377]125   SUBROUTINE zdf_tke( kt, Kbb, Kmm, p_sh2, p_avm, p_avt )
[1239]126      !!----------------------------------------------------------------------
[1531]127      !!                   ***  ROUTINE zdf_tke  ***
[1239]128      !!
129      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
[1492]130      !!              coefficients using a turbulent closure scheme (TKE).
[1239]131      !!
[1492]132      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
133      !!              is computed from a prognostic equation :
134      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
135      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
136      !!                  + avt N^2                      ! stratif. destruc.
137      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
[1239]138      !!      with the boundary conditions:
[1695]139      !!         surface: en = max( rn_emin0, rn_ebb * taum )
[1239]140      !!         bottom : en = rn_emin
[14072]141      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
[1492]142      !!
[14072]143      !!        The now Turbulent kinetic energy is computed using the following
[1492]144      !!      time stepping: implicit for vertical diffusion term, linearized semi
[14072]145      !!      implicit for kolmogoroff dissipation term, and explicit forward for
146      !!      both buoyancy and shear production terms. Therefore a tridiagonal
[1492]147      !!      linear system is solved. Note that buoyancy and shear terms are
148      !!      discretized in a energy conserving form (Bruchard 2002).
149      !!
150      !!        The dissipative and mixing length scale are computed from en and
151      !!      the stratification (see tke_avn)
152      !!
153      !!        The now vertical eddy vicosity and diffusivity coefficients are
[14072]154      !!      given by:
[1492]155      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
[14072]156      !!              avt = max( avmb, pdl * avm                 )
[1239]157      !!              eav = max( avmb, avm )
[1492]158      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
[14072]159      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
[1239]160      !!
161      !! ** Action  :   compute en (now turbulent kinetic energy)
[9019]162      !!                update avt, avm (before vertical eddy coef.)
[1239]163      !!
164      !! References : Gaspar et al., JGR, 1990,
165      !!              Blanke and Delecluse, JPO, 1991
166      !!              Mellor and Blumberg, JPO 2004
167      !!              Axell, JGR, 2002
[1492]168      !!              Bruchard OM 2002
[1239]169      !!----------------------------------------------------------------------
[9019]170      INTEGER                   , INTENT(in   ) ::   kt             ! ocean time step
[12377]171      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]172      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_sh2          ! shear production term
173      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   p_avm, p_avt   !  momentum and tracer Kz (w-points)
[1492]174      !!----------------------------------------------------------------------
[1481]175      !
[12377]176      CALL tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )   ! now tke (en)
[5656]177      !
[12377]178      CALL tke_avn( Kbb, Kmm,        p_avm, p_avt )   ! now avt, avm, dissl
[3632]179      !
[5656]180  END SUBROUTINE zdf_tke
[1239]181
[1492]182
[12377]183   SUBROUTINE tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )
[1239]184      !!----------------------------------------------------------------------
[1492]185      !!                   ***  ROUTINE tke_tke  ***
186      !!
187      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
188      !!
189      !! ** Method  : - TKE surface boundary condition
[2528]190      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
[9019]191      !!              - source term due to shear (= Kz dz[Ub] * dz[Un] )
[1492]192      !!              - Now TKE : resolution of the TKE equation by inverting
193      !!                a tridiagonal linear system by a "methode de chasse"
194      !!              - increase TKE due to surface and internal wave breaking
[14072]195      !!             NB: when sea-ice is present, both LC parameterization
196      !!                 and TKE penetration are turned off when the ice fraction
197      !!                 is smaller than 0.25
[1492]198      !!
199      !! ** Action  : - en : now turbulent kinetic energy)
[1239]200      !! ---------------------------------------------------------------------
[9019]201      USE zdf_oce , ONLY : en   ! ocean vertical physics
202      !!
[12377]203      INTEGER                    , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]204      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_sh2          ! shear production term
205      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
206      !
[13472]207      INTEGER ::   ji, jj, jk                  ! dummy loop arguments
[9019]208      REAL(wp) ::   zetop, zebot, zmsku, zmskv ! local scalars
209      REAL(wp) ::   zrhoa  = 1.22              ! Air density kg/m3
210      REAL(wp) ::   zcdrag = 1.5e-3            ! drag coefficient
[13472]211      REAL(wp) ::   zbbrau, zbbirau, zri       ! local scalars
212      REAL(wp) ::   zfact1, zfact2, zfact3     !   -      -
213      REAL(wp) ::   ztx2  , zty2  , zcof       !   -      -
214      REAL(wp) ::   ztau  , zdif               !   -      -
215      REAL(wp) ::   zus   , zwlc  , zind       !   -      -
216      REAL(wp) ::   zzd_up, zzd_lw             !   -      -
[14007]217      REAL(wp) ::   ztaui, ztauj, z1_norm
[9019]218      INTEGER , DIMENSION(jpi,jpj)     ::   imlc
[14007]219      REAL(wp), DIMENSION(jpi,jpj)     ::   zice_fra, zhlc, zus3, zWlc2
[9019]220      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
[1239]221      !!--------------------------------------------------------------------
[1492]222      !
[13472]223      zbbrau  = rn_ebb / rho0       ! Local constant initialisation
224      zbbirau = 3.75_wp / rho0
[14072]225      zfact1  = -.5_wp * rn_Dt
[13472]226      zfact2  = 1.5_wp * rn_Dt * rn_ediss
227      zfact3  = 0.5_wp         * rn_ediss
[1492]228      !
[14007]229      zpelc(:,:,:) = 0._wp ! need to be initialised in case ln_lc is not used
230      !
[13472]231      ! ice fraction considered for attenuation of langmuir & wave breaking
232      SELECT CASE ( nn_eice )
233      CASE( 0 )   ;   zice_fra(:,:) = 0._wp
234      CASE( 1 )   ;   zice_fra(:,:) =        TANH( fr_i(:,:) * 10._wp )
235      CASE( 2 )   ;   zice_fra(:,:) =              fr_i(:,:)
236      CASE( 3 )   ;   zice_fra(:,:) = MIN( 4._wp * fr_i(:,:) , 1._wp )
237      END SELECT
238      !
[1492]239      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]240      !                     !  Surface/top/bottom boundary condition on tke
[1492]241      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[13472]242      !
[14601]243      ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
244      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14057]245         en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) )
[14007]246         zdiag(ji,jj,1) = 1._wp/en(ji,jj,1)
[14072]247         zd_lw(ji,jj,1) = 1._wp
[14007]248         zd_up(ji,jj,1) = 0._wp
[12377]249      END_2D
[9019]250      !
[1492]251      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
252      !                     !  Bottom boundary condition on tke
253      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[1719]254      !
[12489]255      !   en(bot)   = (ebb0/rho0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
[9019]256      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
257      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
[1492]258      !
[13461]259      IF( .NOT.ln_drg_OFF ) THEN    !== friction used as top/bottom boundary condition on TKE
[9019]260         !
[14601]261         ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )        ! bottom friction
262         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )        ! bottom friction
[12377]263            zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
264            zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
[12489]265            !                       ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
[12377]266            zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mbkt(ji,jj),Kbb)+uu(ji-1,jj,mbkt(ji,jj),Kbb) ) )**2  &
267               &                                           + ( zmskv*( vv(ji,jj,mbkt(ji,jj),Kbb)+vv(ji,jj-1,mbkt(ji,jj),Kbb) ) )**2  )
268            en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
269         END_2D
[13461]270         IF( ln_isfcav ) THEN
[14601]271            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )     ! top friction
272            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )     ! top friction
[12377]273               zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
274               zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
[12489]275               !                             ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
[12377]276               zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mikt(ji,jj),Kbb)+uu(ji-1,jj,mikt(ji,jj),Kbb) ) )**2  &
277                  &                                           + ( zmskv*( vv(ji,jj,mikt(ji,jj),Kbb)+vv(ji,jj-1,mikt(ji,jj),Kbb) ) )**2  )
[12702]278               ! (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj) = 1 where ice shelves are present
279               en(ji,jj,mikt(ji,jj)) = en(ji,jj,1)           * tmask(ji,jj,1) &
280                  &                  + MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj)
[12377]281            END_2D
[9019]282         ENDIF
283         !
284      ENDIF
285      !
[1492]286      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[14007]287      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke (Axell JGR 2002)
[1492]288         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
[1239]289         !
[14007]290         !                       !* Langmuir velocity scale
291         !
292         IF ( cpl_sdrftx )  THEN       ! Surface Stokes Drift available
293            !                                ! Craik-Leibovich velocity scale Wlc = ( u* u_s )^1/2    with u* = (taum/rho0)^1/2
294            !                                ! associated kinetic energy : 1/2 (Wlc)^2 = u* u_s
295            !                                ! more precisely, it is the dot product that must be used :
296            !                                !     1/2  (W_lc)^2 = MAX( u* u_s + v* v_s , 0 )   only the positive part
297!!gm  ! PS: currently we don't have neither the 2 stress components at t-point !nor the angle between u* and u_s
298!!gm  ! so we will overestimate the LC velocity....   !!gm I will do the work if !LC have an effect !
[14601]299            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
300            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]301!!XC                  zWlc2(ji,jj) = 0.5_wp * SQRT( taum(ji,jj) * r1_rho0 * ( ut0sd(ji,jj)**2 +vt0sd(ji,jj)**2 )  )
302                  zWlc2(ji,jj) = 0.5_wp *  ( ut0sd(ji,jj)**2 +vt0sd(ji,jj)**2 )
303            END_2D
304!
305!  Projection of Stokes drift in the wind stress direction
306!
[14601]307            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
308            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]309                  ztaui   = 0.5_wp * ( utau(ji,jj) + utau(ji-1,jj) )
310                  ztauj   = 0.5_wp * ( vtau(ji,jj) + vtau(ji,jj-1) )
311                  z1_norm = 1._wp / MAX( SQRT(ztaui*ztaui+ztauj*ztauj), 1.e-12 ) * tmask(ji,jj,1)
312                  zWlc2(ji,jj) = 0.5_wp * z1_norm * ( MAX( ut0sd(ji,jj)*ztaui + vt0sd(ji,jj)*ztauj, 0._wp ) )**2
313            END_2D
[14601]314            ! [comm_cleanup]
315            IF (nn_hls.eq.1) CALL lbc_lnk      ( 'zdftke', zWlc2, 'T', 1. )
[14007]316!
317         ELSE                          ! Surface Stokes drift deduced from surface stress
318            !                                ! Wlc = u_s   with u_s = 0.016*U_10m, the surface stokes drift  (Axell 2002, Eq.44)
319            !                                ! using |tau| = rho_air Cd |U_10m|^2 , it comes:
320            !                                ! Wlc = 0.016 * [|tau|/(rho_air Cdrag) ]^1/2   and thus:
321            !                                ! 1/2 Wlc^2 = 0.5 * 0.016 * 0.016 |tau| /( rho_air Cdrag )
322            zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )      ! to convert stress in 10m wind using a constant drag
[14601]323            ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
324            DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[14007]325               zWlc2(ji,jj) = zcof * taum(ji,jj)
326            END_2D
327            !
328         ENDIF
329         !
330         !                       !* Depth of the LC circulation  (Axell 2002, Eq.47)
331         !                             !- LHS of Eq.47
[12377]332         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * gdepw(:,:,1,Kmm) * e3w(:,:,1,Kmm)
[1239]333         DO jk = 2, jpk
[13237]334            zpelc(:,:,jk)  = zpelc(:,:,jk-1) +   &
335               &        MAX( rn2b(:,:,jk), 0._wp ) * gdepw(:,:,jk,Kmm) * e3w(:,:,jk,Kmm)
[1239]336         END DO
[14007]337         !
338         !                             !- compare LHS to RHS of Eq.47
[7753]339         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
[14601]340         ! [comm_cleanup] ! DO_3DS( 1, 1, 1, 1, jpkm1, 2, -1 )
341         DO_3DS( nn_hls, nn_hls, nn_hls, nn_hls, jpkm1, 2, -1 )
[14007]342            IF( zpelc(ji,jj,jk) > zWlc2(ji,jj) )   imlc(ji,jj) = jk
[12377]343         END_3D
[1492]344         !                               ! finite LC depth
[14601]345         ! [comm_cleanup] ! DO_2D( 1, 1, 1, 1 )
346         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[12377]347            zhlc(ji,jj) = gdepw(ji,jj,imlc(ji,jj),Kmm)
348         END_2D
[14007]349         !
[1705]350         zcof = 0.016 / SQRT( zrhoa * zcdrag )
[14601]351         ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
352         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]353            zus = SQRT( 2. * zWlc2(ji,jj) )             ! Stokes drift
[13472]354            zus3(ji,jj) = MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok
[12377]355         END_2D
[14601]356         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )                  !* TKE Langmuir circulation source term added to en
357         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )                  !* TKE Langmuir circulation source term added to en
[14072]358            IF ( zus3(ji,jj) /= 0._wp ) THEN
[12377]359               IF ( gdepw(ji,jj,jk,Kmm) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN
360                  !                                           ! vertical velocity due to LC
[13472]361                  zwlc = rn_lc * SIN( rpi * gdepw(ji,jj,jk,Kmm) / zhlc(ji,jj) )
[12377]362                  !                                           ! TKE Langmuir circulation source term
[13472]363                  en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * zus3(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
[12377]364               ENDIF
365            ENDIF
366         END_3D
[1239]367         !
368      ENDIF
[1492]369      !
370      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
371      !                     !  Now Turbulent kinetic energy (output in en)
372      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
373      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
374      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
375      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
376      !
[13497]377      IF( nn_pdl == 1 ) THEN          !* Prandtl number = F( Ri )
[14601]378         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
379         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]380            !                             ! local Richardson number
[13226]381            IF (rn2b(ji,jj,jk) <= 0.0_wp) then
382                zri = 0.0_wp
383            ELSE
384                zri = rn2b(ji,jj,jk) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
385            ENDIF
[12377]386            !                             ! inverse of Prandtl number
387            apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
388         END_3D
[5656]389      ENDIF
[14072]390      !
[14601]391      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )   !* Matrix and right hand side in en
392      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )   !* Matrix and right hand side in en
[12377]393         zcof   = zfact1 * tmask(ji,jj,jk)
394         !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
395         !                                   ! eddy coefficient (ensure numerical stability)
396         zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
[13472]397            &          /    (  e3t(ji,jj,jk  ,Kmm) * e3w(ji,jj,jk  ,Kmm)  )
[12377]398         zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
[13472]399            &          /    (  e3t(ji,jj,jk-1,Kmm) * e3w(ji,jj,jk  ,Kmm)  )
[12377]400         !
401         zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
402         zd_lw(ji,jj,jk) = zzd_lw
403         zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
404         !
405         !                                   ! right hand side in en
[12489]406         en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * (  p_sh2(ji,jj,jk)                        &   ! shear
[12377]407            &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
408            &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
409            &                                ) * wmask(ji,jj,jk)
410      END_3D
[14007]411      !
412      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
413      !                     !  Surface boundary condition on tke if
414      !                     !  coupling with waves
415      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
416      !
417      IF ( cpl_phioc .and. ln_phioc )  THEN
[14072]418         SELECT CASE (nn_bc_surf) ! Boundary Condition using surface TKE flux from waves
[14007]419
420         CASE ( 0 ) ! Dirichlet BC
[14601]421            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )    ! en(1)   = rn_ebb taum / rho0  (min value rn_emin0)
422            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )    ! en(1)   = rn_ebb taum / rho0  (min value rn_emin0)
[14007]423               IF ( phioc(ji,jj) < 0 )  phioc(ji,jj) = 0._wp
424               en(ji,jj,1) = MAX( rn_emin0, .5 * ( 15.8 * phioc(ji,jj) / rho0 )**(2./3.) )  * tmask(ji,jj,1)
425               zdiag(ji,jj,1) = 1._wp/en(ji,jj,1)  ! choose to keep coherence with former estimation of
426            END_2D
427
428         CASE ( 1 ) ! Neumann BC
[14601]429            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
430            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]431               IF ( phioc(ji,jj) < 0 )  phioc(ji,jj) = 0._wp
432               en(ji,jj,2)    = en(ji,jj,2) + ( rn_Dt * phioc(ji,jj) / rho0 ) /e3w(ji,jj,2,Kmm)
433               en(ji,jj,1)    = en(ji,jj,2) + (2 * e3t(ji,jj,1,Kmm) * phioc(ji,jj)/rho0) / ( p_avm(ji,jj,1) + p_avm(ji,jj,2) )
434               zdiag(ji,jj,2) = zdiag(ji,jj,2) + zd_lw(ji,jj,2)
435               zdiag(ji,jj,1) = 1._wp
436               zd_lw(ji,jj,2) = 0._wp
437            END_2D
438
439         END SELECT
440
441      ENDIF
442      !
[5120]443      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
[14601]444      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )                ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
445      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )                ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
[12377]446         zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
447      END_3D
[14007]448!XC : commented to allow for neumann boundary condition
449!      DO_2D( 0, 0, 0, 0 )
450!         zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
451!      END_2D
[14601]452      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
453      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]454         zd_lw(ji,jj,jk) = en(ji,jj,jk) - zd_lw(ji,jj,jk) / zdiag(ji,jj,jk-1) *zd_lw(ji,jj,jk-1)
455      END_3D
[14601]456      ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )                          ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
457      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                          ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
[12377]458         en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
459      END_2D
[14601]460      ! [comm_cleanup] ! DO_3DS( 0, 0, 0, 0, jpk-2, 2, -1 )
461      DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpk-2, 2, -1 )
[12377]462         en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
463      END_3D
[14601]464      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )                ! set the minimum value of tke
465      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )                ! set the minimum value of tke
[12377]466         en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
467      END_3D
[9019]468      !
[1492]469      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
470      !                            !  TKE due to surface and internal wave breaking
471      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[6140]472!!gm BUG : in the exp  remove the depth of ssh !!!
[12377]473!!gm       i.e. use gde3w in argument (gdepw(:,:,:,Kmm))
[13530]474      !
475      ! penetration is partly switched off below sea-ice if nn_eice/=0
476      !
[2528]477      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
[14601]478         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
479         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]480            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
[13472]481               &                                 * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[12377]482         END_3D
[2528]483      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
[14601]484         ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
485         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]486            jk = nmln(ji,jj)
487            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
[13472]488               &                                 * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[12377]489         END_2D
[2528]490      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
[14601]491         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
492         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]493            ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
494            zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
[14072]495            ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
496            zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
[12377]497            zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
498            en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
[13472]499               &                        * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[12377]500         END_3D
[1239]501      ENDIF
[1492]502      !
[1239]503   END SUBROUTINE tke_tke
504
[1492]505
[12377]506   SUBROUTINE tke_avn( Kbb, Kmm, p_avm, p_avt )
[1239]507      !!----------------------------------------------------------------------
[1492]508      !!                   ***  ROUTINE tke_avn  ***
[1239]509      !!
[1492]510      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
511      !!
[14072]512      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
513      !!              the tke_tke routine). First, the now mixing lenth is
[1492]514      !!      computed from en and the strafification (N^2), then the mixings
515      !!      coefficients are computed.
516      !!              - Mixing length : a first evaluation of the mixing lengh
517      !!      scales is:
[14072]518      !!                      mxl = sqrt(2*en) / N
[1492]519      !!      where N is the brunt-vaisala frequency, with a minimum value set
[2528]520      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
[14072]521      !!        The mixing and dissipative length scale are bound as follow :
[1492]522      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
523      !!                        zmxld = zmxlm = mxl
524      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
[14072]525      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
[1492]526      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
527      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
[14072]528      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
529      !!                    the surface to obtain ldown. the resulting length
[1492]530      !!                    scales are:
[14072]531      !!                        zmxld = sqrt( lup * ldown )
[1492]532      !!                        zmxlm = min ( lup , ldown )
533      !!              - Vertical eddy viscosity and diffusivity:
534      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
[14072]535      !!                      avt = max( avmb, pdlr * avm )
[1492]536      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
537      !!
[9019]538      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
[1239]539      !!----------------------------------------------------------------------
[9019]540      USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d   ! ocean vertical physics
541      !!
[12377]542      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]543      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
544      !
[2715]545      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[9019]546      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
547      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
[13007]548      REAL(wp) ::   zemxl, zemlm, zemlp, zmaxice       !   -      -
[9019]549      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmxlm, zmxld   ! 3D workspace
[1239]550      !!--------------------------------------------------------------------
[3294]551      !
[1492]552      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
553      !                     !  Mixing length
554      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
555      !
556      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
557      !
[5120]558      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
[14072]559      zmxlm(:,:,:)  = rmxl_min
[7753]560      zmxld(:,:,:)  = rmxl_min
[14007]561      !
562      IF(ln_sdw .AND. ln_mxhsw) THEN
563         zmxlm(:,:,1)= vkarmn * MAX ( 1.6 * hsw(:,:) , 0.02 )        ! surface mixing length = F(wave height)
564         ! from terray et al 1999 and mellor and blumberg 2004 it should be 0.85 and not 1.6
565         zcoef       = vkarmn * ( (rn_ediff*rn_ediss)**0.25 ) / rn_ediff
566         zmxlm(:,:,1)= zcoef * MAX ( 1.6 * hsw(:,:) , 0.02 )        ! surface mixing length = F(wave height)
567      ELSE
[14072]568      !
[14007]569         IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rho0*g)
[13007]570         !
[14007]571            zraug = vkarmn * 2.e5_wp / ( rho0 * grav )
[13007]572#if ! defined key_si3 && ! defined key_cice
[14601]573            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )                  ! No sea-ice
574            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )                  ! No sea-ice
[14007]575               zmxlm(ji,jj,1) =  zraug * taum(ji,jj) * tmask(ji,jj,1)
576            END_2D
[13007]577#else
[14007]578            SELECT CASE( nn_mxlice )             ! Type of scaling under sea-ice
[13007]579            !
[14007]580            CASE( 0 )                      ! No scaling under sea-ice
[14601]581               ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
582               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]583                  zmxlm(ji,jj,1) = zraug * taum(ji,jj) * tmask(ji,jj,1)
584               END_2D
585               !
586            CASE( 1 )                      ! scaling with constant sea-ice thickness
[14601]587               ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
588               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]589                  zmxlm(ji,jj,1) =  ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
590                     &                          fr_i(ji,jj)   * rn_mxlice           ) * tmask(ji,jj,1)
591               END_2D
592               !
593            CASE( 2 )                      ! scaling with mean sea-ice thickness
[14601]594               ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
595               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[13007]596#if defined key_si3
[14007]597                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
598                     &                         fr_i(ji,jj)   * hm_i(ji,jj) * 2._wp ) * tmask(ji,jj,1)
[13007]599#elif defined key_cice
[14007]600                  zmaxice = MAXVAL( h_i(ji,jj,:) )
601                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
602                     &                         fr_i(ji,jj)   * zmaxice             ) * tmask(ji,jj,1)
[13007]603#endif
[14007]604               END_2D
605               !
606            CASE( 3 )                      ! scaling with max sea-ice thickness
[14601]607               ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
608               DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]609                  zmaxice = MAXVAL( h_i(ji,jj,:) )
610                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
611                     &                         fr_i(ji,jj)   * zmaxice             ) * tmask(ji,jj,1)
612               END_2D
613               !
614            END SELECT
615#endif
[13007]616            !
[14601]617            ! [comm_cleanup] ! DO_2D( 0, 0, 0, 0 )
618            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[14007]619               zmxlm(ji,jj,1) = MAX( rn_mxl0, zmxlm(ji,jj,1) )
[13007]620            END_2D
621            !
[14007]622         ELSE
623            zmxlm(:,:,1) = rn_mxl0
624         ENDIF
[1239]625      ENDIF
626      !
[14601]627      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
628      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]629         zrn2 = MAX( rn2(ji,jj,jk), rsmall )
630         zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
631      END_3D
[1492]632      !
633      !                     !* Physical limits for the mixing length
634      !
[14072]635      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
[7753]636      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
[1492]637      !
[1239]638      SELECT CASE ( nn_mxl )
639      !
[5836]640 !!gm Not sure of that coding for ISF....
[12377]641      ! where wmask = 0 set zmxlm == e3w(:,:,:,Kmm)
[1239]642      CASE ( 0 )           ! bounded by the distance to surface and bottom
[14601]643         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
644         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]645            zemxl = MIN( gdepw(ji,jj,jk,Kmm) - gdepw(ji,jj,mikt(ji,jj),Kmm), zmxlm(ji,jj,jk),   &
646            &            gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) - gdepw(ji,jj,jk,Kmm) )
647            ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
[13237]648            zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk)   &
649               &            + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
650            zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk)   &
651               &            + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
[12377]652         END_3D
[1239]653         !
654      CASE ( 1 )           ! bounded by the vertical scale factor
[14601]655         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
656         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]657            zemxl = MIN( e3w(ji,jj,jk,Kmm), zmxlm(ji,jj,jk) )
658            zmxlm(ji,jj,jk) = zemxl
659            zmxld(ji,jj,jk) = zemxl
660         END_3D
[1239]661         !
662      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
[14601]663         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )       ! from the surface to the bottom :
664         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )       ! from the surface to the bottom :
[13237]665            zmxlm(ji,jj,jk) =   &
666               &    MIN( zmxlm(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
[12377]667         END_3D
[14601]668         ! [comm_cleanup] ! DO_3DS( 0, 0, 0, 0, jpkm1, 2, -1 )   ! from the bottom to the surface :
669         DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpkm1, 2, -1 )   ! from the bottom to the surface :
[12377]670            zemxl = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
671            zmxlm(ji,jj,jk) = zemxl
672            zmxld(ji,jj,jk) = zemxl
673         END_3D
[1239]674         !
675      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
[14601]676         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )        ! from the surface to the bottom : lup
677         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )        ! from the surface to the bottom : lup
[13237]678            zmxld(ji,jj,jk) =    &
679               &    MIN( zmxld(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
[12377]680         END_3D
[14601]681         ! [comm_cleanup] ! DO_3DS( 0, 0, 0, 0, jpkm1, 2, -1 )   ! from the bottom to the surface : ldown
682         DO_3DS( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, jpkm1, 2, -1 )   ! from the bottom to the surface : ldown
[13237]683            zmxlm(ji,jj,jk) =   &
684               &    MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
[12377]685         END_3D
[14601]686         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
687         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12377]688            zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
689            zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
690            zmxlm(ji,jj,jk) = zemlm
691            zmxld(ji,jj,jk) = zemlp
692         END_3D
[1239]693         !
694      END SELECT
[1492]695      !
696      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]697      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
[1492]698      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[14601]699      ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 1, jpkm1 )   !* vertical eddy viscosity & diffivity at w-points
700      DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 1, jpkm1 )   !* vertical eddy viscosity & diffivity at w-points
[12377]701         zsqen = SQRT( en(ji,jj,jk) )
702         zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
703         p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
704         p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
705         dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
706      END_3D
[1492]707      !
708      !
[13497]709      IF( nn_pdl == 1 ) THEN          !* Prandtl number case: update avt
[14601]710         ! [comm_cleanup] ! DO_3D( 0, 0, 0, 0, 2, jpkm1 )
711         DO_3D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1, 2, jpkm1 )
[12698]712            p_avt(ji,jj,jk)   = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
[12377]713         END_3D
[1239]714      ENDIF
[9019]715      !
[12377]716      IF(sn_cfctl%l_prtctl) THEN
[9440]717         CALL prt_ctl( tab3d_1=en   , clinfo1=' tke  - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk)
718         CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke  - m: ', kdim=jpk )
[1239]719      ENDIF
720      !
[1492]721   END SUBROUTINE tke_avn
[1239]722
[1492]723
[12377]724   SUBROUTINE zdf_tke_init( Kmm )
[1239]725      !!----------------------------------------------------------------------
[2528]726      !!                  ***  ROUTINE zdf_tke_init  ***
[14072]727      !!
728      !! ** Purpose :   Initialization of the vertical eddy diffivity and
[1492]729      !!              viscosity when using a tke turbulent closure scheme
[1239]730      !!
[1601]731      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
[1492]732      !!              called at the first timestep (nit000)
[1239]733      !!
[1601]734      !! ** input   :   Namlist namzdf_tke
[1239]735      !!
736      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
737      !!----------------------------------------------------------------------
[9019]738      USE zdf_oce , ONLY : ln_zdfiwm   ! Internal Wave Mixing flag
739      !!
[12377]740      INTEGER, INTENT(in) ::   Kmm          ! time level index
741      INTEGER             ::   ji, jj, jk   ! dummy loop indices
742      INTEGER             ::   ios
[1239]743      !!
[13007]744      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb   , rn_emin  ,  &
745         &                 rn_emin0, rn_bshear, nn_mxl   , ln_mxl0  ,  &
746         &                 rn_mxl0 , nn_mxlice, rn_mxlice,             &
[13461]747         &                 nn_pdl  , ln_lc    , rn_lc    ,             &
[14072]748         &                 nn_etau , nn_htau  , rn_efr   , nn_eice  ,  &
[14007]749         &                 nn_bc_surf, nn_bc_bot, ln_mxhsw
[1239]750      !!----------------------------------------------------------------------
[2715]751      !
[4147]752      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
[11536]753901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' )
[4147]754
755      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
[11536]756902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' )
[4624]757      IF(lwm) WRITE ( numond, namzdf_tke )
[2715]758      !
[2528]759      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
[2715]760      !
[1492]761      IF(lwp) THEN                    !* Control print
[1239]762         WRITE(numout,*)
[2528]763         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
764         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]765         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
[1705]766         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
767         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
768         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
769         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
770         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
[9019]771         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
[1705]772         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
773         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
[9019]774         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
775         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
[13472]776         IF( ln_mxl0 ) THEN
777            WRITE(numout,*) '      type of scaling under sea-ice               nn_mxlice = ', nn_mxlice
778            IF( nn_mxlice == 1 ) &
779            WRITE(numout,*) '      ice thickness when scaling under sea-ice    rn_mxlice = ', rn_mxlice
780            SELECT CASE( nn_mxlice )             ! Type of scaling under sea-ice
781            CASE( 0 )   ;   WRITE(numout,*) '   ==>>>   No scaling under sea-ice'
782            CASE( 1 )   ;   WRITE(numout,*) '   ==>>>   scaling with constant sea-ice thickness'
783            CASE( 2 )   ;   WRITE(numout,*) '   ==>>>   scaling with mean sea-ice thickness'
784            CASE( 3 )   ;   WRITE(numout,*) '   ==>>>   scaling with max sea-ice thickness'
785            CASE DEFAULT
786               CALL ctl_stop( 'zdf_tke_init: wrong value for nn_mxlice, should be 0,1,2,3 or 4')
787            END SELECT
788         ENDIF
[9019]789         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
790         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
[14007]791         IF ( cpl_phioc .and. ln_phioc )  THEN
792            SELECT CASE( nn_bc_surf)             ! Type of scaling under sea-ice
793            CASE( 0 )   ;   WRITE(numout,*) '  nn_bc_surf=0 ==>>> DIRICHLET SBC using surface TKE flux from waves'
794            CASE( 1 )   ;   WRITE(numout,*) '  nn_bc_surf=1 ==>>> NEUMANN SBC using surface TKE flux from waves'
795            END SELECT
796         ENDIF
[1705]797         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
[9019]798         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
799         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
[13472]800         WRITE(numout,*) '      langmuir & surface wave breaking under ice  nn_eice = ', nn_eice
[14072]801         SELECT CASE( nn_eice )
[13472]802         CASE( 0 )   ;   WRITE(numout,*) '   ==>>>   no impact of ice cover on langmuir & surface wave breaking'
803         CASE( 1 )   ;   WRITE(numout,*) '   ==>>>   weigthed by 1-TANH( fr_i(:,:) * 10 )'
804         CASE( 2 )   ;   WRITE(numout,*) '   ==>>>   weighted by 1-fr_i(:,:)'
805         CASE( 3 )   ;   WRITE(numout,*) '   ==>>>   weighted by 1-MIN( 1, 4 * fr_i(:,:) )'
806         CASE DEFAULT
807            CALL ctl_stop( 'zdf_tke_init: wrong value for nn_eice, should be 0,1,2, or 3')
[14072]808         END SELECT
[9019]809         WRITE(numout,*)
[9190]810         WRITE(numout,*) '   ==>>>   critical Richardson nb with your parameters  ri_cri = ', ri_cri
[9019]811         WRITE(numout,*)
[1239]812      ENDIF
[2715]813      !
[9019]814      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
815         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
816         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
[9190]817         IF(lwp) WRITE(numout,*) '   ==>>>   Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3'
[9019]818      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
819         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
[9190]820         IF(lwp) WRITE(numout,*) '   ==>>>   minimum mixing length with your parameters rmxl_min = ', rmxl_min
[9019]821      ENDIF
822      !
[2715]823      !                              ! allocate tke arrays
824      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
825      !
[1492]826      !                               !* Check of some namelist values
[4990]827      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
828      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
829      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
[5407]830      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
[9019]831      !
[2528]832      IF( ln_mxl0 ) THEN
[9169]833         IF(lwp) WRITE(numout,*)
[9190]834         IF(lwp) WRITE(numout,*) '   ==>>>   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
[2528]835         rn_mxl0 = rmxl_min
836      ENDIF
[1239]837
[14072]838      IF( nn_etau == 2  )   CALL zdf_mxl( nit000, Kmm )      ! Initialization of nmln
839
[1492]840      !                               !* depth of penetration of surface tke
[14072]841      IF( nn_etau /= 0 ) THEN
[1601]842         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
[2528]843         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
[7753]844            htau(:,:) = 10._wp
[2528]845         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
[14072]846            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )
[1492]847         END SELECT
848      ENDIF
[9019]849      !                                !* read or initialize all required files
[14072]850      CALL tke_rst( nit000, 'READ' )      ! (en, avt_k, avm_k, dissl)
[1239]851      !
[2528]852   END SUBROUTINE zdf_tke_init
[1239]853
854
[1531]855   SUBROUTINE tke_rst( kt, cdrw )
[9019]856      !!---------------------------------------------------------------------
857      !!                   ***  ROUTINE tke_rst  ***
[14072]858      !!
[9019]859      !! ** Purpose :   Read or write TKE file (en) in restart file
860      !!
861      !! ** Method  :   use of IOM library
[14072]862      !!                if the restart does not contain TKE, en is either
863      !!                set to rn_emin or recomputed
[9019]864      !!----------------------------------------------------------------------
865      USE zdf_oce , ONLY : en, avt_k, avm_k   ! ocean vertical physics
866      !!
867      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
868      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
869      !
870      INTEGER ::   jit, jk              ! dummy loop indices
871      INTEGER ::   id1, id2, id3, id4   ! local integers
872      !!----------------------------------------------------------------------
873      !
[14072]874      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
[9019]875         !                                   ! ---------------
876         IF( ln_rstart ) THEN                   !* Read the restart file
877            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
878            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
879            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
880            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
881            !
882            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
[13970]883               CALL iom_get( numror, jpdom_auto, 'en'   , en    )
884               CALL iom_get( numror, jpdom_auto, 'avt_k', avt_k )
885               CALL iom_get( numror, jpdom_auto, 'avm_k', avm_k )
886               CALL iom_get( numror, jpdom_auto, 'dissl', dissl )
[9019]887            ELSE                                          ! start TKE from rest
[9169]888               IF(lwp) WRITE(numout,*)
[9190]889               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without TKE scheme, set en to background values'
[9019]890               en   (:,:,:) = rn_emin * wmask(:,:,:)
891               dissl(:,:,:) = 1.e-12_wp
892               ! avt_k, avm_k already set to the background value in zdf_phy_init
893            ENDIF
894         ELSE                                   !* Start from rest
[9169]895            IF(lwp) WRITE(numout,*)
[9190]896            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set en to the background value'
[9019]897            en   (:,:,:) = rn_emin * wmask(:,:,:)
898            dissl(:,:,:) = 1.e-12_wp
899            ! avt_k, avm_k already set to the background value in zdf_phy_init
900         ENDIF
901         !
902      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
903         !                                   ! -------------------
[9169]904         IF(lwp) WRITE(numout,*) '---- tke_rst ----'
[13970]905         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en    )
906         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k )
907         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k )
908         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl )
[9019]909         !
910      ENDIF
911      !
[1531]912   END SUBROUTINE tke_rst
[1239]913
914   !!======================================================================
[1531]915END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.