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

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

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 7910

Last change on this file since 7910 was 7910, checked in by timgraham, 7 years ago

All wrk_alloc removed

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