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_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 8863

Last change on this file since 8863 was 8863, checked in by flavoni, 6 years ago

(ENHANCE-09): fix AGRIF reproducibility problem

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