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

source: branches/2014/dev_CNRS1_10_TEOS10_Ediag/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 4619

Last change on this file since 4619 was 4619, checked in by gm, 10 years ago

#1294 : TEOS-10 and Ediag

  • Property svn:keywords set to Id
File size: 45.0 KB
Line 
1MODULE zdftke
2   !!======================================================================
3   !!                       ***  MODULE  zdftke  ***
4   !! Ocean physics:  vertical mixing coefficient computed from the tke
5   !!                 turbulent closure parameterization
6   !!=====================================================================
7   !! History :  OPA  !  1991-03  (b. blanke)  Original code
8   !!            7.0  !  1991-11  (G. Madec)   bug fix
9   !!            7.1  !  1992-10  (G. Madec)   new mixing length and eav
10   !!            7.2  !  1993-03  (M. Guyon)   symetrical conditions
11   !!            7.3  !  1994-08  (G. Madec, M. Imbard)  nn_pdl flag
12   !!            7.5  !  1996-01  (G. Madec)   s-coordinates
13   !!            8.0  !  1997-07  (G. Madec)   lbc
14   !!            8.1  !  1999-01  (E. Stretta) new option for the mixing length
15   !!  NEMO      1.0  !  2002-06  (G. Madec) add tke_init routine
16   !!             -   !  2004-10  (C. Ethe )  1D configuration
17   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
18   !!            3.0  !  2008-05  (C. Ethe,  G.Madec) : update TKE physics:
19   !!                 !           - tke penetration (wind steering)
20   !!                 !           - suface condition for tke & mixing length
21   !!                 !           - Langmuir cells
22   !!             -   !  2008-05  (J.-M. Molines, G. Madec)  2D form of avtb
23   !!             -   !  2008-06  (G. Madec)  style + DOCTOR name for namelist parameters
24   !!             -   !  2008-12  (G. Reffray) stable discretization of the production term
25   !!            3.2  !  2009-06  (G. Madec, S. Masson) TKE restart compatible with key_cpl
26   !!                 !                                + cleaning of the parameters + bugs correction
27   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
28   !!----------------------------------------------------------------------
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            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
244         END DO
245      END DO
246     
247!!bfr   - start commented area
248      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
249      !                     !  Bottom boundary condition on tke
250      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
251      !
252      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
253      ! Tests to date have found the bottom boundary condition on tke to have very little effect.
254      ! The condition is coded here for completion but commented out until there is proof that the
255      ! computational cost is justified
256      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
257      !                     en(bot)   = (rn_ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
258!CDIR NOVERRCHK
259!!    DO jj = 2, jpjm1
260!CDIR NOVERRCHK
261!!       DO ji = fs_2, fs_jpim1   ! vector opt.
262!!          ztx2 = bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) + &
263!!                 bfrua(ji  ,jj) * ub(ji  ,jj,mbku(ji  ,jj) )
264!!          zty2 = bfrva(ji,jj  ) * vb(ji,jj  ,mbkv(ji,jj  )) + &
265!!                 bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1) )
266!!          zebot = 0.001875_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 )   !  where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.
267!!          en (ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * tmask(ji,jj,1)
268!!       END DO
269!!    END DO
270!!bfr   - end commented area
271      !
272      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
273      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke       (Axell JGR 2002)
274         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
275         !
276         !                        !* total energy produce by LC : cumulative sum over jk
277         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * fsdepw(:,:,1) * fse3w(:,:,1)
278         DO jk = 2, jpk
279            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * fsdepw(:,:,jk) * fse3w(:,:,jk)
280         END DO
281         !                        !* finite Langmuir Circulation depth
282         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
283         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
284         DO jk = jpkm1, 2, -1
285            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
286               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
287                  zus  = zcof * taum(ji,jj)
288                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
289               END DO
290            END DO
291         END DO
292         !                               ! finite LC depth
293         DO jj = 1, jpj 
294            DO ji = 1, jpi
295               zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj))
296            END DO
297         END DO
298         zcof = 0.016 / SQRT( zrhoa * zcdrag )
299         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
300            DO jj = 2, jpjm1
301               DO ji = fs_2, fs_jpim1   ! vector opt.
302                  zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
303                  !                                           ! vertical velocity due to LC
304                  zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) )
305                  zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) )
306                  !                                           ! TKE Langmuir circulation source term
307                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj) * tmask(ji,jj,jk)
308               END DO
309            END DO
310         END DO
311         !
312      ENDIF
313      !
314      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
315      !                     !  Now Turbulent kinetic energy (output in en)
316      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
317      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
318      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
319      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
320      !
321      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form)
322         DO jj = 1, jpj                 ! here avmu, avmv used as workspace
323            DO ji = 1, jpi
324               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
325                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   & 
326                  &           / (  fse3uw_n(ji,jj,jk)         &
327                  &              * fse3uw_b(ji,jj,jk) )
328               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
329                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
330                  &                            / (  fse3vw_n(ji,jj,jk)               &
331                  &                              *  fse3vw_b(ji,jj,jk)  )
332            END DO
333         END DO
334      END DO
335      !
336      DO jk = 2, jpkm1           !* Matrix and right hand side in en
337         DO jj = 2, jpjm1
338            DO ji = fs_2, fs_jpim1   ! vector opt.
339               zcof   = zfact1 * tmask(ji,jj,jk)
340               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal
341                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
342               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal
343                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
344                  !                                                           ! shear prod. at w-point weightened by mask
345               zesh2  =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
346                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )   
347                  !
348               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
349               zd_lw(ji,jj,jk) = zzd_lw
350               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk)
351               !
352               !                                   ! right hand side in en
353               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  zesh2  -   avt(ji,jj,jk) * rn2(ji,jj,jk)    &
354                  &                                 + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk)  ) * tmask(ji,jj,jk)
355            END DO
356         END DO
357      END DO
358      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
359      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
360         DO jj = 2, jpjm1
361            DO ji = fs_2, fs_jpim1    ! vector opt.
362               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
363            END DO
364         END DO
365      END DO
366      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
367         DO ji = fs_2, fs_jpim1    ! vector opt.
368            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
369         END DO
370      END DO
371      DO jk = 3, jpkm1
372         DO jj = 2, jpjm1
373            DO ji = fs_2, fs_jpim1    ! vector opt.
374               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)
375            END DO
376         END DO
377      END DO
378      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
379         DO ji = fs_2, fs_jpim1    ! vector opt.
380            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
381         END DO
382      END DO
383      DO jk = jpk-2, 2, -1
384         DO jj = 2, jpjm1
385            DO ji = fs_2, fs_jpim1    ! vector opt.
386               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
387            END DO
388         END DO
389      END DO
390      DO jk = 2, jpkm1                             ! set the minimum value of tke
391         DO jj = 2, jpjm1
392            DO ji = fs_2, fs_jpim1   ! vector opt.
393               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * tmask(ji,jj,jk)
394            END DO
395         END DO
396      END DO
397
398      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
399      !                            !  TKE due to surface and internal wave breaking
400      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
401      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
402         DO jk = 2, jpkm1
403            DO jj = 2, jpjm1
404               DO ji = fs_2, fs_jpim1   ! vector opt.
405                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
406                     &                                               * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk)
407               END DO
408            END DO
409         END DO
410      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
411         DO jj = 2, jpjm1
412            DO ji = fs_2, fs_jpim1   ! vector opt.
413               jk = nmln(ji,jj)
414               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
415                  &                                               * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk)
416            END DO
417         END DO
418      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
419!CDIR NOVERRCHK
420         DO jk = 2, jpkm1
421!CDIR NOVERRCHK
422            DO jj = 2, jpjm1
423!CDIR NOVERRCHK
424               DO ji = fs_2, fs_jpim1   ! vector opt.
425                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
426                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
427                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 )    ! module of the mean stress
428                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
429                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
430                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
431                     &                                        * ( 1._wp - fr_i(ji,jj) ) * tmask(ji,jj,jk)
432               END DO
433            END DO
434         END DO
435      ENDIF
436      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
437      !
438      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer
439      CALL wrk_dealloc( jpi,jpj, zhlc ) 
440      CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
441      !
442      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
443      !
444   END SUBROUTINE tke_tke
445
446
447   SUBROUTINE tke_avn
448      !!----------------------------------------------------------------------
449      !!                   ***  ROUTINE tke_avn  ***
450      !!
451      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
452      !!
453      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
454      !!              the tke_tke routine). First, the now mixing lenth is
455      !!      computed from en and the strafification (N^2), then the mixings
456      !!      coefficients are computed.
457      !!              - Mixing length : a first evaluation of the mixing lengh
458      !!      scales is:
459      !!                      mxl = sqrt(2*en) / N 
460      !!      where N is the brunt-vaisala frequency, with a minimum value set
461      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
462      !!        The mixing and dissipative length scale are bound as follow :
463      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
464      !!                        zmxld = zmxlm = mxl
465      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
466      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
467      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
468      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
469      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
470      !!                    the surface to obtain ldown. the resulting length
471      !!                    scales are:
472      !!                        zmxld = sqrt( lup * ldown )
473      !!                        zmxlm = min ( lup , ldown )
474      !!              - Vertical eddy viscosity and diffusivity:
475      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
476      !!                      avt = max( avmb, pdlr * avm ) 
477      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
478      !!
479      !! ** Action  : - avt : now vertical eddy diffusivity (w-point)
480      !!              - avmu, avmv : now vertical eddy viscosity at uw- and vw-points
481      !!----------------------------------------------------------------------
482      INTEGER  ::   ji, jj, jk   ! dummy loop indices
483      REAL(wp) ::   zrn2, zraug, zcoef, zav     ! local scalars
484      REAL(wp) ::   zdku, zpdlr, zri, zsqen     !   -      -
485      REAL(wp) ::   zdkv, zemxl, zemlm, zemlp   !   -      -
486      REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld
487      !!--------------------------------------------------------------------
488      !
489      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
490
491      CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
492
493      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
494      !                     !  Mixing length
495      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
496      !
497      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
498      !
499      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
500         zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
501         zmxlm(:,:,1) = MAX(  rn_mxl0,  zraug * taum(:,:)  )
502      ELSE                          ! surface set to the minimum value
503         zmxlm(:,:,1) = rn_mxl0
504      ENDIF
505      zmxlm(:,:,jpk)  = rmxl_min     ! last level set to the interior minium value
506      !
507!CDIR NOVERRCHK
508      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
509!CDIR NOVERRCHK
510         DO jj = 2, jpjm1
511!CDIR NOVERRCHK
512            DO ji = fs_2, fs_jpim1   ! vector opt.
513               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
514               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
515            END DO
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 jk = 2, jpkm1
528            DO jj = 2, jpjm1
529               DO ji = fs_2, fs_jpim1   ! vector opt.
530                  zemxl = MIN( fsdepw(ji,jj,jk), 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 jk = 2, jpkm1
540            DO jj = 2, jpjm1
541               DO ji = fs_2, fs_jpim1   ! vector opt.
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 jk = 2, jpkm1         ! from the surface to the bottom :
551            DO jj = 2, jpjm1
552               DO ji = fs_2, fs_jpim1   ! vector opt.
553                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
554               END DO
555            END DO
556         END DO
557         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
558            DO jj = 2, jpjm1
559               DO ji = fs_2, fs_jpim1   ! vector opt.
560                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
561                  zmxlm(ji,jj,jk) = zemxl
562                  zmxld(ji,jj,jk) = zemxl
563               END DO
564            END DO
565         END DO
566         !
567      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
568         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
569            DO jj = 2, jpjm1
570               DO ji = fs_2, fs_jpim1   ! vector opt.
571                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
572               END DO
573            END DO
574         END DO
575         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
576            DO jj = 2, jpjm1
577               DO ji = fs_2, fs_jpim1   ! vector opt.
578                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
579               END DO
580            END DO
581         END DO
582!CDIR NOVERRCHK
583         DO jk = 2, jpkm1
584!CDIR NOVERRCHK
585            DO jj = 2, jpjm1
586!CDIR NOVERRCHK
587               DO ji = fs_2, fs_jpim1   ! vector opt.
588                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
589                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
590                  zmxlm(ji,jj,jk) = zemlm
591                  zmxld(ji,jj,jk) = zemlp
592               END DO
593            END DO
594         END DO
595         !
596      END SELECT
597      !
598# if defined key_c1d
599      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
600      e_mix(:,:,:) = zmxlm(:,:,:)
601# endif
602
603      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
604      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
605      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
606!CDIR NOVERRCHK
607      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
608!CDIR NOVERRCHK
609         DO jj = 2, jpjm1
610!CDIR NOVERRCHK
611            DO ji = fs_2, fs_jpim1   ! vector opt.
612               zsqen = SQRT( en(ji,jj,jk) )
613               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
614               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * tmask(ji,jj,jk)
615               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
616               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
617            END DO
618         END DO
619      END DO
620      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
621      !
622      DO jk = 2, jpkm1            !* vertical eddy viscosity at u- and v-points
623         DO jj = 2, jpjm1
624            DO ji = fs_2, fs_jpim1   ! vector opt.
625               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * umask(ji,jj,jk)
626               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * vmask(ji,jj,jk)
627            END DO
628         END DO
629      END DO
630      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
631      !
632      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
633         DO jk = 2, jpkm1
634            DO jj = 2, jpjm1
635               DO ji = fs_2, fs_jpim1   ! vector opt.
636                  zcoef = avm(ji,jj,jk) * 2._wp * fse3w(ji,jj,jk) * fse3w(ji,jj,jk)
637                  !                                          ! shear
638                  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) )   &
639                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
640                  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) )   &
641                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
642                  !                                          ! local Richardson number
643                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) * zcoef / (zdku + zdkv + rn_bshear )
644                  zpdlr = MAX(  0.1_wp,  0.2 / MAX( 0.2 , zri )  )
645!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
646!!gm              zpdlr = MAX(  0.1_wp,  ri_crit / MAX( ri_crit , zri )  )
647                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
648# if defined key_c1d
649                  e_pdl(ji,jj,jk) = zpdlr * tmask(ji,jj,jk)    ! c1d configuration : save masked Prandlt number
650                  e_ric(ji,jj,jk) = zri * tmask(ji,jj,jk)                            ! c1d config. : save Ri
651# endif
652              END DO
653            END DO
654         END DO
655      ENDIF
656      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
657
658      IF(ln_ctl) THEN
659         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
660         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
661            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
662      ENDIF
663      !
664      CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
665      !
666      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
667      !
668   END SUBROUTINE tke_avn
669
670
671   SUBROUTINE zdf_tke_init
672      !!----------------------------------------------------------------------
673      !!                  ***  ROUTINE zdf_tke_init  ***
674      !!                     
675      !! ** Purpose :   Initialization of the vertical eddy diffivity and
676      !!              viscosity when using a tke turbulent closure scheme
677      !!
678      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
679      !!              called at the first timestep (nit000)
680      !!
681      !! ** input   :   Namlist namzdf_tke
682      !!
683      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
684      !!----------------------------------------------------------------------
685      INTEGER ::   ji, jj, jk   ! dummy loop indices
686      INTEGER ::   ios
687      !!
688      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
689         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
690         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   &
691         &                 nn_etau , nn_htau  , rn_efr   
692      !!----------------------------------------------------------------------
693      !
694      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
695      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
696901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
697
698      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
699      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
700902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
701      WRITE ( numond, namzdf_tke )
702      !
703      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
704      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
705      !
706      IF(lwp) THEN                    !* Control print
707         WRITE(numout,*)
708         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
709         WRITE(numout,*) '~~~~~~~~~~~~'
710         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
711         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
712         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
713         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
714         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
715         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
716         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
717         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
718         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
719         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
720         WRITE(numout,*) '      surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
721         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
722         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
723         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
724         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
725         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
726         WRITE(numout,*)
727         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
728      ENDIF
729      !
730      !                              ! allocate tke arrays
731      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
732      !
733      !                               !* Check of some namelist values
734      IF( nn_mxl  < 0  .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
735      IF( nn_pdl  < 0  .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
736      IF( nn_htau < 0  .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
737#if ! key_coupled
738      IF( nn_etau == 3 )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
739#endif
740
741      IF( ln_mxl0 ) THEN
742         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
743         rn_mxl0 = rmxl_min
744      ENDIF
745     
746      IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
747
748      !                               !* depth of penetration of surface tke
749      IF( nn_etau /= 0 ) THEN     
750         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
751         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
752            htau(:,:) = 10._wp
753         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
754            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
755         END SELECT
756      ENDIF
757      !                               !* set vertical eddy coef. to the background value
758      DO jk = 1, jpk
759         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
760         avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
761         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
762         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
763      END DO
764      dissl(:,:,:) = 1.e-12_wp
765      !                             
766      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
767      !
768   END SUBROUTINE zdf_tke_init
769
770
771   SUBROUTINE tke_rst( kt, cdrw )
772     !!---------------------------------------------------------------------
773     !!                   ***  ROUTINE tke_rst  ***
774     !!                     
775     !! ** Purpose :   Read or write TKE file (en) in restart file
776     !!
777     !! ** Method  :   use of IOM library
778     !!                if the restart does not contain TKE, en is either
779     !!                set to rn_emin or recomputed
780     !!----------------------------------------------------------------------
781     INTEGER         , INTENT(in) ::   kt     ! ocean time-step
782     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
783     !
784     INTEGER ::   jit, jk   ! dummy loop indices
785     INTEGER ::   id1, id2, id3, id4, id5, id6   ! local integers
786     !!----------------------------------------------------------------------
787     !
788     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
789        !                                   ! ---------------
790        IF( ln_rstart ) THEN                   !* Read the restart file
791           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
792           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
793           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
794           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
795           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
796           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
797           !
798           IF( id1 > 0 ) THEN                       ! 'en' exists
799              CALL iom_get( numror, jpdom_autoglo, 'en', en )
800              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
801                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt   )
802                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm   )
803                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu  )
804                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv  )
805                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
806              ELSE                                                 ! one at least array is missing
807                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
808              ENDIF
809           ELSE                                     ! No TKE array found: initialisation
810              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
811              en (:,:,:) = rn_emin * tmask(:,:,:)
812              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
813              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
814           ENDIF
815        ELSE                                   !* Start from rest
816           en(:,:,:) = rn_emin * tmask(:,:,:)
817           DO jk = 1, jpk                           ! set the Kz to the background value
818              avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
819              avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
820              avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
821              avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
822           END DO
823        ENDIF
824        !
825     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
826        !                                   ! -------------------
827        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
828        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     )
829        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
830        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
831        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
832        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
833        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl  )
834        !
835     ENDIF
836     !
837   END SUBROUTINE tke_rst
838
839#else
840   !!----------------------------------------------------------------------
841   !!   Dummy module :                                        NO TKE scheme
842   !!----------------------------------------------------------------------
843   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
844CONTAINS
845   SUBROUTINE zdf_tke_init           ! Dummy routine
846   END SUBROUTINE zdf_tke_init
847   SUBROUTINE zdf_tke( kt )          ! Dummy routine
848      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
849   END SUBROUTINE zdf_tke
850   SUBROUTINE tke_rst( kt, cdrw )
851     CHARACTER(len=*) ::   cdrw
852     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
853   END SUBROUTINE tke_rst
854#endif
855
856   !!======================================================================
857END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.