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

source: branches/UKMO/dev_r5518_obs_oper_update/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 11200

Last change on this file since 11200 was 9188, checked in by kingr, 6 years ago

Small change to zdftke/zdfgls to avoid merge conflicts with AMM15 branch.

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