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 NEMO/branches/2020/r4.0-HEAD_r12713_clem_dan_fixcpl/src/OCE/ZDF – NEMO

source: NEMO/branches/2020/r4.0-HEAD_r12713_clem_dan_fixcpl/src/OCE/ZDF/zdftke.F90 @ 13309

Last change on this file since 13309 was 13309, checked in by clem, 4 years ago

zdftke: orca1 crashes after 6 months unless implicit drags is used or the input of energy below sea-ice is increased. We chose here to remove the dependency of ice fraction in the surface condition on tke (back to old formulation). The rest of the changes are just there to prevent roundoff errors.

  • Property svn:keywords set to Id
File size: 47.8 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
[5120]28   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
[9019]29   !!            4.0  !  2017-04  (G. Madec)  remove CPP ddm key & avm at t-point only
[13279]30   !!             -   !  2017-05  (G. Madec)  add top/bottom friction as boundary condition
[1239]31   !!----------------------------------------------------------------------
[9019]32
[1239]33   !!----------------------------------------------------------------------
[3625]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
[1239]39   !!----------------------------------------------------------------------
[2528]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
[1492]44   USE sbc_oce        ! surface boundary condition: ocean
[9019]45   USE zdfdrg         ! vertical physics: top/bottom drag coef.
[2528]46   USE zdfmxl         ! vertical physics: mixed layer
[9019]47   !
[13006]48#if defined key_si3
49   USE ice, ONLY: hm_i, h_i
50#endif
51#if defined key_cice
52   USE sbc_ice, ONLY: h_i
53#endif
[1492]54   USE in_out_manager ! I/O manager
55   USE iom            ! I/O manager library
[2715]56   USE lib_mpp        ! MPP library
[9019]57   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
58   USE prtctl         ! Print control
[3625]59   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[1239]60
61   IMPLICIT NONE
62   PRIVATE
63
[2528]64   PUBLIC   zdf_tke        ! routine called in step module
65   PUBLIC   zdf_tke_init   ! routine called in opa module
66   PUBLIC   tke_rst        ! routine called in step module
[1239]67
[4147]68   !                      !!** Namelist  namzdf_tke  **
69   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
[13006]70   INTEGER  ::   nn_mxlice ! type of scaling under sea-ice (=0/1/2/3)
71   REAL(wp) ::   rn_mxlice ! ice thickness value when scaling under sea-ice
[4147]72   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
73   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
74   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
75   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
76   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
77   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
78   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
79   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
80   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
81   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
[9019]82   INTEGER  ::      nn_htau   ! type of tke profile of penetration (=0/1)
83   REAL(wp) ::      rn_efr    ! fraction of TKE surface value which penetrates in the ocean
[4147]84   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
[9019]85   REAL(wp) ::      rn_lc     ! coef to compute vertical velocity of Langmuir cells
[13249]86   INTEGER  ::   nn_eice   ! attenutaion of langmuir & surface wave breaking under ice (=0/1/2/3)   
[1239]87
[4147]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]
[2528]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)
[1239]92
[9019]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   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr   ! now mixing lenght of dissipation
[1492]96
[1239]97   !! * Substitutions
98#  include "vectopt_loop_substitute.h90"
99   !!----------------------------------------------------------------------
[9598]100   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[2528]101   !! $Id$
[10068]102   !! Software governed by the CeCILL license (see ./LICENSE)
[1239]103   !!----------------------------------------------------------------------
104CONTAINS
105
[2715]106   INTEGER FUNCTION zdf_tke_alloc()
107      !!----------------------------------------------------------------------
108      !!                ***  FUNCTION zdf_tke_alloc  ***
109      !!----------------------------------------------------------------------
[9019]110      ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) ,   STAT= zdf_tke_alloc )
111      !
[10425]112      CALL mpp_sum ( 'zdftke', zdf_tke_alloc )
113      IF( zdf_tke_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_alloc: failed to allocate arrays' )
[2715]114      !
115   END FUNCTION zdf_tke_alloc
116
117
[9019]118   SUBROUTINE zdf_tke( kt, p_sh2, p_avm, p_avt )
[1239]119      !!----------------------------------------------------------------------
[1531]120      !!                   ***  ROUTINE zdf_tke  ***
[1239]121      !!
122      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
[1492]123      !!              coefficients using a turbulent closure scheme (TKE).
[1239]124      !!
[1492]125      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
126      !!              is computed from a prognostic equation :
127      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
128      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
129      !!                  + avt N^2                      ! stratif. destruc.
130      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
[1239]131      !!      with the boundary conditions:
[1695]132      !!         surface: en = max( rn_emin0, rn_ebb * taum )
[1239]133      !!         bottom : en = rn_emin
[1492]134      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
135      !!
136      !!        The now Turbulent kinetic energy is computed using the following
137      !!      time stepping: implicit for vertical diffusion term, linearized semi
138      !!      implicit for kolmogoroff dissipation term, and explicit forward for
139      !!      both buoyancy and shear production terms. Therefore a tridiagonal
140      !!      linear system is solved. Note that buoyancy and shear terms are
141      !!      discretized in a energy conserving form (Bruchard 2002).
142      !!
143      !!        The dissipative and mixing length scale are computed from en and
144      !!      the stratification (see tke_avn)
145      !!
146      !!        The now vertical eddy vicosity and diffusivity coefficients are
147      !!      given by:
148      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
149      !!              avt = max( avmb, pdl * avm                 ) 
[1239]150      !!              eav = max( avmb, avm )
[1492]151      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
152      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
[1239]153      !!
154      !! ** Action  :   compute en (now turbulent kinetic energy)
[9019]155      !!                update avt, avm (before vertical eddy coef.)
[1239]156      !!
157      !! References : Gaspar et al., JGR, 1990,
158      !!              Blanke and Delecluse, JPO, 1991
159      !!              Mellor and Blumberg, JPO 2004
160      !!              Axell, JGR, 2002
[1492]161      !!              Bruchard OM 2002
[1239]162      !!----------------------------------------------------------------------
[9019]163      INTEGER                   , INTENT(in   ) ::   kt             ! ocean time step
164      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_sh2          ! shear production term
165      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   p_avm, p_avt   !  momentum and tracer Kz (w-points)
[1492]166      !!----------------------------------------------------------------------
[1481]167      !
[9019]168      CALL tke_tke( gdepw_n, e3t_n, e3w_n, p_sh2, p_avm, p_avt )   ! now tke (en)
[5656]169      !
[9019]170      CALL tke_avn( gdepw_n, e3t_n, e3w_n,        p_avm, p_avt )   ! now avt, avm, dissl
[3632]171      !
[5656]172  END SUBROUTINE zdf_tke
[1239]173
[1492]174
[9019]175   SUBROUTINE tke_tke( pdepw, p_e3t, p_e3w, p_sh2, p_avm, p_avt )
[1239]176      !!----------------------------------------------------------------------
[1492]177      !!                   ***  ROUTINE tke_tke  ***
178      !!
179      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
180      !!
181      !! ** Method  : - TKE surface boundary condition
[2528]182      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
[9019]183      !!              - source term due to shear (= Kz dz[Ub] * dz[Un] )
[1492]184      !!              - Now TKE : resolution of the TKE equation by inverting
185      !!                a tridiagonal linear system by a "methode de chasse"
186      !!              - increase TKE due to surface and internal wave breaking
[9019]187      !!             NB: when sea-ice is present, both LC parameterization
188      !!                 and TKE penetration are turned off when the ice fraction
189      !!                 is smaller than 0.25
[1492]190      !!
191      !! ** Action  : - en : now turbulent kinetic energy)
[1239]192      !! ---------------------------------------------------------------------
[9019]193      USE zdf_oce , ONLY : en   ! ocean vertical physics
194      !!
195      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   pdepw          ! depth of w-points
196      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_e3t, p_e3w   ! level thickness (t- & w-points)
197      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_sh2          ! shear production term
198      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
199      !
[13249]200      INTEGER ::   ji, jj, jk                  ! dummy loop arguments
[9019]201      REAL(wp) ::   zetop, zebot, zmsku, zmskv ! local scalars
202      REAL(wp) ::   zrhoa  = 1.22              ! Air density kg/m3
203      REAL(wp) ::   zcdrag = 1.5e-3            ! drag coefficient
[13249]204      REAL(wp) ::   zbbrau, zbbirau, zri       ! local scalars
205      REAL(wp) ::   zfact1, zfact2, zfact3     !   -      -
206      REAL(wp) ::   ztx2  , zty2  , zcof       !   -      -
207      REAL(wp) ::   ztau  , zdif               !   -      -
208      REAL(wp) ::   zus   , zwlc  , zind       !   -      -
209      REAL(wp) ::   zzd_up, zzd_lw             !   -      -
[9019]210      INTEGER , DIMENSION(jpi,jpj)     ::   imlc
[13249]211      REAL(wp), DIMENSION(jpi,jpj)     ::   zice_fra, zhlc, zus3
[9019]212      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
[1239]213      !!--------------------------------------------------------------------
[1492]214      !
[13249]215      zbbrau  =  rn_ebb / rau0       ! Local constant initialisation
216      zbbirau =  3.75_wp / rau0
217      zfact1  = -0.5_wp * rdt 
218      zfact2  =  1.5_wp * rdt * rn_ediss
219      zfact3  =  0.5_wp       * rn_ediss
[1492]220      !
[13249]221      ! ice fraction considered for attenuation of langmuir & wave breaking
222      SELECT CASE ( nn_eice )
223      CASE( 0 )   ;   zice_fra(:,:) = 0._wp
224      CASE( 1 )   ;   zice_fra(:,:) =        TANH( fr_i(:,:) * 10._wp )
225      CASE( 2 )   ;   zice_fra(:,:) =              fr_i(:,:)
226      CASE( 3 )   ;   zice_fra(:,:) = MIN( 4._wp * fr_i(:,:) , 1._wp )
227      END SELECT
228      !
[1492]229      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]230      !                     !  Surface/top/bottom boundary condition on tke
[1492]231      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[13249]232      !
[9019]233      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0)
234         DO ji = fs_2, fs_jpim1   ! vector opt.
[13309]235!! clem: this should be the right formulation but it makes the model unstable unless drags are calculated implicitly
236!!       one way around would be to increase zbbirau
237!!          en(ji,jj,1) = MAX( rn_emin0, ( ( 1._wp - fr_i(ji,jj) ) * zbbrau + &
238!!             &                                     fr_i(ji,jj)   * zbbirau ) * taum(ji,jj) ) * tmask(ji,jj,1)
239            en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
[9019]240         END DO
241      END DO
242      !
[1492]243      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
244      !                     !  Bottom boundary condition on tke
245      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[1719]246      !
[9019]247      !   en(bot)   = (ebb0/rau0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
248      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
249      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
[1492]250      !
[13279]251      IF( .NOT.ln_drg_OFF ) THEN    !== friction used as top/bottom boundary condition on TKE
[9019]252         !
[13279]253         DO jj = 2, jpjm1              ! bottom friction
[9019]254            DO ji = fs_2, fs_jpim1     ! vector opt.
255               zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
256               zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
257               !                       ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
258               zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mbkt(ji,jj))+ub(ji-1,jj,mbkt(ji,jj)) ) )**2  &
259                  &                                           + ( zmskv*( vb(ji,jj,mbkt(ji,jj))+vb(ji,jj-1,mbkt(ji,jj)) ) )**2  )
260               en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
261            END DO
262         END DO
263         IF( ln_isfcav ) THEN       ! top friction
264            DO jj = 2, jpjm1
265               DO ji = fs_2, fs_jpim1   ! vector opt.
266                  zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
267                  zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
268                  !                             ! where 0.001875 = (rn_ebb0/rau0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
269                  zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( ub(ji,jj,mikt(ji,jj))+ub(ji-1,jj,mikt(ji,jj)) ) )**2  &
270                     &                                           + ( zmskv*( vb(ji,jj,mikt(ji,jj))+vb(ji,jj-1,mikt(ji,jj)) ) )**2  )
[13249]271                  en(ji,jj,mikt(ji,jj)) = en(ji,jj,1)           * tmask(ji,jj,1) &     
[12703]272                     &                  + MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj)
[9019]273               END DO
274            END DO
275         ENDIF
276         !
277      ENDIF
278      !
[1492]279      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]280      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke   !   (Axell JGR 2002)
[1492]281         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
[1239]282         !
[1492]283         !                        !* total energy produce by LC : cumulative sum over jk
[9019]284         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * pdepw(:,:,1) * p_e3w(:,:,1)
[1239]285         DO jk = 2, jpk
[9019]286            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * pdepw(:,:,jk) * p_e3w(:,:,jk)
[1239]287         END DO
[1492]288         !                        !* finite Langmuir Circulation depth
[1705]289         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
[7753]290         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
[1239]291         DO jk = jpkm1, 2, -1
[1492]292            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us
293               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1)
[1705]294                  zus  = zcof * taum(ji,jj)
[1239]295                  IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
296               END DO
297            END DO
298         END DO
[1492]299         !                               ! finite LC depth
300         DO jj = 1, jpj 
[1239]301            DO ji = 1, jpi
[9019]302               zhlc(ji,jj) = pdepw(ji,jj,imlc(ji,jj))
[1239]303            END DO
304         END DO
[1705]305         zcof = 0.016 / SQRT( zrhoa * zcdrag )
[10425]306         DO jj = 2, jpjm1
307            DO ji = fs_2, fs_jpim1   ! vector opt.
308               zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
[13309]309               zus3(ji,jj) = MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok
[10425]310            END DO
311         END DO         
[1492]312         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en
[1239]313            DO jj = 2, jpjm1
314               DO ji = fs_2, fs_jpim1   ! vector opt.
[13309]315                  IF ( zus3(ji,jj) /= 0._wp ) THEN               
[10425]316                     ! vertical velocity due to LC   
317                     IF ( pdepw(ji,jj,jk) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN
318                        !                                           ! vertical velocity due to LC
[13249]319                        zwlc = rn_lc * SIN( rpi * pdepw(ji,jj,jk) / zhlc(ji,jj) )
[10425]320                        !                                           ! TKE Langmuir circulation source term
[13249]321                        en(ji,jj,jk) = en(ji,jj,jk) + rdt * zus3(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
[10425]322                     ENDIF
323                  ENDIF
[1239]324               END DO
325            END DO
326         END DO
327         !
328      ENDIF
[1492]329      !
330      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
331      !                     !  Now Turbulent kinetic energy (output in en)
332      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
333      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
334      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
335      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
336      !
[9019]337      IF( nn_pdl == 1 ) THEN      !* Prandtl number = F( Ri )
[5656]338         DO jk = 2, jpkm1
339            DO jj = 2, jpjm1
[9019]340               DO ji = 2, jpim1
341                  !                             ! local Richardson number
342                  zri = MAX( rn2b(ji,jj,jk), 0._wp ) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
343                  !                             ! inverse of Prandtl number
[5656]344                  apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
345               END DO
346            END DO
347         END DO
348      ENDIF
[5836]349      !         
[5120]350      DO jk = 2, jpkm1           !* Matrix and right hand side in en
351         DO jj = 2, jpjm1
352            DO ji = fs_2, fs_jpim1   ! vector opt.
[1492]353               zcof   = zfact1 * tmask(ji,jj,jk)
[9019]354               !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
355               !                                   ! eddy coefficient (ensure numerical stability)
356               zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
357                  &          /    (  p_e3t(ji,jj,jk  ) * p_e3w(ji,jj,jk  )  )
358               zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
359                  &          /    (  p_e3t(ji,jj,jk-1) * p_e3w(ji,jj,jk  )  )
[5656]360               !
[1492]361               zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
362               zd_lw(ji,jj,jk) = zzd_lw
[9019]363               zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
[1239]364               !
[1492]365               !                                   ! right hand side in en
[9019]366               en(ji,jj,jk) = en(ji,jj,jk) + rdt * (  p_sh2(ji,jj,jk)                          &   ! shear
367                  &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
368                  &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
369                  &                                ) * wmask(ji,jj,jk)
[1239]370            END DO
[5120]371         END DO
372      END DO
373      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
374      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
375         DO jj = 2, jpjm1
376            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]377               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]378            END DO
[5120]379         END DO
380      END DO
[5836]381      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
[5120]382         DO ji = fs_2, fs_jpim1   ! vector opt.
383            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
384         END DO
385      END DO
386      DO jk = 3, jpkm1
387         DO jj = 2, jpjm1
388            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]389               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]390            END DO
[5120]391         END DO
392      END DO
[5836]393      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
[5120]394         DO ji = fs_2, fs_jpim1   ! vector opt.
[1492]395            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
[5120]396         END DO
397      END DO
398      DO jk = jpk-2, 2, -1
399         DO jj = 2, jpjm1
400            DO ji = fs_2, fs_jpim1    ! vector opt.
[1492]401               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
[1239]402            END DO
[5120]403         END DO
404      END DO
405      DO jk = 2, jpkm1                             ! set the minimum value of tke
406         DO jj = 2, jpjm1
407            DO ji = fs_2, fs_jpim1   ! vector opt.
408               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
[1239]409            END DO
410         END DO
411      END DO
[9019]412      !
[1492]413      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
414      !                            !  TKE due to surface and internal wave breaking
415      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[6140]416!!gm BUG : in the exp  remove the depth of ssh !!!
[9019]417!!gm       i.e. use gde3w in argument (pdepw)
[6140]418     
419     
[2528]420      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
[13249]421         DO jk = 2, jpkm1                       ! nn_eice=0 : ON below sea-ice ; nn_eice>0 : partly OFF
[1239]422            DO jj = 2, jpjm1
423               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]424                  en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[13309]425                     &                                 * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1239]426               END DO
427            END DO
[1492]428         END DO
[2528]429      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
[1492]430         DO jj = 2, jpjm1
431            DO ji = fs_2, fs_jpim1   ! vector opt.
432               jk = nmln(ji,jj)
[9019]433               en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[13309]434                  &                                 * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1239]435            END DO
[1492]436         END DO
[2528]437      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
[1705]438         DO jk = 2, jpkm1
439            DO jj = 2, jpjm1
440               DO ji = fs_2, fs_jpim1   ! vector opt.
441                  ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
442                  zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
[4990]443                  ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
[2528]444                  zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
445                  zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
[9019]446                  en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -pdepw(ji,jj,jk) / htau(ji,jj) )   &
[13309]447                     &                                 * MAX( 0._wp, 1._wp - zice_fra(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
[1705]448               END DO
449            END DO
450         END DO
[1239]451      ENDIF
[1492]452      !
[1239]453   END SUBROUTINE tke_tke
454
[1492]455
[9019]456   SUBROUTINE tke_avn( pdepw, p_e3t, p_e3w, p_avm, p_avt )
[1239]457      !!----------------------------------------------------------------------
[1492]458      !!                   ***  ROUTINE tke_avn  ***
[1239]459      !!
[1492]460      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
461      !!
462      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
463      !!              the tke_tke routine). First, the now mixing lenth is
464      !!      computed from en and the strafification (N^2), then the mixings
465      !!      coefficients are computed.
466      !!              - Mixing length : a first evaluation of the mixing lengh
467      !!      scales is:
468      !!                      mxl = sqrt(2*en) / N 
469      !!      where N is the brunt-vaisala frequency, with a minimum value set
[2528]470      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
[1492]471      !!        The mixing and dissipative length scale are bound as follow :
472      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
473      !!                        zmxld = zmxlm = mxl
474      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
475      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
476      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
477      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
478      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
479      !!                    the surface to obtain ldown. the resulting length
480      !!                    scales are:
481      !!                        zmxld = sqrt( lup * ldown )
482      !!                        zmxlm = min ( lup , ldown )
483      !!              - Vertical eddy viscosity and diffusivity:
484      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
485      !!                      avt = max( avmb, pdlr * avm ) 
486      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
487      !!
[9019]488      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
[1239]489      !!----------------------------------------------------------------------
[9019]490      USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d   ! ocean vertical physics
491      !!
492      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdepw          ! depth (w-points)
493      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_e3t, p_e3w   ! level thickness (t- & w-points)
494      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
495      !
[2715]496      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[9019]497      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
498      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
[13006]499      REAL(wp) ::   zemxl, zemlm, zemlp, zmaxice       !   -      -
[9019]500      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmxlm, zmxld   ! 3D workspace
[1239]501      !!--------------------------------------------------------------------
[3294]502      !
[1492]503      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
504      !                     !  Mixing length
505      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
506      !
507      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
508      !
[5120]509      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
[7753]510      zmxlm(:,:,:)  = rmxl_min   
511      zmxld(:,:,:)  = rmxl_min
[13249]512      !
513      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g)
[13006]514         !
[9019]515         zraug = vkarmn * 2.e5_wp / ( rau0 * grav )
[13006]516#if ! defined key_si3 && ! defined key_cice
[13249]517         DO jj = 2, jpjm1                     ! No sea-ice
[4990]518            DO ji = fs_2, fs_jpim1
[13006]519               zmxlm(ji,jj,1) =  zraug * taum(ji,jj) * tmask(ji,jj,1)
[4990]520            END DO
521         END DO
[13006]522#else
[13249]523
[13006]524         SELECT CASE( nn_mxlice )             ! Type of scaling under sea-ice
525         !
526         CASE( 0 )                      ! No scaling under sea-ice
527            DO jj = 2, jpjm1
528               DO ji = fs_2, fs_jpim1
529                  zmxlm(ji,jj,1) = zraug * taum(ji,jj) * tmask(ji,jj,1)
530               END DO
531            END DO
532            !
[13249]533         CASE( 1 )                      ! scaling with constant sea-ice thickness
[13006]534            DO jj = 2, jpjm1
535               DO ji = fs_2, fs_jpim1
[13249]536                  zmxlm(ji,jj,1) =  ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
537                     &                          fr_i(ji,jj)   * rn_mxlice           ) * tmask(ji,jj,1)
[13006]538               END DO
539            END DO
540            !
[13249]541         CASE( 2 )                      ! scaling with mean sea-ice thickness
[13006]542            DO jj = 2, jpjm1
543               DO ji = fs_2, fs_jpim1
544#if defined key_si3
[13249]545                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
546                     &                         fr_i(ji,jj)   * hm_i(ji,jj) * 2._wp ) * tmask(ji,jj,1)
[13006]547#elif defined key_cice
548                  zmaxice = MAXVAL( h_i(ji,jj,:) )
[13249]549                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
550                     &                         fr_i(ji,jj)   * zmaxice             ) * tmask(ji,jj,1)
[13006]551#endif
552               END DO
553            END DO
554            !
[13249]555         CASE( 3 )                      ! scaling with max sea-ice thickness
[13006]556            DO jj = 2, jpjm1
557               DO ji = fs_2, fs_jpim1
558                  zmaxice = MAXVAL( h_i(ji,jj,:) )
[13249]559                  zmxlm(ji,jj,1) = ( ( 1._wp - fr_i(ji,jj) ) * zraug * taum(ji,jj) + &
560                     &                         fr_i(ji,jj)   * zmaxice             ) * tmask(ji,jj,1)
[13006]561               END DO
562            END DO
563            !
564         END SELECT
565#endif
566         !
567         DO jj = 2, jpjm1
568            DO ji = fs_2, fs_jpim1
569               zmxlm(ji,jj,1) = MAX( rn_mxl0, zmxlm(ji,jj,1) )
570            END DO
571         END DO
572         !
573      ELSE
[7753]574         zmxlm(:,:,1) = rn_mxl0
[1239]575      ENDIF
576      !
[5120]577      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2)
578         DO jj = 2, jpjm1
579            DO ji = fs_2, fs_jpim1   ! vector opt.
[1239]580               zrn2 = MAX( rn2(ji,jj,jk), rsmall )
[5836]581               zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
[1239]582            END DO
583         END DO
584      END DO
[1492]585      !
586      !                     !* Physical limits for the mixing length
587      !
[7753]588      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
589      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
[1492]590      !
[1239]591      SELECT CASE ( nn_mxl )
592      !
[5836]593 !!gm Not sure of that coding for ISF....
[9019]594      ! where wmask = 0 set zmxlm == p_e3w
[1239]595      CASE ( 0 )           ! bounded by the distance to surface and bottom
[5120]596         DO jk = 2, jpkm1
597            DO jj = 2, jpjm1
598               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]599                  zemxl = MIN( pdepw(ji,jj,jk) - pdepw(ji,jj,mikt(ji,jj)), zmxlm(ji,jj,jk),   &
600                  &            pdepw(ji,jj,mbkt(ji,jj)+1) - pdepw(ji,jj,jk) )
[5120]601                  ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
[9019]602                  zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
603                  zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk) + MIN( zmxlm(ji,jj,jk) , p_e3w(ji,jj,jk) ) * (1 - wmask(ji,jj,jk))
[1239]604               END DO
605            END DO
606         END DO
607         !
608      CASE ( 1 )           ! bounded by the vertical scale factor
[5120]609         DO jk = 2, jpkm1
610            DO jj = 2, jpjm1
611               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]612                  zemxl = MIN( p_e3w(ji,jj,jk), zmxlm(ji,jj,jk) )
[1239]613                  zmxlm(ji,jj,jk) = zemxl
614                  zmxld(ji,jj,jk) = zemxl
615               END DO
616            END DO
617         END DO
618         !
619      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
[5120]620         DO jk = 2, jpkm1         ! from the surface to the bottom :
621            DO jj = 2, jpjm1
622               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]623                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
[1239]624               END DO
[5120]625            END DO
626         END DO
627         DO jk = jpkm1, 2, -1     ! from the bottom to the surface :
628            DO jj = 2, jpjm1
629               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]630                  zemxl = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
[1239]631                  zmxlm(ji,jj,jk) = zemxl
632                  zmxld(ji,jj,jk) = zemxl
633               END DO
634            END DO
635         END DO
636         !
637      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
[5120]638         DO jk = 2, jpkm1         ! from the surface to the bottom : lup
639            DO jj = 2, jpjm1
640               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]641                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + p_e3t(ji,jj,jk-1), zmxlm(ji,jj,jk) )
[1239]642               END DO
[5120]643            END DO
644         END DO
645         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown
646            DO jj = 2, jpjm1
647               DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]648                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + p_e3t(ji,jj,jk+1), zmxlm(ji,jj,jk) )
[1239]649               END DO
650            END DO
651         END DO
652         DO jk = 2, jpkm1
653            DO jj = 2, jpjm1
654               DO ji = fs_2, fs_jpim1   ! vector opt.
655                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
656                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
657                  zmxlm(ji,jj,jk) = zemlm
658                  zmxld(ji,jj,jk) = zemlp
659               END DO
660            END DO
661         END DO
662         !
663      END SELECT
[1492]664      !
665      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]666      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
[1492]667      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
668      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points
[1239]669         DO jj = 2, jpjm1
670            DO ji = fs_2, fs_jpim1   ! vector opt.
671               zsqen = SQRT( en(ji,jj,jk) )
672               zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
[9019]673               p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
674               p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
[1239]675               dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
676            END DO
677         END DO
678      END DO
[1492]679      !
680      !
681      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
[5120]682         DO jk = 2, jpkm1
683            DO jj = 2, jpjm1
684               DO ji = fs_2, fs_jpim1   ! vector opt.
[12697]685                  p_avt(ji,jj,jk)   = MAX( apdlr(ji,jj,jk) * p_avt(ji,jj,jk), avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
[1239]686              END DO
687            END DO
688         END DO
689      ENDIF
[9019]690      !
[1239]691      IF(ln_ctl) THEN
[9440]692         CALL prt_ctl( tab3d_1=en   , clinfo1=' tke  - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk)
693         CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke  - m: ', kdim=jpk )
[1239]694      ENDIF
695      !
[1492]696   END SUBROUTINE tke_avn
[1239]697
[1492]698
[2528]699   SUBROUTINE zdf_tke_init
[1239]700      !!----------------------------------------------------------------------
[2528]701      !!                  ***  ROUTINE zdf_tke_init  ***
[1239]702      !!                     
703      !! ** Purpose :   Initialization of the vertical eddy diffivity and
[1492]704      !!              viscosity when using a tke turbulent closure scheme
[1239]705      !!
[1601]706      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
[1492]707      !!              called at the first timestep (nit000)
[1239]708      !!
[1601]709      !! ** input   :   Namlist namzdf_tke
[1239]710      !!
711      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
712      !!----------------------------------------------------------------------
[9019]713      USE zdf_oce , ONLY : ln_zdfiwm   ! Internal Wave Mixing flag
714      !!
[1239]715      INTEGER ::   ji, jj, jk   ! dummy loop indices
[4147]716      INTEGER ::   ios
[1239]717      !!
[13006]718      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb   , rn_emin  ,  &
719         &                 rn_emin0, rn_bshear, nn_mxl   , ln_mxl0  ,  &
720         &                 rn_mxl0 , nn_mxlice, rn_mxlice,             &
[13279]721         &                 nn_pdl  , ln_lc    , rn_lc,                 &
[13249]722         &                 nn_etau , nn_htau  , rn_efr   , nn_eice 
[1239]723      !!----------------------------------------------------------------------
[2715]724      !
[4147]725      REWIND( numnam_ref )              ! Namelist namzdf_tke in reference namelist : Turbulent Kinetic Energy
726      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
[11536]727901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' )
[4147]728
729      REWIND( numnam_cfg )              ! Namelist namzdf_tke in configuration namelist : Turbulent Kinetic Energy
730      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
[11536]731902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' )
[4624]732      IF(lwm) WRITE ( numond, namzdf_tke )
[2715]733      !
[2528]734      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
[2715]735      !
[1492]736      IF(lwp) THEN                    !* Control print
[1239]737         WRITE(numout,*)
[2528]738         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
739         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]740         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
[1705]741         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
742         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
743         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
744         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
745         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
[9019]746         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
[1705]747         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
748         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
[9019]749         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
750         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
[13006]751         IF( ln_mxl0 ) THEN
752            WRITE(numout,*) '      type of scaling under sea-ice               nn_mxlice = ', nn_mxlice
753            IF( nn_mxlice == 1 ) &
754            WRITE(numout,*) '      ice thickness when scaling under sea-ice    rn_mxlice = ', rn_mxlice
[13249]755            SELECT CASE( nn_mxlice )             ! Type of scaling under sea-ice
756            CASE( 0 )   ;   WRITE(numout,*) '   ==>>>   No scaling under sea-ice'
757            CASE( 1 )   ;   WRITE(numout,*) '   ==>>>   scaling with constant sea-ice thickness'
758            CASE( 2 )   ;   WRITE(numout,*) '   ==>>>   scaling with mean sea-ice thickness'
759            CASE( 3 )   ;   WRITE(numout,*) '   ==>>>   scaling with max sea-ice thickness'
760            CASE DEFAULT
761               CALL ctl_stop( 'zdf_tke_init: wrong value for nn_mxlice, should be 0,1,2,3 or 4')
762            END SELECT
[13006]763         ENDIF
[9019]764         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
765         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
[1705]766         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
[9019]767         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
768         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
[13249]769         WRITE(numout,*) '      langmuir & surface wave breaking under ice  nn_eice = ', nn_eice
770         SELECT CASE( nn_eice ) 
771         CASE( 0 )   ;   WRITE(numout,*) '   ==>>>   no impact of ice cover on langmuir & surface wave breaking'
772         CASE( 1 )   ;   WRITE(numout,*) '   ==>>>   weigthed by 1-TANH( fr_i(:,:) * 10 )'
773         CASE( 2 )   ;   WRITE(numout,*) '   ==>>>   weighted by 1-fr_i(:,:)'
774         CASE( 3 )   ;   WRITE(numout,*) '   ==>>>   weighted by 1-MIN( 1, 4 * fr_i(:,:) )'
775         CASE DEFAULT
776            CALL ctl_stop( 'zdf_tke_init: wrong value for nn_eice, should be 0,1,2, or 3')
777         END SELECT     
[13279]778         IF( .NOT.ln_drg_OFF ) THEN
[9169]779            WRITE(numout,*)
[9019]780            WRITE(numout,*) '   Namelist namdrg_top/_bot:   used values:'
781            WRITE(numout,*) '      top    ocean cavity roughness (m)          rn_z0(_top)= ', r_z0_top
782            WRITE(numout,*) '      Bottom seafloor     roughness (m)          rn_z0(_bot)= ', r_z0_bot
783         ENDIF
784         WRITE(numout,*)
[9190]785         WRITE(numout,*) '   ==>>>   critical Richardson nb with your parameters  ri_cri = ', ri_cri
[9019]786         WRITE(numout,*)
[1239]787      ENDIF
[2715]788      !
[9019]789      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
790         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
791         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
[9190]792         IF(lwp) WRITE(numout,*) '   ==>>>   Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3'
[9019]793      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
794         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
[9190]795         IF(lwp) WRITE(numout,*) '   ==>>>   minimum mixing length with your parameters rmxl_min = ', rmxl_min
[9019]796      ENDIF
797      !
[2715]798      !                              ! allocate tke arrays
799      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
800      !
[1492]801      !                               !* Check of some namelist values
[4990]802      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
803      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
804      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
[5407]805      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
[9019]806      !
[2528]807      IF( ln_mxl0 ) THEN
[9169]808         IF(lwp) WRITE(numout,*)
[9190]809         IF(lwp) WRITE(numout,*) '   ==>>>   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
[2528]810         rn_mxl0 = rmxl_min
811      ENDIF
812     
[1492]813      IF( nn_etau == 2  )   CALL zdf_mxl( nit000 )      ! Initialization of nmln
[1239]814
[1492]815      !                               !* depth of penetration of surface tke
816      IF( nn_etau /= 0 ) THEN     
[1601]817         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
[2528]818         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
[7753]819            htau(:,:) = 10._wp
[2528]820         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
[7753]821            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
[1492]822         END SELECT
823      ENDIF
[9019]824      !                                !* read or initialize all required files
825      CALL tke_rst( nit000, 'READ' )      ! (en, avt_k, avm_k, dissl)
[1239]826      !
[9367]827      IF( lwxios ) THEN
828         CALL iom_set_rstw_var_active('en')
829         CALL iom_set_rstw_var_active('avt_k')
830         CALL iom_set_rstw_var_active('avm_k')
831         CALL iom_set_rstw_var_active('dissl')
832      ENDIF
[2528]833   END SUBROUTINE zdf_tke_init
[1239]834
835
[1531]836   SUBROUTINE tke_rst( kt, cdrw )
[9019]837      !!---------------------------------------------------------------------
838      !!                   ***  ROUTINE tke_rst  ***
839      !!                     
840      !! ** Purpose :   Read or write TKE file (en) in restart file
841      !!
842      !! ** Method  :   use of IOM library
843      !!                if the restart does not contain TKE, en is either
844      !!                set to rn_emin or recomputed
845      !!----------------------------------------------------------------------
846      USE zdf_oce , ONLY : en, avt_k, avm_k   ! ocean vertical physics
847      !!
848      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
849      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
850      !
851      INTEGER ::   jit, jk              ! dummy loop indices
852      INTEGER ::   id1, id2, id3, id4   ! local integers
853      !!----------------------------------------------------------------------
854      !
855      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
856         !                                   ! ---------------
857         IF( ln_rstart ) THEN                   !* Read the restart file
858            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
859            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
860            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
861            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
862            !
863            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
[9367]864               CALL iom_get( numror, jpdom_autoglo, 'en'   , en   , ldxios = lrxios )
865               CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k, ldxios = lrxios )
866               CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k, ldxios = lrxios )
867               CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, ldxios = lrxios )
[9019]868            ELSE                                          ! start TKE from rest
[9169]869               IF(lwp) WRITE(numout,*)
[9190]870               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without TKE scheme, set en to background values'
[9019]871               en   (:,:,:) = rn_emin * wmask(:,:,:)
872               dissl(:,:,:) = 1.e-12_wp
873               ! avt_k, avm_k already set to the background value in zdf_phy_init
874            ENDIF
875         ELSE                                   !* Start from rest
[9169]876            IF(lwp) WRITE(numout,*)
[9190]877            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set en to the background value'
[9019]878            en   (:,:,:) = rn_emin * wmask(:,:,:)
879            dissl(:,:,:) = 1.e-12_wp
880            ! avt_k, avm_k already set to the background value in zdf_phy_init
881         ENDIF
882         !
883      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
884         !                                   ! -------------------
[9169]885         IF(lwp) WRITE(numout,*) '---- tke_rst ----'
[9367]886         IF( lwxios ) CALL iom_swap(      cwxios_context          ) 
887         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en   , ldxios = lwxios )
888         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k, ldxios = lwxios )
889         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k, ldxios = lwxios )
890         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl, ldxios = lwxios )
891         IF( lwxios ) CALL iom_swap(      cxios_context          )
[9019]892         !
893      ENDIF
894      !
[1531]895   END SUBROUTINE tke_rst
[1239]896
897   !!======================================================================
[1531]898END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.