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

source: branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 8882

Last change on this file since 8882 was 8882, checked in by flavoni, 7 years ago

dev_CNRS_2017 branch: merged dev_r7881_ENHANCE09_RK3 with trunk r8864

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