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

source: branches/UKMO/r6232_tracer_advection/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 9295

Last change on this file since 9295 was 9295, checked in by jcastill, 6 years ago

Remove svn keywords

File size: 46.0 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 * ( 1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) /   &
326                     &   zhlc(ji,jj) * wmask(ji,jj,jk) * tmask(ji,jj,1)
327               END DO
328            END DO
329         END DO
330         !
331      ENDIF
332      !
333      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
334      !                     !  Now Turbulent kinetic energy (output in en)
335      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
336      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
337      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
338      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
339      !
340      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form)
341         DO jj = 1, jpj                 ! here avmu, avmv used as workspace
342            DO ji = 1, jpi
343               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
344                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   & 
345                  &                            / (  fse3uw_n(ji,jj,jk)               &
346                  &                              *  fse3uw_b(ji,jj,jk)  )
347               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
348                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
349                  &                            / (  fse3vw_n(ji,jj,jk)               &
350                  &                              *  fse3vw_b(ji,jj,jk)  )
351            END DO
352         END DO
353      END DO
354      !
355      DO jk = 2, jpkm1           !* Matrix and right hand side in en
356         DO jj = 2, jpjm1
357            DO ji = fs_2, fs_jpim1   ! vector opt.
358               zcof   = zfact1 * tmask(ji,jj,jk)
359               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal
360                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
361               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal
362                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
363                  !                                                           ! shear prod. at w-point weightened by mask
364               zesh2  =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
365                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )   
366                  !
367               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
368               zd_lw(ji,jj,jk) = zzd_lw
369               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk)
370               !
371               !                                   ! right hand side in en
372               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  zesh2  -   avt(ji,jj,jk) * rn2(ji,jj,jk)    &
373                  &                                 + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk)  ) &
374                  &                                 * wmask(ji,jj,jk)
375            END DO
376         END DO
377      END DO
378      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
379      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
380         DO jj = 2, jpjm1
381            DO ji = fs_2, fs_jpim1    ! vector opt.
382               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
383            END DO
384         END DO
385      END DO
386      !
387      ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
388      DO jj = 2, jpjm1
389         DO ji = fs_2, fs_jpim1   ! vector opt.
390            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
391         END DO
392      END DO
393      DO jk = 3, jpkm1
394         DO jj = 2, jpjm1
395            DO ji = fs_2, fs_jpim1    ! vector opt.
396               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)
397            END DO
398         END DO
399      END DO
400      !
401      ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
402      DO jj = 2, jpjm1
403         DO ji = fs_2, fs_jpim1   ! vector opt.
404            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
405         END DO
406      END DO
407      DO jk = jpk-2, 2, -1
408         DO jj = 2, jpjm1
409            DO ji = fs_2, fs_jpim1    ! vector opt.
410               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
411            END DO
412         END DO
413      END DO
414      DO jk = 2, jpkm1                             ! set the minimum value of tke
415         DO jj = 2, jpjm1
416            DO ji = fs_2, fs_jpim1   ! vector opt.
417               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
418            END DO
419         END DO
420      END DO
421
422      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
423      !                            !  TKE due to surface and internal wave breaking
424      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
425      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
426         DO jk = 2, jpkm1
427            DO jj = 2, jpjm1
428               DO ji = fs_2, fs_jpim1   ! vector opt.
429                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
430                     &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
431               END DO
432            END DO
433         END DO
434      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
435         DO jj = 2, jpjm1
436            DO ji = fs_2, fs_jpim1   ! vector opt.
437               jk = nmln(ji,jj)
438               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
439                  &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
440            END DO
441         END DO
442      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
443!CDIR NOVERRCHK
444         DO jk = 2, jpkm1
445!CDIR NOVERRCHK
446            DO jj = 2, jpjm1
447!CDIR NOVERRCHK
448               DO ji = fs_2, fs_jpim1   ! vector opt.
449                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
450                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
451                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
452                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
453                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
454                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
455                     &                        * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
456               END DO
457            END DO
458         END DO
459      ENDIF
460      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
461      !
462      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer
463      CALL wrk_dealloc( jpi,jpj, zhlc ) 
464      CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
465      !
466      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
467      !
468   END SUBROUTINE tke_tke
469
470
471   SUBROUTINE tke_avn
472      !!----------------------------------------------------------------------
473      !!                   ***  ROUTINE tke_avn  ***
474      !!
475      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
476      !!
477      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
478      !!              the tke_tke routine). First, the now mixing lenth is
479      !!      computed from en and the strafification (N^2), then the mixings
480      !!      coefficients are computed.
481      !!              - Mixing length : a first evaluation of the mixing lengh
482      !!      scales is:
483      !!                      mxl = sqrt(2*en) / N 
484      !!      where N is the brunt-vaisala frequency, with a minimum value set
485      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
486      !!        The mixing and dissipative length scale are bound as follow :
487      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
488      !!                        zmxld = zmxlm = mxl
489      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
490      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
491      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
492      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
493      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
494      !!                    the surface to obtain ldown. the resulting length
495      !!                    scales are:
496      !!                        zmxld = sqrt( lup * ldown )
497      !!                        zmxlm = min ( lup , ldown )
498      !!              - Vertical eddy viscosity and diffusivity:
499      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
500      !!                      avt = max( avmb, pdlr * avm ) 
501      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
502      !!
503      !! ** Action  : - avt : now vertical eddy diffusivity (w-point)
504      !!              - avmu, avmv : now vertical eddy viscosity at uw- and vw-points
505      !!----------------------------------------------------------------------
506      INTEGER  ::   ji, jj, jk   ! dummy loop indices
507      REAL(wp) ::   zrn2, zraug, zcoef, zav     ! local scalars
508      REAL(wp) ::   zdku, zpdlr, zri, zsqen     !   -      -
509      REAL(wp) ::   zdkv, zemxl, zemlm, zemlp   !   -      -
510      REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld
511      !!--------------------------------------------------------------------
512      !
513      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
514
515      CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
516
517      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
518      !                     !  Mixing length
519      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
520      !
521      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
522      !
523      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
524      zmxlm(:,:,:)  = rmxl_min   
525      zmxld(:,:,:)  = rmxl_min
526      !
527      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
528         DO jj = 2, jpjm1
529            DO ji = fs_2, fs_jpim1
530               zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
531               zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
532            END DO
533         END DO
534      ELSE
535         zmxlm(:,:,1) = rn_mxl0
536      ENDIF
537      !
538!CDIR NOVERRCHK
539      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
540!CDIR NOVERRCHK
541         DO jj = 2, jpjm1
542!CDIR NOVERRCHK
543            DO ji = fs_2, fs_jpim1   ! vector opt.
544               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
545               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) )
546            END DO
547         END DO
548      END DO
549      !
550      !                     !* Physical limits for the mixing length
551      !
552      zmxld(:,:,1  ) = zmxlm(:,:,1)   ! surface set to the minimum value
553      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
554      !
555      SELECT CASE ( nn_mxl )
556      !
557      ! where wmask = 0 set zmxlm == fse3w
558      CASE ( 0 )           ! bounded by the distance to surface and bottom
559         DO jk = 2, jpkm1
560            DO jj = 2, jpjm1
561               DO ji = fs_2, fs_jpim1   ! vector opt.
562                  zemxl = MIN( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
563                  &            fsdepw(ji,jj,mbkt(ji,jj)+1) - fsdepw(ji,jj,jk) )
564                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
565                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
566                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
567               END DO
568            END DO
569         END DO
570         !
571      CASE ( 1 )           ! bounded by the vertical scale factor
572         DO jk = 2, jpkm1
573            DO jj = 2, jpjm1
574               DO ji = fs_2, fs_jpim1   ! vector opt.
575                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
576                  zmxlm(ji,jj,jk) = zemxl
577                  zmxld(ji,jj,jk) = zemxl
578               END DO
579            END DO
580         END DO
581         !
582      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
583         DO jk = 2, jpkm1         ! from the surface to the bottom :
584            DO jj = 2, jpjm1
585               DO ji = fs_2, fs_jpim1   ! vector opt.
586                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
587               END DO
588            END DO
589         END DO
590         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
591            DO jj = 2, jpjm1
592               DO ji = fs_2, fs_jpim1   ! vector opt.
593                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
594                  zmxlm(ji,jj,jk) = zemxl
595                  zmxld(ji,jj,jk) = zemxl
596               END DO
597            END DO
598         END DO
599         !
600      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
601         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
602            DO jj = 2, jpjm1
603               DO ji = fs_2, fs_jpim1   ! vector opt.
604                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
605               END DO
606            END DO
607         END DO
608         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
609            DO jj = 2, jpjm1
610               DO ji = fs_2, fs_jpim1   ! vector opt.
611                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
612               END DO
613            END DO
614         END DO
615!CDIR NOVERRCHK
616         DO jk = 2, jpkm1
617!CDIR NOVERRCHK
618            DO jj = 2, jpjm1
619!CDIR NOVERRCHK
620               DO ji = fs_2, fs_jpim1   ! vector opt.
621                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
622                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
623                  zmxlm(ji,jj,jk) = zemlm
624                  zmxld(ji,jj,jk) = zemlp
625               END DO
626            END DO
627         END DO
628         !
629      END SELECT
630      !
631# if defined key_c1d
632      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
633      e_mix(:,:,:) = zmxlm(:,:,:)
634# endif
635
636      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
637      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
638      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
639!CDIR NOVERRCHK
640      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
641!CDIR NOVERRCHK
642         DO jj = 2, jpjm1
643!CDIR NOVERRCHK
644            DO ji = fs_2, fs_jpim1   ! vector opt.
645               zsqen = SQRT( en(ji,jj,jk) )
646               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
647               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
648               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
649               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
650            END DO
651         END DO
652      END DO
653      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
654      !
655      DO jk = 2, jpkm1            !* vertical eddy viscosity at wu- and wv-points
656         DO jj = 2, jpjm1
657            DO ji = fs_2, fs_jpim1   ! vector opt.
658               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * wumask(ji,jj,jk)
659               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * wvmask(ji,jj,jk)
660            END DO
661         END DO
662      END DO
663      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
664      !
665      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
666         DO jk = 2, jpkm1
667            DO jj = 2, jpjm1
668               DO ji = fs_2, fs_jpim1   ! vector opt.
669                  zcoef = avm(ji,jj,jk) * 2._wp * fse3w(ji,jj,jk) * fse3w(ji,jj,jk)
670                  !                                          ! shear
671                  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) )   &
672                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
673                  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) )   &
674                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
675                  !                                          ! local Richardson number
676                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) * zcoef / (zdku + zdkv + rn_bshear )
677                  zpdlr = MAX(  0.1_wp,  0.2 / MAX( 0.2 , zri )  )
678!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
679!!gm              zpdlr = MAX(  0.1_wp,  ri_crit / MAX( ri_crit , zri )  )
680                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
681# if defined key_c1d
682                  e_pdl(ji,jj,jk) = zpdlr * wmask(ji,jj,jk)  ! c1d configuration : save masked Prandlt number
683                  e_ric(ji,jj,jk) = zri   * wmask(ji,jj,jk)  ! c1d config. : save Ri
684# endif
685              END DO
686            END DO
687         END DO
688      ENDIF
689      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
690
691      IF(ln_ctl) THEN
692         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
693         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
694            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
695      ENDIF
696      !
697      CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
698      !
699      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
700      !
701   END SUBROUTINE tke_avn
702
703
704   SUBROUTINE zdf_tke_init
705      !!----------------------------------------------------------------------
706      !!                  ***  ROUTINE zdf_tke_init  ***
707      !!                     
708      !! ** Purpose :   Initialization of the vertical eddy diffivity and
709      !!              viscosity when using a tke turbulent closure scheme
710      !!
711      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
712      !!              called at the first timestep (nit000)
713      !!
714      !! ** input   :   Namlist namzdf_tke
715      !!
716      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
717      !!----------------------------------------------------------------------
718      INTEGER ::   ji, jj, jk   ! dummy loop indices
719      INTEGER ::   ios, ierr
720      !!
721      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
722         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
723         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   &
724         &                 nn_etau , nn_htau  , rn_efr   
725      !!----------------------------------------------------------------------
726      !
727      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
728      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
729901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
730
731      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
732      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
733902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
734      IF(lwm) WRITE ( numond, namzdf_tke )
735      !
736      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
737      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
738      !
739      IF(lwp) THEN                    !* Control print
740         WRITE(numout,*)
741         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
742         WRITE(numout,*) '~~~~~~~~~~~~'
743         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
744         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
745         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
746         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
747         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
748         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
749         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
750         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
751         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
752         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
753         WRITE(numout,*) '      surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
754         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
755         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
756         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
757         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
758         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
759         WRITE(numout,*)
760         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
761      ENDIF
762      !
763      !                              ! allocate tke arrays
764      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
765      !
766      !                               !* Check of some namelist values
767      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
768      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
769      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
770      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
771
772      IF( ln_mxl0 ) THEN
773         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
774         rn_mxl0 = rmxl_min
775      ENDIF
776     
777      IF( nn_etau == 2  ) THEN
778          ierr = zdf_mxl_alloc()
779          nmln(:,:) = nlb10           ! Initialization of nmln
780      ENDIF
781
782      !                               !* depth of penetration of surface tke
783      IF( nn_etau /= 0 ) THEN     
784         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
785         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
786            htau(:,:) = 10._wp
787         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
788            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
789         END SELECT
790      ENDIF
791      !                               !* set vertical eddy coef. to the background value
792      DO jk = 1, jpk
793         avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
794         avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
795         avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
796         avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
797      END DO
798      dissl(:,:,:) = 1.e-12_wp
799      !                             
800      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
801      !
802   END SUBROUTINE zdf_tke_init
803
804
805   SUBROUTINE tke_rst( kt, cdrw )
806     !!---------------------------------------------------------------------
807     !!                   ***  ROUTINE tke_rst  ***
808     !!                     
809     !! ** Purpose :   Read or write TKE file (en) in restart file
810     !!
811     !! ** Method  :   use of IOM library
812     !!                if the restart does not contain TKE, en is either
813     !!                set to rn_emin or recomputed
814     !!----------------------------------------------------------------------
815     INTEGER         , INTENT(in) ::   kt     ! ocean time-step
816     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
817     !
818     INTEGER ::   jit, jk   ! dummy loop indices
819     INTEGER ::   id1, id2, id3, id4, id5, id6   ! local integers
820     !!----------------------------------------------------------------------
821     !
822     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
823        !                                   ! ---------------
824        IF( ln_rstart ) THEN                   !* Read the restart file
825           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
826           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
827           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
828           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
829           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
830           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
831           !
832           IF( id1 > 0 ) THEN                       ! 'en' exists
833              CALL iom_get( numror, jpdom_autoglo, 'en', en )
834              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
835                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt   )
836                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm   )
837                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu  )
838                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv  )
839                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
840              ELSE                                                 ! one at least array is missing
841                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
842              ENDIF
843           ELSE                                     ! No TKE array found: initialisation
844              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
845              en (:,:,:) = rn_emin * tmask(:,:,:)
846              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
847              !
848              avt_k (:,:,:) = avt (:,:,:)
849              avm_k (:,:,:) = avm (:,:,:)
850              avmu_k(:,:,:) = avmu(:,:,:)
851              avmv_k(:,:,:) = avmv(:,:,:)
852              !
853              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
854           ENDIF
855        ELSE                                   !* Start from rest
856           en(:,:,:) = rn_emin * tmask(:,:,:)
857           DO jk = 1, jpk                           ! set the Kz to the background value
858              avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
859              avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
860              avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
861              avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
862           END DO
863        ENDIF
864        !
865     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
866        !                                   ! -------------------
867        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
868        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     )
869        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
870        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
871        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
872        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
873        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl  )
874        !
875     ENDIF
876     !
877   END SUBROUTINE tke_rst
878
879#else
880   !!----------------------------------------------------------------------
881   !!   Dummy module :                                        NO TKE scheme
882   !!----------------------------------------------------------------------
883   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
884CONTAINS
885   SUBROUTINE zdf_tke_init           ! Dummy routine
886   END SUBROUTINE zdf_tke_init
887   SUBROUTINE zdf_tke( kt )          ! Dummy routine
888      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
889   END SUBROUTINE zdf_tke
890   SUBROUTINE tke_rst( kt, cdrw )
891     CHARACTER(len=*) ::   cdrw
892     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
893   END SUBROUTINE tke_rst
894#endif
895
896   !!======================================================================
897END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.