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 @ 7813

Last change on this file since 7813 was 7813, checked in by clem, 7 years ago

synchronize trunk with 3.6

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