zdftke.F90 50 KB

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