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

source: branches/2015/nemo_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 9294

Last change on this file since 9294 was 8290, checked in by emanuelaclementi, 7 years ago

#1733 bug fixed for Initialization of av*_k in tke

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