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

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 5845

Last change on this file since 5845 was 5845, checked in by gm, 8 years ago

#1613: vvl by default: suppression of domzgr_substitute.h90

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