source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 2375

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

v3.3beta: zdftke namelist simplification associated with an update of DOC

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