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/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90 @ 3680

Last change on this file since 3680 was 3680, checked in by rblod, 11 years ago

First commit of the final branch for 2012 (future nemo_3_5), see ticket #1028

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