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/2011/dev_r2787_LOCEAN3_TRA_TRP/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2011/dev_r2787_LOCEAN3_TRA_TRP/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 2789

Last change on this file since 2789 was 2789, checked in by cetlod, 13 years ago

Implementation of the merge of TRA/TRP : first guess, see ticket #842

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