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/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 7924

Last change on this file since 7924 was 7924, checked in by andmirek, 7 years ago

first commit with XIOS restart read functionality

File size: 51.4 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   USE iom_def, ONLY : lxios_read
60
61
62
63   IMPLICIT NONE
64   PRIVATE
65
66   PUBLIC   zdf_tke        ! routine called in step module
67   PUBLIC   zdf_tke_init   ! routine called in opa module
68   PUBLIC   tke_rst        ! routine called in step module
69
70   LOGICAL , PUBLIC, PARAMETER ::   lk_zdftke = .TRUE.  !: TKE vertical mixing flag
71
72   !                      !!** Namelist  namzdf_tke  **
73   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
74   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
75   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
76   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
77   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
78   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
79   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
80   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
81   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
82   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
83   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
84   INTEGER  ::   nn_htau   ! type of tke profile of penetration (=0/1)
85   REAL(wp) ::   rn_efr    ! fraction of TKE surface value which penetrates in the ocean
86   REAL(wp) ::   rn_c      ! fraction of TKE added within the mixed layer by nn_etau
87   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
88   REAL(wp) ::   rn_lc     ! coef to compute vertical velocity of Langmuir cells
89
90   REAL(wp) ::   ri_cri    ! critic Richardson number (deduced from rn_ediff and rn_ediss values)
91   REAL(wp) ::   rmxl_min  ! minimum mixing length value (deduced from rn_ediff and rn_emin values)  [m]
92   REAL(wp) ::   rhtau                     ! coefficient to relate MLD to htau when nn_htau == 2
93   REAL(wp) ::   rhftau_add = 1.e-3_wp     ! add offset   applied to HF part of taum  (nn_etau=3)
94   REAL(wp) ::   rhftau_scl = 1.0_wp       ! scale factor applied to HF part of taum  (nn_etau=3)
95
96   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau           ! depth of tke penetration (nn_htau)
97   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_niw          !: TKE budget- near-inertial waves term
98   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   efr            ! surface boundary condition for nn_etau = 4
99   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl          ! now mixing lenght of dissipation
100#if defined key_c1d
101   !                                                                        !!** 1D cfg only  **   ('key_c1d')
102   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_dis, e_mix   !: dissipation and mixing turbulent lengh scales
103   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   e_pdl, e_ric   !: prandl and local Richardson numbers
104#endif
105
106   !! * Substitutions
107#  include "domzgr_substitute.h90"
108#  include "vectopt_loop_substitute.h90"
109   !!----------------------------------------------------------------------
110   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
111   !! $Id$
112   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
113   !!----------------------------------------------------------------------
114CONTAINS
115
116   INTEGER FUNCTION zdf_tke_alloc()
117      !!----------------------------------------------------------------------
118      !!                ***  FUNCTION zdf_tke_alloc  ***
119      !!----------------------------------------------------------------------
120      ALLOCATE(                                                                    &
121         &      efr  (jpi,jpj)     , e_niw(jpi,jpj,jpk) ,                         &     
122#if defined key_c1d
123         &      e_dis(jpi,jpj,jpk) , e_mix(jpi,jpj,jpk) ,                          &
124         &      e_pdl(jpi,jpj,jpk) , e_ric(jpi,jpj,jpk) ,                          &
125#endif
126         &      htau  (jpi,jpj)    , dissl(jpi,jpj,jpk) , STAT= zdf_tke_alloc      )
127         !
128      IF( lk_mpp             )   CALL mpp_sum ( zdf_tke_alloc )
129      IF( zdf_tke_alloc /= 0 )   CALL ctl_warn('zdf_tke_alloc: failed to allocate arrays')
130      !
131   END FUNCTION zdf_tke_alloc
132
133
134   SUBROUTINE zdf_tke( kt )
135      !!----------------------------------------------------------------------
136      !!                   ***  ROUTINE zdf_tke  ***
137      !!
138      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
139      !!              coefficients using a turbulent closure scheme (TKE).
140      !!
141      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
142      !!              is computed from a prognostic equation :
143      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
144      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
145      !!                  + avt N^2                      ! stratif. destruc.
146      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
147      !!      with the boundary conditions:
148      !!         surface: en = max( rn_emin0, rn_ebb * taum )
149      !!         bottom : en = rn_emin
150      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
151      !!
152      !!        The now Turbulent kinetic energy is computed using the following
153      !!      time stepping: implicit for vertical diffusion term, linearized semi
154      !!      implicit for kolmogoroff dissipation term, and explicit forward for
155      !!      both buoyancy and shear production terms. Therefore a tridiagonal
156      !!      linear system is solved. Note that buoyancy and shear terms are
157      !!      discretized in a energy conserving form (Bruchard 2002).
158      !!
159      !!        The dissipative and mixing length scale are computed from en and
160      !!      the stratification (see tke_avn)
161      !!
162      !!        The now vertical eddy vicosity and diffusivity coefficients are
163      !!      given by:
164      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
165      !!              avt = max( avmb, pdl * avm                 ) 
166      !!              eav = max( avmb, avm )
167      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
168      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
169      !!
170      !! ** Action  :   compute en (now turbulent kinetic energy)
171      !!                update avt, avmu, avmv (before vertical eddy coef.)
172      !!
173      !! References : Gaspar et al., JGR, 1990,
174      !!              Blanke and Delecluse, JPO, 1991
175      !!              Mellor and Blumberg, JPO 2004
176      !!              Axell, JGR, 2002
177      !!              Bruchard OM 2002
178      !!----------------------------------------------------------------------
179      INTEGER, INTENT(in) ::   kt   ! ocean time step
180      !!----------------------------------------------------------------------
181      !
182      IF( kt /= nit000 ) THEN   ! restore before value to compute tke
183         avt (:,:,:) = avt_k (:,:,:) 
184         avm (:,:,:) = avm_k (:,:,:) 
185         avmu(:,:,:) = avmu_k(:,:,:) 
186         avmv(:,:,:) = avmv_k(:,:,:) 
187      ENDIF 
188      !
189      CALL tke_tke      ! now tke (en)
190      !
191      CALL tke_avn      ! now avt, avm, avmu, avmv
192      !
193      avt_k (:,:,:) = avt (:,:,:) 
194      avm_k (:,:,:) = avm (:,:,:) 
195      avmu_k(:,:,:) = avmu(:,:,:) 
196      avmv_k(:,:,:) = avmv(:,:,:) 
197      !
198#if defined key_agrif
199      ! Update child grid f => parent grid
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
238      !!--------------------------------------------------------------------
239      !
240      IF( nn_timing == 1 )  CALL timing_start('tke_tke')
241      !
242      CALL wrk_alloc( jpi,jpj, imlc )    ! integer
243      CALL wrk_alloc( jpi,jpj, zhlc ) 
244      CALL wrk_alloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
245      !
246      zbbrau = rn_ebb / rau0       ! Local constant initialisation
247      zfact1 = -.5_wp * rdt 
248      zfact2 = 1.5_wp * rdt * rn_ediss
249      zfact3 = 0.5_wp       * rn_ediss
250      !
251      !
252      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
253      !                     !  Surface boundary condition on tke
254      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
255      IF ( ln_isfcav ) THEN
256         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin
257            DO ji = fs_2, fs_jpim1   ! vector opt.
258               en(ji,jj,mikt(ji,jj))=rn_emin * tmask(ji,jj,1)
259            END DO
260         END DO
261      END IF
262      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
263         DO ji = fs_2, fs_jpim1   ! vector opt.
264            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
265         END DO
266      END DO
267     
268!!bfr   - start commented area
269      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
270      !                     !  Bottom boundary condition on tke
271      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
272      !
273      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
274      ! Tests to date have found the bottom boundary condition on tke to have very little effect.
275      ! The condition is coded here for completion but commented out until there is proof that the
276      ! computational cost is justified
277      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
278      !                     en(bot)   = (rn_ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
279!CDIR NOVERRCHK
280!!    DO jj = 2, jpjm1
281!CDIR NOVERRCHK
282!!       DO ji = fs_2, fs_jpim1   ! vector opt.
283!!          ztx2 = bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) + &
284!!                 bfrua(ji  ,jj) * ub(ji  ,jj,mbku(ji  ,jj) )
285!!          zty2 = bfrva(ji,jj  ) * vb(ji,jj  ,mbkv(ji,jj  )) + &
286!!                 bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1) )
287!!          zebot = 0.001875_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 )   !  where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.
288!!          en (ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * tmask(ji,jj,1)
289!!       END DO
290!!    END DO
291!!bfr   - end commented area
292      !
293      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
294      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke       (Axell JGR 2002)
295         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
296         !
297         !                        !* total energy produce by LC : cumulative sum over jk
298         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * fsdepw(:,:,1) * fse3w(:,:,1)
299         DO jk = 2, jpk
300            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * fsdepw(:,:,jk) * fse3w(:,:,jk)
301         END DO
302         !                        !* finite Langmuir Circulation depth
303         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
304         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
305         DO jk = jpkm1, 2, -1
306            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
307               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
308                  zus  = zcof * taum(ji,jj)
309                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
310               END DO
311            END DO
312         END DO
313         !                               ! finite LC depth
314         DO jj = 1, jpj 
315            DO ji = 1, jpi
316               zhlc(ji,jj) = fsdepw(ji,jj,imlc(ji,jj))
317            END DO
318         END DO
319         zcof = 0.016 / SQRT( zrhoa * zcdrag )
320!CDIR NOVERRCHK
321         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
322!CDIR NOVERRCHK
323            DO jj = 2, jpjm1
324!CDIR NOVERRCHK
325               DO ji = fs_2, fs_jpim1   ! vector opt.
326                  zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
327                  !                                           ! vertical velocity due to LC
328                  zind = 0.5 - SIGN( 0.5, fsdepw(ji,jj,jk) - zhlc(ji,jj) )
329                  zwlc = zind * rn_lc * zus * SIN( rpi * fsdepw(ji,jj,jk) / zhlc(ji,jj) )
330                  !                                           ! TKE Langmuir circulation source term
331                  en(ji,jj,jk) = en(ji,jj,jk) + rdt * ( 1._wp - fr_i(ji,jj) ) * ( zwlc * zwlc * zwlc ) /   &
332                     &   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, jpj                 ! here avmu, avmv used as workspace
348            DO ji = 1, jpi
349               avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) )   &
350                  &                            * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) )   & 
351                  &                            / (  fse3uw_n(ji,jj,jk)               &
352                  &                              *  fse3uw_b(ji,jj,jk)  )
353               avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) )   &
354                  &                            * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) )   &
355                  &                            / (  fse3vw_n(ji,jj,jk)               &
356                  &                              *  fse3vw_b(ji,jj,jk)  )
357            END DO
358         END DO
359      END DO
360      !
361      DO jk = 2, jpkm1           !* Matrix and right hand side in en
362         DO jj = 2, jpjm1
363            DO ji = fs_2, fs_jpim1   ! vector opt.
364               zcof   = zfact1 * tmask(ji,jj,jk)
365# if defined key_zdftmx_new
366               ! key_zdftmx_new: New internal wave-driven param: set a minimum value for Kz on TKE (ensure numerical stability)
367               zzd_up = zcof * ( MAX( avm(ji,jj,jk+1) + avm(ji,jj,jk), 2.e-5_wp ) )   &  ! upper diagonal
368                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
369               zzd_lw = zcof * ( MAX( avm(ji,jj,jk) + avm(ji,jj,jk-1), 2.e-5_wp ) )   &  ! lower diagonal
370                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
371# else
372               zzd_up = zcof * ( avm  (ji,jj,jk+1) + avm  (ji,jj,jk  ) )   &  ! upper diagonal
373                  &          / ( fse3t(ji,jj,jk  ) * fse3w(ji,jj,jk  ) )
374               zzd_lw = zcof * ( avm  (ji,jj,jk  ) + avm  (ji,jj,jk-1) )   &  ! lower diagonal
375                  &          / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk  ) )
376# endif
377                  !                                                           ! shear prod. at w-point weightened by mask
378               zesh2  =  ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1._wp , umask(ji-1,jj,jk) + umask(ji,jj,jk) )   &
379                  &    + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1._wp , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )   
380                  !
381               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
382               zd_lw(ji,jj,jk) = zzd_lw
383               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * tmask(ji,jj,jk)
384               !
385               !                                   ! right hand side in en
386               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  zesh2  -   avt(ji,jj,jk) * rn2(ji,jj,jk)    &
387                  &                                 + zfact3 * dissl(ji,jj,jk) * en (ji,jj,jk)  ) &
388                  &                                 * wmask(ji,jj,jk)
389            END DO
390         END DO
391      END DO
392      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
393      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
394         DO jj = 2, jpjm1
395            DO ji = fs_2, fs_jpim1    ! vector opt.
396               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
397            END DO
398         END DO
399      END DO
400      !
401      ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
402      DO jj = 2, jpjm1
403         DO ji = fs_2, fs_jpim1   ! vector opt.
404            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
405         END DO
406      END DO
407      DO jk = 3, jpkm1
408         DO jj = 2, jpjm1
409            DO ji = fs_2, fs_jpim1    ! vector opt.
410               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)
411            END DO
412         END DO
413      END DO
414      !
415      ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
416      DO jj = 2, jpjm1
417         DO ji = fs_2, fs_jpim1   ! vector opt.
418            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
419         END DO
420      END DO
421      DO jk = jpk-2, 2, -1
422         DO jj = 2, jpjm1
423            DO ji = fs_2, fs_jpim1    ! vector opt.
424               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
425            END DO
426         END DO
427      END DO
428      DO jk = 2, jpkm1                             ! set the minimum value of tke
429         DO jj = 2, jpjm1
430            DO ji = fs_2, fs_jpim1   ! vector opt.
431               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
432            END DO
433         END DO
434      END DO
435
436      !                                 ! Save TKE prior to nn_etau addition 
437      e_niw(:,:,:) = en(:,:,:) 
438     
439      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
440      !                            !  TKE due to surface and internal wave breaking
441      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
442      IF( nn_htau == 2 ) THEN           !* mixed-layer depth dependant length scale
443         DO jj = 2, jpjm1
444            DO ji = fs_2, fs_jpim1   ! vector opt.
445               htau(ji,jj) = rhtau * hmlp(ji,jj)
446            END DO
447         END DO
448      ENDIF
449#if defined key_iomput
450      !
451      CALL iom_put( "htau", htau(:,:) )  ! Check htau (even if constant in time)
452#endif
453      !
454      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
455         DO jk = 2, jpkm1
456            DO jj = 2, jpjm1
457               DO ji = fs_2, fs_jpim1   ! vector opt.
458                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
459                     &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
460               END DO
461            END DO
462         END DO
463      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
464         DO jj = 2, jpjm1
465            DO ji = fs_2, fs_jpim1   ! vector opt.
466               jk = nmln(ji,jj)
467               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
468                  &                                 * ( 1._wp - fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
469            END DO
470         END DO
471      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
472!CDIR NOVERRCHK
473         DO jk = 2, jpkm1
474!CDIR NOVERRCHK
475            DO jj = 2, jpjm1
476!CDIR NOVERRCHK
477               DO ji = fs_2, fs_jpim1   ! vector opt.
478                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
479                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
480                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
481                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
482                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
483                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
484                     &                        * ( 1._wp - fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
485               END DO
486            END DO
487         END DO
488      ELSEIF( nn_etau == 4 ) THEN       !* column integral independant of htau (rn_efr must be scaled up)
489         IF( nn_htau == 2 ) THEN        ! efr dependant on time-varying htau
490            DO jj = 2, jpjm1
491               DO ji = fs_2, fs_jpim1   ! vector opt.
492                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
493               END DO
494            END DO
495         ENDIF
496         DO jk = 2, jpkm1
497            DO jj = 2, jpjm1
498               DO ji = fs_2, fs_jpim1   ! vector opt.
499                  en(ji,jj,jk) = en(ji,jj,jk) + efr(ji,jj) * en(ji,jj,1) * EXP( -fsdepw(ji,jj,jk) / htau(ji,jj) )   &
500                     &                                                   * ( 1._wp - fr_i(ji,jj) )  * tmask(ji,jj,jk)
501               END DO
502            END DO
503         END DO
504      ENDIF
505      CALL lbc_lnk( en, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
506      !
507      DO jk = 2, jpkm1                             ! TKE budget: near-inertial waves term 
508         DO jj = 2, jpjm1 
509            DO ji = fs_2, fs_jpim1   ! vector opt. 
510               e_niw(ji,jj,jk) = en(ji,jj,jk) - e_niw(ji,jj,jk) 
511            END DO 
512         END DO 
513      END DO 
514     
515      CALL lbc_lnk( e_niw, 'W', 1. ) 
516      !
517      CALL wrk_dealloc( jpi,jpj, imlc )    ! integer
518      CALL wrk_dealloc( jpi,jpj, zhlc ) 
519      CALL wrk_dealloc( jpi,jpj,jpk, zpelc, zdiag, zd_up, zd_lw ) 
520      !
521      IF( nn_timing == 1 )  CALL timing_stop('tke_tke')
522      !
523   END SUBROUTINE tke_tke
524
525
526   SUBROUTINE tke_avn
527      !!----------------------------------------------------------------------
528      !!                   ***  ROUTINE tke_avn  ***
529      !!
530      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
531      !!
532      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
533      !!              the tke_tke routine). First, the now mixing lenth is
534      !!      computed from en and the strafification (N^2), then the mixings
535      !!      coefficients are computed.
536      !!              - Mixing length : a first evaluation of the mixing lengh
537      !!      scales is:
538      !!                      mxl = sqrt(2*en) / N 
539      !!      where N is the brunt-vaisala frequency, with a minimum value set
540      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
541      !!        The mixing and dissipative length scale are bound as follow :
542      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
543      !!                        zmxld = zmxlm = mxl
544      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
545      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
546      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
547      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
548      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
549      !!                    the surface to obtain ldown. the resulting length
550      !!                    scales are:
551      !!                        zmxld = sqrt( lup * ldown )
552      !!                        zmxlm = min ( lup , ldown )
553      !!              - Vertical eddy viscosity and diffusivity:
554      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
555      !!                      avt = max( avmb, pdlr * avm ) 
556      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
557      !!
558      !! ** Action  : - avt : now vertical eddy diffusivity (w-point)
559      !!              - avmu, avmv : now vertical eddy viscosity at uw- and vw-points
560      !!----------------------------------------------------------------------
561      INTEGER  ::   ji, jj, jk   ! dummy loop indices
562      REAL(wp) ::   zrn2, zraug, zcoef, zav     ! local scalars
563      REAL(wp) ::   zdku, zpdlr, zri, zsqen     !   -      -
564      REAL(wp) ::   zdkv, zemxl, zemlm, zemlp   !   -      -
565      REAL(wp), POINTER, DIMENSION(:,:,:) :: zmpdl, zmxlm, zmxld
566      !!--------------------------------------------------------------------
567      !
568      IF( nn_timing == 1 )  CALL timing_start('tke_avn')
569
570      CALL wrk_alloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
571
572      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
573      !                     !  Mixing length
574      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
575      !
576      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
577      !
578      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
579      zmxlm(:,:,:)  = rmxl_min   
580      zmxld(:,:,:)  = rmxl_min
581      !
582      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
583         DO jj = 2, jpjm1
584            DO ji = fs_2, fs_jpim1
585               zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
586               zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
587            END DO
588         END DO
589      ELSE
590         zmxlm(:,:,1) = rn_mxl0
591      ENDIF
592      !
593!CDIR NOVERRCHK
594      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
595!CDIR NOVERRCHK
596         DO jj = 2, jpjm1
597!CDIR NOVERRCHK
598            DO ji = fs_2, fs_jpim1   ! vector opt.
599               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
600               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 ) )
601            END DO
602         END DO
603      END DO
604      !
605      !                     !* Physical limits for the mixing length
606      !
607      zmxld(:,:,1  ) = zmxlm(:,:,1)   ! surface set to the minimum value
608      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
609      !
610      SELECT CASE ( nn_mxl )
611      !
612      ! where wmask = 0 set zmxlm == fse3w
613      CASE ( 0 )           ! bounded by the distance to surface and bottom
614         DO jk = 2, jpkm1
615            DO jj = 2, jpjm1
616               DO ji = fs_2, fs_jpim1   ! vector opt.
617                  zemxl = MIN( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
618                  &            fsdepw(ji,jj,mbkt(ji,jj)+1) - fsdepw(ji,jj,jk) )
619                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
620                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
621                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN(zmxlm(ji,jj,jk),fse3w(ji,jj,jk)) * (1 - wmask(ji,jj,jk))
622               END DO
623            END DO
624         END DO
625         !
626      CASE ( 1 )           ! bounded by the vertical scale factor
627         DO jk = 2, jpkm1
628            DO jj = 2, jpjm1
629               DO ji = fs_2, fs_jpim1   ! vector opt.
630                  zemxl = MIN( fse3w(ji,jj,jk), zmxlm(ji,jj,jk) )
631                  zmxlm(ji,jj,jk) = zemxl
632                  zmxld(ji,jj,jk) = zemxl
633               END DO
634            END DO
635         END DO
636         !
637      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
638         DO jk = 2, jpkm1         ! from the surface to the bottom :
639            DO jj = 2, jpjm1
640               DO ji = fs_2, fs_jpim1   ! vector opt.
641                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
642               END DO
643            END DO
644         END DO
645         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
646            DO jj = 2, jpjm1
647               DO ji = fs_2, fs_jpim1   ! vector opt.
648                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
649                  zmxlm(ji,jj,jk) = zemxl
650                  zmxld(ji,jj,jk) = zemxl
651               END DO
652            END DO
653         END DO
654         !
655      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
656         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
657            DO jj = 2, jpjm1
658               DO ji = fs_2, fs_jpim1   ! vector opt.
659                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
660               END DO
661            END DO
662         END DO
663         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
664            DO jj = 2, jpjm1
665               DO ji = fs_2, fs_jpim1   ! vector opt.
666                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
667               END DO
668            END DO
669         END DO
670!CDIR NOVERRCHK
671         DO jk = 2, jpkm1
672!CDIR NOVERRCHK
673            DO jj = 2, jpjm1
674!CDIR NOVERRCHK
675               DO ji = fs_2, fs_jpim1   ! vector opt.
676                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
677                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
678                  zmxlm(ji,jj,jk) = zemlm
679                  zmxld(ji,jj,jk) = zemlp
680               END DO
681            END DO
682         END DO
683         !
684      END SELECT
685      !
686# if defined key_c1d
687      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales
688      e_mix(:,:,:) = zmxlm(:,:,:)
689# endif
690
691      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
692      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt)
693      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
694!CDIR NOVERRCHK
695      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
696!CDIR NOVERRCHK
697         DO jj = 2, jpjm1
698!CDIR NOVERRCHK
699            DO ji = fs_2, fs_jpim1   ! vector opt.
700               zsqen = SQRT( en(ji,jj,jk) )
701               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
702               avm  (ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
703               avt  (ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
704               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
705            END DO
706         END DO
707      END DO
708      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged)
709      !
710      DO jk = 2, jpkm1            !* vertical eddy viscosity at wu- and wv-points
711         DO jj = 2, jpjm1
712            DO ji = fs_2, fs_jpim1   ! vector opt.
713               avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj  ,jk) ) * wumask(ji,jj,jk)
714               avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji  ,jj+1,jk) ) * wvmask(ji,jj,jk)
715            END DO
716         END DO
717      END DO
718      CALL lbc_lnk( avmu, 'U', 1. )   ;   CALL lbc_lnk( avmv, 'V', 1. )      ! Lateral boundary conditions
719      !
720      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
721         DO jk = 2, jpkm1
722            DO jj = 2, jpjm1
723               DO ji = fs_2, fs_jpim1   ! vector opt.
724                  zcoef = avm(ji,jj,jk) * 2._wp * fse3w(ji,jj,jk) * fse3w(ji,jj,jk)
725                  !                                          ! shear
726                  zdku = avmu(ji-1,jj,jk) * ( un(ji-1,jj,jk-1) - un(ji-1,jj,jk) ) * ( ub(ji-1,jj,jk-1) - ub(ji-1,jj,jk) )   &
727                    &  + avmu(ji  ,jj,jk) * ( un(ji  ,jj,jk-1) - un(ji  ,jj,jk) ) * ( ub(ji  ,jj,jk-1) - ub(ji  ,jj,jk) )
728                  zdkv = avmv(ji,jj-1,jk) * ( vn(ji,jj-1,jk-1) - vn(ji,jj-1,jk) ) * ( vb(ji,jj-1,jk-1) - vb(ji,jj-1,jk) )   &
729                    &  + avmv(ji,jj  ,jk) * ( vn(ji,jj  ,jk-1) - vn(ji,jj  ,jk) ) * ( vb(ji,jj  ,jk-1) - vb(ji,jj  ,jk) )
730                  !                                          ! local Richardson number
731                  zri   = MAX( rn2b(ji,jj,jk), 0._wp ) * zcoef / (zdku + zdkv + rn_bshear )
732                  zpdlr = MAX(  0.1_wp,  0.2 / MAX( 0.2 , zri )  )
733!!gm and even better with the use of the "true" ri_crit=0.22222...  (this change the results!)
734!!gm              zpdlr = MAX(  0.1_wp,  ri_crit / MAX( ri_crit , zri )  )
735                  avt(ji,jj,jk)   = MAX( zpdlr * avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
736# if defined key_c1d
737                  e_pdl(ji,jj,jk) = zpdlr * wmask(ji,jj,jk)  ! c1d configuration : save masked Prandlt number
738                  e_ric(ji,jj,jk) = zri   * wmask(ji,jj,jk)  ! c1d config. : save Ri
739# endif
740              END DO
741            END DO
742         END DO
743      ENDIF
744      CALL lbc_lnk( avt, 'W', 1. )                      ! Lateral boundary conditions on avt  (sign unchanged)
745
746      IF(ln_ctl) THEN
747         CALL prt_ctl( tab3d_1=en  , clinfo1=' tke  - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
748         CALL prt_ctl( tab3d_1=avmu, clinfo1=' tke  - u: ', mask1=umask,                   &
749            &          tab3d_2=avmv, clinfo2=       ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
750      ENDIF
751      !
752      CALL wrk_dealloc( jpi,jpj,jpk, zmpdl, zmxlm, zmxld ) 
753      !
754      IF( nn_timing == 1 )  CALL timing_stop('tke_avn')
755      !
756   END SUBROUTINE tke_avn
757
758
759   SUBROUTINE zdf_tke_init
760      !!----------------------------------------------------------------------
761      !!                  ***  ROUTINE zdf_tke_init  ***
762      !!                     
763      !! ** Purpose :   Initialization of the vertical eddy diffivity and
764      !!              viscosity when using a tke turbulent closure scheme
765      !!
766      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
767      !!              called at the first timestep (nit000)
768      !!
769      !! ** input   :   Namlist namzdf_tke
770      !!
771      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
772      !!----------------------------------------------------------------------
773      INTEGER ::   ji, jj, jk   ! dummy loop indices
774      INTEGER ::   ios, ierr
775      !!
776      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,   &
777         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,   &
778         &                 rn_mxl0 , nn_pdl   , ln_lc  , rn_lc    ,   &
779         &                 nn_etau , nn_htau  , rn_efr , rn_c   
780      !!----------------------------------------------------------------------
781
782      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
783      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
784901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist', lwp )
785
786      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
787      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
788902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist', lwp )
789      IF(lwm) WRITE ( numond, namzdf_tke )
790      !
791      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
792# if defined key_zdftmx_new
793      ! key_zdftmx_new: New internal wave-driven param: specified value of rn_emin & rmxl_min are used
794      rn_emin  = 1.e-10_wp
795      rmxl_min = 1.e-03_wp
796      IF(lwp) THEN                  ! Control print
797         WRITE(numout,*)
798         WRITE(numout,*) 'zdf_tke_init :  New tidal mixing case: force rn_emin = 1.e-10 and rmxl_min = 1.e-3 '
799         WRITE(numout,*) '~~~~~~~~~~~~'
800      ENDIF
801# else
802      rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
803# endif
804      !
805      IF(lwp) THEN                    !* Control print
806         WRITE(numout,*)
807         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
808         WRITE(numout,*) '~~~~~~~~~~~~'
809         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
810         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
811         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
812         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
813         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
814         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
815         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
816         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
817         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
818         WRITE(numout,*) '      surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
819         WRITE(numout,*) '      surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
820         WRITE(numout,*) '      flag to take into acc.  Langmuir circ.      ln_lc     = ', ln_lc
821         WRITE(numout,*) '      coef to compute verticla velocity of LC     rn_lc     = ', rn_lc
822         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
823         WRITE(numout,*) '      flag for computation of exp. tke profile    nn_htau   = ', nn_htau
824         WRITE(numout,*) '      fraction of en which pene. the thermocline  rn_efr    = ', rn_efr
825         WRITE(numout,*) '      fraction of TKE added within the mixed layer by nn_etau rn_c    = ', rn_c
826         WRITE(numout,*)
827         WRITE(numout,*) '      critical Richardson nb with your parameters  ri_cri = ', ri_cri
828      ENDIF
829      !
830      !                              ! allocate tke arrays
831      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
832      !
833      !                               !* Check of some namelist values
834      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
835      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
836      IF( nn_htau < 0  .OR.  nn_htau > 5 )   CALL ctl_stop( 'bad flag: nn_htau is 0 to 5    ' )
837      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
838
839      IF( ln_mxl0 ) THEN
840         IF(lwp) WRITE(numout,*) '   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
841         rn_mxl0 = rmxl_min
842      ENDIF
843     
844      IF( nn_etau == 2  ) THEN
845          ierr = zdf_mxl_alloc()
846          nmln(:,:) = nlb10           ! Initialization of nmln
847      ENDIF
848
849      IF( nn_etau /= 0 .and. nn_htau == 2 ) THEN
850          ierr = zdf_mxl_alloc()
851          nmln(:,:) = nlb10           ! Initialization of nmln
852      ENDIF
853
854      !                               !* depth of penetration of surface tke
855      IF( nn_etau /= 0 ) THEN     
856         htau(:,:) = 0._wp
857         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
858         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
859            htau(:,:) = 10._wp
860         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
861            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
862         CASE( 2 )                                 ! fraction of depth-integrated TKE within mixed-layer
863            rhtau = -1._wp / LOG( 1._wp - rn_c )
864         CASE( 3 )                                 ! F(latitude) : 0.5m to 15m poleward of 20 degrees
865            htau(:,:) = MAX(  0.5_wp, MIN( 15._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )
866         CASE( 4 )                                 ! F(latitude) : 0.5m to 10m/30m poleward of 13/40 degrees north/south
867            DO jj = 2, jpjm1
868               DO ji = fs_2, fs_jpim1   ! vector opt.
869                  IF( gphit(ji,jj) <= 0._wp ) THEN
870                     htau(ji,jj) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
871                  ELSE
872                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
873                  ENDIF
874               END DO
875            END DO
876         CASE ( 5 )                                ! F(latitude) : 0.5m to 10m poleward of 13 degrees north/south,
877            DO jj = 2, jpjm1                       !               10m to 30m between 30/45 degrees south
878               DO ji = fs_2, fs_jpim1   ! vector opt.
879                  IF( gphit(ji,jj) <= -30._wp ) THEN
880                     htau(ji,jj) = MAX(  10._wp, MIN( 30._wp, 55._wp* ABS( SIN( rpi/120._wp * ( gphit(ji,jj) + 23._wp ) ) ) )   )
881                  ELSE
882                     htau(ji,jj) = MAX(  0.5_wp, MIN( 10._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(ji,jj) ) ) )   )
883                  ENDIF
884               END DO
885            END DO
886         END SELECT
887         !
888         IF( nn_etau == 4 .AND. nn_htau /= 2 ) THEN            ! efr dependant on constant htau
889            DO jj = 2, jpjm1
890               DO ji = fs_2, fs_jpim1   ! vector opt.
891                  efr(ji,jj) = rn_efr / ( htau(ji,jj) * ( 1._wp - EXP( -bathy(ji,jj) / htau(ji,jj) ) ) )
892               END DO
893            END DO
894         ENDIF
895      ENDIF
896      !                               !* set vertical eddy coef. to the background value
897      DO jk = 1, jpk
898         avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
899         avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
900         avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
901         avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
902      END DO
903      dissl(:,:,:) = 1.e-12_wp
904      !                             
905      CALL tke_rst( nit000, 'READ' )  !* read or initialize all required files
906      !
907   END SUBROUTINE zdf_tke_init
908
909
910   SUBROUTINE tke_rst( kt, cdrw )
911     !!---------------------------------------------------------------------
912     !!                   ***  ROUTINE tke_rst  ***
913     !!                     
914     !! ** Purpose :   Read or write TKE file (en) in restart file
915     !!
916     !! ** Method  :   use of IOM library
917     !!                if the restart does not contain TKE, en is either
918     !!                set to rn_emin or recomputed
919     !!----------------------------------------------------------------------
920     INTEGER         , INTENT(in) ::   kt     ! ocean time-step
921     CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
922     !
923     INTEGER ::   jit, jk   ! dummy loop indices
924     INTEGER ::   id1, id2, id3, id4, id5, id6   ! local integers
925     !!----------------------------------------------------------------------
926     !
927     IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
928        !                                   ! ---------------
929        IF( ln_rstart ) THEN                   !* Read the restart file
930           id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
931           id2 = iom_varid( numror, 'avt'  , ldstop = .FALSE. )
932           id3 = iom_varid( numror, 'avm'  , ldstop = .FALSE. )
933           id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
934           id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
935           id6 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
936           !
937           IF( id1 > 0 ) THEN                       ! 'en' exists
938              CALL iom_get( numror, jpdom_autoglo, 'en', en, lrxios = lxios_read )
939              IF( MIN( id2, id3, id4, id5, id6 ) > 0 ) THEN        ! all required arrays exist
940                 CALL iom_get( numror, jpdom_autoglo, 'avt'  , avt, lrxios = lxios_read   )
941                 CALL iom_get( numror, jpdom_autoglo, 'avm'  , avm, lrxios = lxios_read   )
942                 CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu, lrxios = lxios_read  )
943                 CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv, lrxios = lxios_read  )
944                 CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, lrxios = lxios_read )
945              ELSE                                                 ! one at least array is missing
946                 CALL tke_avn                                          ! compute avt, avm, avmu, avmv and dissl (approximation)
947              ENDIF
948           ELSE                                     ! No TKE array found: initialisation
949              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop'
950              en (:,:,:) = rn_emin * tmask(:,:,:)
951              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation)
952              !
953              avt_k (:,:,:) = avt (:,:,:)
954              avm_k (:,:,:) = avm (:,:,:)
955              avmu_k(:,:,:) = avmu(:,:,:)
956              avmv_k(:,:,:) = avmv(:,:,:)
957              !
958              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO
959           ENDIF
960        ELSE                                   !* Start from rest
961           en(:,:,:) = rn_emin * tmask(:,:,:)
962           DO jk = 1, jpk                           ! set the Kz to the background value
963              avt (:,:,jk) = avtb(jk) * wmask (:,:,jk)
964              avm (:,:,jk) = avmb(jk) * wmask (:,:,jk)
965              avmu(:,:,jk) = avmb(jk) * wumask(:,:,jk)
966              avmv(:,:,jk) = avmb(jk) * wvmask(:,:,jk)
967           END DO
968        ENDIF
969        !
970     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
971        !                                   ! -------------------
972        IF(lwp) WRITE(numout,*) '---- tke-rst ----'
973        CALL iom_rstput( kt, nitrst, numrow, 'en'   , en     )
974        CALL iom_rstput( kt, nitrst, numrow, 'avt'  , avt_k  )
975        CALL iom_rstput( kt, nitrst, numrow, 'avm'  , avm_k  )
976        CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
977        CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
978        CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl  )
979        !
980     ENDIF
981     !
982   END SUBROUTINE tke_rst
983
984#else
985   !!----------------------------------------------------------------------
986   !!   Dummy module :                                        NO TKE scheme
987   !!----------------------------------------------------------------------
988   LOGICAL, PUBLIC, PARAMETER ::   lk_zdftke = .FALSE.   !: TKE flag
989CONTAINS
990   SUBROUTINE zdf_tke_init           ! Dummy routine
991   END SUBROUTINE zdf_tke_init
992   SUBROUTINE zdf_tke( kt )          ! Dummy routine
993      WRITE(*,*) 'zdf_tke: You should not have seen this print! error?', kt
994   END SUBROUTINE zdf_tke
995   SUBROUTINE tke_rst( kt, cdrw )
996     CHARACTER(len=*) ::   cdrw
997     WRITE(*,*) 'tke_rst: You should not have seen this print! error?', kt, cdwr
998   END SUBROUTINE tke_rst
999#endif
1000
1001   !!======================================================================
1002END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.