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

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

source: branches/UKMO/dev_r10171_test_crs_AMM7/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 10207

Last change on this file since 10207 was 10207, checked in by cmao, 5 years ago

remove svn keyword

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