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/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 7990

Last change on this file since 7990 was 7990, checked in by gm, 7 years ago

#1880 (HPC-09): OPA remove avmu, avmv from zdf modules + move CALL tke(gls)_rst & gls_rst in zdftke(gls) + rename zdftmx and zdfqiao

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