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/dev_r12527_Gurvan_ShallowWater/src/OCE/ZDF – NEMO

source: NEMO/branches/2020/dev_r12527_Gurvan_ShallowWater/src/OCE/ZDF/zdftke.F90 @ 13151

Last change on this file since 13151 was 13151, checked in by gm, 4 years ago

result from merge with qco r12983

  • Property svn:keywords set to Id
File size: 39.4 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
30   !!             -   !  2017-05  (G. Madec)  add top/bottom friction as boundary condition (ln_drg)
[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   !
[1492]48   USE in_out_manager ! I/O manager
49   USE iom            ! I/O manager library
[2715]50   USE lib_mpp        ! MPP library
[9019]51   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
52   USE prtctl         ! Print control
[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
[4147]62   !                      !!** Namelist  namzdf_tke  **
63   LOGICAL  ::   ln_mxl0   ! mixing length scale surface value as function of wind stress or not
64   INTEGER  ::   nn_mxl    ! type of mixing length (=0/1/2/3)
65   REAL(wp) ::   rn_mxl0   ! surface  min value of mixing length (kappa*z_o=0.4*0.1 m)  [m]
66   INTEGER  ::   nn_pdl    ! Prandtl number or not (ratio avt/avm) (=0/1)
67   REAL(wp) ::   rn_ediff  ! coefficient for avt: avt=rn_ediff*mxl*sqrt(e)
68   REAL(wp) ::   rn_ediss  ! coefficient of the Kolmogoroff dissipation
69   REAL(wp) ::   rn_ebb    ! coefficient of the surface input of tke
70   REAL(wp) ::   rn_emin   ! minimum value of tke           [m2/s2]
71   REAL(wp) ::   rn_emin0  ! surface minimum value of tke   [m2/s2]
72   REAL(wp) ::   rn_bshear ! background shear (>0) currently a numerical threshold (do not change it)
[9019]73   LOGICAL  ::   ln_drg    ! top/bottom friction forcing flag
[4147]74   INTEGER  ::   nn_etau   ! type of depth penetration of surface tke (=0/1/2/3)
[9019]75   INTEGER  ::      nn_htau   ! type of tke profile of penetration (=0/1)
76   REAL(wp) ::      rn_efr    ! fraction of TKE surface value which penetrates in the ocean
[9546]77   REAL(wp) ::      rn_eice   ! =0 ON below sea-ice, =4 OFF when ice fraction > 1/4   
[4147]78   LOGICAL  ::   ln_lc     ! Langmuir cells (LC) as a source term of TKE or not
[9019]79   REAL(wp) ::      rn_lc     ! coef to compute vertical velocity of Langmuir cells
[1239]80
[4147]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]
[2528]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
[9019]86   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   htau    ! depth of tke penetration (nn_htau)
87   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dissl   ! now mixing lenght of dissipation
88   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   apdlr   ! now mixing lenght of dissipation
[1492]89
[1239]90   !! * Substitutions
[12377]91#  include "do_loop_substitute.h90"
[13151]92#  include "domzgr_substitute.h90"
[1239]93   !!----------------------------------------------------------------------
[9598]94   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[2528]95   !! $Id$
[10068]96   !! Software governed by the CeCILL license (see ./LICENSE)
[1239]97   !!----------------------------------------------------------------------
98CONTAINS
99
[2715]100   INTEGER FUNCTION zdf_tke_alloc()
101      !!----------------------------------------------------------------------
102      !!                ***  FUNCTION zdf_tke_alloc  ***
103      !!----------------------------------------------------------------------
[9019]104      ALLOCATE( htau(jpi,jpj) , dissl(jpi,jpj,jpk) , apdlr(jpi,jpj,jpk) ,   STAT= zdf_tke_alloc )
105      !
[10425]106      CALL mpp_sum ( 'zdftke', zdf_tke_alloc )
107      IF( zdf_tke_alloc /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_alloc: failed to allocate arrays' )
[2715]108      !
109   END FUNCTION zdf_tke_alloc
110
111
[12377]112   SUBROUTINE zdf_tke( kt, Kbb, Kmm, p_sh2, p_avm, p_avt )
[1239]113      !!----------------------------------------------------------------------
[1531]114      !!                   ***  ROUTINE zdf_tke  ***
[1239]115      !!
116      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
[1492]117      !!              coefficients using a turbulent closure scheme (TKE).
[1239]118      !!
[1492]119      !! ** Method  :   The time evolution of the turbulent kinetic energy (tke)
120      !!              is computed from a prognostic equation :
121      !!         d(en)/dt = avm (d(u)/dz)**2             ! shear production
122      !!                  + d( avm d(en)/dz )/dz         ! diffusion of tke
123      !!                  + avt N^2                      ! stratif. destruc.
124      !!                  - rn_ediss / emxl en**(2/3)    ! Kolmogoroff dissipation
[1239]125      !!      with the boundary conditions:
[1695]126      !!         surface: en = max( rn_emin0, rn_ebb * taum )
[1239]127      !!         bottom : en = rn_emin
[1492]128      !!      The associated critical Richardson number is: ri_cri = 2/(2+rn_ediss/rn_ediff)
129      !!
130      !!        The now Turbulent kinetic energy is computed using the following
131      !!      time stepping: implicit for vertical diffusion term, linearized semi
132      !!      implicit for kolmogoroff dissipation term, and explicit forward for
133      !!      both buoyancy and shear production terms. Therefore a tridiagonal
134      !!      linear system is solved. Note that buoyancy and shear terms are
135      !!      discretized in a energy conserving form (Bruchard 2002).
136      !!
137      !!        The dissipative and mixing length scale are computed from en and
138      !!      the stratification (see tke_avn)
139      !!
140      !!        The now vertical eddy vicosity and diffusivity coefficients are
141      !!      given by:
142      !!              avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
143      !!              avt = max( avmb, pdl * avm                 ) 
[1239]144      !!              eav = max( avmb, avm )
[1492]145      !!      where pdl, the inverse of the Prandtl number is 1 if nn_pdl=0 and
146      !!      given by an empirical funtion of the localRichardson number if nn_pdl=1
[1239]147      !!
148      !! ** Action  :   compute en (now turbulent kinetic energy)
[9019]149      !!                update avt, avm (before vertical eddy coef.)
[1239]150      !!
151      !! References : Gaspar et al., JGR, 1990,
152      !!              Blanke and Delecluse, JPO, 1991
153      !!              Mellor and Blumberg, JPO 2004
154      !!              Axell, JGR, 2002
[1492]155      !!              Bruchard OM 2002
[1239]156      !!----------------------------------------------------------------------
[9019]157      INTEGER                   , INTENT(in   ) ::   kt             ! ocean time step
[12377]158      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]159      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   p_sh2          ! shear production term
160      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   p_avm, p_avt   !  momentum and tracer Kz (w-points)
[1492]161      !!----------------------------------------------------------------------
[1481]162      !
[12377]163      CALL tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )   ! now tke (en)
[5656]164      !
[12377]165      CALL tke_avn( Kbb, Kmm,        p_avm, p_avt )   ! now avt, avm, dissl
[3632]166      !
[5656]167  END SUBROUTINE zdf_tke
[1239]168
[1492]169
[12377]170   SUBROUTINE tke_tke( Kbb, Kmm, p_sh2, p_avm, p_avt )
[1239]171      !!----------------------------------------------------------------------
[1492]172      !!                   ***  ROUTINE tke_tke  ***
173      !!
174      !! ** Purpose :   Compute the now Turbulente Kinetic Energy (TKE)
175      !!
176      !! ** Method  : - TKE surface boundary condition
[2528]177      !!              - source term due to Langmuir cells (Axell JGR 2002) (ln_lc=T)
[9019]178      !!              - source term due to shear (= Kz dz[Ub] * dz[Un] )
[1492]179      !!              - Now TKE : resolution of the TKE equation by inverting
180      !!                a tridiagonal linear system by a "methode de chasse"
181      !!              - increase TKE due to surface and internal wave breaking
[9019]182      !!             NB: when sea-ice is present, both LC parameterization
183      !!                 and TKE penetration are turned off when the ice fraction
184      !!                 is smaller than 0.25
[1492]185      !!
186      !! ** Action  : - en : now turbulent kinetic energy)
[1239]187      !! ---------------------------------------------------------------------
[9019]188      USE zdf_oce , ONLY : en   ! ocean vertical physics
189      !!
[12377]190      INTEGER                    , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]191      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_sh2          ! shear production term
192      REAL(wp), DIMENSION(:,:,:) , INTENT(in   ) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
193      !
194      INTEGER ::   ji, jj, jk              ! dummy loop arguments
195      REAL(wp) ::   zetop, zebot, zmsku, zmskv ! local scalars
196      REAL(wp) ::   zrhoa  = 1.22              ! Air density kg/m3
197      REAL(wp) ::   zcdrag = 1.5e-3            ! drag coefficient
198      REAL(wp) ::   zbbrau, zri                ! local scalars
199      REAL(wp) ::   zfact1, zfact2, zfact3     !   -         -
200      REAL(wp) ::   ztx2  , zty2  , zcof       !   -         -
201      REAL(wp) ::   ztau  , zdif               !   -         -
202      REAL(wp) ::   zus   , zwlc  , zind       !   -         -
203      REAL(wp) ::   zzd_up, zzd_lw             !   -         -
204      INTEGER , DIMENSION(jpi,jpj)     ::   imlc
[10425]205      REAL(wp), DIMENSION(jpi,jpj)     ::   zhlc, zfr_i
[9019]206      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zpelc, zdiag, zd_up, zd_lw
[1239]207      !!--------------------------------------------------------------------
[1492]208      !
[12489]209      zbbrau = rn_ebb / rho0       ! Local constant initialisation
210      zfact1 = -.5_wp * rn_Dt 
211      zfact2 = 1.5_wp * rn_Dt * rn_ediss
[2528]212      zfact3 = 0.5_wp       * rn_ediss
[1492]213      !
214      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]215      !                     !  Surface/top/bottom boundary condition on tke
[1492]216      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[13151]217      !
[12377]218      DO_2D_00_00
219         en(ji,jj,1) = MAX( rn_emin0, zbbrau * taum(ji,jj) ) * tmask(ji,jj,1)
220      END_2D
[9019]221      !
[1492]222      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
223      !                     !  Bottom boundary condition on tke
224      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[1719]225      !
[12489]226      !   en(bot)   = (ebb0/rho0)*0.5*sqrt(u_botfr^2+v_botfr^2) (min value rn_emin)
[9019]227      ! where ebb0 does not includes surface wave enhancement (i.e. ebb0=3.75)
228      ! Note that stress averaged is done using an wet-only calculation of u and v at t-point like in zdfsh2
[1492]229      !
[9019]230      IF( ln_drg ) THEN       !== friction used as top/bottom boundary condition on TKE
231         !
[12377]232         DO_2D_00_00
233            zmsku = ( 2. - umask(ji-1,jj,mbkt(ji,jj)) * umask(ji,jj,mbkt(ji,jj)) )
234            zmskv = ( 2. - vmask(ji,jj-1,mbkt(ji,jj)) * vmask(ji,jj,mbkt(ji,jj)) )
[12489]235            !                       ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000. (CAUTION CdU<0)
[12377]236            zebot = - 0.001875_wp * rCdU_bot(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mbkt(ji,jj),Kbb)+uu(ji-1,jj,mbkt(ji,jj),Kbb) ) )**2  &
237               &                                           + ( zmskv*( vv(ji,jj,mbkt(ji,jj),Kbb)+vv(ji,jj-1,mbkt(ji,jj),Kbb) ) )**2  )
238            en(ji,jj,mbkt(ji,jj)+1) = MAX( zebot, rn_emin ) * ssmask(ji,jj)
239         END_2D
[9019]240         IF( ln_isfcav ) THEN       ! top friction
[12377]241            DO_2D_00_00
242               zmsku = ( 2. - umask(ji-1,jj,mikt(ji,jj)) * umask(ji,jj,mikt(ji,jj)) )
243               zmskv = ( 2. - vmask(ji,jj-1,mikt(ji,jj)) * vmask(ji,jj,mikt(ji,jj)) )
[12489]244               !                             ! where 0.001875 = (rn_ebb0/rho0) * 0.5 = 3.75*0.5/1000.  (CAUTION CdU<0)
[12377]245               zetop = - 0.001875_wp * rCdU_top(ji,jj) * SQRT(  ( zmsku*( uu(ji,jj,mikt(ji,jj),Kbb)+uu(ji-1,jj,mikt(ji,jj),Kbb) ) )**2  &
246                  &                                           + ( zmskv*( vv(ji,jj,mikt(ji,jj),Kbb)+vv(ji,jj-1,mikt(ji,jj),Kbb) ) )**2  )
[13151]247               en(ji,jj,mikt(ji,jj)) = en(ji,jj,1) * tmask(ji,jj,1) + MAX( zetop, rn_emin ) * (1._wp - tmask(ji,jj,1)) * ssmask(ji,jj)   ! masked at ocean surface
[12377]248            END_2D
[9019]249         ENDIF
250         !
251      ENDIF
252      !
[1492]253      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]254      IF( ln_lc ) THEN      !  Langmuir circulation source term added to tke   !   (Axell JGR 2002)
[1492]255         !                  !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
[1239]256         !
[1492]257         !                        !* total energy produce by LC : cumulative sum over jk
[13151]258         zpelc(:,:,1) = MAX( rn2b(:,:,1), 0._wp ) * gdepw(:,:,1,Kmm) * e3w(:,:,1,Kmm)
[1239]259         DO jk = 2, jpk
[13151]260            zpelc(:,:,jk)  = zpelc(:,:,jk-1) +   &
261               &        MAX( rn2b(:,:,jk), 0._wp ) * gdepw(:,:,jk,Kmm) * e3w(:,:,jk,Kmm)
[1239]262         END DO
[1492]263         !                        !* finite Langmuir Circulation depth
[1705]264         zcof = 0.5 * 0.016 * 0.016 / ( zrhoa * zcdrag )
[7753]265         imlc(:,:) = mbkt(:,:) + 1       ! Initialization to the number of w ocean point (=2 over land)
[12377]266         DO_3DS_11_11( jpkm1, 2, -1 )
267            zus  = zcof * taum(ji,jj)
268            IF( zpelc(ji,jj,jk) > zus )   imlc(ji,jj) = jk
269         END_3D
[1492]270         !                               ! finite LC depth
[12377]271         DO_2D_11_11
272            zhlc(ji,jj) = gdepw(ji,jj,imlc(ji,jj),Kmm)
273         END_2D
[1705]274         zcof = 0.016 / SQRT( zrhoa * zcdrag )
[12377]275         DO_2D_00_00
276            zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift
277            zfr_i(ji,jj) = ( 1._wp - 4._wp * fr_i(ji,jj) ) * zus * zus * zus * tmask(ji,jj,1) ! zus > 0. ok
278            IF (zfr_i(ji,jj) < 0. ) zfr_i(ji,jj) = 0.
279         END_2D
280         DO_3D_00_00( 2, jpkm1 )
281            IF ( zfr_i(ji,jj) /= 0. ) THEN               
282               ! vertical velocity due to LC   
283               IF ( gdepw(ji,jj,jk,Kmm) - zhlc(ji,jj) < 0 .AND. wmask(ji,jj,jk) /= 0. ) THEN
284                  !                                           ! vertical velocity due to LC
285                  zwlc = rn_lc * SIN( rpi * gdepw(ji,jj,jk,Kmm) / zhlc(ji,jj) )   ! warning: optimization: zus^3 is in zfr_i
286                  !                                           ! TKE Langmuir circulation source term
[12489]287                  en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * zfr_i(ji,jj) * ( zwlc * zwlc * zwlc ) / zhlc(ji,jj)
[12377]288               ENDIF
289            ENDIF
290         END_3D
[1239]291         !
292      ENDIF
[1492]293      !
294      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
295      !                     !  Now Turbulent kinetic energy (output in en)
296      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
297      !                     ! Resolution of a tridiagonal linear system by a "methode de chasse"
298      !                     ! computation from level 2 to jpkm1  (e(1) already computed and e(jpk)=0 ).
299      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal
300      !
[9019]301      IF( nn_pdl == 1 ) THEN      !* Prandtl number = F( Ri )
[12377]302         DO_3D_00_00( 2, jpkm1 )
303            !                             ! local Richardson number
304            zri = MAX( rn2b(ji,jj,jk), 0._wp ) * p_avm(ji,jj,jk) / ( p_sh2(ji,jj,jk) + rn_bshear )
305            !                             ! inverse of Prandtl number
306            apdlr(ji,jj,jk) = MAX(  0.1_wp,  ri_cri / MAX( ri_cri , zri )  )
307         END_3D
[5656]308      ENDIF
[5836]309      !         
[12377]310      DO_3D_00_00( 2, jpkm1 )
311         zcof   = zfact1 * tmask(ji,jj,jk)
312         !                                   ! A minimum of 2.e-5 m2/s is imposed on TKE vertical
313         !                                   ! eddy coefficient (ensure numerical stability)
314         zzd_up = zcof * MAX(  p_avm(ji,jj,jk+1) + p_avm(ji,jj,jk  ) , 2.e-5_wp  )   &  ! upper diagonal
[13151]315            &          /    (  e3t(ji,jj,jk  ,Kmm)   &
316            &                * e3w(ji,jj,jk  ,Kmm)  )
[12377]317         zzd_lw = zcof * MAX(  p_avm(ji,jj,jk  ) + p_avm(ji,jj,jk-1) , 2.e-5_wp  )   &  ! lower diagonal
[13151]318            &          /    (  e3t(ji,jj,jk-1,Kmm)   &
319            &                * e3w(ji,jj,jk  ,Kmm)  )
[12377]320         !
321         zd_up(ji,jj,jk) = zzd_up            ! Matrix (zdiag, zd_up, zd_lw)
322         zd_lw(ji,jj,jk) = zzd_lw
323         zdiag(ji,jj,jk) = 1._wp - zzd_lw - zzd_up + zfact2 * dissl(ji,jj,jk) * wmask(ji,jj,jk)
324         !
325         !                                   ! right hand side in en
[12489]326         en(ji,jj,jk) = en(ji,jj,jk) + rn_Dt * (  p_sh2(ji,jj,jk)                        &   ! shear
[12377]327            &                                 - p_avt(ji,jj,jk) * rn2(ji,jj,jk)          &   ! stratification
328            &                                 + zfact3 * dissl(ji,jj,jk) * en(ji,jj,jk)  &   ! dissipation
329            &                                ) * wmask(ji,jj,jk)
330      END_3D
[5120]331      !                          !* Matrix inversion from level 2 (tke prescribed at level 1)
[12377]332      DO_3D_00_00( 3, jpkm1 )
333         zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1)
334      END_3D
335      DO_2D_00_00
336         zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke
337      END_2D
338      DO_3D_00_00( 3, jpkm1 )
339         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)
340      END_3D
341      DO_2D_00_00
342         en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1)
343      END_2D
344      DO_3DS_00_00( jpk-2, 2, -1 )
345         en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk)
346      END_3D
347      DO_3D_00_00( 2, jpkm1 )
348         en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * wmask(ji,jj,jk)
349      END_3D
[9019]350      !
[1492]351      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
352      !                            !  TKE due to surface and internal wave breaking
353      !                            !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[6140]354!!gm BUG : in the exp  remove the depth of ssh !!!
[12377]355!!gm       i.e. use gde3w in argument (gdepw(:,:,:,Kmm))
[6140]356     
357     
[2528]358      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction)
[12377]359         DO_3D_00_00( 2, jpkm1 )
360            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
361               &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
362         END_3D
[2528]363      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction)
[12377]364         DO_2D_00_00
365            jk = nmln(ji,jj)
366            en(ji,jj,jk) = en(ji,jj,jk) + rn_efr * en(ji,jj,1) * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
367               &                                 * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) )  * wmask(ji,jj,jk) * tmask(ji,jj,1)
368         END_2D
[2528]369      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability)
[12377]370         DO_3D_00_00( 2, jpkm1 )
371            ztx2 = utau(ji-1,jj  ) + utau(ji,jj)
372            zty2 = vtau(ji  ,jj-1) + vtau(ji,jj)
373            ztau = 0.5_wp * SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)    ! module of the mean stress
374            zdif = taum(ji,jj) - ztau                            ! mean of modulus - modulus of the mean
375            zdif = rhftau_scl * MAX( 0._wp, zdif + rhftau_add )  ! apply some modifications...
376            en(ji,jj,jk) = en(ji,jj,jk) + zbbrau * zdif * EXP( -gdepw(ji,jj,jk,Kmm) / htau(ji,jj) )   &
377               &                        * MAX(0.,1._wp - rn_eice *fr_i(ji,jj) ) * wmask(ji,jj,jk) * tmask(ji,jj,1)
378         END_3D
[1239]379      ENDIF
[1492]380      !
[1239]381   END SUBROUTINE tke_tke
382
[1492]383
[12377]384   SUBROUTINE tke_avn( Kbb, Kmm, p_avm, p_avt )
[1239]385      !!----------------------------------------------------------------------
[1492]386      !!                   ***  ROUTINE tke_avn  ***
[1239]387      !!
[1492]388      !! ** Purpose :   Compute the vertical eddy viscosity and diffusivity
389      !!
390      !! ** Method  :   At this stage, en, the now TKE, is known (computed in
391      !!              the tke_tke routine). First, the now mixing lenth is
392      !!      computed from en and the strafification (N^2), then the mixings
393      !!      coefficients are computed.
394      !!              - Mixing length : a first evaluation of the mixing lengh
395      !!      scales is:
396      !!                      mxl = sqrt(2*en) / N 
397      !!      where N is the brunt-vaisala frequency, with a minimum value set
[2528]398      !!      to rmxl_min (rn_mxl0) in the interior (surface) ocean.
[1492]399      !!        The mixing and dissipative length scale are bound as follow :
400      !!         nn_mxl=0 : mxl bounded by the distance to surface and bottom.
401      !!                        zmxld = zmxlm = mxl
402      !!         nn_mxl=1 : mxl bounded by the e3w and zmxld = zmxlm = mxl
403      !!         nn_mxl=2 : mxl bounded such that the vertical derivative of mxl is
404      !!                    less than 1 (|d/dz(mxl)|<1) and zmxld = zmxlm = mxl
405      !!         nn_mxl=3 : mxl is bounded from the surface to the bottom usings
406      !!                    |d/dz(xml)|<1 to obtain lup, and from the bottom to
407      !!                    the surface to obtain ldown. the resulting length
408      !!                    scales are:
409      !!                        zmxld = sqrt( lup * ldown )
410      !!                        zmxlm = min ( lup , ldown )
411      !!              - Vertical eddy viscosity and diffusivity:
412      !!                      avm = max( avtb, rn_ediff * zmxlm * en^1/2 )
413      !!                      avt = max( avmb, pdlr * avm ) 
414      !!      with pdlr=1 if nn_pdl=0, pdlr=1/pdl=F(Ri) otherwise.
415      !!
[9019]416      !! ** Action  : - avt, avm : now vertical eddy diffusivity and viscosity (w-point)
[1239]417      !!----------------------------------------------------------------------
[9019]418      USE zdf_oce , ONLY : en, avtb, avmb, avtb_2d   ! ocean vertical physics
419      !!
[12377]420      INTEGER                   , INTENT(in   ) ::   Kbb, Kmm       ! ocean time level indices
[9019]421      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_avm, p_avt   ! vertical eddy viscosity & diffusivity (w-points)
422      !
[2715]423      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[9019]424      REAL(wp) ::   zrn2, zraug, zcoef, zav   ! local scalars
425      REAL(wp) ::   zdku,   zdkv, zsqen       !   -      -
426      REAL(wp) ::   zemxl, zemlm, zemlp       !   -      -
427      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmxlm, zmxld   ! 3D workspace
[1239]428      !!--------------------------------------------------------------------
[3294]429      !
[1492]430      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
431      !                     !  Mixing length
432      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
433      !
434      !                     !* Buoyancy length scale: l=sqrt(2*e/n**2)
435      !
[5120]436      ! initialisation of interior minimum value (avoid a 2d loop with mikt)
[7753]437      zmxlm(:,:,:)  = rmxl_min   
438      zmxld(:,:,:)  = rmxl_min
[5120]439      !
[12489]440      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rho0*g)
441         zraug = vkarmn * 2.e5_wp / ( rho0 * grav )
[12377]442         DO_2D_00_00
443            zmxlm(ji,jj,1) = MAX( rn_mxl0, zraug * taum(ji,jj) * tmask(ji,jj,1) )
444         END_2D
[4990]445      ELSE
[7753]446         zmxlm(:,:,1) = rn_mxl0
[1239]447      ENDIF
448      !
[12377]449      DO_3D_00_00( 2, jpkm1 )
450         zrn2 = MAX( rn2(ji,jj,jk), rsmall )
451         zmxlm(ji,jj,jk) = MAX(  rmxl_min,  SQRT( 2._wp * en(ji,jj,jk) / zrn2 )  )
452      END_3D
[1492]453      !
454      !                     !* Physical limits for the mixing length
455      !
[7753]456      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value
457      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value
[1492]458      !
[1239]459      SELECT CASE ( nn_mxl )
460      !
[5836]461 !!gm Not sure of that coding for ISF....
[12377]462      ! where wmask = 0 set zmxlm == e3w(:,:,:,Kmm)
[1239]463      CASE ( 0 )           ! bounded by the distance to surface and bottom
[12377]464         DO_3D_00_00( 2, jpkm1 )
465            zemxl = MIN( gdepw(ji,jj,jk,Kmm) - gdepw(ji,jj,mikt(ji,jj),Kmm), zmxlm(ji,jj,jk),   &
466            &            gdepw(ji,jj,mbkt(ji,jj)+1,Kmm) - gdepw(ji,jj,jk,Kmm) )
467            ! wmask prevent zmxlm = 0 if jk = mikt(ji,jj)
[13151]468            zmxlm(ji,jj,jk) = zemxl * wmask(ji,jj,jk)   &
469               &            + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
470            zmxld(ji,jj,jk) = zemxl * wmask(ji,jj,jk)   &
471               &            + MIN( zmxlm(ji,jj,jk) , e3w(ji,jj,jk,Kmm) ) * (1 - wmask(ji,jj,jk))
[12377]472         END_3D
[1239]473         !
474      CASE ( 1 )           ! bounded by the vertical scale factor
[12377]475         DO_3D_00_00( 2, jpkm1 )
476            zemxl = MIN( e3w(ji,jj,jk,Kmm), zmxlm(ji,jj,jk) )
477            zmxlm(ji,jj,jk) = zemxl
478            zmxld(ji,jj,jk) = zemxl
479         END_3D
[1239]480         !
481      CASE ( 2 )           ! |dk[xml]| bounded by e3t :
[12377]482         DO_3D_00_00( 2, jpkm1 )
[13151]483            zmxlm(ji,jj,jk) =   &
484               &    MIN( zmxlm(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
[12377]485         END_3D
486         DO_3DS_00_00( jpkm1, 2, -1 )
487            zemxl = MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
488            zmxlm(ji,jj,jk) = zemxl
489            zmxld(ji,jj,jk) = zemxl
490         END_3D
[1239]491         !
492      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t :
[12377]493         DO_3D_00_00( 2, jpkm1 )
[13151]494            zmxld(ji,jj,jk) =    &
495               &    MIN( zmxld(ji,jj,jk-1) + e3t(ji,jj,jk-1,Kmm), zmxlm(ji,jj,jk) )
[12377]496         END_3D
497         DO_3DS_00_00( jpkm1, 2, -1 )
[13151]498            zmxlm(ji,jj,jk) =   &
499               &    MIN( zmxlm(ji,jj,jk+1) + e3t(ji,jj,jk+1,Kmm), zmxlm(ji,jj,jk) )
[12377]500         END_3D
501         DO_3D_00_00( 2, jpkm1 )
502            zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) )
503            zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) )
504            zmxlm(ji,jj,jk) = zemlm
505            zmxld(ji,jj,jk) = zemlp
506         END_3D
[1239]507         !
508      END SELECT
[1492]509      !
510      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[9019]511      !                     !  Vertical eddy viscosity and diffusivity  (avm and avt)
[1492]512      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
[12377]513      DO_3D_00_00( 1, jpkm1 )
514         zsqen = SQRT( en(ji,jj,jk) )
515         zav   = rn_ediff * zmxlm(ji,jj,jk) * zsqen
516         p_avm(ji,jj,jk) = MAX( zav,                  avmb(jk) ) * wmask(ji,jj,jk)
517         p_avt(ji,jj,jk) = MAX( zav, avtb_2d(ji,jj) * avtb(jk) ) * wmask(ji,jj,jk)
518         dissl(ji,jj,jk) = zsqen / zmxld(ji,jj,jk)
519      END_3D
[1492]520      !
521      !
522      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt
[12377]523         DO_3D_00_00( 2, jpkm1 )
[13151]524            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)
[12377]525         END_3D
[1239]526      ENDIF
[9019]527      !
[12377]528      IF(sn_cfctl%l_prtctl) THEN
[9440]529         CALL prt_ctl( tab3d_1=en   , clinfo1=' tke  - e: ', tab3d_2=p_avt, clinfo2=' t: ', kdim=jpk)
530         CALL prt_ctl( tab3d_1=p_avm, clinfo1=' tke  - m: ', kdim=jpk )
[1239]531      ENDIF
532      !
[1492]533   END SUBROUTINE tke_avn
[1239]534
[1492]535
[12377]536   SUBROUTINE zdf_tke_init( Kmm )
[1239]537      !!----------------------------------------------------------------------
[2528]538      !!                  ***  ROUTINE zdf_tke_init  ***
[1239]539      !!                     
540      !! ** Purpose :   Initialization of the vertical eddy diffivity and
[1492]541      !!              viscosity when using a tke turbulent closure scheme
[1239]542      !!
[1601]543      !! ** Method  :   Read the namzdf_tke namelist and check the parameters
[1492]544      !!              called at the first timestep (nit000)
[1239]545      !!
[1601]546      !! ** input   :   Namlist namzdf_tke
[1239]547      !!
548      !! ** Action  :   Increase by 1 the nstop flag is setting problem encounter
549      !!----------------------------------------------------------------------
[9019]550      USE zdf_oce , ONLY : ln_zdfiwm   ! Internal Wave Mixing flag
551      !!
[12377]552      INTEGER, INTENT(in) ::   Kmm          ! time level index
553      INTEGER             ::   ji, jj, jk   ! dummy loop indices
554      INTEGER             ::   ios
[1239]555      !!
[9019]556      NAMELIST/namzdf_tke/ rn_ediff, rn_ediss , rn_ebb , rn_emin  ,          &
557         &                 rn_emin0, rn_bshear, nn_mxl , ln_mxl0  ,          &
558         &                 rn_mxl0 , nn_pdl   , ln_drg , ln_lc    , rn_lc,   &
[9558]559         &                 nn_etau , nn_htau  , rn_efr , rn_eice 
[1239]560      !!----------------------------------------------------------------------
[2715]561      !
[4147]562      READ  ( numnam_ref, namzdf_tke, IOSTAT = ios, ERR = 901)
[11536]563901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_tke in reference namelist' )
[4147]564
565      READ  ( numnam_cfg, namzdf_tke, IOSTAT = ios, ERR = 902 )
[11536]566902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namzdf_tke in configuration namelist' )
[4624]567      IF(lwm) WRITE ( numond, namzdf_tke )
[2715]568      !
[2528]569      ri_cri   = 2._wp    / ( 2._wp + rn_ediss / rn_ediff )   ! resulting critical Richardson number
[2715]570      !
[1492]571      IF(lwp) THEN                    !* Control print
[1239]572         WRITE(numout,*)
[2528]573         WRITE(numout,*) 'zdf_tke_init : tke turbulent closure scheme - initialisation'
574         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]575         WRITE(numout,*) '   Namelist namzdf_tke : set tke mixing parameters'
[1705]576         WRITE(numout,*) '      coef. to compute avt                        rn_ediff  = ', rn_ediff
577         WRITE(numout,*) '      Kolmogoroff dissipation coef.               rn_ediss  = ', rn_ediss
578         WRITE(numout,*) '      tke surface input coef.                     rn_ebb    = ', rn_ebb
579         WRITE(numout,*) '      minimum value of tke                        rn_emin   = ', rn_emin
580         WRITE(numout,*) '      surface minimum value of tke                rn_emin0  = ', rn_emin0
[9019]581         WRITE(numout,*) '      prandl number flag                          nn_pdl    = ', nn_pdl
[1705]582         WRITE(numout,*) '      background shear (>0)                       rn_bshear = ', rn_bshear
583         WRITE(numout,*) '      mixing length type                          nn_mxl    = ', nn_mxl
[9019]584         WRITE(numout,*) '         surface mixing length = F(stress) or not    ln_mxl0   = ', ln_mxl0
585         WRITE(numout,*) '         surface  mixing length minimum value        rn_mxl0   = ', rn_mxl0
586         WRITE(numout,*) '      top/bottom friction forcing flag            ln_drg    = ', ln_drg
587         WRITE(numout,*) '      Langmuir cells parametrization              ln_lc     = ', ln_lc
588         WRITE(numout,*) '         coef to compute vertical velocity of LC     rn_lc  = ', rn_lc
[1705]589         WRITE(numout,*) '      test param. to add tke induced by wind      nn_etau   = ', nn_etau
[9019]590         WRITE(numout,*) '          type of tke penetration profile            nn_htau   = ', nn_htau
591         WRITE(numout,*) '          fraction of TKE that penetrates            rn_efr    = ', rn_efr
[9546]592         WRITE(numout,*) '          below sea-ice:  =0 ON                      rn_eice   = ', rn_eice
593         WRITE(numout,*) '          =4 OFF when ice fraction > 1/4   '
[9019]594         IF( ln_drg ) THEN
[9169]595            WRITE(numout,*)
[9019]596            WRITE(numout,*) '   Namelist namdrg_top/_bot:   used values:'
597            WRITE(numout,*) '      top    ocean cavity roughness (m)          rn_z0(_top)= ', r_z0_top
598            WRITE(numout,*) '      Bottom seafloor     roughness (m)          rn_z0(_bot)= ', r_z0_bot
599         ENDIF
600         WRITE(numout,*)
[9190]601         WRITE(numout,*) '   ==>>>   critical Richardson nb with your parameters  ri_cri = ', ri_cri
[9019]602         WRITE(numout,*)
[1239]603      ENDIF
[2715]604      !
[9019]605      IF( ln_zdfiwm ) THEN          ! Internal wave-driven mixing
606         rn_emin  = 1.e-10_wp             ! specific values of rn_emin & rmxl_min are used
607         rmxl_min = 1.e-03_wp             ! associated avt minimum = molecular salt diffusivity (10^-9 m2/s)
[9190]608         IF(lwp) WRITE(numout,*) '   ==>>>   Internal wave-driven mixing case:   force   rn_emin = 1.e-10 and rmxl_min = 1.e-3'
[9019]609      ELSE                          ! standard case : associated avt minimum = molecular viscosity (10^-6 m2/s)
610         rmxl_min = 1.e-6_wp / ( rn_ediff * SQRT( rn_emin ) )    ! resulting minimum length to recover molecular viscosity
[9190]611         IF(lwp) WRITE(numout,*) '   ==>>>   minimum mixing length with your parameters rmxl_min = ', rmxl_min
[9019]612      ENDIF
613      !
[2715]614      !                              ! allocate tke arrays
615      IF( zdf_tke_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_tke_init : unable to allocate arrays' )
616      !
[1492]617      !                               !* Check of some namelist values
[4990]618      IF( nn_mxl  < 0   .OR.  nn_mxl  > 3 )   CALL ctl_stop( 'bad flag: nn_mxl is  0, 1 or 2 ' )
619      IF( nn_pdl  < 0   .OR.  nn_pdl  > 1 )   CALL ctl_stop( 'bad flag: nn_pdl is  0 or 1    ' )
620      IF( nn_htau < 0   .OR.  nn_htau > 1 )   CALL ctl_stop( 'bad flag: nn_htau is 0, 1 or 2 ' )
[5407]621      IF( nn_etau == 3 .AND. .NOT. ln_cpl )   CALL ctl_stop( 'nn_etau == 3 : HF taum only known in coupled mode' )
[9019]622      !
[2528]623      IF( ln_mxl0 ) THEN
[9169]624         IF(lwp) WRITE(numout,*)
[9190]625         IF(lwp) WRITE(numout,*) '   ==>>>   use a surface mixing length = F(stress) :   set rn_mxl0 = rmxl_min'
[2528]626         rn_mxl0 = rmxl_min
627      ENDIF
628     
[12377]629      IF( nn_etau == 2  )   CALL zdf_mxl( nit000, Kmm )      ! Initialization of nmln
[1239]630
[1492]631      !                               !* depth of penetration of surface tke
632      IF( nn_etau /= 0 ) THEN     
[1601]633         SELECT CASE( nn_htau )             ! Choice of the depth of penetration
[2528]634         CASE( 0 )                                 ! constant depth penetration (here 10 meters)
[7753]635            htau(:,:) = 10._wp
[2528]636         CASE( 1 )                                 ! F(latitude) : 0.5m to 30m poleward of 40 degrees
[7753]637            htau(:,:) = MAX(  0.5_wp, MIN( 30._wp, 45._wp* ABS( SIN( rpi/180._wp * gphit(:,:) ) ) )   )           
[1492]638         END SELECT
639      ENDIF
[9019]640      !                                !* read or initialize all required files
641      CALL tke_rst( nit000, 'READ' )      ! (en, avt_k, avm_k, dissl)
[1239]642      !
[9367]643      IF( lwxios ) THEN
644         CALL iom_set_rstw_var_active('en')
645         CALL iom_set_rstw_var_active('avt_k')
646         CALL iom_set_rstw_var_active('avm_k')
647         CALL iom_set_rstw_var_active('dissl')
648      ENDIF
[2528]649   END SUBROUTINE zdf_tke_init
[1239]650
651
[1531]652   SUBROUTINE tke_rst( kt, cdrw )
[9019]653      !!---------------------------------------------------------------------
654      !!                   ***  ROUTINE tke_rst  ***
655      !!                     
656      !! ** Purpose :   Read or write TKE file (en) in restart file
657      !!
658      !! ** Method  :   use of IOM library
659      !!                if the restart does not contain TKE, en is either
660      !!                set to rn_emin or recomputed
661      !!----------------------------------------------------------------------
662      USE zdf_oce , ONLY : en, avt_k, avm_k   ! ocean vertical physics
663      !!
664      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
665      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
666      !
667      INTEGER ::   jit, jk              ! dummy loop indices
668      INTEGER ::   id1, id2, id3, id4   ! local integers
669      !!----------------------------------------------------------------------
670      !
671      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
672         !                                   ! ---------------
673         IF( ln_rstart ) THEN                   !* Read the restart file
674            id1 = iom_varid( numror, 'en'   , ldstop = .FALSE. )
675            id2 = iom_varid( numror, 'avt_k', ldstop = .FALSE. )
676            id3 = iom_varid( numror, 'avm_k', ldstop = .FALSE. )
677            id4 = iom_varid( numror, 'dissl', ldstop = .FALSE. )
678            !
679            IF( MIN( id1, id2, id3, id4 ) > 0 ) THEN      ! fields exist
[9367]680               CALL iom_get( numror, jpdom_autoglo, 'en'   , en   , ldxios = lrxios )
681               CALL iom_get( numror, jpdom_autoglo, 'avt_k', avt_k, ldxios = lrxios )
682               CALL iom_get( numror, jpdom_autoglo, 'avm_k', avm_k, ldxios = lrxios )
683               CALL iom_get( numror, jpdom_autoglo, 'dissl', dissl, ldxios = lrxios )
[9019]684            ELSE                                          ! start TKE from rest
[9169]685               IF(lwp) WRITE(numout,*)
[9190]686               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without TKE scheme, set en to background values'
[9019]687               en   (:,:,:) = rn_emin * wmask(:,:,:)
688               dissl(:,:,:) = 1.e-12_wp
689               ! avt_k, avm_k already set to the background value in zdf_phy_init
690            ENDIF
691         ELSE                                   !* Start from rest
[9169]692            IF(lwp) WRITE(numout,*)
[9190]693            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set en to the background value'
[9019]694            en   (:,:,:) = rn_emin * wmask(:,:,:)
695            dissl(:,:,:) = 1.e-12_wp
696            ! avt_k, avm_k already set to the background value in zdf_phy_init
697         ENDIF
698         !
699      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
700         !                                   ! -------------------
[9169]701         IF(lwp) WRITE(numout,*) '---- tke_rst ----'
[9367]702         IF( lwxios ) CALL iom_swap(      cwxios_context          ) 
703         CALL iom_rstput( kt, nitrst, numrow, 'en'   , en   , ldxios = lwxios )
704         CALL iom_rstput( kt, nitrst, numrow, 'avt_k', avt_k, ldxios = lwxios )
705         CALL iom_rstput( kt, nitrst, numrow, 'avm_k', avm_k, ldxios = lwxios )
706         CALL iom_rstput( kt, nitrst, numrow, 'dissl', dissl, ldxios = lwxios )
707         IF( lwxios ) CALL iom_swap(      cxios_context          )
[9019]708         !
709      ENDIF
710      !
[1531]711   END SUBROUTINE tke_rst
[1239]712
713   !!======================================================================
[1531]714END MODULE zdftke
Note: See TracBrowser for help on using the repository browser.