zdfgls.F90 58 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243
  1. MODULE zdfgls
  2. !!======================================================================
  3. !! *** MODULE zdfgls ***
  4. !! Ocean physics: vertical mixing coefficient computed from the gls
  5. !! turbulent closure parameterization
  6. !!======================================================================
  7. !! History : 3.0 ! 2009-09 (G. Reffray) Original code
  8. !! 3.3 ! 2010-10 (C. Bricaud) Add in the reference
  9. !!----------------------------------------------------------------------
  10. #if defined key_zdfgls || defined key_esopa
  11. !!----------------------------------------------------------------------
  12. !! 'key_zdfgls' Generic Length Scale vertical physics
  13. !!----------------------------------------------------------------------
  14. !! zdf_gls : update momentum and tracer Kz from a gls scheme
  15. !! zdf_gls_init : initialization, namelist read, and parameters control
  16. !! gls_rst : read/write gls restart in ocean restart file
  17. !!----------------------------------------------------------------------
  18. USE oce ! ocean dynamics and active tracers
  19. USE dom_oce ! ocean space and time domain
  20. USE domvvl ! ocean space and time domain : variable volume layer
  21. USE zdf_oce ! ocean vertical physics
  22. USE zdfbfr ! bottom friction (only for rn_bfrz0)
  23. USE sbc_oce ! surface boundary condition: ocean
  24. USE phycst ! physical constants
  25. USE zdfmxl ! mixed layer
  26. USE lbclnk ! ocean lateral boundary conditions (or mpp link)
  27. USE lib_mpp ! MPP manager
  28. USE wrk_nemo ! work arrays
  29. USE prtctl ! Print control
  30. USE in_out_manager ! I/O manager
  31. USE iom ! I/O manager library
  32. USE timing ! Timing
  33. USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
  34. IMPLICIT NONE
  35. PRIVATE
  36. PUBLIC zdf_gls ! routine called in step module
  37. PUBLIC zdf_gls_init ! routine called in opa module
  38. PUBLIC gls_rst ! routine called in step module
  39. LOGICAL , PUBLIC, PARAMETER :: lk_zdfgls = .TRUE. !: TKE vertical mixing flag
  40. !
  41. REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: mxln !: now mixing length
  42. REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: zwall !: wall function
  43. REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ustars2 !: Squared surface velocity scale at T-points
  44. REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ustarb2 !: Squared bottom velocity scale at T-points
  45. ! !! ** Namelist namzdf_gls **
  46. LOGICAL :: ln_length_lim ! use limit on the dissipation rate under stable stratification (Galperin et al. 1988)
  47. LOGICAL :: ln_sigpsi ! Activate Burchard (2003) modification for k-eps closure & wave breaking mixing
  48. INTEGER :: nn_bc_surf ! surface boundary condition (=0/1)
  49. INTEGER :: nn_bc_bot ! bottom boundary condition (=0/1)
  50. INTEGER :: nn_z0_met ! Method for surface roughness computation
  51. INTEGER :: nn_stab_func ! stability functions G88, KC or Canuto (=0/1/2)
  52. INTEGER :: nn_clos ! closure 0/1/2/3 MY82/k-eps/k-w/gen
  53. REAL(wp) :: rn_clim_galp ! Holt 2008 value for k-eps: 0.267
  54. REAL(wp) :: rn_epsmin ! minimum value of dissipation (m2/s3)
  55. REAL(wp) :: rn_emin ! minimum value of TKE (m2/s2)
  56. REAL(wp) :: rn_charn ! Charnock constant for surface breaking waves mixing : 1400. (standard) or 2.e5 (Stacey value)
  57. REAL(wp) :: rn_crban ! Craig and Banner constant for surface breaking waves mixing
  58. REAL(wp) :: rn_hsro ! Minimum surface roughness
  59. REAL(wp) :: rn_frac_hs ! Fraction of wave height as surface roughness (if nn_z0_met > 1)
  60. REAL(wp) :: rcm_sf = 0.73_wp ! Shear free turbulence parameters
  61. REAL(wp) :: ra_sf = -2.0_wp ! Must be negative -2 < ra_sf < -1
  62. REAL(wp) :: rl_sf = 0.2_wp ! 0 <rl_sf<vkarmn
  63. REAL(wp) :: rghmin = -0.28_wp
  64. REAL(wp) :: rgh0 = 0.0329_wp
  65. REAL(wp) :: rghcri = 0.03_wp
  66. REAL(wp) :: ra1 = 0.92_wp
  67. REAL(wp) :: ra2 = 0.74_wp
  68. REAL(wp) :: rb1 = 16.60_wp
  69. REAL(wp) :: rb2 = 10.10_wp
  70. REAL(wp) :: re2 = 1.33_wp
  71. REAL(wp) :: rl1 = 0.107_wp
  72. REAL(wp) :: rl2 = 0.0032_wp
  73. REAL(wp) :: rl3 = 0.0864_wp
  74. REAL(wp) :: rl4 = 0.12_wp
  75. REAL(wp) :: rl5 = 11.9_wp
  76. REAL(wp) :: rl6 = 0.4_wp
  77. REAL(wp) :: rl7 = 0.0_wp
  78. REAL(wp) :: rl8 = 0.48_wp
  79. REAL(wp) :: rm1 = 0.127_wp
  80. REAL(wp) :: rm2 = 0.00336_wp
  81. REAL(wp) :: rm3 = 0.0906_wp
  82. REAL(wp) :: rm4 = 0.101_wp
  83. REAL(wp) :: rm5 = 11.2_wp
  84. REAL(wp) :: rm6 = 0.4_wp
  85. REAL(wp) :: rm7 = 0.0_wp
  86. REAL(wp) :: rm8 = 0.318_wp
  87. REAL(wp) :: rtrans = 0.1_wp
  88. REAL(wp) :: rc02, rc02r, rc03, rc04 ! coefficients deduced from above parameters
  89. REAL(wp) :: rsbc_tke1, rsbc_tke2, rfact_tke ! - - - -
  90. REAL(wp) :: rsbc_psi1, rsbc_psi2, rfact_psi ! - - - -
  91. REAL(wp) :: rsbc_zs1, rsbc_zs2 ! - - - -
  92. REAL(wp) :: rc0, rc2, rc3, rf6, rcff, rc_diff ! - - - -
  93. REAL(wp) :: rs0, rs1, rs2, rs4, rs5, rs6 ! - - - -
  94. REAL(wp) :: rd0, rd1, rd2, rd3, rd4, rd5 ! - - - -
  95. REAL(wp) :: rsc_tke, rsc_psi, rpsi1, rpsi2, rpsi3, rsc_psi0 ! - - - -
  96. REAL(wp) :: rpsi3m, rpsi3p, rpp, rmm, rnn ! - - - -
  97. !! * Substitutions
  98. # include "domzgr_substitute.h90"
  99. # include "vectopt_loop_substitute.h90"
  100. !!----------------------------------------------------------------------
  101. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  102. !! $Id: zdfgls.F90 5610 2015-07-17 17:42:27Z cetlod $
  103. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  104. !!----------------------------------------------------------------------
  105. CONTAINS
  106. INTEGER FUNCTION zdf_gls_alloc()
  107. !!----------------------------------------------------------------------
  108. !! *** FUNCTION zdf_gls_alloc ***
  109. !!----------------------------------------------------------------------
  110. ALLOCATE( mxln(jpi,jpj,jpk), zwall(jpi,jpj,jpk) , &
  111. & ustars2(jpi,jpj) , ustarb2(jpi,jpj) , STAT= zdf_gls_alloc )
  112. !
  113. IF( lk_mpp ) CALL mpp_sum ( zdf_gls_alloc )
  114. IF( zdf_gls_alloc /= 0 ) CALL ctl_warn('zdf_gls_alloc: failed to allocate arrays')
  115. END FUNCTION zdf_gls_alloc
  116. SUBROUTINE zdf_gls( kt )
  117. !!----------------------------------------------------------------------
  118. !! *** ROUTINE zdf_gls ***
  119. !!
  120. !! ** Purpose : Compute the vertical eddy viscosity and diffusivity
  121. !! coefficients using the GLS turbulent closure scheme.
  122. !!----------------------------------------------------------------------
  123. INTEGER, INTENT(in) :: kt ! ocean time step
  124. INTEGER :: ji, jj, jk, ibot, ibotm1, dir ! dummy loop arguments
  125. REAL(wp) :: zesh2, zsigpsi, zcoef, zex1, zex2 ! local scalars
  126. REAL(wp) :: ztx2, zty2, zup, zdown, zcof ! - -
  127. REAL(wp) :: zratio, zrn2, zflxb, sh ! - -
  128. REAL(wp) :: prod, buoy, diss, zdiss, sm ! - -
  129. REAL(wp) :: gh, gm, shr, dif, zsqen, zav ! - -
  130. REAL(wp), POINTER, DIMENSION(:,: ) :: zdep
  131. REAL(wp), POINTER, DIMENSION(:,: ) :: zkar
  132. REAL(wp), POINTER, DIMENSION(:,: ) :: zflxs ! Turbulence fluxed induced by internal waves
  133. REAL(wp), POINTER, DIMENSION(:,: ) :: zhsro ! Surface roughness (surface waves)
  134. REAL(wp), POINTER, DIMENSION(:,:,:) :: eb ! tke at time before
  135. REAL(wp), POINTER, DIMENSION(:,:,:) :: mxlb ! mixing length at time before
  136. REAL(wp), POINTER, DIMENSION(:,:,:) :: shear ! vertical shear
  137. REAL(wp), POINTER, DIMENSION(:,:,:) :: eps ! dissipation rate
  138. REAL(wp), POINTER, DIMENSION(:,:,:) :: zwall_psi ! Wall function use in the wb case (ln_sigpsi)
  139. REAL(wp), POINTER, DIMENSION(:,:,:) :: psi ! psi at time now
  140. REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_a ! element of the first matrix diagonal
  141. REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_b ! element of the second matrix diagonal
  142. REAL(wp), POINTER, DIMENSION(:,:,:) :: z_elem_c ! element of the third matrix diagonal
  143. !!--------------------------------------------------------------------
  144. !
  145. IF( nn_timing == 1 ) CALL timing_start('zdf_gls')
  146. !
  147. CALL wrk_alloc( jpi,jpj, zdep, zkar, zflxs, zhsro )
  148. CALL wrk_alloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi )
  149. ! Preliminary computing
  150. ustars2 = 0._wp ; ustarb2 = 0._wp ; psi = 0._wp ; zwall_psi = 0._wp
  151. IF( kt /= nit000 ) THEN ! restore before value to compute tke
  152. avt (:,:,:) = avt_k (:,:,:)
  153. avm (:,:,:) = avm_k (:,:,:)
  154. avmu(:,:,:) = avmu_k(:,:,:)
  155. avmv(:,:,:) = avmv_k(:,:,:)
  156. ENDIF
  157. ! Compute surface and bottom friction at T-points
  158. !CDIR NOVERRCHK
  159. DO jj = 2, jpjm1
  160. !CDIR NOVERRCHK
  161. DO ji = fs_2, fs_jpim1 ! vector opt.
  162. !
  163. ! surface friction
  164. ustars2(ji,jj) = r1_rau0 * taum(ji,jj) * tmask(ji,jj,1)
  165. !
  166. ! bottom friction (explicit before friction)
  167. ! Note that we chose here not to bound the friction as in dynbfr)
  168. ztx2 = ( bfrua(ji,jj) * ub(ji,jj,mbku(ji,jj)) + bfrua(ji-1,jj) * ub(ji-1,jj,mbku(ji-1,jj)) ) &
  169. & * ( 1._wp - 0.5_wp * umask(ji,jj,1) * umask(ji-1,jj,1) )
  170. zty2 = ( bfrva(ji,jj) * vb(ji,jj,mbkv(ji,jj)) + bfrva(ji,jj-1) * vb(ji,jj-1,mbkv(ji,jj-1)) ) &
  171. & * ( 1._wp - 0.5_wp * vmask(ji,jj,1) * vmask(ji,jj-1,1) )
  172. ustarb2(ji,jj) = SQRT( ztx2 * ztx2 + zty2 * zty2 ) * tmask(ji,jj,1)
  173. END DO
  174. END DO
  175. ! Set surface roughness length
  176. SELECT CASE ( nn_z0_met )
  177. !
  178. CASE ( 0 ) ! Constant roughness
  179. zhsro(:,:) = rn_hsro
  180. CASE ( 1 ) ! Standard Charnock formula
  181. zhsro(:,:) = MAX(rsbc_zs1 * ustars2(:,:), rn_hsro)
  182. CASE ( 2 ) ! Roughness formulae according to Rascle et al., Ocean Modelling (2008)
  183. zdep(:,:) = 30.*TANH(2.*0.3/(28.*SQRT(MAX(ustars2(:,:),rsmall)))) ! Wave age (eq. 10)
  184. zhsro(:,:) = MAX(rsbc_zs2 * ustars2(:,:) * zdep(:,:)**1.5, rn_hsro) ! zhsro = rn_frac_hs * Hsw (eq. 11)
  185. !
  186. END SELECT
  187. ! Compute shear and dissipation rate
  188. DO jk = 2, jpkm1
  189. DO jj = 2, jpjm1
  190. DO ji = fs_2, fs_jpim1 ! vector opt.
  191. avmu(ji,jj,jk) = avmu(ji,jj,jk) * ( un(ji,jj,jk-1) - un(ji,jj,jk) ) &
  192. & * ( ub(ji,jj,jk-1) - ub(ji,jj,jk) ) &
  193. & / ( fse3uw_n(ji,jj,jk) &
  194. & * fse3uw_b(ji,jj,jk) )
  195. avmv(ji,jj,jk) = avmv(ji,jj,jk) * ( vn(ji,jj,jk-1) - vn(ji,jj,jk) ) &
  196. & * ( vb(ji,jj,jk-1) - vb(ji,jj,jk) ) &
  197. & / ( fse3vw_n(ji,jj,jk) &
  198. & * fse3vw_b(ji,jj,jk) )
  199. eps(ji,jj,jk) = rc03 * en(ji,jj,jk) * SQRT(en(ji,jj,jk)) / mxln(ji,jj,jk)
  200. END DO
  201. END DO
  202. END DO
  203. !
  204. ! Lateral boundary conditions (avmu,avmv) (sign unchanged)
  205. CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. )
  206. ! Save tke at before time step
  207. eb (:,:,:) = en (:,:,:)
  208. mxlb(:,:,:) = mxln(:,:,:)
  209. IF( nn_clos == 0 ) THEN ! Mellor-Yamada
  210. DO jk = 2, jpkm1
  211. DO jj = 2, jpjm1
  212. DO ji = fs_2, fs_jpim1 ! vector opt.
  213. zup = mxln(ji,jj,jk) * fsdepw(ji,jj,mbkt(ji,jj)+1)
  214. zdown = vkarmn * fsdepw(ji,jj,jk) * ( -fsdepw(ji,jj,jk) + fsdepw(ji,jj,mbkt(ji,jj)+1) )
  215. zcoef = ( zup / MAX( zdown, rsmall ) )
  216. zwall (ji,jj,jk) = ( 1._wp + re2 * zcoef*zcoef ) * tmask(ji,jj,jk)
  217. END DO
  218. END DO
  219. END DO
  220. ENDIF
  221. !!---------------------------------!!
  222. !! Equation to prognostic k !!
  223. !!---------------------------------!!
  224. !
  225. ! Now Turbulent kinetic energy (output in en)
  226. ! -------------------------------
  227. ! Resolution of a tridiagonal linear system by a "methode de chasse"
  228. ! computation from level 2 to jpkm1 (e(1) computed after and e(jpk)=0 ).
  229. ! The surface boundary condition are set after
  230. ! The bottom boundary condition are also set after. In standard e(bottom)=0.
  231. ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal
  232. ! Warning : after this step, en : right hand side of the matrix
  233. DO jk = 2, jpkm1
  234. DO jj = 2, jpjm1
  235. DO ji = fs_2, fs_jpim1 ! vector opt.
  236. !
  237. ! shear prod. at w-point weightened by mask
  238. shear(ji,jj,jk) = ( avmu(ji-1,jj,jk) + avmu(ji,jj,jk) ) / MAX( 1.e0 , umask(ji-1,jj,jk) + umask(ji,jj,jk) ) &
  239. & + ( avmv(ji,jj-1,jk) + avmv(ji,jj,jk) ) / MAX( 1.e0 , vmask(ji,jj-1,jk) + vmask(ji,jj,jk) )
  240. !
  241. ! stratif. destruction
  242. buoy = - avt(ji,jj,jk) * rn2(ji,jj,jk)
  243. !
  244. ! shear prod. - stratif. destruction
  245. diss = eps(ji,jj,jk)
  246. !
  247. dir = 0.5_wp + SIGN( 0.5_wp, shear(ji,jj,jk) + buoy ) ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0)
  248. !
  249. zesh2 = dir*(shear(ji,jj,jk)+buoy)+(1._wp-dir)*shear(ji,jj,jk) ! production term
  250. zdiss = dir*(diss/en(ji,jj,jk)) +(1._wp-dir)*(diss-buoy)/en(ji,jj,jk) ! dissipation term
  251. !
  252. ! Compute a wall function from 1. to rsc_psi*zwall/rsc_psi0
  253. ! Note that as long that Dirichlet boundary conditions are NOT set at the first and last levels (GOTM style)
  254. ! there is no need to set a boundary condition for zwall_psi at the top and bottom boundaries.
  255. ! Otherwise, this should be rsc_psi/rsc_psi0
  256. IF( ln_sigpsi ) THEN
  257. zsigpsi = MIN( 1._wp, zesh2 / eps(ji,jj,jk) ) ! 0. <= zsigpsi <= 1.
  258. zwall_psi(ji,jj,jk) = rsc_psi / &
  259. & ( zsigpsi * rsc_psi + (1._wp-zsigpsi) * rsc_psi0 / MAX( zwall(ji,jj,jk), 1._wp ) )
  260. ELSE
  261. zwall_psi(ji,jj,jk) = 1._wp
  262. ENDIF
  263. !
  264. ! building the matrix
  265. zcof = rfact_tke * tmask(ji,jj,jk)
  266. !
  267. ! lower diagonal
  268. z_elem_a(ji,jj,jk) = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) &
  269. & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) )
  270. !
  271. ! upper diagonal
  272. z_elem_c(ji,jj,jk) = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) &
  273. & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) )
  274. !
  275. ! diagonal
  276. z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) &
  277. & + rdt * zdiss * tmask(ji,jj,jk)
  278. !
  279. ! right hand side in en
  280. en(ji,jj,jk) = en(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk)
  281. END DO
  282. END DO
  283. END DO
  284. !
  285. z_elem_b(:,:,jpk) = 1._wp
  286. !
  287. ! Set surface condition on zwall_psi (1 at the bottom)
  288. zwall_psi(:,:,1) = zwall_psi(:,:,2)
  289. zwall_psi(:,:,jpk) = 1.
  290. !
  291. ! Surface boundary condition on tke
  292. ! ---------------------------------
  293. !
  294. SELECT CASE ( nn_bc_surf )
  295. !
  296. CASE ( 0 ) ! Dirichlet case
  297. ! First level
  298. en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1)**(2._wp/3._wp)
  299. en(:,:,1) = MAX(en(:,:,1), rn_emin)
  300. z_elem_a(:,:,1) = en(:,:,1)
  301. z_elem_c(:,:,1) = 0._wp
  302. z_elem_b(:,:,1) = 1._wp
  303. !
  304. ! One level below
  305. en(:,:,2) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1 * ((zhsro(:,:)+fsdepw(:,:,2)) &
  306. & / zhsro(:,:) )**(1.5_wp*ra_sf))**(2._wp/3._wp)
  307. en(:,:,2) = MAX(en(:,:,2), rn_emin )
  308. z_elem_a(:,:,2) = 0._wp
  309. z_elem_c(:,:,2) = 0._wp
  310. z_elem_b(:,:,2) = 1._wp
  311. !
  312. !
  313. CASE ( 1 ) ! Neumann boundary condition on d(e)/dz
  314. !
  315. ! Dirichlet conditions at k=1
  316. en(:,:,1) = rc02r * ustars2(:,:) * (1._wp + rsbc_tke1)**(2._wp/3._wp)
  317. en(:,:,1) = MAX(en(:,:,1), rn_emin)
  318. z_elem_a(:,:,1) = en(:,:,1)
  319. z_elem_c(:,:,1) = 0._wp
  320. z_elem_b(:,:,1) = 1._wp
  321. !
  322. ! at k=2, set de/dz=Fw
  323. !cbr
  324. z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
  325. z_elem_a(:,:,2) = 0._wp
  326. zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1.-exp(-rtrans*fsdept(:,:,1)/zhsro(:,:)) ))
  327. zflxs(:,:) = rsbc_tke2 * ustars2(:,:)**1.5_wp * zkar(:,:) &
  328. & * ((zhsro(:,:)+fsdept(:,:,1))/zhsro(:,:) )**(1.5_wp*ra_sf)
  329. en(:,:,2) = en(:,:,2) + zflxs(:,:)/fse3w(:,:,2)
  330. !
  331. !
  332. END SELECT
  333. ! Bottom boundary condition on tke
  334. ! --------------------------------
  335. !
  336. SELECT CASE ( nn_bc_bot )
  337. !
  338. CASE ( 0 ) ! Dirichlet
  339. ! ! en(ibot) = u*^2 / Co2 and mxln(ibot) = rn_lmin
  340. ! ! Balance between the production and the dissipation terms
  341. !CDIR NOVERRCHK
  342. DO jj = 2, jpjm1
  343. !CDIR NOVERRCHK
  344. DO ji = fs_2, fs_jpim1 ! vector opt.
  345. ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point
  346. ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1
  347. !
  348. ! Bottom level Dirichlet condition:
  349. z_elem_a(ji,jj,ibot ) = 0._wp
  350. z_elem_c(ji,jj,ibot ) = 0._wp
  351. z_elem_b(ji,jj,ibot ) = 1._wp
  352. en(ji,jj,ibot ) = MAX( rc02r * ustarb2(ji,jj), rn_emin )
  353. !
  354. ! Just above last level, Dirichlet condition again
  355. z_elem_a(ji,jj,ibotm1) = 0._wp
  356. z_elem_c(ji,jj,ibotm1) = 0._wp
  357. z_elem_b(ji,jj,ibotm1) = 1._wp
  358. en(ji,jj,ibotm1) = MAX( rc02r * ustarb2(ji,jj), rn_emin )
  359. END DO
  360. END DO
  361. !
  362. CASE ( 1 ) ! Neumman boundary condition
  363. !
  364. !CDIR NOVERRCHK
  365. DO jj = 2, jpjm1
  366. !CDIR NOVERRCHK
  367. DO ji = fs_2, fs_jpim1 ! vector opt.
  368. ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point
  369. ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1
  370. !
  371. ! Bottom level Dirichlet condition:
  372. z_elem_a(ji,jj,ibot) = 0._wp
  373. z_elem_c(ji,jj,ibot) = 0._wp
  374. z_elem_b(ji,jj,ibot) = 1._wp
  375. en(ji,jj,ibot) = MAX( rc02r * ustarb2(ji,jj), rn_emin )
  376. !
  377. ! Just above last level: Neumann condition
  378. z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1) ! Remove z_elem_c from z_elem_b
  379. z_elem_c(ji,jj,ibotm1) = 0._wp
  380. END DO
  381. END DO
  382. !
  383. END SELECT
  384. ! Matrix inversion (en prescribed at surface and the bottom)
  385. ! ----------------------------------------------------------
  386. !
  387. DO jk = 2, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
  388. DO jj = 2, jpjm1
  389. DO ji = fs_2, fs_jpim1 ! vector opt.
  390. z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1)
  391. END DO
  392. END DO
  393. END DO
  394. DO jk = 2, jpk ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
  395. DO jj = 2, jpjm1
  396. DO ji = fs_2, fs_jpim1 ! vector opt.
  397. z_elem_a(ji,jj,jk) = en(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1)
  398. END DO
  399. END DO
  400. END DO
  401. DO jk = jpk-1, 2, -1 ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
  402. DO jj = 2, jpjm1
  403. DO ji = fs_2, fs_jpim1 ! vector opt.
  404. en(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * en(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk)
  405. END DO
  406. END DO
  407. END DO
  408. ! ! set the minimum value of tke
  409. en(:,:,:) = MAX( en(:,:,:), rn_emin )
  410. !!----------------------------------------!!
  411. !! Solve prognostic equation for psi !!
  412. !!----------------------------------------!!
  413. ! Set psi to previous time step value
  414. !
  415. SELECT CASE ( nn_clos )
  416. !
  417. CASE( 0 ) ! k-kl (Mellor-Yamada)
  418. DO jk = 2, jpkm1
  419. DO jj = 2, jpjm1
  420. DO ji = fs_2, fs_jpim1 ! vector opt.
  421. psi(ji,jj,jk) = eb(ji,jj,jk) * mxlb(ji,jj,jk)
  422. END DO
  423. END DO
  424. END DO
  425. !
  426. CASE( 1 ) ! k-eps
  427. DO jk = 2, jpkm1
  428. DO jj = 2, jpjm1
  429. DO ji = fs_2, fs_jpim1 ! vector opt.
  430. psi(ji,jj,jk) = eps(ji,jj,jk)
  431. END DO
  432. END DO
  433. END DO
  434. !
  435. CASE( 2 ) ! k-w
  436. DO jk = 2, jpkm1
  437. DO jj = 2, jpjm1
  438. DO ji = fs_2, fs_jpim1 ! vector opt.
  439. psi(ji,jj,jk) = SQRT( eb(ji,jj,jk) ) / ( rc0 * mxlb(ji,jj,jk) )
  440. END DO
  441. END DO
  442. END DO
  443. !
  444. CASE( 3 ) ! generic
  445. DO jk = 2, jpkm1
  446. DO jj = 2, jpjm1
  447. DO ji = fs_2, fs_jpim1 ! vector opt.
  448. psi(ji,jj,jk) = rc02 * eb(ji,jj,jk) * mxlb(ji,jj,jk)**rnn
  449. END DO
  450. END DO
  451. END DO
  452. !
  453. END SELECT
  454. !
  455. ! Now gls (output in psi)
  456. ! -------------------------------
  457. ! Resolution of a tridiagonal linear system by a "methode de chasse"
  458. ! computation from level 2 to jpkm1 (e(1) already computed and e(jpk)=0 ).
  459. ! z_elem_b : diagonal z_elem_c : upper diagonal z_elem_a : lower diagonal
  460. ! Warning : after this step, en : right hand side of the matrix
  461. DO jk = 2, jpkm1
  462. DO jj = 2, jpjm1
  463. DO ji = fs_2, fs_jpim1 ! vector opt.
  464. !
  465. ! psi / k
  466. zratio = psi(ji,jj,jk) / eb(ji,jj,jk)
  467. !
  468. ! psi3+ : stable : B=-KhN²<0 => N²>0 if rn2>0 dir = 1 (stable) otherwise dir = 0 (unstable)
  469. dir = 0.5_wp + SIGN( 0.5_wp, rn2(ji,jj,jk) )
  470. !
  471. rpsi3 = dir * rpsi3m + ( 1._wp - dir ) * rpsi3p
  472. !
  473. ! shear prod. - stratif. destruction
  474. prod = rpsi1 * zratio * shear(ji,jj,jk)
  475. !
  476. ! stratif. destruction
  477. buoy = rpsi3 * zratio * (- avt(ji,jj,jk) * rn2(ji,jj,jk) )
  478. !
  479. ! shear prod. - stratif. destruction
  480. diss = rpsi2 * zratio * zwall(ji,jj,jk) * eps(ji,jj,jk)
  481. !
  482. dir = 0.5_wp + SIGN( 0.5_wp, prod + buoy ) ! dir =1(=0) if shear(ji,jj,jk)+buoy >0(<0)
  483. !
  484. zesh2 = dir * ( prod + buoy ) + (1._wp - dir ) * prod ! production term
  485. zdiss = dir * ( diss / psi(ji,jj,jk) ) + (1._wp - dir ) * (diss-buoy) / psi(ji,jj,jk) ! dissipation term
  486. !
  487. ! building the matrix
  488. zcof = rfact_psi * zwall_psi(ji,jj,jk) * tmask(ji,jj,jk)
  489. ! lower diagonal
  490. z_elem_a(ji,jj,jk) = zcof * ( avm (ji,jj,jk ) + avm (ji,jj,jk-1) ) &
  491. & / ( fse3t(ji,jj,jk-1) * fse3w(ji,jj,jk ) )
  492. ! upper diagonal
  493. z_elem_c(ji,jj,jk) = zcof * ( avm (ji,jj,jk+1) + avm (ji,jj,jk ) ) &
  494. & / ( fse3t(ji,jj,jk ) * fse3w(ji,jj,jk) )
  495. ! diagonal
  496. z_elem_b(ji,jj,jk) = 1._wp - z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) &
  497. & + rdt * zdiss * tmask(ji,jj,jk)
  498. !
  499. ! right hand side in psi
  500. psi(ji,jj,jk) = psi(ji,jj,jk) + rdt * zesh2 * tmask(ji,jj,jk)
  501. END DO
  502. END DO
  503. END DO
  504. !
  505. z_elem_b(:,:,jpk) = 1._wp
  506. ! Surface boundary condition on psi
  507. ! ---------------------------------
  508. !
  509. SELECT CASE ( nn_bc_surf )
  510. !
  511. CASE ( 0 ) ! Dirichlet boundary conditions
  512. !
  513. ! Surface value
  514. zdep(:,:) = zhsro(:,:) * rl_sf ! Cosmetic
  515. psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1)
  516. z_elem_a(:,:,1) = psi(:,:,1)
  517. z_elem_c(:,:,1) = 0._wp
  518. z_elem_b(:,:,1) = 1._wp
  519. !
  520. ! One level below
  521. zkar(:,:) = (rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans*fsdepw(:,:,2)/zhsro(:,:) )))
  522. zdep(:,:) = (zhsro(:,:) + fsdepw(:,:,2)) * zkar(:,:)
  523. psi (:,:,2) = rc0**rpp * en(:,:,2)**rmm * zdep(:,:)**rnn * tmask(:,:,1)
  524. z_elem_a(:,:,2) = 0._wp
  525. z_elem_c(:,:,2) = 0._wp
  526. z_elem_b(:,:,2) = 1._wp
  527. !
  528. !
  529. CASE ( 1 ) ! Neumann boundary condition on d(psi)/dz
  530. !
  531. ! Surface value: Dirichlet
  532. zdep(:,:) = zhsro(:,:) * rl_sf
  533. psi (:,:,1) = rc0**rpp * en(:,:,1)**rmm * zdep(:,:)**rnn * tmask(:,:,1)
  534. z_elem_a(:,:,1) = psi(:,:,1)
  535. z_elem_c(:,:,1) = 0._wp
  536. z_elem_b(:,:,1) = 1._wp
  537. !
  538. ! Neumann condition at k=2
  539. z_elem_b(:,:,2) = z_elem_b(:,:,2) + z_elem_a(:,:,2) ! Remove z_elem_a from z_elem_b
  540. z_elem_a(:,:,2) = 0._wp
  541. !
  542. ! Set psi vertical flux at the surface:
  543. zkar(:,:) = rl_sf + (vkarmn-rl_sf)*(1._wp-exp(-rtrans*fsdept(:,:,1)/zhsro(:,:) )) ! Lengh scale slope
  544. zdep(:,:) = ((zhsro(:,:) + fsdept(:,:,1)) / zhsro(:,:))**(rmm*ra_sf)
  545. zflxs(:,:) = (rnn + rsbc_tke1 * (rnn + rmm*ra_sf) * zdep(:,:))*(1._wp + rsbc_tke1*zdep(:,:))**(2._wp*rmm/3._wp-1_wp)
  546. zdep(:,:) = rsbc_psi1 * (zwall_psi(:,:,1)*avm(:,:,1)+zwall_psi(:,:,2)*avm(:,:,2)) * &
  547. & ustars2(:,:)**rmm * zkar(:,:)**rnn * (zhsro(:,:) + fsdept(:,:,1))**(rnn-1.)
  548. zflxs(:,:) = zdep(:,:) * zflxs(:,:)
  549. psi(:,:,2) = psi(:,:,2) + zflxs(:,:) / fse3w(:,:,2)
  550. !
  551. !
  552. END SELECT
  553. ! Bottom boundary condition on psi
  554. ! --------------------------------
  555. !
  556. SELECT CASE ( nn_bc_bot )
  557. !
  558. !
  559. CASE ( 0 ) ! Dirichlet
  560. ! ! en(ibot) = u*^2 / Co2 and mxln(ibot) = vkarmn * rn_bfrz0
  561. ! ! Balance between the production and the dissipation terms
  562. !CDIR NOVERRCHK
  563. DO jj = 2, jpjm1
  564. !CDIR NOVERRCHK
  565. DO ji = fs_2, fs_jpim1 ! vector opt.
  566. ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point
  567. ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1
  568. zdep(ji,jj) = vkarmn * rn_bfrz0
  569. psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn
  570. z_elem_a(ji,jj,ibot) = 0._wp
  571. z_elem_c(ji,jj,ibot) = 0._wp
  572. z_elem_b(ji,jj,ibot) = 1._wp
  573. !
  574. ! Just above last level, Dirichlet condition again (GOTM like)
  575. zdep(ji,jj) = vkarmn * ( rn_bfrz0 + fse3t(ji,jj,ibotm1) )
  576. psi (ji,jj,ibotm1) = rc0**rpp * en(ji,jj,ibot )**rmm * zdep(ji,jj)**rnn
  577. z_elem_a(ji,jj,ibotm1) = 0._wp
  578. z_elem_c(ji,jj,ibotm1) = 0._wp
  579. z_elem_b(ji,jj,ibotm1) = 1._wp
  580. END DO
  581. END DO
  582. !
  583. CASE ( 1 ) ! Neumman boundary condition
  584. !
  585. !CDIR NOVERRCHK
  586. DO jj = 2, jpjm1
  587. !CDIR NOVERRCHK
  588. DO ji = fs_2, fs_jpim1 ! vector opt.
  589. ibot = mbkt(ji,jj) + 1 ! k bottom level of w-point
  590. ibotm1 = mbkt(ji,jj) ! k-1 bottom level of w-point but >=1
  591. !
  592. ! Bottom level Dirichlet condition:
  593. zdep(ji,jj) = vkarmn * rn_bfrz0
  594. psi (ji,jj,ibot) = rc0**rpp * en(ji,jj,ibot)**rmm * zdep(ji,jj)**rnn
  595. !
  596. z_elem_a(ji,jj,ibot) = 0._wp
  597. z_elem_c(ji,jj,ibot) = 0._wp
  598. z_elem_b(ji,jj,ibot) = 1._wp
  599. !
  600. ! Just above last level: Neumann condition with flux injection
  601. z_elem_b(ji,jj,ibotm1) = z_elem_b(ji,jj,ibotm1) + z_elem_c(ji,jj,ibotm1) ! Remove z_elem_c from z_elem_b
  602. z_elem_c(ji,jj,ibotm1) = 0.
  603. !
  604. ! Set psi vertical flux at the bottom:
  605. zdep(ji,jj) = rn_bfrz0 + 0.5_wp*fse3t(ji,jj,ibotm1)
  606. zflxb = rsbc_psi2 * ( avm(ji,jj,ibot) + avm(ji,jj,ibotm1) ) &
  607. & * (0.5_wp*(en(ji,jj,ibot)+en(ji,jj,ibotm1)))**rmm * zdep(ji,jj)**(rnn-1._wp)
  608. psi(ji,jj,ibotm1) = psi(ji,jj,ibotm1) + zflxb / fse3w(ji,jj,ibotm1)
  609. END DO
  610. END DO
  611. !
  612. END SELECT
  613. ! Matrix inversion
  614. ! ----------------
  615. !
  616. DO jk = 2, jpkm1 ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1
  617. DO jj = 2, jpjm1
  618. DO ji = fs_2, fs_jpim1 ! vector opt.
  619. z_elem_b(ji,jj,jk) = z_elem_b(ji,jj,jk) - z_elem_a(ji,jj,jk) * z_elem_c(ji,jj,jk-1) / z_elem_b(ji,jj,jk-1)
  620. END DO
  621. END DO
  622. END DO
  623. DO jk = 2, jpk ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1
  624. DO jj = 2, jpjm1
  625. DO ji = fs_2, fs_jpim1 ! vector opt.
  626. z_elem_a(ji,jj,jk) = psi(ji,jj,jk) - z_elem_a(ji,jj,jk) / z_elem_b(ji,jj,jk-1) * z_elem_a(ji,jj,jk-1)
  627. END DO
  628. END DO
  629. END DO
  630. DO jk = jpk-1, 2, -1 ! Third recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk
  631. DO jj = 2, jpjm1
  632. DO ji = fs_2, fs_jpim1 ! vector opt.
  633. psi(ji,jj,jk) = ( z_elem_a(ji,jj,jk) - z_elem_c(ji,jj,jk) * psi(ji,jj,jk+1) ) / z_elem_b(ji,jj,jk)
  634. END DO
  635. END DO
  636. END DO
  637. ! Set dissipation
  638. !----------------
  639. SELECT CASE ( nn_clos )
  640. !
  641. CASE( 0 ) ! k-kl (Mellor-Yamada)
  642. DO jk = 1, jpkm1
  643. DO jj = 2, jpjm1
  644. DO ji = fs_2, fs_jpim1 ! vector opt.
  645. eps(ji,jj,jk) = rc03 * en(ji,jj,jk) * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / MAX( psi(ji,jj,jk), rn_epsmin)
  646. END DO
  647. END DO
  648. END DO
  649. !
  650. CASE( 1 ) ! k-eps
  651. DO jk = 1, jpkm1
  652. DO jj = 2, jpjm1
  653. DO ji = fs_2, fs_jpim1 ! vector opt.
  654. eps(ji,jj,jk) = psi(ji,jj,jk)
  655. END DO
  656. END DO
  657. END DO
  658. !
  659. CASE( 2 ) ! k-w
  660. DO jk = 1, jpkm1
  661. DO jj = 2, jpjm1
  662. DO ji = fs_2, fs_jpim1 ! vector opt.
  663. eps(ji,jj,jk) = rc04 * en(ji,jj,jk) * psi(ji,jj,jk)
  664. END DO
  665. END DO
  666. END DO
  667. !
  668. CASE( 3 ) ! generic
  669. zcoef = rc0**( 3._wp + rpp/rnn )
  670. zex1 = ( 1.5_wp + rmm/rnn )
  671. zex2 = -1._wp / rnn
  672. DO jk = 1, jpkm1
  673. DO jj = 2, jpjm1
  674. DO ji = fs_2, fs_jpim1 ! vector opt.
  675. eps(ji,jj,jk) = zcoef * en(ji,jj,jk)**zex1 * psi(ji,jj,jk)**zex2
  676. END DO
  677. END DO
  678. END DO
  679. !
  680. END SELECT
  681. ! Limit dissipation rate under stable stratification
  682. ! --------------------------------------------------
  683. DO jk = 1, jpkm1 ! Note that this set boundary conditions on mxln at the same time
  684. DO jj = 2, jpjm1
  685. DO ji = fs_2, fs_jpim1 ! vector opt.
  686. ! limitation
  687. eps(ji,jj,jk) = MAX( eps(ji,jj,jk), rn_epsmin )
  688. mxln(ji,jj,jk) = rc03 * en(ji,jj,jk) * SQRT( en(ji,jj,jk) ) / eps(ji,jj,jk)
  689. ! Galperin criterium (NOTE : Not required if the proper value of C3 in stable cases is calculated)
  690. zrn2 = MAX( rn2(ji,jj,jk), rsmall )
  691. IF (ln_length_lim) mxln(ji,jj,jk) = MIN( rn_clim_galp * SQRT( 2._wp * en(ji,jj,jk) / zrn2 ), mxln(ji,jj,jk) )
  692. END DO
  693. END DO
  694. END DO
  695. !
  696. ! Stability function and vertical viscosity and diffusivity
  697. ! ---------------------------------------------------------
  698. !
  699. SELECT CASE ( nn_stab_func )
  700. !
  701. CASE ( 0 , 1 ) ! Galperin or Kantha-Clayson stability functions
  702. DO jk = 2, jpkm1
  703. DO jj = 2, jpjm1
  704. DO ji = fs_2, fs_jpim1 ! vector opt.
  705. ! zcof = l²/q²
  706. zcof = mxlb(ji,jj,jk) * mxlb(ji,jj,jk) / ( 2._wp*eb(ji,jj,jk) )
  707. ! Gh = -N²l²/q²
  708. gh = - rn2(ji,jj,jk) * zcof
  709. gh = MIN( gh, rgh0 )
  710. gh = MAX( gh, rghmin )
  711. ! Stability functions from Kantha and Clayson (if C2=C3=0 => Galperin)
  712. sh = ra2*( 1._wp-6._wp*ra1/rb1 ) / ( 1.-3.*ra2*gh*(6.*ra1+rb2*( 1._wp-rc3 ) ) )
  713. sm = ( rb1**(-1._wp/3._wp) + ( 18._wp*ra1*ra1 + 9._wp*ra1*ra2*(1._wp-rc2) )*sh*gh ) / (1._wp-9._wp*ra1*ra2*gh)
  714. !
  715. ! Store stability function in avmu and avmv
  716. avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk)
  717. avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk)
  718. END DO
  719. END DO
  720. END DO
  721. !
  722. CASE ( 2, 3 ) ! Canuto stability functions
  723. DO jk = 2, jpkm1
  724. DO jj = 2, jpjm1
  725. DO ji = fs_2, fs_jpim1 ! vector opt.
  726. ! zcof = l²/q²
  727. zcof = mxlb(ji,jj,jk)*mxlb(ji,jj,jk) / ( 2._wp * eb(ji,jj,jk) )
  728. ! Gh = -N²l²/q²
  729. gh = - rn2(ji,jj,jk) * zcof
  730. gh = MIN( gh, rgh0 )
  731. gh = MAX( gh, rghmin )
  732. gh = gh * rf6
  733. ! Gm = M²l²/q² Shear number
  734. shr = shear(ji,jj,jk) / MAX( avm(ji,jj,jk), rsmall )
  735. gm = MAX( shr * zcof , 1.e-10 )
  736. gm = gm * rf6
  737. gm = MIN ( (rd0 - rd1*gh + rd3*gh*gh) / (rd2-rd4*gh) , gm )
  738. ! Stability functions from Canuto
  739. rcff = rd0 - rd1*gh +rd2*gm + rd3*gh*gh - rd4*gh*gm + rd5*gm*gm
  740. sm = (rs0 - rs1*gh + rs2*gm) / rcff
  741. sh = (rs4 - rs5*gh + rs6*gm) / rcff
  742. !
  743. ! Store stability function in avmu and avmv
  744. avmu(ji,jj,jk) = rc_diff * sh * tmask(ji,jj,jk)
  745. avmv(ji,jj,jk) = rc_diff * sm * tmask(ji,jj,jk)
  746. END DO
  747. END DO
  748. END DO
  749. !
  750. END SELECT
  751. ! Boundary conditions on stability functions for momentum (Neumann):
  752. ! Lines below are useless if GOTM style Dirichlet conditions are used
  753. avmv(:,:,1) = avmv(:,:,2)
  754. DO jj = 2, jpjm1
  755. DO ji = fs_2, fs_jpim1 ! vector opt.
  756. avmv(ji,jj,mbkt(ji,jj)+1) = avmv(ji,jj,mbkt(ji,jj))
  757. END DO
  758. END DO
  759. ! Compute diffusivities/viscosities
  760. ! The computation below could be restrained to jk=2 to jpkm1 if GOTM style Dirichlet conditions are used
  761. DO jk = 1, jpk
  762. DO jj = 2, jpjm1
  763. DO ji = fs_2, fs_jpim1 ! vector opt.
  764. zsqen = SQRT( 2._wp * en(ji,jj,jk) ) * mxln(ji,jj,jk)
  765. zav = zsqen * avmu(ji,jj,jk)
  766. avt(ji,jj,jk) = MAX( zav, avtb(jk) )*tmask(ji,jj,jk) ! apply mask for zdfmxl routine
  767. zav = zsqen * avmv(ji,jj,jk)
  768. avm(ji,jj,jk) = MAX( zav, avmb(jk) ) ! Note that avm is not masked at the surface and the bottom
  769. END DO
  770. END DO
  771. END DO
  772. !
  773. ! Lateral boundary conditions (sign unchanged)
  774. avt(:,:,1) = 0._wp
  775. CALL lbc_lnk( avm, 'W', 1. ) ; CALL lbc_lnk( avt, 'W', 1. )
  776. DO jk = 2, jpkm1 !* vertical eddy viscosity at u- and v-points
  777. DO jj = 2, jpjm1
  778. DO ji = fs_2, fs_jpim1 ! vector opt.
  779. avmu(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji+1,jj ,jk) ) * umask(ji,jj,jk)
  780. avmv(ji,jj,jk) = 0.5 * ( avm(ji,jj,jk) + avm(ji ,jj+1,jk) ) * vmask(ji,jj,jk)
  781. END DO
  782. END DO
  783. END DO
  784. avmu(:,:,1) = 0._wp ; avmv(:,:,1) = 0._wp ! set surface to zero
  785. CALL lbc_lnk( avmu, 'U', 1. ) ; CALL lbc_lnk( avmv, 'V', 1. ) ! Lateral boundary conditions
  786. IF(ln_ctl) THEN
  787. CALL prt_ctl( tab3d_1=en , clinfo1=' gls - e: ', tab3d_2=avt, clinfo2=' t: ', ovlap=1, kdim=jpk)
  788. CALL prt_ctl( tab3d_1=avmu, clinfo1=' gls - u: ', mask1=umask, &
  789. & tab3d_2=avmv, clinfo2= ' v: ', mask2=vmask, ovlap=1, kdim=jpk )
  790. ENDIF
  791. !
  792. avt_k (:,:,:) = avt (:,:,:)
  793. avm_k (:,:,:) = avm (:,:,:)
  794. avmu_k(:,:,:) = avmu(:,:,:)
  795. avmv_k(:,:,:) = avmv(:,:,:)
  796. !
  797. CALL wrk_dealloc( jpi,jpj, zdep, zkar, zflxs, zhsro )
  798. CALL wrk_dealloc( jpi,jpj,jpk, eb, mxlb, shear, eps, zwall_psi, z_elem_a, z_elem_b, z_elem_c, psi )
  799. !
  800. IF( nn_timing == 1 ) CALL timing_stop('zdf_gls')
  801. !
  802. !
  803. END SUBROUTINE zdf_gls
  804. SUBROUTINE zdf_gls_init
  805. !!----------------------------------------------------------------------
  806. !! *** ROUTINE zdf_gls_init ***
  807. !!
  808. !! ** Purpose : Initialization of the vertical eddy diffivity and
  809. !! viscosity when using a gls turbulent closure scheme
  810. !!
  811. !! ** Method : Read the namzdf_gls namelist and check the parameters
  812. !! called at the first timestep (nit000)
  813. !!
  814. !! ** input : Namlist namzdf_gls
  815. !!
  816. !! ** Action : Increase by 1 the nstop flag is setting problem encounter
  817. !!
  818. !!----------------------------------------------------------------------
  819. USE dynzdf_exp
  820. USE trazdf_exp
  821. !
  822. INTEGER :: jk ! dummy loop indices
  823. INTEGER :: ios ! Local integer output status for namelist read
  824. REAL(wp):: zcr ! local scalar
  825. !!
  826. NAMELIST/namzdf_gls/rn_emin, rn_epsmin, ln_length_lim, &
  827. & rn_clim_galp, ln_sigpsi, rn_hsro, &
  828. & rn_crban, rn_charn, rn_frac_hs, &
  829. & nn_bc_surf, nn_bc_bot, nn_z0_met, &
  830. & nn_stab_func, nn_clos
  831. !!----------------------------------------------------------
  832. !
  833. IF( nn_timing == 1 ) CALL timing_start('zdf_gls_init')
  834. !
  835. REWIND( numnam_ref ) ! Namelist namzdf_gls in reference namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme
  836. READ ( numnam_ref, namzdf_gls, IOSTAT = ios, ERR = 901)
  837. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in reference namelist', lwp )
  838. REWIND( numnam_cfg ) ! Namelist namzdf_gls in configuration namelist : Vertical eddy diffivity and viscosity using gls turbulent closure scheme
  839. READ ( numnam_cfg, namzdf_gls, IOSTAT = ios, ERR = 902 )
  840. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_gls in configuration namelist', lwp )
  841. IF(lwm) WRITE ( numond, namzdf_gls )
  842. IF(lwp) THEN !* Control print
  843. WRITE(numout,*)
  844. WRITE(numout,*) 'zdf_gls_init : gls turbulent closure scheme'
  845. WRITE(numout,*) '~~~~~~~~~~~~'
  846. WRITE(numout,*) ' Namelist namzdf_gls : set gls mixing parameters'
  847. WRITE(numout,*) ' minimum value of en rn_emin = ', rn_emin
  848. WRITE(numout,*) ' minimum value of eps rn_epsmin = ', rn_epsmin
  849. WRITE(numout,*) ' Limit dissipation rate under stable stratif. ln_length_lim = ', ln_length_lim
  850. WRITE(numout,*) ' Galperin limit (Standard: 0.53, Holt: 0.26) rn_clim_galp = ', rn_clim_galp
  851. WRITE(numout,*) ' TKE Surface boundary condition nn_bc_surf = ', nn_bc_surf
  852. WRITE(numout,*) ' TKE Bottom boundary condition nn_bc_bot = ', nn_bc_bot
  853. WRITE(numout,*) ' Modify psi Schmidt number (wb case) ln_sigpsi = ', ln_sigpsi
  854. WRITE(numout,*) ' Craig and Banner coefficient rn_crban = ', rn_crban
  855. WRITE(numout,*) ' Charnock coefficient rn_charn = ', rn_charn
  856. WRITE(numout,*) ' Surface roughness formula nn_z0_met = ', nn_z0_met
  857. WRITE(numout,*) ' Wave height frac. (used if nn_z0_met=2) rn_frac_hs = ', rn_frac_hs
  858. WRITE(numout,*) ' Stability functions nn_stab_func = ', nn_stab_func
  859. WRITE(numout,*) ' Type of closure nn_clos = ', nn_clos
  860. WRITE(numout,*) ' Surface roughness (m) rn_hsro = ', rn_hsro
  861. WRITE(numout,*) ' Bottom roughness (m) (nambfr namelist) rn_bfrz0 = ', rn_bfrz0
  862. ENDIF
  863. ! !* allocate gls arrays
  864. IF( zdf_gls_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_gls_init : unable to allocate arrays' )
  865. ! !* Check of some namelist values
  866. IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'bad flag: nn_bc_surf is 0 or 1' )
  867. IF( nn_bc_surf < 0 .OR. nn_bc_surf > 1 ) CALL ctl_stop( 'bad flag: nn_bc_surf is 0 or 1' )
  868. IF( nn_z0_met < 0 .OR. nn_z0_met > 2 ) CALL ctl_stop( 'bad flag: nn_z0_met is 0, 1 or 2' )
  869. IF( nn_stab_func < 0 .OR. nn_stab_func > 3 ) CALL ctl_stop( 'bad flag: nn_stab_func is 0, 1, 2 and 3' )
  870. IF( nn_clos < 0 .OR. nn_clos > 3 ) CALL ctl_stop( 'bad flag: nn_clos is 0, 1, 2 or 3' )
  871. SELECT CASE ( nn_clos ) !* set the parameters for the chosen closure
  872. !
  873. CASE( 0 ) ! k-kl (Mellor-Yamada)
  874. !
  875. IF(lwp) WRITE(numout,*) 'The choosen closure is k-kl closed to the classical Mellor-Yamada'
  876. rpp = 0._wp
  877. rmm = 1._wp
  878. rnn = 1._wp
  879. rsc_tke = 1.96_wp
  880. rsc_psi = 1.96_wp
  881. rpsi1 = 0.9_wp
  882. rpsi3p = 1._wp
  883. rpsi2 = 0.5_wp
  884. !
  885. SELECT CASE ( nn_stab_func )
  886. CASE( 0, 1 ) ; rpsi3m = 2.53_wp ! G88 or KC stability functions
  887. CASE( 2 ) ; rpsi3m = 2.62_wp ! Canuto A stability functions
  888. CASE( 3 ) ; rpsi3m = 2.38 ! Canuto B stability functions (caution : constant not identified)
  889. END SELECT
  890. !
  891. CASE( 1 ) ! k-eps
  892. !
  893. IF(lwp) WRITE(numout,*) 'The choosen closure is k-eps'
  894. rpp = 3._wp
  895. rmm = 1.5_wp
  896. rnn = -1._wp
  897. rsc_tke = 1._wp
  898. rsc_psi = 1.2_wp ! Schmidt number for psi
  899. rpsi1 = 1.44_wp
  900. rpsi3p = 1._wp
  901. rpsi2 = 1.92_wp
  902. !
  903. SELECT CASE ( nn_stab_func )
  904. CASE( 0, 1 ) ; rpsi3m = -0.52_wp ! G88 or KC stability functions
  905. CASE( 2 ) ; rpsi3m = -0.629_wp ! Canuto A stability functions
  906. CASE( 3 ) ; rpsi3m = -0.566 ! Canuto B stability functions
  907. END SELECT
  908. !
  909. CASE( 2 ) ! k-omega
  910. !
  911. IF(lwp) WRITE(numout,*) 'The choosen closure is k-omega'
  912. rpp = -1._wp
  913. rmm = 0.5_wp
  914. rnn = -1._wp
  915. rsc_tke = 2._wp
  916. rsc_psi = 2._wp
  917. rpsi1 = 0.555_wp
  918. rpsi3p = 1._wp
  919. rpsi2 = 0.833_wp
  920. !
  921. SELECT CASE ( nn_stab_func )
  922. CASE( 0, 1 ) ; rpsi3m = -0.58_wp ! G88 or KC stability functions
  923. CASE( 2 ) ; rpsi3m = -0.64_wp ! Canuto A stability functions
  924. CASE( 3 ) ; rpsi3m = -0.64_wp ! Canuto B stability functions caution : constant not identified)
  925. END SELECT
  926. !
  927. CASE( 3 ) ! generic
  928. !
  929. IF(lwp) WRITE(numout,*) 'The choosen closure is generic'
  930. rpp = 2._wp
  931. rmm = 1._wp
  932. rnn = -0.67_wp
  933. rsc_tke = 0.8_wp
  934. rsc_psi = 1.07_wp
  935. rpsi1 = 1._wp
  936. rpsi3p = 1._wp
  937. rpsi2 = 1.22_wp
  938. !
  939. SELECT CASE ( nn_stab_func )
  940. CASE( 0, 1 ) ; rpsi3m = 0.1_wp ! G88 or KC stability functions
  941. CASE( 2 ) ; rpsi3m = 0.05_wp ! Canuto A stability functions
  942. CASE( 3 ) ; rpsi3m = 0.05_wp ! Canuto B stability functions caution : constant not identified)
  943. END SELECT
  944. !
  945. END SELECT
  946. !
  947. SELECT CASE ( nn_stab_func ) !* set the parameters of the stability functions
  948. !
  949. CASE ( 0 ) ! Galperin stability functions
  950. !
  951. IF(lwp) WRITE(numout,*) 'Stability functions from Galperin'
  952. rc2 = 0._wp
  953. rc3 = 0._wp
  954. rc_diff = 1._wp
  955. rc0 = 0.5544_wp
  956. rcm_sf = 0.9884_wp
  957. rghmin = -0.28_wp
  958. rgh0 = 0.0233_wp
  959. rghcri = 0.02_wp
  960. !
  961. CASE ( 1 ) ! Kantha-Clayson stability functions
  962. !
  963. IF(lwp) WRITE(numout,*) 'Stability functions from Kantha-Clayson'
  964. rc2 = 0.7_wp
  965. rc3 = 0.2_wp
  966. rc_diff = 1._wp
  967. rc0 = 0.5544_wp
  968. rcm_sf = 0.9884_wp
  969. rghmin = -0.28_wp
  970. rgh0 = 0.0233_wp
  971. rghcri = 0.02_wp
  972. !
  973. CASE ( 2 ) ! Canuto A stability functions
  974. !
  975. IF(lwp) WRITE(numout,*) 'Stability functions from Canuto A'
  976. rs0 = 1.5_wp * rl1 * rl5*rl5
  977. rs1 = -rl4*(rl6+rl7) + 2._wp*rl4*rl5*(rl1-(1._wp/3._wp)*rl2-rl3) + 1.5_wp*rl1*rl5*rl8
  978. rs2 = -(3._wp/8._wp) * rl1*(rl6*rl6-rl7*rl7)
  979. rs4 = 2._wp * rl5
  980. rs5 = 2._wp * rl4
  981. rs6 = (2._wp/3._wp) * rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.5_wp * rl5*rl1 * (3._wp*rl3-rl2) &
  982. & + 0.75_wp * rl1 * ( rl6 - rl7 )
  983. rd0 = 3._wp * rl5*rl5
  984. rd1 = rl5 * ( 7._wp*rl4 + 3._wp*rl8 )
  985. rd2 = rl5*rl5 * ( 3._wp*rl3*rl3 - rl2*rl2 ) - 0.75_wp*(rl6*rl6 - rl7*rl7 )
  986. rd3 = rl4 * ( 4._wp*rl4 + 3._wp*rl8)
  987. rd4 = rl4 * ( rl2 * rl6 - 3._wp*rl3*rl7 - rl5*(rl2*rl2 - rl3*rl3 ) ) + rl5*rl8 * ( 3._wp*rl3*rl3 - rl2*rl2 )
  988. rd5 = 0.25_wp * ( rl2*rl2 - 3._wp *rl3*rl3 ) * ( rl6*rl6 - rl7*rl7 )
  989. rc0 = 0.5268_wp
  990. rf6 = 8._wp / (rc0**6._wp)
  991. rc_diff = SQRT(2._wp) / (rc0**3._wp)
  992. rcm_sf = 0.7310_wp
  993. rghmin = -0.28_wp
  994. rgh0 = 0.0329_wp
  995. rghcri = 0.03_wp
  996. !
  997. CASE ( 3 ) ! Canuto B stability functions
  998. !
  999. IF(lwp) WRITE(numout,*) 'Stability functions from Canuto B'
  1000. rs0 = 1.5_wp * rm1 * rm5*rm5
  1001. rs1 = -rm4 * (rm6+rm7) + 2._wp * rm4*rm5*(rm1-(1._wp/3._wp)*rm2-rm3) + 1.5_wp * rm1*rm5*rm8
  1002. rs2 = -(3._wp/8._wp) * rm1 * (rm6*rm6-rm7*rm7 )
  1003. rs4 = 2._wp * rm5
  1004. rs5 = 2._wp * rm4
  1005. rs6 = (2._wp/3._wp) * rm5 * (3._wp*rm3*rm3-rm2*rm2) - 0.5_wp * rm5*rm1*(3._wp*rm3-rm2) + 0.75_wp * rm1*(rm6-rm7)
  1006. rd0 = 3._wp * rm5*rm5
  1007. rd1 = rm5 * (7._wp*rm4 + 3._wp*rm8)
  1008. rd2 = rm5*rm5 * (3._wp*rm3*rm3 - rm2*rm2) - 0.75_wp * (rm6*rm6 - rm7*rm7)
  1009. rd3 = rm4 * ( 4._wp*rm4 + 3._wp*rm8 )
  1010. rd4 = rm4 * ( rm2*rm6 -3._wp*rm3*rm7 - rm5*(rm2*rm2 - rm3*rm3) ) + rm5 * rm8 * ( 3._wp*rm3*rm3 - rm2*rm2 )
  1011. rd5 = 0.25_wp * ( rm2*rm2 - 3._wp*rm3*rm3 ) * ( rm6*rm6 - rm7*rm7 )
  1012. rc0 = 0.5268_wp !! rc0 = 0.5540_wp (Warner ...) to verify !
  1013. rf6 = 8._wp / ( rc0**6._wp )
  1014. rc_diff = SQRT(2._wp)/(rc0**3.)
  1015. rcm_sf = 0.7470_wp
  1016. rghmin = -0.28_wp
  1017. rgh0 = 0.0444_wp
  1018. rghcri = 0.0414_wp
  1019. !
  1020. END SELECT
  1021. ! !* Set Schmidt number for psi diffusion in the wave breaking case
  1022. ! ! See Eq. (13) of Carniel et al, OM, 30, 225-239, 2009
  1023. ! ! or Eq. (17) of Burchard, JPO, 31, 3133-3145, 2001
  1024. IF( ln_sigpsi ) THEN
  1025. ra_sf = -1.5 ! Set kinetic energy slope, then deduce rsc_psi and rl_sf
  1026. ! Verification: retrieve Burchard (2001) results by uncomenting the line below:
  1027. ! Note that the results depend on the value of rn_cm_sf which is constant (=rc0) in his work
  1028. ! ra_sf = -SQRT(2./3.*rc0**3./rn_cm_sf*rn_sc_tke)/vkarmn
  1029. rsc_psi0 = rsc_tke/(24.*rpsi2)*(-1.+(4.*rnn + ra_sf*(1.+4.*rmm))**2./(ra_sf**2.))
  1030. ELSE
  1031. rsc_psi0 = rsc_psi
  1032. ENDIF
  1033. ! !* Shear free turbulence parameters
  1034. !
  1035. ra_sf = -4._wp*rnn*SQRT(rsc_tke) / ( (1._wp+4._wp*rmm)*SQRT(rsc_tke) &
  1036. & - SQRT(rsc_tke + 24._wp*rsc_psi0*rpsi2 ) )
  1037. IF ( rn_crban==0._wp ) THEN
  1038. rl_sf = vkarmn
  1039. ELSE
  1040. rl_sf = rc0 * SQRT(rc0/rcm_sf) * SQRT( ( (1._wp + 4._wp*rmm + 8._wp*rmm**2_wp)*rsc_tke &
  1041. & + 12._wp * rsc_psi0*rpsi2 - (1._wp + 4._wp*rmm) &
  1042. & *SQRT(rsc_tke*(rsc_tke &
  1043. & + 24._wp*rsc_psi0*rpsi2)) ) &
  1044. & /(12._wp*rnn**2.) &
  1045. & )
  1046. ENDIF
  1047. !
  1048. IF(lwp) THEN !* Control print
  1049. WRITE(numout,*)
  1050. WRITE(numout,*) 'Limit values'
  1051. WRITE(numout,*) '~~~~~~~~~~~~'
  1052. WRITE(numout,*) 'Parameter m = ',rmm
  1053. WRITE(numout,*) 'Parameter n = ',rnn
  1054. WRITE(numout,*) 'Parameter p = ',rpp
  1055. WRITE(numout,*) 'rpsi1 = ',rpsi1
  1056. WRITE(numout,*) 'rpsi2 = ',rpsi2
  1057. WRITE(numout,*) 'rpsi3m = ',rpsi3m
  1058. WRITE(numout,*) 'rpsi3p = ',rpsi3p
  1059. WRITE(numout,*) 'rsc_tke = ',rsc_tke
  1060. WRITE(numout,*) 'rsc_psi = ',rsc_psi
  1061. WRITE(numout,*) 'rsc_psi0 = ',rsc_psi0
  1062. WRITE(numout,*) 'rc0 = ',rc0
  1063. WRITE(numout,*)
  1064. WRITE(numout,*) 'Shear free turbulence parameters:'
  1065. WRITE(numout,*) 'rcm_sf = ',rcm_sf
  1066. WRITE(numout,*) 'ra_sf = ',ra_sf
  1067. WRITE(numout,*) 'rl_sf = ',rl_sf
  1068. WRITE(numout,*)
  1069. ENDIF
  1070. ! !* Constants initialization
  1071. rc02 = rc0 * rc0 ; rc02r = 1. / rc02
  1072. rc03 = rc02 * rc0
  1073. rc04 = rc03 * rc0
  1074. rsbc_tke1 = -3._wp/2._wp*rn_crban*ra_sf*rl_sf ! Dirichlet + Wave breaking
  1075. rsbc_tke2 = rdt * rn_crban / rl_sf ! Neumann + Wave breaking
  1076. zcr = MAX(rsmall, rsbc_tke1**(1./(-ra_sf*3._wp/2._wp))-1._wp )
  1077. rtrans = 0.2_wp / zcr ! Ad. inverse transition length between log and wave layer
  1078. rsbc_zs1 = rn_charn/grav ! Charnock formula for surface roughness
  1079. rsbc_zs2 = rn_frac_hs / 0.85_wp / grav * 665._wp ! Rascle formula for surface roughness
  1080. rsbc_psi1 = -0.5_wp * rdt * rc0**(rpp-2._wp*rmm) / rsc_psi
  1081. rsbc_psi2 = -0.5_wp * rdt * rc0**rpp * rnn * vkarmn**rnn / rsc_psi ! Neumann + NO Wave breaking
  1082. rfact_tke = -0.5_wp / rsc_tke * rdt ! Cst used for the Diffusion term of tke
  1083. rfact_psi = -0.5_wp / rsc_psi * rdt ! Cst used for the Diffusion term of tke
  1084. ! !* Wall proximity function
  1085. zwall (:,:,:) = 1._wp * tmask(:,:,:)
  1086. ! !* set vertical eddy coef. to the background value
  1087. DO jk = 1, jpk
  1088. avt (:,:,jk) = avtb(jk) * tmask(:,:,jk)
  1089. avm (:,:,jk) = avmb(jk) * tmask(:,:,jk)
  1090. avmu(:,:,jk) = avmb(jk) * umask(:,:,jk)
  1091. avmv(:,:,jk) = avmb(jk) * vmask(:,:,jk)
  1092. END DO
  1093. !
  1094. CALL gls_rst( nit000, 'READ' ) !* read or initialize all required files
  1095. !
  1096. IF( nn_timing == 1 ) CALL timing_stop('zdf_gls_init')
  1097. !
  1098. END SUBROUTINE zdf_gls_init
  1099. SUBROUTINE gls_rst( kt, cdrw )
  1100. !!---------------------------------------------------------------------
  1101. !! *** ROUTINE ts_rst ***
  1102. !!
  1103. !! ** Purpose : Read or write TKE file (en) in restart file
  1104. !!
  1105. !! ** Method : use of IOM library
  1106. !! if the restart does not contain TKE, en is either
  1107. !! set to rn_emin or recomputed (nn_igls/=0)
  1108. !!----------------------------------------------------------------------
  1109. INTEGER , INTENT(in) :: kt ! ocean time-step
  1110. CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag
  1111. !
  1112. INTEGER :: jit, jk ! dummy loop indices
  1113. INTEGER :: id1, id2, id3, id4, id5, id6
  1114. INTEGER :: ji, jj, ikbu, ikbv
  1115. REAL(wp):: cbx, cby
  1116. !!----------------------------------------------------------------------
  1117. !
  1118. IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise
  1119. ! ! ---------------
  1120. IF( ln_rstart ) THEN !* Read the restart file
  1121. id1 = iom_varid( numror, 'en' , ldstop = .FALSE. )
  1122. id2 = iom_varid( numror, 'avt' , ldstop = .FALSE. )
  1123. id3 = iom_varid( numror, 'avm' , ldstop = .FALSE. )
  1124. id4 = iom_varid( numror, 'avmu' , ldstop = .FALSE. )
  1125. id5 = iom_varid( numror, 'avmv' , ldstop = .FALSE. )
  1126. id6 = iom_varid( numror, 'mxln' , ldstop = .FALSE. )
  1127. !
  1128. IF( MIN( id1, id2, id3, id4, id5, id6 ) > 0 ) THEN ! all required arrays exist
  1129. CALL iom_get( numror, jpdom_autoglo, 'en' , en )
  1130. CALL iom_get( numror, jpdom_autoglo, 'avt' , avt )
  1131. CALL iom_get( numror, jpdom_autoglo, 'avm' , avm )
  1132. CALL iom_get( numror, jpdom_autoglo, 'avmu' , avmu )
  1133. CALL iom_get( numror, jpdom_autoglo, 'avmv' , avmv )
  1134. CALL iom_get( numror, jpdom_autoglo, 'mxln' , mxln )
  1135. ELSE
  1136. IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without gls scheme, en and mxln computed by iterative loop'
  1137. en (:,:,:) = rn_emin
  1138. mxln(:,:,:) = 0.05
  1139. avt_k (:,:,:) = avt (:,:,:)
  1140. avm_k (:,:,:) = avm (:,:,:)
  1141. avmu_k(:,:,:) = avmu(:,:,:)
  1142. avmv_k(:,:,:) = avmv(:,:,:)
  1143. DO jit = nit000 + 1, nit000 + 10 ; CALL zdf_gls( jit ) ; END DO
  1144. ENDIF
  1145. ELSE !* Start from rest
  1146. IF(lwp) WRITE(numout,*) ' ===>>>> : Initialisation of en and mxln by background values'
  1147. en (:,:,:) = rn_emin
  1148. mxln(:,:,:) = 0.05
  1149. ENDIF
  1150. !
  1151. ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file
  1152. ! ! -------------------
  1153. IF(lwp) WRITE(numout,*) '---- gls-rst ----'
  1154. CALL iom_rstput( kt, nitrst, numrow, 'en' , en )
  1155. CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt_k )
  1156. CALL iom_rstput( kt, nitrst, numrow, 'avm' , avm_k )
  1157. CALL iom_rstput( kt, nitrst, numrow, 'avmu' , avmu_k )
  1158. CALL iom_rstput( kt, nitrst, numrow, 'avmv' , avmv_k )
  1159. CALL iom_rstput( kt, nitrst, numrow, 'mxln' , mxln )
  1160. !
  1161. ENDIF
  1162. !
  1163. END SUBROUTINE gls_rst
  1164. #else
  1165. !!----------------------------------------------------------------------
  1166. !! Dummy module : NO TKE scheme
  1167. !!----------------------------------------------------------------------
  1168. LOGICAL, PUBLIC, PARAMETER :: lk_zdfgls = .FALSE. !: TKE flag
  1169. CONTAINS
  1170. SUBROUTINE zdf_gls_init ! Empty routine
  1171. WRITE(*,*) 'zdf_gls_init: You should not have seen this print! error?'
  1172. END SUBROUTINE zdf_gls_init
  1173. SUBROUTINE zdf_gls( kt ) ! Empty routine
  1174. WRITE(*,*) 'zdf_gls: You should not have seen this print! error?', kt
  1175. END SUBROUTINE zdf_gls
  1176. SUBROUTINE gls_rst( kt, cdrw ) ! Empty routine
  1177. INTEGER , INTENT(in) :: kt ! ocean time-step
  1178. CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag
  1179. WRITE(*,*) 'gls_rst: You should not have seen this print! error?', kt, cdrw
  1180. END SUBROUTINE gls_rst
  1181. #endif
  1182. !!======================================================================
  1183. END MODULE zdfgls