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

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

source: branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 5782

Last change on this file since 5782 was 5782, checked in by cetlod, 9 years ago

Final step of improvements/simplifications of ADV & LDF momentum trends. The branch is now phased with the trunk at revision 5721 and is ready to be used

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