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

source: branches/2015/dev_r5776_UKMO2_OBS_efficiency_improvs/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 6041

Last change on this file since 6041 was 6041, checked in by timgraham, 8 years ago

Merged head of trunk into branch

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