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

source: branches/UKMO/dev_r5518_GO6_under_ice_relax_dr_hook/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 11738

Last change on this file since 11738 was 11738, checked in by marc, 5 years ago

The Dr Hook changes from my perl code.

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