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 branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 8143

Last change on this file since 8143 was 8143, checked in by gm, 7 years ago

#1880 (HPC-09) - step-7: top/bottom drag computed at T-points, zdfbfr.F90 replaced by zdfdrg.F90 + changes in namelist

  • Property svn:keywords set to Id
File size: 42.8 KB
RevLine 
[1531]1MODULE zdftke
[1239]2   !!======================================================================
[1531]3   !!                       ***  MODULE  zdftke  ***
[1239]4   !! Ocean physics:  vertical mixing coefficient computed from the tke
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
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
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
[7990]29   !!            4.0  !  2017-04  (G. Madec)  remove CPP ddm key & avm at t-point only
[8143]30   !!             -   !  2017-05  (G. Madec)  add top/bottom friction as boundary condition (ln_drg)
[1239]31   !!----------------------------------------------------------------------
[7953]32
[1239]33   !!----------------------------------------------------------------------
[3625]34   !!   zdf_tke       : update momentum and tracer Kz from a tke scheme
35   !!   tke_tke       : tke time stepping: update tke at now time step (en)
36   !!   tke_avn       : compute mixing length scale and deduce avm and avt
37   !!   zdf_tke_init  : initialization, namelist read, and parameters control
38   !!   tke_rst       : read/write tke restart in ocean restart file
[1239]39   !!----------------------------------------------------------------------
[2528]40   USE oce            ! ocean: dynamics and active tracers variables
41   USE phycst         ! physical constants
42   USE dom_oce        ! domain: ocean
43   USE domvvl         ! domain: variable volume layer
[1492]44   USE sbc_oce        ! surface boundary condition: ocean
[2528]45   USE zdf_oce        ! vertical physics: ocean variables
[8093]46   USE zdfdrg         ! vertical physics: top/bottom drag coef.
[2528]47   USE zdfmxl         ! vertical physics: mixed layer
[8093]48#if defined key_agrif
49   USE agrif_opa_interp
50   USE agrif_opa_update
51#endif
52   !
[1492]53   USE in_out_manager ! I/O manager
54   USE iom            ! I/O manager library
[2715]55   USE lib_mpp        ! MPP library
[8093]56   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
57   USE prtctl         ! Print control
[3294]58   USE timing         ! Timing
[3625]59   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[1239]60
61   IMPLICIT NONE
62   PRIVATE
63
[2528]64   PUBLIC   zdf_tke        ! routine called in step module
65   PUBLIC   zdf_tke_init   ! routine called in opa module
66   PUBLIC   tke_rst        ! routine called in step module
[1239]67
[4147]68   !                      !!** Namelist  namzdf_tke  **
69   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
70   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
71   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
72   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
73   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
74   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
75   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
76   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
77   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
78   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
[8143]79   LOGICAL  ::   ln_drg    ! top/bottom friction forcing flag
[4147]80   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
[8143]81   INTEGER  ::      nn_htau   ! type of tke profile of penetration (=0/1)
82   REAL(wp) ::      rn_efr    ! fraction of TKE surface value which penetrates in the ocean
[4147]83   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
[8143]84   REAL(wp) ::      rn_lc     ! coef to compute vertical velocity of Langmuir cells
[1239]85
[4147]86   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values)
87   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m]
[2528]88   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3)
89   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3)
[1239]90
[8093]91   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau    ! depth of tke penetration (nn_htau)
92   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl   ! now mixing lenght of dissipation
93   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr   ! now mixing lenght of dissipation
[1492]94
[1239]95   !! * Substitutions
96#  include "vectopt_loop_substitute.h90"
97   !!----------------------------------------------------------------------
[5836]98   !! NEMO/OPA 3.7 , NEMO Consortium (2015)
[2528]99   !! $Id$
100   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[1239]101   !!----------------------------------------------------------------------
102CONTAINS
103
[2715]104   INTEGER FUNCTION zdf_tke_alloc()
105      !!----------------------------------------------------------------------
106      !!                ***  FUNCTION zdf_tke_alloc  ***
107      !!----------------------------------------------------------------------
[7991]108      ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) ,   STAT= zdf_tke_alloc )
[2715]109         !
110      IF( lk_mpp             )   CALL mpp_sum ( zdf_tke_alloc )
111      IF( zdf_tke_alloc /= 0 )   CALL ctl_warn('zdf_tke_alloc: failed to allocate arrays')
112      !
113   END FUNCTION zdf_tke_alloc
114
115
[8093]116   SUBROUTINE zdf_tke( kt, p_sh2, p_avm, p_avt )
[1239]117      !!----------------------------------------------------------------------
[1531]118      !!                   ***  ROUTINE zdf_tke  ***
[1239]119      !!
120      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
[1492]121      !!              coefficients using a turbulent closure scheme (TKE).
[1239]122      !!
[1492]123      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
124      !!              is computed from a prognostic equation :
125      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
126      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
127      !!                  + avt N^2                      ! stratif. destruc.
128      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
[1239]129      !!      with the boundary conditions:
[1695]130      !!         surface: en = max( rn_emin0, rn_ebb * taum )
[1239]131      !!         bottom : en = rn_emin
[1492]132      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
133      !!
134      !!        The now Turbulent kinetic energy is computed using the following
135      !!      time stepping: implicit for vertical diffusion term, linearized semi
136      !!      implicit for kolmogoroff dissipation term, and explicit forward for
137      !!      both buoyancy and shear production terms. Therefore a tridiagonal
138      !!      linear system is solved. Note that buoyancy and shear terms are
139      !!      discretized in a energy conserving form (Bruchard 2002).
140      !!
141      !!        The dissipative and mixing length scale are computed from en and
142      !!      the stratification (see tke_avn)
143      !!
144      !!        The now vertical eddy vicosity and diffusivity coefficients are
145      !!      given by:
146      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
147      !!              avt = max( avmb, pdl * avm                 ) 
[1239]148      !!              eav = max( avmb, avm )
[1492]149      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
150      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
[1239]151      !!
152      !! ** Action  :   compute en (now turbulent kinetic energy)
[7990]153      !!                update avt, avm (before vertical eddy coef.)
[1239]154      !!
155      !! References : Gaspar et al., JGR, 1990,
156      !!              Blanke and Delecluse, JPO, 1991
157      !!              Mellor and Blumberg, JPO 2004
158      !!              Axell, JGR, 2002
[1492]159      !!              Bruchard OM 2002
[1239]160      !!----------------------------------------------------------------------
[8093]161      INTEGER                    , INTENT(in   ) ::   kt             ! ocean time step
162      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_sh2          ! shear production term
163      REAL(wp), DIMENSION(:,:,:) , INTENT(inout) ::   p_avm, p_avt   !  momentum and tracer Kz (w-points)
[1492]164      !!----------------------------------------------------------------------
[1481]165      !
[5656]166#if defined key_agrif 
167      ! interpolation parent grid => child grid for avm_k ( ex : at west border: update column 1 and 2)
168      IF( .NOT.Agrif_Root() )   CALL Agrif_Tke
169#endif
170      !
[8093]171      CALL tke_tke( gdepw_n, e3t_n, e3w_n, p_sh2, p_avm, p_avt )   ! now tke (en)
[3632]172      !
[8093]173      CALL tke_avn( gdepw_n, e3t_n, e3w_n,        p_avm, p_avt )   ! now avt, avm, dissl
[1492]174      !
[5656]175#if defined key_agrif
176      ! Update child grid f => parent grid
177      IF( .NOT.Agrif_Root() )   CALL Agrif_Update_Tke( kt )      ! children only
178#endif     
[7990]179      !
[5656]180  END SUBROUTINE zdf_tke
[1239]181
[1492]182
[8093]183   SUBROUTINE tke_tke( pdepw, p_e3t, p_e3w, p_sh2    &
184      &                            , p_avm, p_avt )
[1239]185      !!----------------------------------------------------------------------
[1492]186      !!                   ***  ROUTINE tke_tke  ***
187      !!
188      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
189      !!
190      !! ** Method  : - TKE surface boundary condition
[2528]191      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
[7990]192      !!              - source term due to shear (= Kz dz[Ub] * dz[Un] )
[1492]193      !!              - Now TKE : resolution of the TKE equation by inverting
194      !!                a tridiagonal linear system by a "methode de chasse"
195      !!              - increase TKE due to surface and internal wave breaking
196      !!
197      !! ** Action  : - en : now turbulent kinetic energy)
[1239]198      !! ---------------------------------------------------------------------
[8093]199      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   pdepw          ! depth of w-points
200      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_e3t, p_e3w   ! level thickness (t- & w-points)
201      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_sh2          ! shear production term
202      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
203      !
204      INTEGER ::   ji, jj, jk              ! dummy loop arguments
[8143]205      REAL(wp) ::   zetop, zebot, zmsku, zmskv ! local scalars
206      REAL(wp) ::   zrhoa  = 1.22              ! Air density kg/m3
207      REAL(wp) ::   zcdrag = 1.5e-3            ! drag coefficient
208      REAL(wp) ::   zbbrau, zri                ! local scalars
209      REAL(wp) ::   zfact1, zfact2, zfact3     !   -         -
210      REAL(wp) ::   ztx2  , zty2  , zcof       !   -         -
211      REAL(wp) ::   ztau  , zdif               !   -         -
212      REAL(wp) ::   zus   , zwlc  , zind       !   -         -
213      REAL(wp) ::   zzd_up, zzd_lw             !   -         -
[8093]214      INTEGER , DIMENSION(jpi,jpj)     ::   imlc
215      REAL(wp), DIMENSION(jpi,jpj)     ::   zhlc
216      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
[1239]217      !!--------------------------------------------------------------------
[1492]218      !
[3294]219      IF( nn_timing == 1 )  CALL timing_start('tke_tke')
220      !
[1695]221      zbbrau = rn_ebb / rau0       ! Local constant initialisation
[2528]222      zfact1 = -.5_wp * rdt 
223      zfact2 = 1.5_wp * rdt * rn_ediss
224      zfact3 = 0.5_wp       * rn_ediss
[1492]225      !
[5120]226      !
[1492]227      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[8143]228      !                     !  Surface/top/bottom boundary condition on tke
[1492]229      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[8143]230     
231      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
232         DO ji = fs_2, fs_jpim1   ! vector opt.
233            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
234         END DO
235      END DO
[5120]236      IF ( ln_isfcav ) THEN
237         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin
238            DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]239               en(ji,jj,mikt(ji,jj)) = rn_emin * tmask(ji,jj,1)
[5120]240            END DO
241         END DO
[8143]242      ENDIF
[2528]243     
[1492]244      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
245      !                     !  Bottom boundary condition on tke
246      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[1719]247      !
[8143]248      !   en(bot)   = (ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
249      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
250      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
[1492]251      !
[8143]252      IF( ln_drg ) THEN       !== friction used as top/bottom boundary condition on TKE
253         !
254         DO jj = 2, jpjm1           ! bottom friction
255            DO ji = fs_2, fs_jpim1     ! vector opt.
256               zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
257               zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
258               !                       ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
259               zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mbkt(ji,jj))+ub(ji-1,jj,mbkt(ji,jj)) ) )**2  &
260                  &                                           + ( zmskv*( vb(ji,jj,mbkt(ji,jj))+vb(ji,jj-1,mbkt(ji,jj)) ) )**2  )
261               en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
262            END DO
263         END DO
264         IF( ln_isfcav ) THEN       ! top friction
265            DO jj = 2, jpjm1
266               DO ji = fs_2, fs_jpim1   ! vector opt.
267                  zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
268                  zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
269                  !                             ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
270                  zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mikt(ji,jj))+ub(ji-1,jj,mikt(ji,jj)) ) )**2  &
271                     &                                           + ( zmskv*( vb(ji,jj,mikt(ji,jj))+vb(ji,jj-1,mikt(ji,jj)) ) )**2  )
272                  en(ji,jj,mikt(ji,jj)) = MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1))   ! masked at ocean surface
273               END DO
274            END DO
275         ENDIF
276         !
277      ENDIF
278      !
[1492]279      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[8093]280      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke   !   (Axell JGR 2002)
[1492]281         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
[1239]282         !
[1492]283         !                        !* total energy produce by LC : cumulative sum over jk
[8093]284         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * pdepw(:,:,1) * p_e3w(:,:,1)
[1239]285         DO jk = 2, jpk
[8093]286            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * pdepw(:,:,jk) * p_e3w(:,:,jk)
[1239]287         END DO
[1492]288         !                        !* finite Langmuir Circulation depth
[1705]289         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
[7753]290         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
[1239]291         DO jk = jpkm1, 2, -1
[1492]292            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
293               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
[1705]294                  zus  = zcof * taum(ji,jj)
[1239]295                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
296               END DO
297            END DO
298         END DO
[1492]299         !                               ! finite LC depth
300         DO jj = 1, jpj 
[1239]301            DO ji = 1, jpi
[8093]302               zhlc(ji,jj) = pdepw(ji,jj,imlc(ji,jj))
[1239]303            END DO
304         END DO
[1705]305         zcof = 0.016 / SQRT( zrhoa * zcdrag )
[1492]306         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
[1239]307            DO jj = 2, jpjm1
308               DO ji = fs_2, fs_jpim1   ! vector opt.
[1705]309                  zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
[1492]310                  !                                           ! vertical velocity due to LC
[8093]311                  zind = 0.5 - SIGN( 0.5, pdepw(ji,jj,jk) - zhlc(ji,jj) )
312                  zwlc = zind * rn_lc * zus * SIN( rpi * pdepw(ji,jj,jk) / zhlc(ji,jj) )
[1492]313                  !                                           ! TKE Langmuir circulation source term
[7813]314                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * MAX(0.,1._wp - 2.*fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc )   &
[6497]315                     &                              / zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1239]316               END DO
317            END DO
318         END DO
319         !
320      ENDIF
[1492]321      !
322      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
323      !                     !  Now Turbulent kinetic energy (output in en)
324      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
325      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
326      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
327      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
328      !
[7991]329      IF( nn_pdl == 1 ) THEN      !* Prandtl number = F( Ri )
[5656]330         DO jk = 2, jpkm1
331            DO jj = 2, jpjm1
[7991]332               DO ji = 2, jpim1
333                  !                             ! local Richardson number
[8093]334                  zri = MAX( rn2b(ji,jj,jk), 0._wp ) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
[7991]335                  !                             ! inverse of Prandtl number
[5656]336                  apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
337               END DO
338            END DO
339         END DO
340      ENDIF
[5836]341      !         
[5120]342      DO jk = 2, jpkm1           !* Matrix and right hand side in en
343         DO jj = 2, jpjm1
344            DO ji = fs_2, fs_jpim1   ! vector opt.
[1492]345               zcof   = zfact1 * tmask(ji,jj,jk)
[7953]346               !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
347               !                                   ! eddy coefficient (ensure numerical stability)
[8093]348               zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
349                  &          /    (  p_e3t(ji,jj,jk  ) * p_e3w(ji,jj,jk  )  )
350               zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
351                  &          /    (  p_e3t(ji,jj,jk-1) * p_e3w(ji,jj,jk  )  )
[7953]352               !
[1492]353               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
354               zd_lw(ji,jj,jk) = zzd_lw
[7991]355               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
[1239]356               !
[1492]357               !                                   ! right hand side in en
[8093]358               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  p_sh2(ji,jj,jk)                          &   ! shear
359                  &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
[7991]360                  &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
361                  &                                ) * wmask(ji,jj,jk)
[1239]362            END DO
[5120]363         END DO
364      END DO
365      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
366      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
367         DO jj = 2, jpjm1
368            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]369               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
[1239]370            END DO
[5120]371         END DO
372      END DO
[5836]373      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
[5120]374         DO ji = fs_2, fs_jpim1   ! vector opt.
375            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
376         END DO
377      END DO
378      DO jk = 3, jpkm1
379         DO jj = 2, jpjm1
380            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]381               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)
[1239]382            END DO
[5120]383         END DO
384      END DO
[5836]385      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
[5120]386         DO ji = fs_2, fs_jpim1   ! vector opt.
[1492]387            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
[5120]388         END DO
389      END DO
390      DO jk = jpk-2, 2, -1
391         DO jj = 2, jpjm1
392            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]393               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
[1239]394            END DO
[5120]395         END DO
396      END DO
397      DO jk = 2, jpkm1                             ! set the minimum value of tke
398         DO jj = 2, jpjm1
399            DO ji = fs_2, fs_jpim1   ! vector opt.
400               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
[1239]401            END DO
402         END DO
403      END DO
404
[1492]405      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
406      !                            !  TKE due to surface and internal wave breaking
407      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[6140]408!!gm BUG : in the exp  remove the depth of ssh !!!
[8143]409!!gm       i.e. use gde3w in argument (pdepw)
[6140]410     
411     
[2528]412      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
[1492]413         DO jk = 2, jpkm1
[1239]414            DO jj = 2, jpjm1
415               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]416                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[7813]417                     &                                 * MAX(0.,1._wp - 2.*fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1239]418               END DO
419            END DO
[1492]420         END DO
[2528]421      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
[1492]422         DO jj = 2, jpjm1
423            DO ji = fs_2, fs_jpim1   ! vector opt.
424               jk = nmln(ji,jj)
[8093]425               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[7813]426                  &                                 * MAX(0.,1._wp - 2.*fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1239]427            END DO
[1492]428         END DO
[2528]429      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
[1705]430         DO jk = 2, jpkm1
431            DO jj = 2, jpjm1
432               DO ji = fs_2, fs_jpim1   ! vector opt.
433                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
434                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
[4990]435                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
[2528]436                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
437                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
[8093]438                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[7813]439                     &                        * MAX(0.,1._wp - 2.*fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1705]440               END DO
441            END DO
442         END DO
[1239]443      ENDIF
[1492]444      !
[3294]445      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
446      !
[1239]447   END SUBROUTINE tke_tke
448
[1492]449
[8093]450   SUBROUTINE tke_avn( pdepw, p_e3t, p_e3w, p_avm, p_avt )
[1239]451      !!----------------------------------------------------------------------
[1492]452      !!                   ***  ROUTINE tke_avn  ***
[1239]453      !!
[1492]454      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
455      !!
456      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
457      !!              the tke_tke routine). First, the now mixing lenth is
458      !!      computed from en and the strafification (N^2), then the mixings
459      !!      coefficients are computed.
460      !!              - Mixing length : a first evaluation of the mixing lengh
461      !!      scales is:
462      !!                      mxl = sqrt(2*en) / N 
463      !!      where N is the brunt-vaisala frequency, with a minimum value set
[2528]464      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
[1492]465      !!        The mixing and dissipative length scale are bound as follow :
466      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
467      !!                        zmxld = zmxlm = mxl
468      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
469      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
470      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
471      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
472      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
473      !!                    the surface to obtain ldown. the resulting length
474      !!                    scales are:
475      !!                        zmxld = sqrt( lup * ldown )
476      !!                        zmxlm = min ( lup , ldown )
477      !!              - Vertical eddy viscosity and diffusivity:
478      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
479      !!                      avt = max( avmb, pdlr * avm ) 
480      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
481      !!
[7990]482      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
[1239]483      !!----------------------------------------------------------------------
[8093]484      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdepw          ! depth (w-points)
485      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_e3t, p_e3w   ! level thickness (t- & w-points)
486      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
487      !
[2715]488      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[7990]489      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
490      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
491      REAL(wp) ::   zemxl, zemlm, zemlp       !   -      -
[8093]492      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmxlm, zmxld
[1239]493      !!--------------------------------------------------------------------
[3294]494      !
495      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
[1239]496
[1492]497      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
498      !                     !  Mixing length
499      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
500      !
501      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
502      !
[5120]503      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
[7753]504      zmxlm(:,:,:)  = rmxl_min   
505      zmxld(:,:,:)  = rmxl_min
[5120]506      !
[2528]507      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
[8093]508         zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
[4990]509         DO jj = 2, jpjm1
510            DO ji = fs_2, fs_jpim1
[5120]511               zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
[4990]512            END DO
513         END DO
514      ELSE
[7753]515         zmxlm(:,:,1) = rn_mxl0
[1239]516      ENDIF
517      !
[5120]518      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
519         DO jj = 2, jpjm1
520            DO ji = fs_2, fs_jpim1   ! vector opt.
[1239]521               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
[5836]522               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
[1239]523            END DO
524         END DO
525      END DO
[1492]526      !
527      !                     !* Physical limits for the mixing length
528      !
[7753]529      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
530      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
[1492]531      !
[1239]532      SELECT CASE ( nn_mxl )
533      !
[5836]534 !!gm Not sure of that coding for ISF....
[8093]535      ! where wmask = 0 set zmxlm == p_e3w
[1239]536      CASE ( 0 )           ! bounded by the distance to surface and bottom
[5120]537         DO jk = 2, jpkm1
538            DO jj = 2, jpjm1
539               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]540                  zemxl = MIN( pdepw(ji,jj,jk) - pdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
541                  &            pdepw(ji,jj,mbkt(ji,jj)+1) - pdepw(ji,jj,jk) )
[5120]542                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
[8093]543                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
544                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
[1239]545               END DO
546            END DO
547         END DO
548         !
549      CASE ( 1 )           ! bounded by the vertical scale factor
[5120]550         DO jk = 2, jpkm1
551            DO jj = 2, jpjm1
552               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]553                  zemxl = MIN( p_e3w(ji,jj,jk), zmxlm(ji,jj,jk) )
[1239]554                  zmxlm(ji,jj,jk) = zemxl
555                  zmxld(ji,jj,jk) = zemxl
556               END DO
557            END DO
558         END DO
559         !
560      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
[5120]561         DO jk = 2, jpkm1         ! from the surface to the bottom :
562            DO jj = 2, jpjm1
563               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]564                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
[1239]565               END DO
[5120]566            END DO
567         END DO
568         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
569            DO jj = 2, jpjm1
570               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]571                  zemxl = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
[1239]572                  zmxlm(ji,jj,jk) = zemxl
573                  zmxld(ji,jj,jk) = zemxl
574               END DO
575            END DO
576         END DO
577         !
578      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
[5120]579         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
580            DO jj = 2, jpjm1
581               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]582                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
[1239]583               END DO
[5120]584            END DO
585         END DO
586         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
587            DO jj = 2, jpjm1
588               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]589                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
[1239]590               END DO
591            END DO
592         END DO
593         DO jk = 2, jpkm1
594            DO jj = 2, jpjm1
595               DO ji = fs_2, fs_jpim1   ! vector opt.
596                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
597                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
598                  zmxlm(ji,jj,jk) = zemlm
599                  zmxld(ji,jj,jk) = zemlp
600               END DO
601            END DO
602         END DO
603         !
604      END SELECT
[1492]605      !
[1239]606
[1492]607      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[7990]608      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
[1492]609      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
610      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
[1239]611         DO jj = 2, jpjm1
612            DO ji = fs_2, fs_jpim1   ! vector opt.
613               zsqen = SQRT( en(ji,jj,jk) )
614               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
[8093]615               p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
616               p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
[1239]617               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
618            END DO
619         END DO
620      END DO
[1492]621      !
622      !
623      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
[5120]624         DO jk = 2, jpkm1
625            DO jj = 2, jpjm1
626               DO ji = fs_2, fs_jpim1   ! vector opt.
[8093]627                  p_avt(ji,jj,jk)   = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
[1239]628              END DO
629            END DO
630         END DO
631      ENDIF
632
633      IF(ln_ctl) THEN
[7990]634         CALL prt_ctl( tab3d_1=en , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
635         CALL prt_ctl( tab3d_1=avm, clinfo1=' tke  - m: ', ovlap=1, kdim=jpk )
[1239]636      ENDIF
637      !
[3294]638      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
639      !
[1492]640   END SUBROUTINE tke_avn
[1239]641
[1492]642
[2528]643   SUBROUTINE zdf_tke_init
[1239]644      !!----------------------------------------------------------------------
[2528]645      !!                  ***  ROUTINE zdf_tke_init  ***
[1239]646      !!                     
647      !! ** Purpose :   Initialization of the vertical eddy diffivity and
[1492]648      !!              viscosity when using a tke turbulent closure scheme
[1239]649      !!
[1601]650      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
[1492]651      !!              called at the first timestep (nit000)
[1239]652      !!
[1601]653      !! ** input   :   Namlist namzdf_tke
[1239]654      !!
655      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
656      !!----------------------------------------------------------------------
657      INTEGER ::   ji, jj, jk   ! dummy loop indices
[4147]658      INTEGER ::   ios
[1239]659      !!
[2528]660      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
661         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
[8143]662         &                 rn_mxl0 , nn_pdl   , ln_drg , ln_lc    , rn_lc    ,   &
[2528]663         &                 nn_etau , nn_htau  , rn_efr   
[1239]664      !!----------------------------------------------------------------------
[2715]665      !
[4147]666      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
667      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
668901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
669
670      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
671      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
672902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
[4624]673      IF(lwm) WRITE ( numond, namzdf_tke )
[2715]674      !
[2528]675      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
[2715]676      !
[1492]677      IF(lwp) THEN                    !* Control print
[1239]678         WRITE(numout,*)
[2528]679         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
680         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]681         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
[1705]682         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
683         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
684         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
685         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
686         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
[8143]687         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
[1705]688         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
689         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
[8143]690         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
691         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
692         WRITE(numout,*) '      top/bottom friction forcing flag            ln_drg    = ', ln_drg
693         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
694         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
[1705]695         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
[8143]696         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
697         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
[1239]698         WRITE(numout,*)
[8143]699         IF( ln_drg ) THEN
700            WRITE(numout,*) '   Namelist namdrg_top/_bot:   used values:'
701            WRITE(numout,*) '      top    ocean cavity roughness (m)          rn_z0(_top)= ', r_z0_top
702            WRITE(numout,*) '      Bottom seafloor     roughness (m)          rn_z0(_bot)= ', r_z0_bot
703         ENDIF
[7953]704         WRITE(numout,*)
[8143]705         WRITE(numout,*)
706         WRITE(numout,*) '   ==>> critical Richardson nb with your parameters  ri_cri = ', ri_cri
707         WRITE(numout,*)
[1239]708      ENDIF
[2715]709      !
[7990]710      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
711         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
712         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
[7953]713         IF(lwp) WRITE(numout,*) '      Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3 '
[7990]714      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
[7953]715         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
716         IF(lwp) WRITE(numout,*) '      minimum mixing length with your parameters rmxl_min = ', rmxl_min
717      ENDIF
718      !
[2715]719      !                              ! allocate tke arrays
720      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
721      !
[1492]722      !                               !* Check of some namelist values
[4990]723      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
724      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
725      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
[5407]726      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
[1239]727
[2528]728      IF( ln_mxl0 ) THEN
729         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
730         rn_mxl0 = rmxl_min
731      ENDIF
732     
[1492]733      IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
[1239]734
[1492]735      !                               !* depth of penetration of surface tke
736      IF( nn_etau /= 0 ) THEN     
[1601]737         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
[2528]738         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
[7753]739            htau(:,:) = 10._wp
[2528]740         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
[7753]741            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
[1492]742         END SELECT
743      ENDIF
744      !                               !* set vertical eddy coef. to the background value
[1239]745      DO jk = 1, jpk
[8093]746         avt(:,:,jk) = avtb(jk) * wmask(:,:,jk)
747         avm(:,:,jk) = avmb(jk) * wmask(:,:,jk)
[1239]748      END DO
[7753]749      dissl(:,:,:) = 1.e-12_wp
[2715]750      !                             
751      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
[1239]752      !
[2528]753   END SUBROUTINE zdf_tke_init
[1239]754
755
[1531]756   SUBROUTINE tke_rst( kt, cdrw )
[8093]757      !!---------------------------------------------------------------------
758      !!                   ***  ROUTINE tke_rst  ***
759      !!                     
760      !! ** Purpose :   Read or write TKE file (en) in restart file
761      !!
762      !! ** Method  :   use of IOM library
763      !!                if the restart does not contain TKE, en is either
764      !!                set to rn_emin or recomputed
765      !!----------------------------------------------------------------------
766      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
767      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
768      !
769      INTEGER ::   jit, jk              ! dummy loop indices
770      INTEGER ::   id1, id2, id3, id4   ! local integers
771      !!----------------------------------------------------------------------
772      !
773      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
774         !                                   ! ---------------
775         IF( ln_rstart ) THEN                   !* Read the restart file
776            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
777            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
778            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
779            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
780            !
781            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
782               CALL iom_get( numror, jpdom_autoglo, 'en', en )
783               CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k )
784               CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k )
785               CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
786            ELSE                                          ! start TKE from rest
787               IF(lwp) WRITE(numout,*) '   ==>>   previous run without TKE scheme, set en to background values'
788               en(:,:,:) = rn_emin * wmask(:,:,:)
789               ! avt_k, avm_k already set to the background value in zdf_phy_init
790            ENDIF
791         ELSE                                   !* Start from rest
792            IF(lwp) WRITE(numout,*) '   ==>>   start from rest: set en to the background value'
793            en(:,:,:) = rn_emin * wmask(:,:,:)
794            ! avt_k, avm_k already set to the background value in zdf_phy_init
795         ENDIF
796         !
797      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
798         !                                   ! -------------------
799         IF(lwp) WRITE(numout,*) '---- tke-rst ----'
800         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en    )
801         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k )
802         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k )
803         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl )
804         !
805      ENDIF
806      !
[1531]807   END SUBROUTINE tke_rst
[1239]808
809   !!======================================================================
[1531]810END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.