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

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

source: branches/2015/dev_r5918_nemo_v3_6_STABLE_XIOS2/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 6156

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

Merged in dev_agrif_v3_6_STABLE

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