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/2015/dev_r5072_UKMO2_OBS_simplification/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2015/dev_r5072_UKMO2_OBS_simplification/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 5682

Last change on this file since 5682 was 5682, checked in by mattmartin, 9 years ago

OBS simplification changes committed to branch after running SETTE tests to make sure we get the same results as the trunk for ORCA2_LIM_OBS.

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