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

source: branches/DEV_r2006_merge_TRA_TRC/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 2083

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

Reorganisation of the initialisation phase, see ticket:695

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