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

source: trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 5798

Last change on this file since 5798 was 5798, checked in by timgraham, 9 years ago

Minor change so that results are preserved if nn_pdl=0 after changes from AGRIF branch

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