New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdftke.F90 in branches/DEV_r2191_3partymerge2010/NEMO/OPA_SRC/ZDF – NEMO

source: branches/DEV_r2191_3partymerge2010/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 2207

Last change on this file since 2207 was 2207, checked in by acc, 13 years ago

#733 DEV_r2191_3partymerge2010. Merged in changes from devukmo2010 branch

File size: 43.8 KB
Line 
1MODULE zdftke
2   !!======================================================================
3   !!                       ***  MODULE  zdftke  ***
4   !! Ocean physics:  vertical mixing coefficient computed from the tke
5   !!                 turbulent closure parameterization
6   !!=====================================================================
7   !! History :  OPA  !  1991-03  (b. blanke)  Original code
8   !!            7.0  !  1991-11  (G. Madec)   bug fix
9   !!            7.1  !  1992-10  (G. Madec)   new mixing length and eav
10   !!            7.2  !  1993-03  (M. Guyon)   symetrical conditions
11   !!            7.3  !  1994-08  (G. Madec, M. Imbard)  nn_pdl flag
12   !!            7.5  !  1996-01  (G. Madec)   s-coordinates
13   !!            8.0  !  1997-07  (G. Madec)   lbc
14   !!            8.1  !  1999-01  (E. Stretta) new option for the mixing length
15   !!  NEMO      1.0  !  2002-06  (G. Madec) add tke_init routine
16   !!             -   !  2004-10  (C. Ethe )  1D configuration
17   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
18   !!            3.0  !  2008-05  (C. Ethe,  G.Madec) : update TKE physics:
19   !!                 !           - tke penetration (wind steering)
20   !!                 !           - suface condition for tke & mixing length
21   !!                 !           - Langmuir cells
22   !!             -   !  2008-05  (J.-M. Molines, G. Madec)  2D form of avtb
23   !!             -   !  2008-06  (G. Madec)  style + DOCTOR name for namelist parameters
24   !!             -   !  2008-12  (G. Reffray) stable discretization of the production term
25   !!            3.2  !  2009-06  (G. Madec, S. Masson) TKE restart compatible with key_cpl
26   !!                 !                                + cleaning of the parameters + bugs correction
27   !!----------------------------------------------------------------------
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), PUBLIC, 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      !
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
248         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
249         !
250         !                        !* total energy produce by LC : cumulative sum over jk
251         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0. ) * fsdepw(:,:,1) * fse3w(:,:,1)
252         DO jk = 2, jpk
253            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0. ) * 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.e0 , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
331                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1.e0 , 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.e0 - 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 throughout the water column
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.e0 - 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)
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.e0 - fr_i(ji,jj) )  * tmask(ji,jj,jk)
401            END DO
402         END DO
403      ELSEIF( nn_etau == 3 ) THEN       !* penetration throughout the water column
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 * SQRT( ztx2 * ztx2 + zty2 * zty2 )      ! module of the mean stress
413                  zdif = taum(ji,jj) - ztau                           ! mean of the module - module of the mean
414                  zdif = rn_sclhft * MAX( 0.e0, zdif + rn_addhft )    ! apply some modifications...
415                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
416                     &                                        * ( 1.e0 - 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 rn_lmin (rn_lmin0) 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!!gm  this should be useless
483         zmxlm(:,:,1) = 0.e0
484!!gm end
485         zraug = vkarmn * 2.e5 / ( rau0 * grav )
486         DO jj = 2, jpjm1
487            DO ji = fs_2, fs_jpim1   ! vector opt.
488               zmxlm(ji,jj,1) = MAX(  rn_lmin0,  zraug * taum(ji,jj)  )
489            END DO
490         END DO
491      ELSE                       ! surface set to the minimum value
492         zmxlm(:,:,1) = rn_lmin0
493      ENDIF
494      zmxlm(:,:,jpk) = rn_lmin   ! bottom set to the interior minium value
495      !
496!CDIR NOVERRCHK
497      DO jk = 2, jpkm1           ! interior value : l=sqrt(2*e/n**2)
498!CDIR NOVERRCHK
499         DO jj = 2, jpjm1
500!CDIR NOVERRCHK
501            DO ji = fs_2, fs_jpim1   ! vector opt.
502               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
503               zmxlm(ji,jj,jk) = MAX(  rn_lmin,  SQRT( 2. * en(ji,jj,jk) / zrn2 )  )
504            END DO
505         END DO
506      END DO
507      !
508!!gm  CAUTION: to be added here the bottom boundary condition on zmxlm
509      !
510      !                     !* Physical limits for the mixing length
511      !
512      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
513      zmxld(:,:,jpk) = rn_lmin        ! bottom  set to the minimum value
514      !
515      SELECT CASE ( nn_mxl )
516      !
517      CASE ( 0 )           ! bounded by the distance to surface and bottom
518         DO jk = 2, jpkm1
519            DO jj = 2, jpjm1
520               DO ji = fs_2, fs_jpim1   ! vector opt.
521                  zemxl = MIN( fsdepw(ji,jj,jk), zmxlm(ji,jj,jk),   &
522                  &            fsdepw(ji,jj,mbathy(ji,jj)) - fsdepw(ji,jj,jk) )
523                  zmxlm(ji,jj,jk) = zemxl
524                  zmxld(ji,jj,jk) = zemxl
525               END DO
526            END DO
527         END DO
528         !
529      CASE ( 1 )           ! bounded by the vertical scale factor
530         DO jk = 2, jpkm1
531            DO jj = 2, jpjm1
532               DO ji = fs_2, fs_jpim1   ! vector opt.
533                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
534                  zmxlm(ji,jj,jk) = zemxl
535                  zmxld(ji,jj,jk) = zemxl
536               END DO
537            END DO
538         END DO
539         !
540      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
541         DO jk = 2, jpkm1         ! from the surface to the bottom :
542            DO jj = 2, jpjm1
543               DO ji = fs_2, fs_jpim1   ! vector opt.
544                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
545               END DO
546            END DO
547         END DO
548         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
549            DO jj = 2, jpjm1
550               DO ji = fs_2, fs_jpim1   ! vector opt.
551                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
552                  zmxlm(ji,jj,jk) = zemxl
553                  zmxld(ji,jj,jk) = zemxl
554               END DO
555            END DO
556         END DO
557         !
558      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
559         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
560            DO jj = 2, jpjm1
561               DO ji = fs_2, fs_jpim1   ! vector opt.
562                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
563               END DO
564            END DO
565         END DO
566         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
567            DO jj = 2, jpjm1
568               DO ji = fs_2, fs_jpim1   ! vector opt.
569                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
570               END DO
571            END DO
572         END DO
573!CDIR NOVERRCHK
574         DO jk = 2, jpkm1
575!CDIR NOVERRCHK
576            DO jj = 2, jpjm1
577!CDIR NOVERRCHK
578               DO ji = fs_2, fs_jpim1   ! vector opt.
579                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
580                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
581                  zmxlm(ji,jj,jk) = zemlm
582                  zmxld(ji,jj,jk) = zemlp
583               END DO
584            END DO
585         END DO
586         !
587      END SELECT
588      !
589# if defined key_c1d
590      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
591      e_mix(:,:,:) = zmxlm(:,:,:)
592# endif
593
594      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
595      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
596      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
597!CDIR NOVERRCHK
598      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
599!CDIR NOVERRCHK
600         DO jj = 2, jpjm1
601!CDIR NOVERRCHK
602            DO ji = fs_2, fs_jpim1   ! vector opt.
603               zsqen = SQRT( en(ji,jj,jk) )
604               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
605               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * tmask(ji,jj,jk)
606               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
607               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
608            END DO
609         END DO
610      END DO
611      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
612      !
613      DO jk = 2, jpkm1            !* vertical eddy viscosity at u- and v-points
614         DO jj = 2, jpjm1
615            DO ji = fs_2, fs_jpim1   ! vector opt.
616               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * umask(ji,jj,jk)
617               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * vmask(ji,jj,jk)
618            END DO
619         END DO
620      END DO
621      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
622      !
623      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
624         DO jk = 2, jpkm1
625            DO jj = 2, jpjm1
626               DO ji = fs_2, fs_jpim1   ! vector opt.
627                  zcoef = 0.5 / ( fse3w(ji,jj,jk) * fse3w(ji,jj,jk) )
628                  !                                          ! shear
629                  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) )   &
630                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
631                  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) )   &
632                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
633                  !                                          ! local Richardson number
634                  zri   = MAX( rn2b(ji,jj,jk), 0. ) * avm(ji,jj,jk) / (zcoef * (zdku + zdkv + rn_bshear ) )
635                  zpdlr = MAX(  0.1,  0.2 / MAX( 0.2 , zri )  )
636!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
637!!gm              zpdlr = MAX(  0.1,  ri_crit / MAX( ri_crit , zri )  )
638                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * tmask(ji,jj,jk)
639# if defined key_c1d
640                  e_pdl(ji,jj,jk) = zpdlr * tmask(ji,jj,jk)    ! c1d configuration : save masked Prandlt number
641                  e_ric(ji,jj,jk) = zri * tmask(ji,jj,jk)                            ! c1d config. : save Ri
642# endif
643              END DO
644            END DO
645         END DO
646      ENDIF
647      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
648
649      IF(ln_ctl) THEN
650         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
651         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
652            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
653      ENDIF
654      !
655   END SUBROUTINE tke_avn
656
657
658   SUBROUTINE tke_init
659      !!----------------------------------------------------------------------
660      !!                  ***  ROUTINE tke_init  ***
661      !!                     
662      !! ** Purpose :   Initialization of the vertical eddy diffivity and
663      !!              viscosity when using a tke turbulent closure scheme
664      !!
665      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
666      !!              called at the first timestep (nit000)
667      !!
668      !! ** input   :   Namlist namzdf_tke
669      !!
670      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
671      !!----------------------------------------------------------------------
672      INTEGER ::   ji, jj, jk   ! dummy loop indices
673      !!
674      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb   , rn_emin  ,   &
675         &                 rn_emin0, rn_bshear, nn_mxl   , ln_mxl0  ,   &
676         &                 rn_lmin , rn_lmin0 , nn_pdl   , nn_etau  ,   &
677         &                 nn_htau , rn_efr   , rn_addhft, rn_sclhft,   &
678         &                 ln_lc   , rn_lc 
679      !!----------------------------------------------------------------------
680
681      REWIND ( numnam )               !* Read Namelist namzdf_tke : Turbulente Kinetic Energy
682      READ   ( numnam, namzdf_tke )
683     
684      ri_cri = 2. / ( 2. + rn_ediss / rn_ediff )      ! resulting critical Richardson number
685
686      IF(lwp) THEN                    !* Control print
687         WRITE(numout,*)
688         WRITE(numout,*) 'zdf_tke : tke turbulent closure scheme - initialisation'
689         WRITE(numout,*) '~~~~~~~~'
690         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
691         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
692         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
693         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
694         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
695         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
696         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
697         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
698         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
699         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
700         WRITE(numout,*) '      surface  mixing length minimum value        rn_lmin0  = ', rn_lmin0
701         WRITE(numout,*) '      interior mixing length minimum value        rn_lmin0  = ', rn_lmin0
702         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
703         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
704         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
705         WRITE(numout,*) '      add offset   applied to HF tau              rn_addhft = ', rn_addhft
706         WRITE(numout,*) '      scale factor applied to HF tau              rn_sclhft = ', rn_sclhft
707         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
708         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
709         WRITE(numout,*)
710         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
711      ENDIF
712
713      !                               !* Check of some namelist values
714      IF( nn_mxl  < 0     .OR.  nn_mxl  > 3   )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
715      IF( nn_pdl  < 0     .OR.  nn_pdl  > 1   )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
716      IF( nn_htau < 0     .OR.  nn_htau > 2   )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
717      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 ' )
718      IF( rn_sclhft == 0. .AND. nn_etau == 3  )   CALL ctl_stop( 'force null HF tau to penetrate the thermocline...' )
719      IF( .NOT. lhftau    .AND. nn_etau == 3  )   CALL ctl_stop( 'bad flag: nn_etau == 3 must be used with HF tau' )
720
721      IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
722
723      !                               !* depth of penetration of surface tke
724      IF( nn_etau /= 0 ) THEN     
725         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
726         CASE( 0 )                                    ! constant depth penetration (here 10 meters)
727            htau(:,:) = 10.e0
728         CASE( 1 )                                    ! F(latitude) : 0.5m to 30m at high lat.
729            DO jj = 1, jpj
730               DO ji = 1, jpi 
731                  htau(ji,jj) = MAX( 0.5, 3./4. * MIN( 40., 60.*ABS( SIN( rpi/180. * gphit(ji,jj) ) ) )   )
732               END DO
733            END DO
734         CASE( 2 )                                    ! constant depth penetration (here 30 meters)
735            htau(:,:) = 30.e0
736         END SELECT
737      ENDIF
738
739      !                               !* set vertical eddy coef. to the background value
740      DO jk = 1, jpk
741         avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
742         avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
743         avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
744         avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
745      END DO
746      dissl(:,:,:) = 1.e-12
747      !                               !* read or initialize all required files
748      CALL tke_rst( nit000, 'READ' )
749      !
750   END SUBROUTINE tke_init
751
752
753   SUBROUTINE tke_rst( kt, cdrw )
754     !!---------------------------------------------------------------------
755     !!                   ***  ROUTINE tke_rst  ***
756     !!                     
757     !! ** Purpose :   Read or write TKE file (en) in restart file
758     !!
759     !! ** Method  :   use of IOM library
760     !!                if the restart does not contain TKE, en is either
761     !!                set to rn_emin or recomputed
762     !!----------------------------------------------------------------------
763     INTEGER         , INTENT(in) ::   kt         ! ocean time-step
764     CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
765     !
766     INTEGER ::   jit, jk   ! dummy loop indices
767     INTEGER ::   id1, id2, id3, id4, id5, id6
768     !!----------------------------------------------------------------------
769     !
770     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
771        !                                   ! ---------------
772        IF( ln_rstart ) THEN                   !* Read the restart file
773           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
774           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
775           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
776           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
777           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
778           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
779           !
780           IF( id1 > 0 ) THEN                       ! 'en' exists
781              CALL iom_get( numror, jpdom_autoglo, 'en', en )
782              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
783                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt   )
784                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm   )
785                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu  )
786                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv  )
787                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl )
788              ELSE                                                 ! one at least array is missing
789                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
790              ENDIF
791           ELSE                                     ! No TKE array found: initialisation
792              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
793              en (:,:,:) = rn_emin * tmask(:,:,:)
794              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
795              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
796           ENDIF
797        ELSE                                   !* Start from rest
798           en(:,:,:) = rn_emin * tmask(:,:,:)
799           DO jk = 1, jpk                           ! set the Kz to the background value
800              avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
801              avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
802              avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
803              avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
804           END DO
805        ENDIF
806        !
807     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
808        !                                   ! -------------------
809        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
810        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en    )
811        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt   )
812        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm   )
813        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu  )
814        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv  )
815        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl )
816        !
817     ENDIF
818     !
819   END SUBROUTINE tke_rst
820
821#else
822   !!----------------------------------------------------------------------
823   !!   Dummy module :                                        NO TKE scheme
824   !!----------------------------------------------------------------------
825   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
826CONTAINS
827   SUBROUTINE zdf_tke( kt )          ! Empty routine
828      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
829   END SUBROUTINE zdf_tke
830   SUBROUTINE tke_rst( kt, cdrw )
831     CHARACTER(len=*) ::   cdrw
832     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
833   END SUBROUTINE tke_rst
834#endif
835
836   !!======================================================================
837END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.