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

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 4616

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

#1260 : see the associated wiki page for explanation

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