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

source: trunk/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 1708

Last change on this file since 1708 was 1705, checked in by smasson, 15 years ago

impact of HF winds in TKE, see ticket:585

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