New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdftke.F90 in branches/2017/dev_r8600_xios_read_write/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r8600_xios_read_write/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 8793

Last change on this file since 8793 was 8793, checked in by andmirek, 6 years ago

#1953 and #1962 change lxios_read to lrxios to be consistent with write branch

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