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

source: branches/UKMO/dev_r6501_GO6_package_trunk/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 6825

Last change on this file since 6825 was 6507, checked in by timgraham, 8 years ago

First attempt at merging in science changes from GO6 package branch at v3.6 stable (Note-namelists not yet dealt with)

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