New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
zdftke.F90 in branches/2017/dev_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2017/dev_r7963_nemo_v3_6_AGRIF-3_AGRIFVVL/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 8010

Last change on this file since 8010 was 7988, checked in by jchanut, 7 years ago

Add AGRIF proper AGRIF bcs to GLS and TKE + vvl update

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