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_old.F90 in branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/ZDF – NEMO

source: branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/ZDF/zdftke_old.F90 @ 2236

Last change on this file since 2236 was 2236, checked in by cetlod, 14 years ago

First guess of NEMO_v3.3

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 46.9 KB
Line 
1MODULE zdftke_old
2   !!======================================================================
3   !!                       ***  MODULE  zdftke_old  ***
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 zdf_tke_init routine
16   !!              -   !  2002-08  (G. Madec)  rn_cri and Free form, F90
17   !!              -   !  2004-10  (C. Ethe )  1D configuration
18   !!             2.0  !  2006-07  (S. Masson)  distributed restart using iom
19   !!             3.0  !  2008-05  (C. Ethe,  G.Madec) : update TKE physics:
20   !!                              - tke penetration (wind steering)
21   !!                              - suface condition for tke & mixing length
22   !!                              - Langmuir cells
23   !!              -   !  2008-05  (J.-M. Molines, G. Madec)  2D form of avtb
24   !!              -   !  2008-06  (G. Madec)  style + DOCTOR name for namelist parameters
25   !!----------------------------------------------------------------------
26#if defined key_zdftke_old   ||   defined key_esopa
27   !!----------------------------------------------------------------------
28   !!   'key_zdftke_old'                               TKE vertical physics
29   !!----------------------------------------------------------------------
30   !!----------------------------------------------------------------------
31   !!   zdf_tke_old  : update momentum and tracer Kz from a tke scheme
32   !!   zdf_tke_init : initialization, namelist read, and parameters control
33   !!   tke_rst      : read/write tke restart in ocean restart file
34   !!----------------------------------------------------------------------
35   USE oce             ! ocean dynamics and active tracers
36   USE dom_oce         ! ocean space and time domain
37   USE zdf_oce         ! ocean vertical physics
38   USE sbc_oce         ! surface boundary condition: ocean
39   USE phycst          ! physical constants
40   USE zdfmxl          ! mixed layer
41   USE restart         ! only for lrst_oce
42   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
43   USE prtctl          ! Print control
44   USE in_out_manager  ! I/O manager
45   USE iom             ! I/O manager library
46
47   IMPLICIT NONE
48   PRIVATE
49
50   PUBLIC   zdf_tke_old     ! routine called in step module
51   PUBLIC   zdf_tke_init_o  ! routine called in opa module
52
53   LOGICAL , PUBLIC, PARAMETER              ::   lk_zdftke_old = .TRUE.  !: TKE vertical mixing flag
54   REAL(wp), PUBLIC                         ::   eboost              !: multiplicative coeff of the shear product.
55   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   en                  !: now turbulent kinetic energy
56# if defined key_vectopt_memory
57   !                                                                !!! key_vectopt_memory
58   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   etmean              !: coefficient used for horizontal smoothing
59   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   eumean, evmean      !: at t-, u- and v-points
60# endif
61#if defined key_c1d
62   !                                                                !!! 1D cfg only
63   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   e_dis, e_mix        !: dissipation and mixing turbulent lengh scales
64   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   e_pdl, e_ric        !: prandl and local Richardson numbers
65   REAL(wp), PUBLIC, DIMENSION(jpi,jpj)     ::   hlc                 !: save finite Langmuir Circulation depth
66#endif
67
68   !                                       !!! ** Namelist  namzdf_tke  **
69   LOGICAL  ::   ln_rstke = .FALSE.         ! =T restart with tke from a run without tke
70   LOGICAL  ::   ln_mxl0  = .FALSE.         ! mixing length scale surface value as function of wind stress or not
71   LOGICAL  ::   ln_lc    = .FALSE.         ! Langmuir cells (LC) as a source term of TKE or not
72   INTEGER  ::   nn_itke  = 50              ! number of restart iterative loops
73   INTEGER  ::   nn_mxl   =  2              ! type of mixing length (=0/1/2/3)
74   INTEGER  ::   nn_pdl   =  1              ! Prandtl number or not (ratio avt/avm) (=0/1)
75   INTEGER  ::   nn_ave   =  1              ! horizontal average or not on avt, avmu, avmv (=0/1)
76   REAL(wp) ::   rn_ediff = 0.1_wp          ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
77   REAL(wp) ::   rn_ediss = 0.7_wp          ! coefficient of the Kolmogoroff dissipation
78   REAL(wp) ::   rn_ebb   = 3.75_wp         ! coefficient of the surface input of tke
79   REAL(wp) ::   rn_efave = 1._wp           ! coefficient for ave : ave=rn_efave*avm
80   REAL(wp) ::   rn_emin  = 0.7071e-6_wp    ! minimum value of tke (m2/s2)
81   REAL(wp) ::   rn_emin0 = 1.e-4_wp        ! surface minimum value of tke (m2/s2)
82   REAL(wp) ::   rn_cri   = 2._wp / 9._wp   ! critic Richardson number
83   INTEGER  ::   nn_etau  = 0               ! type of depth penetration of surface tke (=0/1/2)
84   INTEGER  ::   nn_htau  = 0               ! type of tke profile of penetration (=0/1)
85   REAL(wp) ::   rn_lmin0 = 0.4_wp          ! surface  min value of mixing length
86   REAL(wp) ::   rn_lmin  = 0.1_wp          ! interior min value of mixing length
87   REAL(wp) ::   rn_efr   = 1.0_wp          ! fraction of TKE surface value which penetrates in the ocean
88   REAL(wp) ::   rn_lc    = 0.15_wp         ! coef to compute vertical velocity of Langmuir cells
89
90   !! * Substitutions
91#  include "domzgr_substitute.h90"
92#  include "vectopt_loop_substitute.h90"
93   !!----------------------------------------------------------------------
94   !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008)
95   !! $Id$
96   !! Software governed by the CeCILL licence  (NEMOGCM/License_CeCILL.txt)
97   !!----------------------------------------------------------------------
98
99CONTAINS
100
101   SUBROUTINE zdf_tke_old( kt )
102      !!----------------------------------------------------------------------
103      !!                   ***  ROUTINE zdf_tke_old  ***
104      !!
105      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
106      !!      coefficients using a 1.5 turbulent closure scheme.
107      !!
108      !! ** Method  :   The time evolution of the turbulent kinetic energy
109      !!      (tke) is computed from a prognostic equation :
110      !!         d(en)/dt = eboost eav (d(u)/dz)**2         ! shear production
111      !!                  + d( rn_efave eav d(en)/dz )/dz   ! diffusion of tke
112      !!                  + grav/rau0 pdl eav d(rau)/dz     ! stratif. destruc.
113      !!                  - rn_ediss / emxl en**(2/3)       ! dissipation
114      !!      with the boundary conditions:
115      !!         surface: en = max( rn_emin0, rn_ebb sqrt(utau^2 + vtau^2) )
116      !!         bottom : en = rn_emin
117      !!      -1- The dissipation and mixing turbulent lengh scales are computed
118      !!         from the usual diagnostic buoyancy length scale: 
119      !!         mxl= sqrt(2*en)/N  where N is the brunt-vaisala frequency
120      !!         with mxl = rn_lmin at the bottom minimum value of 0.4
121      !!      Four cases :
122      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
123      !!                  zmxld = zmxlm = mxl
124      !!         nn_mxl=1 : mxl bounded by the vertical scale factor.
125      !!                  zmxld = zmxlm = mxl
126      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl
127      !!                  is less than 1 (|d/dz(xml)|<1).
128      !!                  zmxld = zmxlm = mxl
129      !!         nn_mxl=3 : lup = mxl bounded using |d/dz(xml)|<1 from the surface
130      !!                        to the bottom
131      !!                  ldown = mxl bounded using |d/dz(xml)|<1 from the bottom
132      !!                        to the surface
133      !!                  zmxld = sqrt (lup*ldown) ; zmxlm = min(lup,ldown)
134      !!      -2- Compute the now Turbulent kinetic energy. The time differencing
135      !!      is implicit for vertical diffusion term, linearized for kolmo-
136      !!      goroff dissipation term, and explicit forward for both buoyancy
137      !!      and dynamic production terms. Thus a tridiagonal linear system is
138      !!      solved.
139      !!         Note that - the shear production is multiplied by eboost in order
140      !!      to set the critic richardson number to rn_cri (namelist parameter)
141      !!                   - the destruction by stratification term is multiplied
142      !!      by the Prandtl number (defined by an empirical funtion of the local
143      !!      Richardson number) if nn_pdl=1 (namelist parameter)
144      !!      coefficient (zesh2):
145      !!      -3- Compute the now vertical eddy vicosity and diffusivity
146      !!      coefficients from en (before the time stepping) and zmxlm:
147      !!              avm = max( avtb, rn_ediff*zmxlm*en^1/2 )
148      !!              avt = max( avmb, pdl*avm )  (pdl=1 if nn_pdl=0)
149      !!              eav = max( avmb, avm )
150      !!      avt and avm are horizontally averaged to avoid numerical insta-
151      !!      bilities.
152      !!        N.B. The computation is done from jk=2 to jpkm1 except for
153      !!      en. Surface value of avt avmu avmv are set once a time to
154      !!      their background value in routine zdf_tke_init.
155      !!
156      !! ** Action  :   compute en (now turbulent kinetic energy)
157      !!                update avt, avmu, avmv (before vertical eddy coef.)
158      !!
159      !! References : Gaspar et al., JGR, 1990,
160      !!              Blanke and Delecluse, JPO, 1991
161      !!              Mellor and Blumberg, JPO 2004
162      !!              Axell, JGR, 2002
163      !!----------------------------------------------------------------------
164      USE oce,     zwd    =>   ua   ! use ua as workspace
165      USE oce,     zmxlm  =>   va   ! use va as workspace
166      USE oce,     zmxld  =>   ta   ! use ta as workspace
167      USE oce,     ztkelc =>   sa   ! use sa as workspace
168      !
169      INTEGER, INTENT(in) ::   kt   ! ocean time step
170      !
171      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
172      REAL(wp) ::   zbbrau, zrn2, zesurf            ! temporary scalars
173      REAL(wp) ::   zfact1, ztx2, zdku              !    -         -
174      REAL(wp) ::   zfact2, zty2, zdkv              !    -         -
175      REAL(wp) ::   zfact3, zcoef, zcof, zav        !    -         -
176      REAL(wp) ::   zsh2, zpdl, zri, zsqen, zesh2   !    -         -
177      REAL(wp) ::   zemxl, zemlm, zemlp             !    -         -
178      REAL(wp) ::   zraug, zus, zwlc, zind          !    -         -
179      INTEGER , DIMENSION(jpi,jpj)     ::   imlc    ! 2D workspace
180      REAL(wp), DIMENSION(jpi,jpj)     ::   zhtau   !  -      -
181      REAL(wp), DIMENSION(jpi,jpj)     ::   zhlc    !  -      -
182      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc   ! 3D workspace
183      !!--------------------------------------------------------------------
184
185      !                                            ! Local constant initialization
186      zbbrau =  .5 * rn_ebb / rau0
187      zfact1 = -.5 * rdt * rn_efave
188      zfact2 = 1.5 * rdt * rn_ediss
189      zfact3 = 0.5 * rdt * rn_ediss
190
191      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
192      ! I.  Mixing length
193      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
194
195      ! Buoyancy length scale: l=sqrt(2*e/n**2)
196      ! ---------------------
197      IF( ln_mxl0 ) THEN         ! surface mixing length = F(stress) : l=vkarmn*2.e5*sqrt(utau^2 + vtau^2)/(rau0*g)
198!!gm  this should be useless
199         zmxlm(:,:,1) = 0.e0
200!!gm end
201         zraug = 0.5 * vkarmn * 2.e5 / ( rau0 * grav )
202         DO jj = 2, jpjm1
203!CDIR NOVERRCHK
204            DO ji = fs_2, fs_jpim1   ! vector opt.
205               ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
206               zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
207               zmxlm(ji,jj,1) = MAX(  rn_lmin0,  zraug * SQRT( ztx2 * ztx2 + zty2 * zty2 )  )
208            END DO
209         END DO
210      ELSE                       ! surface set to the minimum value
211         zmxlm(:,:,1) = rn_lmin0
212      ENDIF
213      zmxlm(:,:,jpk) = rn_lmin   ! bottom set to the interior minium value
214      !
215!CDIR NOVERRCHK
216      DO jk = 2, jpkm1           ! interior value : l=sqrt(2*e/n**2)
217!CDIR NOVERRCHK
218         DO jj = 2, jpjm1
219!CDIR NOVERRCHK
220            DO ji = fs_2, fs_jpim1   ! vector opt.
221               zrn2 = MAX( rn2b(ji,jj,jk), rsmall )
222               zmxlm(ji,jj,jk) = MAX(  rn_lmin,  SQRT( 2. * en(ji,jj,jk) / zrn2 )  )
223            END DO
224         END DO
225      END DO
226
227      ! Physical limits for the mixing length
228      ! -------------------------------------
229      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
230      zmxld(:,:,jpk) = rn_lmin        ! bottom  set to the minimum value
231
232      SELECT CASE ( nn_mxl )
233      !
234      CASE ( 0 )           ! bounded by the distance to surface and bottom
235         DO jk = 2, jpkm1
236            DO jj = 2, jpjm1
237               DO ji = fs_2, fs_jpim1   ! vector opt.
238                  zemxl = MIN( fsdepw(ji,jj,jk), zmxlm(ji,jj,jk),   &
239                  &            fsdepw(ji,jj,mbathy(ji,jj)) - fsdepw(ji,jj,jk) )
240                  zmxlm(ji,jj,jk) = zemxl
241                  zmxld(ji,jj,jk) = zemxl
242               END DO
243            END DO
244         END DO
245         !
246      CASE ( 1 )           ! bounded by the vertical scale factor
247         DO jk = 2, jpkm1
248            DO jj = 2, jpjm1
249               DO ji = fs_2, fs_jpim1   ! vector opt.
250                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
251                  zmxlm(ji,jj,jk) = zemxl
252                  zmxld(ji,jj,jk) = zemxl
253               END DO
254            END DO
255         END DO
256         !
257      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
258         DO jk = 2, jpkm1         ! from the surface to the bottom :
259            DO jj = 2, jpjm1
260               DO ji = fs_2, fs_jpim1   ! vector opt.
261                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
262               END DO
263            END DO
264         END DO
265         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
266            DO jj = 2, jpjm1
267               DO ji = fs_2, fs_jpim1   ! vector opt.
268                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
269                  zmxlm(ji,jj,jk) = zemxl
270                  zmxld(ji,jj,jk) = zemxl
271               END DO
272            END DO
273         END DO
274         !
275      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
276         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
277            DO jj = 2, jpjm1
278               DO ji = fs_2, fs_jpim1   ! vector opt.
279                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
280               END DO
281            END DO
282         END DO
283         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
284            DO jj = 2, jpjm1
285               DO ji = fs_2, fs_jpim1   ! vector opt.
286                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
287               END DO
288            END DO
289         END DO
290!CDIR NOVERRCHK
291         DO jk = 2, jpkm1
292!CDIR NOVERRCHK
293            DO jj = 2, jpjm1
294!CDIR NOVERRCHK
295               DO ji = fs_2, fs_jpim1   ! vector opt.
296                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
297                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
298                  zmxlm(ji,jj,jk) = zemlm
299                  zmxld(ji,jj,jk) = zemlp
300               END DO
301            END DO
302         END DO
303         !
304      END SELECT
305
306# if defined key_c1d
307      ! c1d configuration : save mixing and dissipation turbulent length scales
308      e_dis(:,:,:) = zmxld(:,:,:)
309      e_mix(:,:,:) = zmxlm(:,:,:)
310# endif
311
312      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
313      ! II  TKE Langmuir circulation source term
314      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
315      IF( ln_lc ) THEN
316         !
317         ! Computation of total energy produce by LC : cumulative sum over jk
318         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0. ) * fsdepw(:,:,1) * fse3w(:,:,1)
319         DO jk = 2, jpk
320            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0. ) * fsdepw(:,:,jk) * fse3w(:,:,jk)
321         END DO
322         !
323         ! Computation of finite Langmuir Circulation depth
324         ! Initialization to the number of w ocean point mbathy
325         imlc(:,:) = mbathy(:,:)
326         DO jk = jpkm1, 2, -1
327            DO jj = 1, jpj
328               DO ji = 1, jpi
329                  ! Last w-level at which zpelc>=0.5*us*us with us=0.016*wind(starting from jpk-1)
330                  zus  = 0.000128 * wndm(ji,jj) * wndm(ji,jj)
331                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
332               END DO
333            END DO
334         END DO
335         !
336         ! finite LC depth
337         DO jj = 1, jpj
338            DO ji = 1, jpi
339               zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj))
340            END DO
341         END DO
342         !
343# if defined key_c1d
344         hlc(:,:) = zhlc(:,:) * tmask(:,:,1)      ! c1d configuration: save finite Langmuir Circulation depth
345# endif
346         !
347         ! TKE Langmuir circulation source term
348!CDIR NOVERRCHK
349         DO jk = 2, jpkm1
350!CDIR NOVERRCHK
351            DO jj = 2, jpjm1
352!CDIR NOVERRCHK
353               DO ji = fs_2, fs_jpim1   ! vector opt.
354                  ! Stokes drift
355                  zus  = 0.016 * wndm(ji,jj)
356                  ! computation of vertical velocity due to LC
357                  zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) )
358                  zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) )
359                  ! TKE Langmuir circulation source term
360                  ztkelc(ji,jj,jk) = ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
361               END DO
362            END DO
363         END DO
364         !
365      ENDIF
366
367      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
368      ! III  Tubulent kinetic energy time stepping
369      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
370
371      ! 1. Vertical eddy viscosity on tke (put in zmxlm) and first estimate of avt
372      ! ---------------------------------------------------------------------
373!CDIR NOVERRCHK
374      DO jk = 2, jpkm1
375!CDIR NOVERRCHK
376         DO jj = 2, jpjm1
377!CDIR NOVERRCHK
378            DO ji = fs_2, fs_jpim1   ! vector opt.
379               zsqen = SQRT( en(ji,jj,jk) )
380               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
381               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
382               zmxlm(ji,jj,jk) = MAX( zav, avmb(jk) ) * tmask(ji,jj,jk)
383               zmxld(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
384            END DO
385         END DO
386      END DO
387
388      ! 2. Surface boundary condition on tke and its eddy viscosity (zmxlm)
389      ! -------------------------------------------------
390      ! en(1)   = rn_ebb sqrt(utau^2+vtau^2) / rau0  (min value rn_emin0)
391      ! zmxlm(1) = avmb(1) and zmxlm(jpk) = 0.
392!CDIR NOVERRCHK
393      DO jj = 2, jpjm1
394!CDIR NOVERRCHK
395         DO ji = fs_2, fs_jpim1   ! vector opt.
396            ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
397            zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
398            zesurf = zbbrau * SQRT( ztx2 * ztx2 + zty2 * zty2 )
399            en (ji,jj,1) = MAX( zesurf, rn_emin0 ) * tmask(ji,jj,1)
400            zav  =  rn_ediff * zmxlm(ji,jj,1) * SQRT( en(ji,jj,1) )
401            zmxlm(ji,jj,1  ) = MAX( zav, avmb(1)                  ) * tmask(ji,jj,1)
402            avt  (ji,jj,1  ) = MAX( zav, avtb(1) * avtb_2d(ji,jj) ) * tmask(ji,jj,1)
403            zmxlm(ji,jj,jpk) = 0.e0
404         END DO
405      END DO
406
407      ! 3. Now Turbulent kinetic energy (output in en)
408      ! -------------------------------
409      ! Resolution of a tridiagonal linear system by a "methode de chasse"
410      ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
411
412      SELECT CASE ( nn_pdl )
413      !
414      CASE ( 0 )           ! No Prandtl number
415         DO jk = 2, jpkm1
416            DO jj = 2, jpjm1
417               DO ji = fs_2, fs_jpim1   ! vector opt.
418                  !                                                             ! shear prod. - stratif. destruction
419                  zcoef = 0.5 / fse3w(ji,jj,jk)
420                  zdku = zcoef * (  ub(ji-1, jj ,jk-1) + ub(ji,jj,jk-1)   &          ! shear
421                     &            - ub(ji-1, jj ,jk  ) - ub(ji,jj,jk  )  )
422                  zdkv = zcoef * (  vb( ji ,jj-1,jk-1) + vb(ji,jj,jk-1)   &
423                     &            - vb( ji ,jj-1,jk  ) - vb(ji,jj,jk  )  )
424                  zesh2 =  eboost * ( zdku*zdku + zdkv*zdkv ) - rn2b(ji,jj,jk)        ! coefficient (zesh2)
425                  !
426                  !                                                             ! Matrix
427                  zcof = zfact1 * tmask(ji,jj,jk)
428                  !                                                                  ! lower diagonal
429                  avmv(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk  ) + zmxlm(ji,jj,jk-1) )   &
430                     &                  / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
431                  !                                                                  ! upper diagonal
432                  avmu(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk+1) + zmxlm(ji,jj,jk  ) )   &
433                     &                  / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk) )
434                  !                                                                  ! diagonal
435                  zwd(ji,jj,jk) = 1. - avmv(ji,jj,jk) - avmu(ji,jj,jk) + zfact2 * zmxld(ji,jj,jk)
436                  !
437                  !                                                             ! right hand side in en
438                  IF( .NOT. ln_lc ) THEN                                             ! No Langmuir cells
439                     en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en(ji,jj,jk)   &
440                        &                        +   rdt  * zmxlm (ji,jj,jk) * zesh2
441                  ELSE                                                               ! Langmuir cell source term
442                     en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en(ji,jj,jk)   &
443                        &                        +   rdt  * zmxlm (ji,jj,jk) * zesh2          &
444                        &                        +   rdt  * ztkelc(ji,jj,jk)
445                  ENDIF
446               END DO
447            END DO
448         END DO
449         !
450      CASE ( 1 )           ! Prandtl number
451         DO jk = 2, jpkm1
452            DO jj = 2, jpjm1
453               DO ji = fs_2, fs_jpim1   ! vector opt.
454                  !                                                             ! shear prod. - stratif. destruction
455                  zcoef = 0.5 / fse3w(ji,jj,jk)                                     
456                  zdku = zcoef * (  ub(ji-1,jj  ,jk-1) + ub(ji,jj,jk-1)   &          ! shear
457                  &               - ub(ji-1,jj  ,jk  ) - ub(ji,jj,jk  )  )
458                  zdkv = zcoef * (  vb(ji  ,jj-1,jk-1) + vb(ji,jj,jk-1)   &
459                  &               - vb(ji  ,jj-1,jk  ) - vb(ji,jj,jk  )  )
460                  zsh2 = zdku * zdku + zdkv * zdkv                                   ! square of shear
461                  zri  = MAX( rn2b(ji,jj,jk), 0. ) / ( zsh2 + 1.e-20 )                ! local Richardson number
462# if defined key_c1d
463                  e_ric(ji,jj,jk) = zri * tmask(ji,jj,jk)                            ! c1d config. : save Ri
464# endif
465                  zpdl = 1.0                                                         ! Prandtl number
466                  IF( zri >= 0.2 )   zpdl = 0.2 / zri
467                  zpdl = MAX( 0.1, zpdl )
468                  zesh2 = eboost * zsh2 - zpdl * rn2b(ji,jj,jk)                       ! coefficient (esh2)
469                  !
470                  !                                                             ! Matrix
471                  zcof = zfact1 * tmask(ji,jj,jk)
472                  !                                                                 ! lower diagonal
473                  avmv(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk  ) + zmxlm(ji,jj,jk-1) )   &
474                  &                     / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
475                  !                                                                 ! upper diagonal
476                  avmu(ji,jj,jk) = zcof * ( zmxlm(ji,jj,jk+1) + zmxlm(ji,jj,jk  ) )   &
477                  &                     / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk) )
478                  !                                                                 ! diagonal
479                  zwd(ji,jj,jk) = 1. - avmv(ji,jj,jk) - avmu(ji,jj,jk) + zfact2 * zmxld(ji,jj,jk)
480                  !
481                  !                                                             ! right hand side in en
482                  IF( .NOT. ln_lc ) THEN                                             ! No Langmuir cells
483                     en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en   (ji,jj,jk)   &
484                        &                        +   rdt  * zmxlm (ji,jj,jk) * zesh2
485                  ELSE                                                               ! Langmuir cell source term
486                     en(ji,jj,jk) = en(ji,jj,jk) + zfact3 * zmxld (ji,jj,jk) * en   (ji,jj,jk)   &
487                        &                        +   rdt  * zmxlm (ji,jj,jk) * zesh2             &
488                        &                        +   rdt  * ztkelc(ji,jj,jk)
489                  ENDIF
490                  zmxld(ji,jj,jk) = zpdl * tmask(ji,jj,jk)                      ! store masked Prandlt number in zmxld array
491               END DO
492            END DO
493         END DO
494         !
495      END SELECT
496
497# if defined key_c1d
498      e_pdl(:,:,2:jpkm1) = zmxld(:,:,2:jpkm1)      ! c1d configuration : save masked Prandlt number
499      e_pdl(:,:,      1) = e_pdl(:,:,      2)
500      e_pdl(:,:,    jpk) = e_pdl(:,:,  jpkm1)     
501# endif
502
503      ! 4. Matrix inversion from level 2 (tke prescribed at level 1)
504      !!--------------------------------
505      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
506         DO jj = 2, jpjm1
507            DO ji = fs_2, fs_jpim1    ! vector opt.
508               zwd(ji,jj,jk) = zwd(ji,jj,jk) - avmv(ji,jj,jk) * avmu(ji,jj,jk-1) / zwd(ji,jj,jk-1)
509            END DO
510         END DO
511      END DO
512      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
513         DO ji = fs_2, fs_jpim1    ! vector opt.
514            avmv(ji,jj,2) = en(ji,jj,2) - avmv(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
515         END DO
516      END DO
517      DO jk = 3, jpkm1
518         DO jj = 2, jpjm1
519            DO ji = fs_2, fs_jpim1    ! vector opt.
520               avmv(ji,jj,jk) = en(ji,jj,jk) - avmv(ji,jj,jk) / zwd(ji,jj,jk-1) *avmv(ji,jj,jk-1)
521            END DO
522         END DO
523      END DO
524      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
525         DO ji = fs_2, fs_jpim1    ! vector opt.
526            en(ji,jj,jpkm1) = avmv(ji,jj,jpkm1) / zwd(ji,jj,jpkm1)
527         END DO
528      END DO
529      DO jk = jpk-2, 2, -1
530         DO jj = 2, jpjm1
531            DO ji = fs_2, fs_jpim1    ! vector opt.
532               en(ji,jj,jk) = ( avmv(ji,jj,jk) - avmu(ji,jj,jk) * en(ji,jj,jk+1) ) / zwd(ji,jj,jk)
533            END DO
534         END DO
535      END DO
536      DO jk = 2, jpkm1                             ! set the minimum value of tke
537         DO jj = 2, jpjm1
538            DO ji = fs_2, fs_jpim1   ! vector opt.
539               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * tmask(ji,jj,jk)
540            END DO
541         END DO
542      END DO
543     
544      ! 5. Add extra TKE due to surface and internal wave breaking (nn_etau /= 0)
545      !!----------------------------------------------------------
546      IF( nn_etau /= 0 ) THEN        ! extra tke : en = en + rn_efr * en(1) * exp( -z/zhtau )
547         !
548         SELECT CASE( nn_htau )           ! Choice of the depth of penetration
549         CASE( 0 )                                    ! constant depth penetration (here 10 meters)
550            DO jj = 2, jpjm1
551               DO ji = fs_2, fs_jpim1   ! vector opt.
552                  zhtau(ji,jj) = 10.
553               END DO
554            END DO
555         CASE( 1 )                                    ! meridional profile  1
556            DO jj = 2, jpjm1                          ! ( 0.5m in the tropics to a maximum of 30 m at high lat.)
557               DO ji = fs_2, fs_jpim1   ! vector opt.
558                  zhtau(ji,jj) = MAX( 0.5, 3./4. * MIN( 40., 60.*ABS( SIN( rpi/180. * gphit(ji,jj) ) ) )   )
559               END DO
560            END DO
561         END SELECT
562         !
563         IF( nn_etau == 1 ) THEN          ! extra term throughout the water column
564            DO jk = 2, jpkm1
565               DO jj = 2, jpjm1
566                  DO ji = fs_2, fs_jpim1   ! vector opt.
567                     en(ji,jj,jk) = en(ji,jj,jk)   &
568                        &         + rn_efr * en(ji,jj,1)*EXP( -fsdepw(ji,jj,jk) / zhtau(ji,jj) ) &
569                        &                  * ( 1.e0 - fr_i(ji,jj) )  * tmask(ji,jj,jk)
570                  END DO
571               END DO
572            END DO
573         ELSEIF( nn_etau == 2 ) THEN      ! extra term only at the base of the mixed layer
574            DO jj = 2, jpjm1
575               DO ji = fs_2, fs_jpim1   ! vector opt.
576                  jk = nmln(ji,jj)
577                  en(ji,jj,jk) = en(ji,jj,jk)    &
578                     &         + rn_efr * en(ji,jj,1)*EXP( -fsdepw(ji,jj,jk) / zhtau(ji,jj) ) &
579                     &                   * ( 1.e0 - fr_i(ji,jj) )  * tmask(ji,jj,jk)
580               END DO
581            END DO
582         ENDIF
583         !
584      ENDIF
585
586
587      ! Lateral boundary conditions (sign unchanged)
588      CALL lbc_lnk( en, 'W', 1. )   ;   CALL lbc_lnk( avt, 'W', 1. )
589
590
591      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
592      ! IV.  Before vertical eddy vicosity and diffusivity coefficients
593      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
594      !
595      SELECT CASE ( nn_ave )
596      CASE ( 0 )                                     ! no horizontal average
597         DO jk = 2, jpkm1                                 ! only vertical eddy viscosity
598            DO jj = 2, jpjm1
599               DO ji = fs_2, fs_jpim1   ! vector opt.
600                  avmu(ji,jj,jk) = ( avt  (ji,jj,jk) + avt  (ji+1,jj  ,jk) ) * umask(ji,jj,jk)   &
601                  &     / MAX( 1.,   tmask(ji,jj,jk) + tmask(ji+1,jj  ,jk) )
602                  avmv(ji,jj,jk) = ( avt  (ji,jj,jk) + avt  (ji  ,jj+1,jk) ) * vmask(ji,jj,jk)   &
603                  &     / MAX( 1.,   tmask(ji,jj,jk) + tmask(ji  ,jj+1,jk) )
604               END DO
605            END DO
606         END DO
607         CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
608         !
609      CASE ( 1 )                                     ! horizontal average                         ( 1/2  1/2 )
610         !                                                ! Vertical eddy viscosity    avmu = 1/4 ( 1    1   )
611         !                                                !                                       ( 1/2  1/2 )
612         !                                                !
613         !                                                !                                       ( 1/2  1   1/2 )           
614         !                                                !                            avmv = 1/4 ( 1/2  1   1/2 )     
615         DO jk = 2, jpkm1                               
616            DO jj = 2, jpjm1
617               DO ji = fs_2, fs_jpim1   ! vector opt.
618# if defined key_vectopt_memory
619                  !                                       ! caution : vectopt_memory change the solution
620                  !                                       !           (last digit of the solver stat)
621                  avmu(ji,jj,jk) = (      avt(ji,jj  ,jk) + avt(ji+1,jj  ,jk)   &
622                     &              +.5*( avt(ji,jj-1,jk) + avt(ji+1,jj-1,jk)   &
623                     &                   +avt(ji,jj+1,jk) + avt(ji+1,jj+1,jk) ) ) * eumean(ji,jj,jk)
624
625                  avmv(ji,jj,jk) = (      avt(ji  ,jj,jk) + avt(ji  ,jj+1,jk)   &
626                     &              +.5*( avt(ji-1,jj,jk) + avt(ji-1,jj+1,jk)   &
627                     &                   +avt(ji+1,jj,jk) + avt(ji+1,jj+1,jk) ) ) * evmean(ji,jj,jk)
628# else
629                  avmu(ji,jj,jk) = (   avt  (ji,jj  ,jk) + avt  (ji+1,jj  ,jk)   &
630                     &           +.5*( avt  (ji,jj-1,jk) + avt  (ji+1,jj-1,jk)   &
631                     &                +avt  (ji,jj+1,jk) + avt  (ji+1,jj+1,jk) ) ) * umask(ji,jj,jk)   &
632                     &    / MAX( 1.,   tmask(ji,jj  ,jk) + tmask(ji+1,jj  ,jk)   &
633                     &           +.5*( tmask(ji,jj-1,jk) + tmask(ji+1,jj-1,jk)   &
634                     &                +tmask(ji,jj+1,jk) + tmask(ji+1,jj+1,jk) )  )
635
636                  avmv(ji,jj,jk) = (   avt  (ji  ,jj,jk) + avt  (ji  ,jj+1,jk)   &
637                     &           +.5*( avt  (ji-1,jj,jk) + avt  (ji-1,jj+1,jk)   &
638                     &                +avt  (ji+1,jj,jk) + avt  (ji+1,jj+1,jk) ) ) * vmask(ji,jj,jk)   &
639                     &   /  MAX( 1.,   tmask(ji  ,jj,jk) + tmask(ji  ,jj+1,jk)   &
640                     &           +.5*( tmask(ji-1,jj,jk) + tmask(ji-1,jj+1,jk)   &
641                     &                +tmask(ji+1,jj,jk) + tmask(ji+1,jj+1,jk) )  )
642# endif
643               END DO
644            END DO
645         END DO
646         CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
647         !
648         !                                                ! Vertical eddy diffusivity             (1 2 1)
649         !                                                !                            avt = 1/16 (2 4 2)
650         !                                                !                                       (1 2 1)
651         DO jk = 2, jpkm1         
652            DO jj = 2, jpjm1
653               DO ji = fs_2, fs_jpim1   ! vector opt.
654#  if defined key_vectopt_memory
655                  avt(ji,jj,jk) = ( avmu(ji,jj,jk) + avmu(ji-1,jj  ,jk)                          &
656                     &            + avmv(ji,jj,jk) + avmv(ji  ,jj-1,jk)  ) * etmean(ji,jj,jk)
657#  else
658                  avt(ji,jj,jk) = ( avmu (ji,jj,jk) + avmu (ji-1,jj  ,jk)                        &
659                     &            + avmv (ji,jj,jk) + avmv (ji  ,jj-1,jk)  ) * tmask(ji,jj,jk)   &
660                     &  / MAX( 1.,  umask(ji,jj,jk) + umask(ji-1,jj  ,jk)                        &
661                     &            + vmask(ji,jj,jk) + vmask(ji  ,jj-1,jk)  )
662#  endif
663               END DO
664            END DO
665         END DO
666         !
667      END SELECT
668      !
669      IF( nn_pdl == 1 ) THEN                            ! Ponderation by the Prandtl number (nn_pdl=1)
670         DO jk = 2, jpkm1 
671            DO jj = 2, jpjm1
672               DO ji = fs_2, fs_jpim1   ! vector opt.
673                  zpdl = zmxld(ji,jj,jk)
674                  avt(ji,jj,jk) = MAX( zpdl * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
675               END DO
676            END DO
677         END DO
678      ENDIF
679      !
680      DO jk = 2, jpkm1                                  ! Minimum value on the eddy viscosity
681         avmu(:,:,jk) = MAX( avmu(:,:,jk), avmb(jk) ) * umask(:,:,jk)
682         avmv(:,:,jk) = MAX( avmv(:,:,jk), avmb(jk) ) * vmask(:,:,jk)
683      END DO
684
685      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
686
687      IF( lrst_oce )   CALL tke_rst( kt, 'WRITE' )      ! write en in restart file
688
689      IF(ln_ctl) THEN
690         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
691         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
692            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
693      ENDIF
694      !
695   END SUBROUTINE zdf_tke_old
696
697
698   SUBROUTINE zdf_tke_init_o
699      !!----------------------------------------------------------------------
700      !!                  ***  ROUTINE zdf_tke_init_o  ***
701      !!                     
702      !! ** Purpose :   Initialization of the vertical eddy diffivity and
703      !!      viscosity when using a tke turbulent closure scheme
704      !!
705      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
706      !!      called at the first timestep (nit000)
707      !!
708      !! ** input   :   Namlist namzdf_tke
709      !!
710      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
711      !!
712      !!----------------------------------------------------------------------
713      USE dynzdf_exp
714      USE trazdf_exp
715      !
716# if defined key_vectopt_memory
717      INTEGER ::   ji, jj, jk   ! dummy loop indices
718# else
719      INTEGER ::           jk   ! dummy loop indices
720# endif
721      !!
722      NAMELIST/namzdf_tke/ ln_rstke, rn_ediff, rn_ediss, rn_ebb  , rn_efave, rn_emin,   &
723         &                 rn_emin0, rn_cri  , nn_itke , nn_mxl  , nn_pdl  , nn_ave ,   &
724         &                 ln_mxl0 , rn_lmin , rn_lmin0, nn_etau,   &
725         &                 nn_htau , rn_efr  , ln_lc   , rn_lc 
726      !!----------------------------------------------------------------------
727
728      ! Read Namelist namzdf_tke : Turbulente Kinetic Energy
729      ! --------------------
730      REWIND ( numnam )
731      READ   ( numnam, namzdf_tke )
732
733      ! Compute boost associated with the Richardson critic
734      !     (control values: rn_cri = 0.3   ==> eboost=1.25 for nn_pdl=1)
735      !     (                rn_cri = 0.222 ==> eboost=1.               )
736      eboost = rn_cri * ( 2. + rn_ediss / rn_ediff ) / 2.
737
738
739
740      ! Parameter control and print
741      ! ---------------------------
742      ! Control print
743      IF(lwp) THEN
744         WRITE(numout,*)
745         WRITE(numout,*) 'zdf_tke_init_o : tke turbulent closure scheme (old scheme)'
746         WRITE(numout,*) '~~~~~~~~~~~~~~'
747         WRITE(numout,*) '          Namelist namzdf_tke : set tke mixing parameters'
748         WRITE(numout,*) '             restart with tke from no tke              ln_rstke = ', ln_rstke
749         WRITE(numout,*) '             coef. to compute avt                      rn_ediff = ', rn_ediff
750         WRITE(numout,*) '             Kolmogoroff dissipation coef.             rn_ediss = ', rn_ediss
751         WRITE(numout,*) '             tke surface input coef.                   rn_ebb   = ', rn_ebb
752         WRITE(numout,*) '             tke diffusion coef.                       rn_efave = ', rn_efave
753         WRITE(numout,*) '             minimum value of tke                      rn_emin  = ', rn_emin
754         WRITE(numout,*) '             surface minimum value of tke              rn_emin0 = ', rn_emin0
755         WRITE(numout,*) '             number of restart iter loops              nn_itke  = ', nn_itke
756         WRITE(numout,*) '             mixing length type                        nn_mxl   = ', nn_mxl
757         WRITE(numout,*) '             prandl number flag                        nn_pdl   = ', nn_pdl
758         WRITE(numout,*) '             horizontal average flag                   nn_ave   = ', nn_ave
759         WRITE(numout,*) '             critic Richardson nb                      rn_cri   = ', rn_cri
760         WRITE(numout,*) '                and its associated coeff.              eboost   = ', eboost
761         WRITE(numout,*) '             surface mixing length = F(stress) or not  ln_mxl0  = ', ln_mxl0
762         WRITE(numout,*) '             surface  mixing length minimum value      rn_lmin0 = ', rn_lmin0
763         WRITE(numout,*) '             interior mixing length minimum value      rn_lmin0 = ', rn_lmin0
764         WRITE(numout,*) '             test param. to add tke induced by wind    nn_etau  = ', nn_etau
765         WRITE(numout,*) '             flag for computation of exp. tke profile  nn_htau  = ', nn_htau
766         WRITE(numout,*) '             % of rn_emin0 which pene. the thermocline rn_efr   = ', rn_efr
767         WRITE(numout,*) '             flag to take into acc.  Langmuir circ.    ln_lc    = ', ln_lc
768         WRITE(numout,*) '             coef to compute verticla velocity of LC   rn_lc    = ', rn_lc
769         WRITE(numout,*)
770      ENDIF
771
772      ! Check of some namelist values
773      IF( nn_mxl  < 0    .OR. nn_mxl  > 3   )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
774      IF( nn_pdl  < 0    .OR. nn_pdl  > 1   )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
775      IF( nn_ave  < 0    .OR. nn_ave  > 1   )   CALL ctl_stop( 'bad flag: nn_ave is  0 or 1    ' )
776      IF( nn_htau < 0    .OR. nn_htau > 1   )   CALL ctl_stop( 'bad flag: nn_htau is 0 or 1    ' )
777      IF( rn_lc   < 0.15 .OR. rn_lc   > 0.2 )   CALL ctl_stop( 'bad value: rn_lc must be between 0.15 and 0.2 ' )
778
779      IF( nn_etau == 2 )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
780
781      ! Horizontal average : initialization of weighting arrays
782      ! -------------------
783      !
784      SELECT CASE ( nn_ave )
785      ! Notice: mean arrays have not to by defined at local domain boundaries due to the vertical nature of TKE
786      !
787      CASE ( 0 )                ! no horizontal average
788         IF(lwp) WRITE(numout,*) '          no horizontal average on avt, avmu, avmv'
789         IF(lwp) WRITE(numout,*) '          only in very high horizontal resolution !'
790# if defined key_vectopt_memory
791         ! caution vectopt_memory change the solution (last digit of the solver stat)
792         ! weighting mean arrays etmean, eumean and evmean
793         !           ( 1  1 )                                          ( 1 )
794         ! avt = 1/4 ( 1  1 )     avmu = 1/2 ( 1  1 )       avmv=  1/2 ( 1 )
795         !                         
796         etmean(:,:,:) = 0.e0
797         eumean(:,:,:) = 0.e0
798         evmean(:,:,:) = 0.e0
799         !
800         DO jk = 1, jpkm1
801            DO jj = 2, jpjm1
802               DO ji = fs_2, fs_jpim1   ! vector opt.
803                  etmean(ji,jj,jk) = tmask(ji,jj,jk) / MAX( 1.,  umask(ji,jj,jk) + umask(ji-1,jj  ,jk)   &
804                  &                                            + vmask(ji,jj,jk) + vmask(ji  ,jj-1,jk) )
805                  eumean(ji,jj,jk) = umask(ji,jj,jk) / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji+1,jj  ,jk) )
806                  evmean(ji,jj,jk) = vmask(ji,jj,jk) / MAX( 1.,  tmask(ji,jj,jk) + tmask(ji  ,jj+1,jk) )
807               END DO
808            END DO
809         END DO
810# endif
811         !
812      CASE ( 1 )                ! horizontal average
813         IF(lwp) WRITE(numout,*) '          horizontal average on avt, avmu, avmv'
814# if defined key_vectopt_memory
815         ! caution vectopt_memory change the solution (last digit of the solver stat)
816         ! weighting mean arrays etmean, eumean and evmean
817         !           ( 1  1 )              ( 1/2  1/2 )             ( 1/2  1  1/2 )
818         ! avt = 1/4 ( 1  1 )   avmu = 1/4 ( 1    1   )   avmv= 1/4 ( 1/2  1  1/2 )
819         !                                 ( 1/2  1/2 )
820         etmean(:,:,:) = 0.e0
821         eumean(:,:,:) = 0.e0
822         evmean(:,:,:) = 0.e0
823         !
824         DO jk = 1, jpkm1
825            DO jj = 2, jpjm1
826               DO ji = fs_2, fs_jpim1   ! vector opt.
827                  etmean(ji,jj,jk) = tmask(ji,jj,jk) / MAX( 1.,   umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk)   &
828                  &                                       +       vmask(ji  ,jj-1,jk) + vmask(ji  ,jj  ,jk) )
829                  eumean(ji,jj,jk) = umask(ji,jj,jk) / MAX( 1.,   tmask(ji  ,jj  ,jk) + tmask(ji+1,jj  ,jk)   &
830                  &                                       +.5 * ( tmask(ji  ,jj-1,jk) + tmask(ji+1,jj-1,jk)   &
831                  &                                             + tmask(ji  ,jj+1,jk) + tmask(ji+1,jj+1,jk) )  )
832                  evmean(ji,jj,jk) = vmask(ji,jj,jk) / MAX( 1.,   tmask(ji  ,jj  ,jk) + tmask(ji  ,jj+1,jk)   &
833                  &                                       +.5 * ( tmask(ji-1,jj  ,jk) + tmask(ji-1,jj+1,jk)   &
834                  &                                             + tmask(ji+1,jj  ,jk) + tmask(ji+1,jj+1,jk) )  )
835               END DO
836            END DO
837         END DO
838# endif
839         !
840      END SELECT
841
842      ! Initialization of vertical eddy coef. to the background value
843      ! -------------------------------------------------------------
844      DO jk = 1, jpk
845         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
846         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
847         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
848      END DO
849
850
851      ! read or initialize turbulent kinetic energy ( en )
852      ! -------------------------------------------------
853      CALL tke_rst( nit000, 'READ' )
854      !
855   END SUBROUTINE zdf_tke_init_o
856
857
858   SUBROUTINE tke_rst( kt, cdrw )
859     !!---------------------------------------------------------------------
860     !!                   ***  ROUTINE ts_rst  ***
861     !!                     
862     !! ** Purpose :   Read or write TKE file (en) in restart file
863     !!
864     !! ** Method  :   use of IOM library
865     !!                if the restart does not contain TKE, en is either
866     !!                set to rn_emin or recomputed (nn_itke/=0)
867     !!----------------------------------------------------------------------
868     INTEGER         , INTENT(in) ::   kt         ! ocean time-step
869     CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
870     !
871     INTEGER ::   jit   ! dummy loop indices
872     !!----------------------------------------------------------------------
873     !
874     IF( TRIM(cdrw) == 'READ' ) THEN
875        IF( ln_rstart ) THEN
876           IF( iom_varid( numror, 'en', ldstop = .FALSE. ) > 0 .AND. .NOT.(ln_rstke) ) THEN
877              CALL iom_get( numror, jpdom_autoglo, 'en', en )
878           ELSE
879              IF( lwp .AND. iom_varid( numror, 'en', ldstop = .FALSE. ) > 0 )   &
880                 &                       WRITE(numout,*) ' ===>>>> : previous run without tke scheme'
881              IF( lwp .AND. ln_rstke )   WRITE(numout,*) ' ===>>>> : We do not use en from the restart file'
882              IF( lwp                )   WRITE(numout,*) ' ===>>>> : en is set by iterative loop'
883              en (:,:,:) = rn_emin * tmask(:,:,:)
884              DO jit = 2, nn_itke + 1
885                 CALL zdf_tke_old( jit )
886              END DO
887           ENDIF
888        ELSE
889           en(:,:,:) = rn_emin * tmask(:,:,:)      ! no restart: en set to its minimum value
890        ENDIF
891     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
892        CALL iom_rstput( kt, nitrst, numrow, 'en', en )
893     ENDIF
894     !
895   END SUBROUTINE tke_rst
896
897#else
898   !!----------------------------------------------------------------------
899   !!   Dummy module :                                        NO TKE scheme
900   !!----------------------------------------------------------------------
901   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke_old = .FALSE.   !: TKE flag
902CONTAINS
903   SUBROUTINE zdf_tke_init_o             ! Dummy routine
904   END SUBROUTINE zdf_tke_init_o
905   SUBROUTINE zdf_tke_old( kt )          ! Dummy routine
906      WRITE(*,*) 'zdf_tke_old: You should not have seen this print! error?', kt
907   END SUBROUTINE zdf_tke_old
908#endif
909
910   !!======================================================================
911END MODULE zdftke_old
Note: See TracBrowser for help on using the repository browser.