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

source: branches/UKMO/dev_r5107_xios_initialize_toyoce/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 5277

Last change on this file since 5277 was 5277, checked in by dancopsey, 9 years ago

Removed SVN keywords.

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