diffusion.F90 56 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584
  1. !#################################################################
  2. !BOP
  3. !
  4. ! !MODULE: DIFFUSION
  5. !
  6. ! !DESCRIPTION:
  7. !
  8. ! Diffusion routines.
  9. !
  10. !
  11. ! !REVISION HISTORY: -----------------------------
  12. !
  13. ! ????-?? Bram Bregman.
  14. ! Adjusted for TM5 by
  15. !
  16. ! 2005-04 ??
  17. ! Adjoint version
  18. !
  19. ! 2010-08 Arjo Segers
  20. ! Debugged computation of u* over sea:
  21. ! abs. surface stress should be devided by air density.
  22. ! Set minimum value for u* to trap division by zero (MK).
  23. !
  24. ! 2012-01 Philippe Le Sager
  25. ! Adapted for Lat-Lon MPI domain decomposition
  26. !
  27. ! 2019-06 Andy Jacobson
  28. ! Changed DKG file format to CF-compliant netCDF.
  29. ! Also added diffusion_filename function (Sourish Basu)
  30. ! to place DKG files in a meteo file hierarchy.
  31. !
  32. !EOP
  33. !#################################################################
  34. !
  35. ! global input fields:
  36. ! phlb(:,:,lm+1) ! half level pressure [Pa]
  37. ! gph(:,:,lm+1) ! height of half-levels [m] !CMK note gph starts NOW at 1!
  38. ! t ! temperature [K]
  39. ! q ! specific humidity [kg/kg]
  40. ! uwind ! mean wind W-E [m/s]
  41. ! vwind ! mean wind S-N [m/s]
  42. ! sshf ! surface sensible heat flux [W / m^2]
  43. ! slhf ! surface latent heat flux [W / m^2]
  44. ! ustar (ustar_loc) ! velocity scale
  45. !
  46. ! global output fields:
  47. ! kvh ! eddy diffusivity for heat at top [m^2/s]
  48. ! pblh ! boundary layer height[m]
  49. ! dkg ! vertical diffusion coefficient [ks s-1]
  50. !
  51. !### macro's #####################################################
  52. !
  53. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  54. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  55. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  56. !
  57. #include "tm5.inc"
  58. !
  59. !#################################################################
  60. MODULE DIFFUSION
  61. use GO, only : gol, goPr, goErr
  62. use dims, only : nregions
  63. use global_types, only : emis_data
  64. implicit none
  65. private
  66. ! methods
  67. public :: Diffusion_Init, Diffusion_Done
  68. public :: Calc_Kzz, Read_Diffusion, Write_Diffusion
  69. ! flag
  70. public :: diffusion_write
  71. ! --- const --------------------------------------
  72. character(len=*), parameter :: mname = 'diffusion'
  73. logical :: diffusion_write = .false.
  74. character(len=2000) :: diffusion_dir = ''
  75. character(len=20) :: mlevs
  76. ! --- local ------------------------------------
  77. type(emis_data), dimension(nregions), target :: ustar_loc, sr_mix
  78. CONTAINS
  79. ! ================================================================
  80. subroutine Diffusion_Init( status )
  81. use global_data, only : rcfile
  82. use MeteoData, only : Set
  83. use MeteoData, only : pu_dat, pv_dat
  84. use MeteoData, only : temper_dat, humid_dat, gph_dat
  85. use MeteoData, only : oro_dat, lsmask_dat
  86. use MeteoData, only : sr_ecm_dat, sr_ols_dat
  87. use MeteoData, only : wspd_dat
  88. use MeteoData, only : slhf_dat, sshf_dat
  89. use MeteoData, only : ewss_dat, nsss_dat
  90. use GO, only : TrcFile, Init, Done, ReadRc
  91. ! --- in/out --------------------------------
  92. integer, intent(out) :: status
  93. ! --- const ------------------------------
  94. character(len=*), parameter :: rname = mname//'/Diffusion_Init'
  95. ! --- local -------------------------------------
  96. integer :: region
  97. type(TrcFile) :: rcF
  98. ! --- begin --------------------------------
  99. ! loop over all (zoom) regions:
  100. do region = 1, nregions
  101. ! meteo used by diffusion
  102. call Set( pu_dat(region), status, used=.true. )
  103. call Set( pv_dat(region), status, used=.true. )
  104. call Set( temper_dat(region), status, used=.true. )
  105. call Set( humid_dat(region), status, used=.true. )
  106. call Set( gph_dat(region), status, used=.true. )
  107. call Set( oro_dat(region), status, used=.true. )
  108. call Set( lsmask_dat(region), status, used=.true. )
  109. call Set( sr_ecm_dat(region), status, used=.true. )
  110. call Set( sr_ols_dat(region), status, used=.true. )
  111. call Set( wspd_dat(region), status, used=.true. )
  112. call Set( sshf_dat(region), status, used=.true. )
  113. call Set( slhf_dat(region), status, used=.true. )
  114. call Set( ewss_dat(region), status, used=.true. )
  115. call Set( nsss_dat(region), status, used=.true. )
  116. end do
  117. call Init( rcF, rcfile , status)
  118. IF_NOTOK_RETURN(status=1)
  119. call ReadRc(rcF, 'diffusion.write', diffusion_write, status, default=.false.)
  120. IF_ERROR_RETURN(status=1)
  121. if (diffusion_write) then
  122. call ReadRc(rcF, 'diffusion.dir', diffusion_dir, status)
  123. IF_NOTOK_RETURN(status=1)
  124. endif
  125. call ReadRc(rcF, 'my.levs', mlevs, status)
  126. IF_NOTOK_RETURN(status=1)
  127. call Done(rcF, status)
  128. IF_NOTOK_RETURN(status=1)
  129. ! ok
  130. status = 0
  131. end subroutine Diffusion_Init
  132. ! ***
  133. subroutine Diffusion_Done( status )
  134. ! --- in/out --------------------------------
  135. integer, intent(out) :: status
  136. ! --- const ------------------------------
  137. character(len=*), parameter :: rname = mname//'/Diffusion_Done'
  138. ! --- begin --------------------------------
  139. ! nothing to be done
  140. ! ok
  141. status = 0
  142. end subroutine Diffusion_Done
  143. ! ================================================================
  144. !------------------------------------
  145. !
  146. ! Purpose:
  147. ! -------
  148. ! this subroutine reads and prepares all fields needed for the calculation of kzz
  149. !
  150. ! External
  151. ! --------
  152. ! dd_get_3_hourly_surface_fields
  153. ! dd_calc_ustar_raero_rb
  154. ! dd_coarsen_fields
  155. !
  156. ! Reference
  157. ! ---------
  158. ! Ganzeveld and Lelieveld (1996) and references therein.
  159. ! Adjusted for TM5, Bram Bregman, August 2003
  160. !
  161. !
  162. subroutine calc_kzz( status )
  163. use binas, only : grav
  164. use dims, only : im, jm, lmax_conv, idate, lm, okdebug
  165. use dims, only : nread, ndyn, tref, at, bt
  166. use dims, only : revert
  167. use global_data, only : mass_dat, wind_dat, region_dat
  168. #ifndef without_convection
  169. #ifndef without_diffusion
  170. use global_data, only : conv_dat
  171. #endif
  172. #endif
  173. use meteoData, only : spm_dat, pu_dat, pv_dat, m_dat
  174. use meteoData, only : temper_dat, humid_dat, gph_dat
  175. use tm5_distgrid, only : dgrid, Get_DistGrid, update_halo, dist_arr_stat
  176. integer, intent(out) :: status
  177. ! --- local ------------------------------
  178. character(len=10) :: c_time
  179. real,dimension(:,:,:),pointer :: phlb,gph,t,q,m,pu,pv
  180. real,dimension(:,:,:),pointer :: am,bm,cm ! fluxes in kg/timestep(region)
  181. #ifndef without_convection
  182. #ifndef without_diffusion
  183. real,dimension(:,:,:),pointer :: dkg !vert diff (ks/s)
  184. #endif
  185. #endif
  186. real,dimension(:),pointer :: dxyp
  187. integer,parameter :: avg_field=1, nlon360=360, nlat180=180
  188. integer :: i,region, nsend, ntimes
  189. real, dimension(:,:,:),allocatable,target :: m_adj, phlb_adj
  190. real, dimension(:,:),allocatable :: p_adj
  191. integer :: imr, jmr, lmr, j,l
  192. real, allocatable, target :: phlb_t(:,:,:)
  193. real, allocatable, target :: m_t(:,:,:)
  194. integer :: i1, i2, j1, j2
  195. character(len=*), parameter :: rname = mname//'/Calc_kzz'
  196. ! --- begin ----------------------------------
  197. if (okdebug) write(*,*) 'start of calc_kzz'
  198. ! Allocate work arrays, and fill ghost cells
  199. DO region = 1, nregions
  200. CALL GET_DISTGRID( DGRID(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  201. ALLOCATE(ustar_loc(region)%surf( i1:i2, j1:j2))
  202. ALLOCATE( sr_mix(region)%surf( i1:i2, j1:j2))
  203. CALL UPDATE_HALO( dgrid(region), wind_dat(region)%am, wind_dat(region)%halo, status)
  204. IF_NOTOK_RETURN(status=1)
  205. CALL UPDATE_HALO( dgrid(region), wind_dat(region)%bm, wind_dat(region)%halo, status)
  206. IF_NOTOK_RETURN(status=1)
  207. CALL UPDATE_HALO( dgrid(region), wind_dat(region)%cm, wind_dat(region)%halo, status)
  208. IF_NOTOK_RETURN(status=1)
  209. CALL UPDATE_HALO( dgrid(region), pu_dat(region)%data, pu_dat(region)%halo, status)
  210. IF_NOTOK_RETURN(status=1)
  211. CALL UPDATE_HALO( dgrid(region), pv_dat(region)%data, pv_dat(region)%halo, status)
  212. IF_NOTOK_RETURN(status=1)
  213. ENDDO
  214. ! -- for output of time header on debug fields
  215. write(c_time,'(i4,3(i2))') idate(1:4) !CMK corrected
  216. do i=1,10
  217. if (c_time(i:i)==' ') c_time(i:i)='0'
  218. enddo
  219. !
  220. ! Surface data sets may have the following characteristics
  221. ! 1. instanteneous. The time written in the file is valid from t-1.5 to t+1.5
  222. ! 2. Accumulated. The time written in the file is valid from t to t+3
  223. ! 3. Daily averaged. The time on the file is always 00 hours, and valid until t+24
  224. !
  225. ! *** It is essential to understand this timing when reading and applying these sets.
  226. !
  227. !
  228. ! fd2mk it has to be decided to 'save' fields like vgrat_low and daily average surface fields
  229. ! for later use, or to read it every 3 hours time step
  230. REG: do region = 1,nregions
  231. ! local grid sizes
  232. lmr=lm(region)
  233. call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  234. if (revert == -1) then
  235. !mkadj
  236. !mass and pressure are at t=t+dt in adjoint
  237. !these should be brought back to t to get
  238. !exactly the same diffusion as in forward run
  239. allocate ( m_adj (i1:i2, j1:j2, lmr ) )
  240. allocate ( phlb_adj(i1:i2, j1:j2, lmr+1) )
  241. allocate ( p_adj (i1:i2, j1:j2 ) )
  242. am => wind_dat(region)%am
  243. bm => wind_dat(region)%bm
  244. cm => wind_dat(region)%cm
  245. dxyp => region_dat(region)%dxyp
  246. ntimes = (nread/(ndyn))*tref(region)
  247. m_adj = m_dat(region)%data(i1:i2, j1:j2, 1:lmr)
  248. do i=1,ntimes
  249. m_adj = m_adj + revert*am(i1-1:i2-1, j1 :j2 , 1:lmr ) & !east
  250. - revert*am(i1 :i2 , j1 :j2 , 1:lmr ) & !west
  251. + revert*bm(i1 :i2 , j1 :j2 , 1:lmr ) & !south
  252. - revert*bm(i1 :i2 , j1+1:j2+1, 1:lmr ) & !north
  253. + revert*cm(i1 :i2 , j1 :j2 , 0:lmr-1) & !lower
  254. - revert*cm(i1 :i2 , j1 :j2 , 1:lmr ) !upper
  255. enddo
  256. p_adj(:,:) = 0.0
  257. do l=1,lmr
  258. ! note that on the EDGES masses are coarse>
  259. ! diffusion iis applied only on core zoom>
  260. do j = j1, j2
  261. do i = i1, i2
  262. p_adj(i,j) = p_adj(i,j) + m_adj(i,j,l)*grav/dxyp(j)
  263. end do
  264. end do
  265. end do
  266. do l=1,lmr+1
  267. do j = j1, j2
  268. do i = i1, i2
  269. phlb_adj(i,j,l) = at(l)+bt(l)*p_adj(i,j)
  270. end do
  271. end do
  272. end do
  273. deallocate(p_adj)
  274. nullify(am)
  275. nullify(bm)
  276. nullify(cm)
  277. nullify(dxyp)
  278. phlb => phlb_adj
  279. m => m_adj
  280. else
  281. ! Bug fix 2006-09-21 (AJS)
  282. ! The two half level pressure arrays phbl_t and phlb_k
  283. ! are sometimes not updated together when running in parallel
  284. ! with zoom regions. To avoid problems, pressure and mass
  285. ! are computed here given the surface pressure in the middle
  286. ! of the current time interval, stored in spm_dat.
  287. !--
  288. !phlb => mass_dat(region)%phlb_t
  289. !m => mass_dat(region)%m_t
  290. !--
  291. allocate( phlb_t(i1:i2, j1:j2, 1:lmr+1) )
  292. do l = 1, lmr+1
  293. phlb_t(:,:,l) = at(l) + bt(l)*spm_dat(region)%data(i1:i2, j1:j2, 1)
  294. end do
  295. allocate( m_t(i1-2:i2+2, j1-2:j2+2, 1:lmr) )
  296. m_t = 0.0
  297. do l = 1, lmr
  298. do j = j1, j2
  299. m_t(i1:i2,j,l) = ( phlb_t(:,j,l) - phlb_t(:,j,l+1) )*region_dat(region)%dxyp(j)/grav
  300. end do
  301. end do
  302. phlb => phlb_t
  303. m => m_t
  304. !--
  305. endif
  306. pu => pu_dat(region)%data
  307. pv => pv_dat(region)%data
  308. gph => gph_dat(region)%data
  309. t => temper_dat(region)%data
  310. q => humid_dat(region)%data
  311. #ifndef without_convection
  312. #ifndef without_diffusion
  313. dkg => conv_dat(region)%dkg
  314. #endif
  315. #endif
  316. ! ustar_loc
  317. call dd_calc_ustar(region, i1, i2, j1, j2)
  318. if (okdebug) then
  319. ! call dd_field_statistics(region, 'ustar_loc', ustar_loc(region)%surf, i1, j1, status)
  320. call dist_arr_stat(dgrid(region), 'ustar_loc', ustar_loc(region)%surf, 0, status)
  321. end if
  322. #ifndef without_convection
  323. #ifndef without_diffusion
  324. ! calculate kzz
  325. call bldiff( region, phlb, m, i1, i2, j1, j2 )
  326. #endif
  327. #endif
  328. nullify(phlb)
  329. nullify(m)
  330. !--
  331. deallocate( phlb_t )
  332. deallocate( m_t )
  333. !--
  334. nullify(pu)
  335. nullify(pv)
  336. #ifndef without_convection
  337. #ifndef without_diffusion
  338. nullify(dkg)
  339. #endif
  340. #endif
  341. nullify(t)
  342. nullify(q)
  343. nullify(gph)
  344. if(revert == -1) then
  345. deallocate(m_adj)
  346. deallocate(phlb_adj)
  347. endif
  348. end do REG
  349. ! Done
  350. do region = 1,nregions
  351. deallocate(sr_mix(region)%surf)
  352. deallocate(ustar_loc(region)%surf)
  353. end do
  354. if (okdebug) write(*,*) 'end of calc_kzz'
  355. end subroutine calc_kzz
  356. !------------------------------------------------------------------------------
  357. ! TM5 !
  358. !------------------------------------------------------------------------------
  359. !BOP
  360. !
  361. ! !IROUTINE: dd_calc_ustar
  362. !
  363. ! !DESCRIPTION:
  364. !
  365. ! Calculate friction velocity (ustar)
  366. !
  367. ! Uses the log normal wind profile information stored by ECMWF in 10 meter wind
  368. ! to estimate a local ustar over land
  369. ! uses Charnock equation over sea to estimate ustar
  370. ! aerodynamic resistance is calculated from heat fluxes and ustar using a constant reference height
  371. ! sub laminar resistance from ustar
  372. !
  373. ! Reference:
  374. ! o M.Z. Jacobson, "Fundamentals of Atmospheric Modeling", sections 8.3-8.4
  375. !
  376. ! !REVISION HISTORY:
  377. !
  378. ! 2003 ??; Ge Verver (personal communication, 2003)
  379. !
  380. ! 2010-08 Arjo Segers
  381. ! Debugged computation of u* over sea:
  382. ! abs. surface stress should be devided by air density.
  383. ! Set minum value for u* to trap division by zero (MK).
  384. !
  385. ! 2012-01 P. Le Sager
  386. ! adapted for lon-lat MPI domain decomposition
  387. !------------------------
  388. subroutine dd_calc_ustar(region, i1, i2, j1, j2)
  389. use binas, only : grav, vKarman
  390. use dims, only : okdebug
  391. use meteoData, only : lsmask_dat, sr_ecm_dat, sr_ols_dat
  392. use meteoData, only : wspd_dat
  393. use meteoData, only : slhf_dat, sshf_dat
  394. use meteoData, only : ewss_dat, nsss_dat
  395. !--- in/out ------------------------------------
  396. integer,intent(in) :: region, i1, i2, j1, j2
  397. ! --- const -------------------------------------
  398. ! Garret relation:
  399. real,parameter :: alfa_garret = 0.11
  400. ! air density at seal level and 293 K :
  401. real, parameter :: rho_air = 1.2 ! kg/m3
  402. ! Charnock relation:
  403. real,parameter :: alfa_charnock = 0.018
  404. ! kinematic viscosity of air at about 300 K :
  405. real,parameter :: nu_air = 1.5e-5 ! m2/s
  406. real, parameter :: Href=30. ! constant reference height for calculations of raero
  407. ! some constants specific to calculation of ustar and raero
  408. real,parameter :: rhoCp = 1231.0
  409. real,parameter :: rhoLv = 3013000.0
  410. real, parameter :: rz0 =2.0
  411. !--- local -------------------------------------
  412. integer :: i,j
  413. real :: tau_z_abs
  414. real :: ustar_sea,ustar_land,xland,sr_sea,sr_land,sr_help,wind10m
  415. ! --- begin ----------------------------------
  416. do j=j1, j2
  417. do i=i1, i2
  418. xland=lsmask_dat(region)%data(i,j,1)/100. ! 0-1
  419. ! SEA:
  420. !
  421. ! From M.Z. Jacobson, "Fundamentals of Atmospheric Modeling",
  422. ! section "8.3 Friction velocity" :
  423. ! The friction velocity :
  424. ! u* = sqrt( |tau_z|/rho_a ) [m/s]
  425. ! where:
  426. ! tau_z = (ewss,nsss) surface stress [N/m2=kg/m/s2]
  427. ! rho_air air density [kg/m3] ; about 1.2 at sea level and 293 K
  428. !
  429. tau_z_abs = sqrt(ewss_dat(region)%data(i,j,1)**2+&
  430. nsss_dat(region)%data(i,j,1)**2 ) ! N/m2
  431. ustar_sea = sqrt(tau_z_abs/rho_air) ! m/s
  432. !
  433. ! limitation to avoid division by zero:
  434. ustar_sea = max( 0.001, ustar_sea ) ! m/s ; minium of 1 mm/s
  435. !
  436. ! From M.Z. Jacobson, "Fundamentals of Atmospheric Modeling",
  437. ! section "8.4 Surface roughness lengths" :
  438. ! "For smooth surfaces, such as over a smooth ocean with low wind speeds (Garret):
  439. ! z0m = 0.11 nu_a / u*
  440. ! where nu_a is the kinematic viscosity of air.
  441. ! Over a rough ocean with high wind speeds:
  442. ! z0m = alpha_c (u*)**2 / grav
  443. ! which is the 'Charnock relation',
  444. ! where alpha_c ~ 0.016 is the 'Charnock constant'."
  445. !
  446. ! Here the sum of these two parameterizations is used,
  447. ! where the first part is dominant for u*<0.1 and the second for u*>0.1 ;
  448. ! minimum value for u* prevents division by zero:
  449. !
  450. sr_sea = alfa_garret * nu_air / ustar_sea + &
  451. alfa_charnock * ustar_sea**2 / grav
  452. !
  453. ! LAND
  454. ! calculate the 'local' surface roughness for momentum that is consistent with 10 m wind
  455. ! and the Olsson data base, we assume that the windspeed at 75 m is independent of
  456. ! surface roughness
  457. !
  458. wind10m = wspd_dat(region)%data(i,j,1)
  459. if ( xland > 0.0 ) then
  460. !>>> TvN
  461. ! Over land, the friction velocity ustar
  462. ! is calculated from the 10 m wind speed, u10.
  463. ! In IFS u10 is a diagnostic variable, calculated to be
  464. ! compatible with "open-terrain" wind speed observations.
  465. ! Over rough or inhomogeneous terrain (z0M > 0.03 m),
  466. ! it is calculated using an aerodynamic roughness
  467. ! length typical for open terrain with grassland (0.03 m)
  468. ! (see IFS cy31r1 documentation, p. 46).
  469. ! In the expressions applied in TM5,
  470. ! the stability profile functions Psi_M have disappeared,
  471. ! and only the logarithmic terms are kept.
  472. ! Apart from this, the expression for ustar was incorrect:
  473. ! 10./sr_land should be 75./sr_land.
  474. ! This has now been corrected.
  475. ! Moreover, over islands and near coast lines,
  476. ! zeros can occur in sr_ols_dat,
  477. ! which have to be removed.
  478. ! However, a lower bound of 1e-2 m = 1 cm
  479. ! seems too high for smooth surfaces,
  480. ! like deserts and ice caps.
  481. ! The minimum positive value in sr_ols_dat
  482. ! is 2.5e-4 m = 0.025 cm,
  483. ! which is obtained in parts of Antarctica.
  484. ! This is likely the result of regridding of
  485. ! a field with minimum value of 0.1 cm
  486. ! from 0.5x0.5 to 1x1 degrees.
  487. ! Please note that with the introduction of CY31R1,
  488. ! the orographic contribution to the aerodynamic roughness length
  489. ! in IFS has been replaced by an explicit specification of stress
  490. ! on model levels due to turbulent orographic form drag,
  491. ! and the climatological variable previously used
  492. ! has been replaced by a prognostic variable.
  493. ! In TM5, the climatological variable is still used,
  494. ! but it would be better to use the prognostic variable instead.
  495. ! Note also that the same calculation
  496. ! is done in dry_deposition.F90.
  497. !sr_land=max(1e-2,sr_ols_dat(region)%data(i,j,1)) ! occurs at Islands, etc.
  498. sr_land=max(1e-3,sr_ols_dat(region)%data(i,j,1))
  499. sr_help=min(sr_ecm_dat(region)%data(i,j,1),0.03)
  500. !fd ustar_land=vKarman*wind10m(i,j)/alog(10./sr_help) !ustar consistent with 'clipped' large scale roughness
  501. !ustar_land=vKarman*wind10m/&
  502. ! alog(10./sr_help)*alog(75./sr_help)/alog(10./sr_land)
  503. ustar_land=vKarman*wind10m/&
  504. alog(10./sr_help)*alog(75./sr_help)/alog(75./sr_land)
  505. ! <<< TvN
  506. else
  507. sr_land = 0.0
  508. ustar_land = 0.0
  509. endif
  510. ! interpolate between land and sea for mixed pixels:
  511. ustar_loc(region)%surf(i,j) = xland*ustar_land+(1-xland)*ustar_sea
  512. sr_mix (region)%surf(i,j) = xland* sr_land+(1-xland)* sr_sea
  513. end do !i
  514. end do !j
  515. if (okdebug) write(*,*) 'end calc_ustar'
  516. end subroutine dd_calc_ustar
  517. ! ***
  518. #ifndef without_convection
  519. #ifndef without_diffusion
  520. !------------------------------------------------------------------
  521. ! Procedure:
  522. ! (1) calc_Kv: Diffusion coefficients full atmosphere from shear
  523. ! (2) pblhght: Determine the boundary layer height
  524. ! (3) difcoef: Replace diffusion coefficients in the boundary layer
  525. !------------------------------------------------------------------
  526. ! subroutine adjusted for TM5, Bram Bregman, August 2003.
  527. !-------------------------IO---------------------------------------
  528. subroutine bldiff( region, phlb, m, i1, i2, j1, j2 )
  529. use binas, only: ae, cp_air, Rgas, grav, Lvap
  530. use binas, only: xmair
  531. use binas, only: vkarman
  532. use dims, only: im, jm, lm, lmax_conv
  533. use dims, only: okdebug
  534. use dims, only: gtor, dx, dy, ybeg, xref, yref
  535. use dims, only: isr, ier, jsr, jer
  536. use toolbox, only: dumpfield
  537. use global_data, only : conv_dat
  538. use meteodata , only : slhf_dat, sshf_dat, pu_dat, pv_dat, gph_dat, temper_dat, humid_dat
  539. ! --- in/out ----------------------------------------------
  540. integer,intent(in) :: region, i1, i2, j1, j2
  541. real,dimension(:,:,:),pointer :: phlb, gph, t, q, m, pu, pv, dkg
  542. ! --- const ----------------------------------
  543. ! minimum value diffusion coefficient
  544. real, parameter :: ckmin = 1.e-15
  545. real, parameter :: sffrac= 0.1
  546. real, parameter :: onet = 1./3.
  547. real, parameter :: rair = Rgas*1.e3/xmair
  548. real, parameter :: ccpq = 1869.46/cp_air -1.
  549. real, parameter :: epsilo = rair/461.51
  550. real, parameter :: apzero = 100000. ! reference pressure (usually 1000 hPa)
  551. real, parameter :: ccon = sffrac*vkarman
  552. real, parameter :: betam=15.,betah=15.
  553. real, parameter :: binm = betam*sffrac
  554. real, parameter :: binh = betah*sffrac
  555. real, parameter :: fr_excess = 1.
  556. real, parameter :: fak = fr_excess * 8.5
  557. real, parameter :: fakn = fr_excess * 7.2
  558. real, parameter :: zkappa = rair / cp_air
  559. real, parameter :: zcrdq = 1.0/epsilo - 1.0
  560. real, parameter :: tiny = 1.e-9 ! bound wind**2 to avoid dividing by 0
  561. real, parameter :: zacb = 100. !factor in ustr-shearproduction term
  562. real, parameter :: ricr = 0.3 ! critical richardson number
  563. real, parameter :: ssfrac=0.1
  564. ! --- local ------------------------------------
  565. real, dimension(:,:,:), allocatable:: wkvf,& ! outer layer Kv (at the top of each layer)
  566. zdup2,& ! vertical shear squared
  567. zrinub,& ! Richardson number
  568. wthm,& ! potential temperature
  569. zthvk,& ! virtual potential temperature
  570. pf,& ! full level pressure
  571. whm, & ! height of layer centers
  572. uwind,vwind,& ! wind velocities [m s-1]
  573. kvh !
  574. real, dimension(:,:), allocatable :: ustr,& ! velocity scale
  575. wheat,& ! surface sensible heat flux [K m/s]
  576. wqflx,& ! surface latent heat flux [kg m/ (kg s)]
  577. wobkl,& ! Monin Obukhov length
  578. wheatv,& ! surface virtual heat flux [K m/s]
  579. wtseff ! theta_v at lowest level [K]
  580. integer :: i, j, l
  581. real :: intfac !cmk BUG was integer!
  582. ! variables for the calculation of free atmosphere vertical diffusion coefficients wkvf
  583. real :: cml2,zlambdac,arg,zdz,zthva,zfunst, zfstabh,zsstab,zkvn,dyy
  584. real, dimension(j1:j2) :: dxx, lat
  585. integer :: jq! jqif function
  586. ! variabales for the calculation of the boundary layer height
  587. real :: wrino ! Richardson number
  588. real :: thvparc ! parcel virtual potential temperature
  589. real :: fmt,wsc,tkv,zrino,vvk,dtkv,zexcess
  590. integer :: jwork
  591. real,dimension(:,:),pointer :: pblh ! planetary boundary layer height [m]
  592. ! variables to calculate difusion coefficient in the boundary layer
  593. !real :: wsc,cml2
  594. real :: z,zwstr,zkvh,kvhmin,zm,zp,zslask,zsl1, &
  595. zh,zzh,zl,pr,zpblk,zfstabh1,zfstabh2,&
  596. term1,term2,term3,obukov
  597. integer :: jq1,jq2,jq3,jq4,jq5,jck,jqq, jq6, jq7
  598. real :: yres
  599. real :: xres
  600. real :: thvgrad, kvhentr
  601. ! --- begin ------------------------------------------
  602. if (okdebug) write(*,*) 'start of bldiff'
  603. allocate(wkvf ( i1:i2, j1:j2, lm(region)))
  604. allocate(zdup2 ( i1:i2, j1:j2, lm(region)))
  605. allocate(zrinub( i1:i2, j1:j2, lm(region)))
  606. allocate(wthm ( i1:i2, j1:j2, lm(region)))
  607. allocate(zthvk ( i1:i2, j1:j2, lm(region)))
  608. allocate(pf ( i1:i2, j1:j2, lm(region)))
  609. allocate(whm ( i1:i2, j1:j2, lm(region)))
  610. allocate(uwind ( i1:i2, j1:j2, lm(region)))
  611. allocate(vwind ( i1:i2, j1:j2, lm(region)))
  612. allocate(kvh ( i1:i2, j1:j2, lm(region)))
  613. allocate(ustr ( i1:i2, j1:j2 ))
  614. allocate(wheat ( i1:i2, j1:j2 ))
  615. allocate(wqflx ( i1:i2, j1:j2 ))
  616. allocate(wobkl ( i1:i2, j1:j2 ))
  617. allocate(wheatv( i1:i2, j1:j2 ))
  618. allocate(wtseff( i1:i2, j1:j2 ))
  619. pblh => conv_dat(region)%blh
  620. pu => pu_dat(region)%data
  621. pv => pv_dat(region)%data
  622. gph => gph_dat(region)%data
  623. t => temper_dat(region)%data
  624. q => humid_dat(region)%data
  625. dkg => conv_dat(region)%dkg
  626. ! mid-layer pressure levels
  627. ! potential temperature field
  628. ! virtual potential temperature field
  629. !
  630. !PLS $OMP PARALLEL &
  631. !PLS $OMP default (none) &
  632. !PLS $OMP shared (lm, jsr, jer, isr, ier, region, pf, phlb, t, q, wthm, &
  633. !PLS $OMP zthvk, whm, gph) &
  634. !PLS $OMP private (i, j, l, jss, jee, intfac)
  635. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  636. do l=1,lm(region)
  637. do j=j1,j2
  638. do i=i1,i2
  639. pf(i,j,l) = (phlb(i,j,l) + phlb(i,j,l+1) )*0.5
  640. wthm(i,j,l) = t(i,j,l) * ( apzero/pf(i,j,l) )**zkappa
  641. zthvk(i,j,l) = wthm(i,j,l) * (1.0 + zcrdq*q(i,j,l))
  642. ! mid-layer heights (note: gph array is for levels 1:lm+1)
  643. intfac = (phlb(i,j,l)-pf(i,j,l)) / (phlb(i,j,l)-phlb(i,j,l+1))
  644. ! cmk: intfac is now real!
  645. whm(i,j,l) = (gph(i,j,l)-gph(i,j,1)) * (1.0-intfac) + (gph(i,j,l+1) -gph(i,j,1)) * intfac
  646. ! cmk bug corrected:
  647. ! cmk old: whm(i,j,l) = (gph(i,j,l+1)-gph(i,j,1)) * (1.0-intfac) + (gph(i,j,l+1) -gph(i,j,1)) * intfac
  648. end do
  649. end do
  650. end do
  651. !PLS $OMP END PARALLEL
  652. ! convert units of sensible and latent heat fluxes
  653. ! sshf: W/m2 / (J/kg/K) / (kg/m3) = J/s/m2 /J kg K /kg m3 = K m/s
  654. ! slhf: W/m2 / (J/kg) / (kg/m3) = J/s/m2 /J kg /kg m3 = m/s
  655. ! ustr is local velocity scale, NOT from ECMWF surface stresses
  656. !
  657. !PLS $OMP PARALLEL &
  658. !PLS $OMP default (none) &
  659. !PLS $OMP shared (jm, im, region, wheat, sshf_dat, gph, phlb, slhf_dat, &
  660. !PLS $OMP wqflx, ustr, ustar_loc) &
  661. !PLS $OMP private (i, j, js, je)
  662. ! call my_omp_domain (1, jm(region), js, je)
  663. do j = j1,j2
  664. do i = i1,i2
  665. wheat(i,j) = -sshf_dat(region)%data(i,j,1) * (gph(i,j,2) - &
  666. gph(i,j,1)) * grav / (phlb(i,j,1) - &
  667. phlb(i,j,2)) / cp_air
  668. wqflx(i,j) = -slhf_dat(region)%data(i,j,1) * (gph(i,j,2) - &
  669. gph(i,j,1)) * grav / (phlb(i,j,1) - &
  670. phlb(i,j,2)) / lvap
  671. ustr(i,j) = max (ustar_loc(region)%surf(i,j), 0.01)
  672. end do
  673. end do
  674. !PLS $OMP END PARALLEL
  675. !
  676. !-----------------------------------------------------------------
  677. ! compute the free atmosphere vertical diffusion coefficients wkvf
  678. ! height-dependent asymptotic mixing length implemented
  679. ! Holtslag and Boville, 1993, J. Climate, 6, 1825-1842.
  680. !------------------------------------------------------------------
  681. ! first calculate wind velocities [m s-1] from the mass fluxes
  682. yres = dy/yref(region)
  683. xres = dx/xref(region)
  684. do j = j1,j2
  685. lat(j) = ybeg(region) + 0.5 * yres + yres * (j-1)
  686. enddo
  687. dxx(j1:j2) = ae * xres * gtor * cos(lat(j1:j2)*gtor)
  688. dyy = ae * yres * gtor
  689. !PLS $OMP PARALLEL &
  690. !PLS $OMP default (none) &
  691. !PLS $OMP shared (lm, jsr, jer, isr, ier, region, dxx, dyy, pu, pv, m, &
  692. !PLS $OMP uwind, vwind) &
  693. !PLS $OMP private (i, j, l, jss, jee)
  694. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  695. do l=1,lm(region)
  696. do j=j1,j2
  697. do i=i1,i2
  698. uwind(i,j,l) = dxx(j)*(pu(i,j,l) + pu(i-1,j,l))*0.5 / m(i,j,l)
  699. vwind(i,j,l) = dyy *(pv(i,j,l) + pv(i,j+1,l))*0.5 / m(i,j,l)
  700. enddo
  701. enddo
  702. enddo
  703. !PLS $OMP END PARALLEL
  704. !
  705. ! Initialize gradient Richardson number and free tropospheric k values
  706. !
  707. if (okdebug) write(*,*) ' Richardson number'
  708. !PLS $OMP PARALLEL &
  709. !PLS $OMP default (none) &
  710. !PLS $OMP shared (lm, jsr, jer, isr, ier, region, gph, zdup2, uwind, &
  711. !PLS $OMP vwind, whm, phlb, pf, zthvk, zrinub, wkvf) &
  712. !PLS $OMP private (i, j, l, jss, jee, jq, zlambdac, cml2, zdz, arg, &
  713. !PLS $OMP zthva, zsstab, zfunst, zfstabh, zkvn)
  714. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  715. do l=1,lm(region)-1
  716. do j=j1,j2
  717. do i=i1,i2
  718. !
  719. ! if (delta gph < 1000) then jq=1
  720. jq = jqif( 1000.0, (gph(i,j,l+1)-gph(i,j,1)) )
  721. zlambdac=jq*300.+(1-jq)*(30. + 270. * exp (1.-(gph(i,j,l+1)-gph(i,j,1))/1000.))
  722. cml2=(1./( 1./zlambdac + 1./(vkarman*(gph(i,j,l+1)-gph(i,j,1))) ))**2.
  723. !
  724. ! vertical shear squared,
  725. ! min value of (delta v)**2 prevents zero shear
  726. !
  727. zdup2(i,j,l) = (uwind(i,j,l+1)-uwind(i,j,l))**2 + &
  728. (vwind(i,j,l+1)-vwind(i,j,l))**2
  729. zdup2(i,j,l) = max(zdup2(i,j,l),1.e-10)
  730. zdz = whm(i,j,l+1) - whm(i,j,l)
  731. zdup2(i,j,l) = zdup2(i,j,l)/(zdz**2)
  732. !
  733. ! static stability (use virtual potential temperature)
  734. ! dv now calculated at interface height (lineair in pressure)
  735. !
  736. arg=(log(phlb(i,j,l+1))-log(pf(i,j,l)))/(log(pf(i,j,l+1))-log(pf(i,j,l)))
  737. zthva = zthvk(i,j,l) + arg * (zthvk(i,j,l+1)-zthvk(i,j,l))
  738. zsstab = grav/zthva*(zthvk(i,j,l+1) - zthvk(i,j,l))/zdz
  739. !
  740. ! gradient Richardson number
  741. !
  742. zrinub(i,j,l) = zsstab/zdup2(i,j,l)
  743. !
  744. ! stability functions
  745. !
  746. zfunst = max(1. - 18.*zrinub(i,j,l),0.)
  747. zfunst = sqrt(zfunst)
  748. !
  749. zfstabh = 1.0 / (1.0 + 10.0*zrinub(i,j,l)*(1.0+8.0*zrinub(i,j,l)))
  750. !
  751. ! neutral diffusion coefficient
  752. !
  753. zkvn = cml2 * sqrt(zdup2(i,j,l))
  754. !
  755. ! correction with stability functions
  756. !
  757. jq = jqif( zrinub(i,j,l),0.0 ) ! thus: if (zrinub>0) then jq=1
  758. wkvf(i,j,l)=jq*(max(ckmin,zkvn*zfstabh))+(1-jq)*(max(ckmin,zkvn*zfunst))
  759. end do
  760. end do
  761. end do
  762. !PLS $OMP END PARALLEL
  763. !
  764. ! Determine the boundary layer height
  765. !--------------------------------------------------------------------------------
  766. ! based on: Vogelezang and Holtslag, 1996, Bound. Layer Meteorol., 81, 245-269.
  767. !--------------------------------------------------------------------------------
  768. if (okdebug) write(*,*) ' boundary layer height'
  769. !PLS $OMP PARALLEL &
  770. !PLS $OMP default (none) &
  771. !PLS $OMP shared (lm, jsr, jer, isr, ier, region, wtseff, wthm, q, &
  772. !PLS $OMP wheatv, wheat, wqflx, wobkl, ustr, pblh, uwind, vwind, &
  773. !PLS $OMP whm ) &
  774. !PLS $OMP private (i, j, l, jss, jee, jwork, wrino, vvk, tkv, dtkv, &
  775. !PLS $OMP zrino, jq, fmt, wsc, zexcess, thvparc )
  776. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  777. do j=j1,j2
  778. do i=i1,i2
  779. ! compute bottom level virtual temperature and heat flux and Monin-Obukhov Length
  780. wtseff(i,j) = wthm(i,j,1)*(1. + zcrdq*q(i,j,1))
  781. wheatv(i,j) = wheat(i,j) + zcrdq*wthm(i,j,1)*wqflx(i,j)
  782. wobkl(i,j) = -wtseff(i,j)*ustr(i,j)**3 / (grav*vkarman*(wheatv(i,j)+sign(1.e-10,wheatv(i,j))) )
  783. ! initialize
  784. pblh(i,j) = 0.0
  785. jwork = 1
  786. wrino = 0.0
  787. do l=1,lm(region)
  788. vvk = (uwind(i,j,l)-uwind(i,j,1))**2 + (vwind(i,j,l)-vwind(i,j,1))**2 + &
  789. zacb*ustr(i,j)*ustr(i,j)
  790. vvk = max(vvk,tiny)
  791. tkv = wthm(i,j,l)*(1. + zcrdq*q(i,j,l))
  792. dtkv = tkv-wtseff(i,j)
  793. !
  794. ! Bulk Richardson number
  795. !
  796. zrino = grav*dtkv * (whm(i,j,l)-whm(i,j,1)) / (wtseff(i,j)*vvk)
  797. zrino = zrino+sign(1.e-10,zrino) ! prevent zrino becoming zero
  798. zrino = jwork*zrino
  799. jq = jqif(ricr,zrino) ! thus: if (zrino < ricr) then jq=1
  800. !
  801. ! calculate pblh from linear interpolation in Ri(bulk)
  802. !
  803. if ( l == 1 ) then
  804. pblh(i,j) = jq*pblh(i,j) + (1-jq) * &
  805. ( (ricr-wrino)/(zrino-wrino)*(whm(i,j,l) ) )
  806. else
  807. pblh(i,j) = jq*pblh(i,j) + (1-jq) * &
  808. ( whm(i,j,l-1) + (ricr-wrino)/(zrino-wrino)*(whm(i,j,l)-whm(i,j,l-1)) )
  809. end if
  810. !
  811. ! first time zrino > ricr we set jwork to zero, thus all further times zrino=0
  812. ! and pblh does not change anymore
  813. !
  814. jwork = jq*jwork
  815. !
  816. ! if pblh is found already avoid dividing by zero next level by a faked value
  817. ! of wrino (jwork=0 already)
  818. !
  819. wrino = zrino+(1-jwork)*0.1
  820. end do
  821. pblh(i,j) = (1-jwork)*pblh(i,j) + whm(i,j,lm(region))*jwork
  822. !
  823. ! second calculation of pblh including excess surface temperature using
  824. ! a convective velocity scale wsc (wm, not equal to w*!!!), using result
  825. ! of the first calculation of pblh
  826. !
  827. jq = jqif(wheatv(i,j),0.) ! thus: if (wheatv(i,j) > 0.) then jq=1
  828. jwork = jq
  829. fmt = ( jq*(1.0 - binm*pblh(i,j)/wobkl(i,j)) + (1.0-jq) )**onet
  830. wsc = ustr(i,j)*fmt
  831. zexcess = wheatv(i,j)*fak/wsc
  832. !
  833. ! improve pblh estimate under convective conditions using
  834. ! convective temperature excess (zexcess)
  835. !
  836. thvparc = wtseff(i,j)+jq*zexcess
  837. !
  838. jwork = 1
  839. wrino = 0.0
  840. !
  841. do l=2,lm(region)
  842. !
  843. ! Bulk Richardson number:
  844. ! if stable corrected with extra shear term
  845. ! if unstable now also corrected for temperature excess
  846. !
  847. vvk = (uwind(i,j,l)-uwind(i,j,1))**2 + (vwind(i,j,l)-vwind(i,j,1))**2 + &
  848. zacb*ustr(i,j)*ustr(i,j)
  849. vvk = max(vvk,tiny)
  850. tkv = wthm(i,j,l)*(1. + zcrdq*q(i,j,l))
  851. zrino = grav*(tkv-thvparc) * (whm(i,j,l)-whm(i,j,1))/(vvk*wtseff(i,j))
  852. zrino = zrino+sign(1.e-10,zrino) ! prevent zrino becoming zero
  853. zrino = zrino*jwork
  854. jq = jqif(ricr,zrino) ! thus: if (zrino < ricr) then jq=0
  855. !
  856. ! calculate pblh from linear interpolation in Ri(bulk)
  857. !
  858. pblh(i,j) = jq*pblh(i,j)+(1-jq)* (whm(i,j,l-1)+(ricr-wrino)/ &
  859. (zrino-wrino)*(whm(i,j,l)-whm(i,j,l-1)))
  860. !
  861. ! first time zrino > ricr we set jwork to zero, thus all further times zrino=0
  862. ! and pblh does not change anymore
  863. !
  864. jwork = jq*jwork
  865. !
  866. ! if pblh is found already avoid dividing by zero next level by a faked value
  867. ! of wrino (jwork=0 already)
  868. !
  869. wrino = zrino+(1-jwork)*0.1
  870. end do
  871. pblh(i,j) = (1-jwork)*pblh(i,j) + jwork*whm(i,j,lm(region))
  872. pblh(i,j) = max(pblh(i,j), 100.0) ! set minimum value for pblh
  873. end do
  874. end do
  875. !PLS $OMP END PARALLEL
  876. !
  877. ! Diffusion coefficients in the surface layer
  878. !
  879. !--------------------------------------------------------------------------------
  880. ! the revised LTG scheme (Beljaars and Viterbo, 1998)
  881. ! Modification - April 1, 1999: prevent zpblk to become too small in very stable conditions
  882. !--------------------------------------------------------------------------------
  883. if (okdebug) write(*,*) ' diffusion coeff'
  884. ! initialize output array with minimum value
  885. kvhmin = 0.1
  886. !
  887. ! Loop to calculate the diffusivity
  888. !
  889. !PLS $OMP PARALLEL &
  890. !PLS $OMP default (none) &
  891. !PLS $OMP shared (lm, jsr, jer, isr, ier, region, pblh, gph, wobkl, &
  892. !PLS $OMP wheatv, whm, zrinub, zdup2, kvhmin, wkvf, ustr, &
  893. !PLS $OMP wtseff, kvh, zthvk) &
  894. !PLS $OMP private (i, j, l, jss, jee, jq, z, zh, zl, zzh, jq1, jq2, jq3, &
  895. !PLS $OMP jq4, jq5, jq6, jq7, jqq, cml2, zfstabh, zpblk, zkvh, &
  896. !PLS $OMP zslask, fmt, wsc, zwstr, obukov, term1, term2, term3, &
  897. !PLS $OMP pr, thvgrad, kvhentr)
  898. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  899. do l=1,lm(region)-1
  900. do j=j1,j2
  901. do i=i1,i2
  902. jq = jqif(pblh(i,j),(gph(i,j,l+1)-gph(i,j,1))) ! if top of level hh < pblh jq=1 inside BL
  903. z = gph(i,j,l+1)-gph(i,j,1)
  904. zh = jq*z/pblh(i,j) ! only defined in BL (jq=1)
  905. ! zl must not be zero
  906. zl = jq*z/wobkl(i,j)+(1-jq) ! z/L outside BL = 1
  907. zzh = jq*(1.-zh)**2 ! (1-z/h)^2 only in BL
  908. jq1 = jqif(wheatv(i,j),0.) ! jq1 is 1 in unstable case
  909. jq2 = jq*(1-jq1) ! jq2 is 1 in stable BL
  910. jq3 = jq*jq1 ! jq3 is 1 in unstable BL
  911. jq6 = jqif(pblh(i,j), whm(i,j,l))
  912. jq7 = jqif(whm(i,j,l+1), pblh(i,j))
  913. jq6 = jq1*jq6*jq7 ! jq6 is 1 at at the kvh-level that is located between the
  914. ! mid-levels that surround pblh
  915. jq1 = jqif(sffrac,zh) ! jq1 is 1 in surface layer
  916. jq4 = jq3*jq1 ! jq4 is 1 in unstable surface layer
  917. jq5 = jq3*(1-jq1) ! jq5 is 1 outside unstable surface layer
  918. !
  919. ! calculate coefficients for momentum, the values for heat are obtained by dividing by Prantl number
  920. !
  921. ! stable and neutral:
  922. !
  923. cml2 = (1./( 1./(vkarman*z) + 1./450. ))**2.
  924. jqq = jqif(zrinub(i,j,l),0.)
  925. !
  926. ! stability function for momentum
  927. !
  928. zfstabh = jqq/(1.+10.*zrinub(i,j,l)*sqrt(1.+jqq*zrinub(i,j,l))) + (1-jqq)
  929. zpblk = cml2 * sqrt(zdup2(i,j,l))*zfstabh
  930. zkvh = jq2*max(zpblk,kvhmin)+(1-jq2)*wkvf(i,j,l) ! take free troposp. values outside BL
  931. !
  932. ! unstable case in surface layer
  933. !
  934. zslask = 1.-betah*zl
  935. zslask = jq4*zslask+(1-jq4)
  936. zpblk = (1-jq4)*zkvh+jq4*(ustr(i,j)*vkarman*z*zzh*zslask**(1./3.))
  937. !
  938. ! unstable case above surface layer
  939. !
  940. ! evaluate stability function for momentum at height z = 0.1*pblh (top of surface layer)
  941. !
  942. zslask = 1.-ssfrac*betah*pblh(i,j)/wobkl(i,j)
  943. zslask = jq5*zslask+(1-jq5)
  944. fmt = zslask**(1./3.)
  945. !
  946. ! use surface layer definition for w_m
  947. wsc = ustr(i,j)*fmt
  948. zpblk = (1-jq5)*zpblk + jq5*wsc*vkarman*z*zzh
  949. !
  950. ! Determine Prandt Number
  951. !
  952. ! Pr-number is assumed to be constant throughout the CBL, i.e. in surface and mixed layer
  953. ! NOTE: this is different from the original formulation
  954. ! zwstr not used if wheatv<0.
  955. !
  956. zwstr = (abs(wheatv(i,j))*grav*pblh(i,j)/wtseff(i,j))**(1./3.) ! bh
  957. obukov = jq3*wobkl(i,j)-(1-jq3)
  958. term1 = (1.-ssfrac*betah*pblh(i,j)/obukov)**(1./3.)
  959. term2 = sqrt(1.-ssfrac*betah*pblh(i,j)/obukov)
  960. term3 = ccon*fakn*zwstr/(ustr(i,j)*(1.-ssfrac*betah*pblh(i,j)/obukov)**(1./3.))
  961. pr = jq3*(term1/term2+term3)+(1-jq3)
  962. !
  963. ! NOTE that kvh applies to the top of the level
  964. !
  965. kvh(i,j,l) = jq3*max(zpblk/pr,kvhmin)+(1-jq3)*zkvh
  966. ! if in the entrainment zone, override with prescribed entrainment
  967. thvgrad = (zthvk(i,j,l+1)-zthvk(i,j,l))/(whm(i,j,l+1)-whm(i,j,l))
  968. kvhentr = 0.2*wheatv(i,j)/thvgrad
  969. kvh(i,j,l) = jq6*kvhentr + (1-jq6)*kvh(i,j,l)
  970. end do
  971. end do
  972. end do
  973. !PLS $OMP END PARALLEl
  974. ! calculate vertical diffusion
  975. !PLS $OMP PARALLEL &
  976. !PLS $OMP default (none) &
  977. !PLS $OMP shared (jsr, jer, isr, ier, region, dkg, kvh, m, gph) &
  978. !PLS $OMP private (i, j, l, jss, jee)
  979. ! call my_omp_domain (jsr(region), jer(region), jss, jee)
  980. do l=1,lmax_conv-1
  981. do j=j1,j2
  982. do i=i1, i2
  983. dkg(i,j,l)=max(0.,kvh(i,j,l))*2.*(m(i,j,l+1)+m(i,j,l))/(gph(i,j,l+2)-gph(i,j,l))**2
  984. ! CMK dec 2002 ----> gph changed to 1--->lm+1 boundaries
  985. enddo
  986. enddo
  987. enddo
  988. dkg(i1:i2,j1:j2,lmax_conv) = 0.0
  989. !PLS $OMP END PARALLEL
  990. deallocate(wkvf)
  991. deallocate(zdup2)
  992. deallocate(zrinub)
  993. deallocate(wthm)
  994. deallocate(zthvk)
  995. deallocate(pf)
  996. deallocate(whm)
  997. deallocate(ustr)
  998. deallocate(wheat)
  999. deallocate(wqflx)
  1000. deallocate(wobkl)
  1001. deallocate(wheatv)
  1002. deallocate(wtseff)
  1003. deallocate(uwind)
  1004. deallocate(vwind)
  1005. deallocate(kvh)
  1006. nullify(pblh)
  1007. if (okdebug) write(*,*) 'end of bldiff'
  1008. end subroutine bldiff
  1009. #endif /* diffusion */
  1010. #endif /* convection */
  1011. !
  1012. ! define vectorizable 'if' function
  1013. ! if xxxz>yyyz jqif=1 else jqif=0
  1014. !
  1015. INTEGER FUNCTION JQIF( xxxz, yyyz )
  1016. ! --- in/out ---------------------------
  1017. real, intent(in) :: xxxz, yyyz
  1018. ! --- begin ----------------------------
  1019. !
  1020. !jqif = nint( 0.5 - sign( 0.5, -dim(xxxz,yyyz)) )
  1021. !
  1022. jqif = nint( 0.5 + sign( 0.5, xxxz-yyyz ) )
  1023. end function jqif
  1024. !=====================================================================================================
  1025. !=====================================================================================================
  1026. subroutine diffusion_filename(tau, region, fname)
  1027. use datetime, only : tau2date
  1028. use Go, only : pathsep
  1029. use dims, only : region_name, revert
  1030. integer(kind=8), intent(in) :: tau
  1031. integer, intent(in) :: region
  1032. character(len=1000), intent(out) :: fname
  1033. character(len=*), parameter :: rname = mname//'/diffusion_filename'
  1034. integer, dimension(6) :: idater
  1035. if ( revert == 1 ) then
  1036. call tau2date( tau, idater )
  1037. else
  1038. ! read field from 3 hours before
  1039. call tau2date( tau - 3600*3, idater )
  1040. end if
  1041. write(fname, '(6a, i4.4, a1, i2.2, a1, i2.2, a1, "dkg_", i4, 4i2.2, ".nc")') &
  1042. trim(diffusion_dir), pathsep, trim(mlevs), pathsep, region_name(region), pathsep, &
  1043. idater(1), pathsep, idater(2), pathsep, idater(3), pathsep, idater(1:5)
  1044. end subroutine diffusion_filename
  1045. !=====================================================================================================
  1046. !=====================================================================================================
  1047. subroutine read_diffusion( status )
  1048. !
  1049. ! Read vertical diffusion from file and store in conv_dat%dkg
  1050. !
  1051. ! Output:
  1052. ! o status 0 = ok
  1053. ! -1 = file not available for at least one region
  1054. ! 1 = read error
  1055. !
  1056. ! 19 Dec 2012 - P. Le Sager - TM5-MP version, using MDF
  1057. ! --- modules ------------------------------
  1058. use dims, only : itau, revert, region_name, lm, okdebug
  1059. use global_data, only : conv_dat
  1060. use datetime, only : tau2date
  1061. use MDF, only : MDF_OPEN, MDF_NETCDF, MDF_Inq_VarID, MDF_Get_Var, MDF_Close, MDF_READ
  1062. use Partools, only : isRoot, par_broadcast
  1063. use meteodata, only : global_lli
  1064. use tm5_distgrid, only : dgrid, scatter
  1065. ! --- in/out ----------------------------------------------
  1066. integer, intent(out) :: status
  1067. ! --- const -----------------------------------------------
  1068. character(len=*), parameter :: rname = mname//', read_diffusion'
  1069. ! --- local -----------------------------------------------
  1070. character(len=1000) :: fname
  1071. integer :: region, fid, varid, imr, jmr, sz(3)
  1072. integer, dimension(6) :: idater
  1073. logical :: exist
  1074. real, allocatable :: global_arr(:,:,:)
  1075. ! --- begin -----------------------------------------------
  1076. if ( revert == 1 ) then
  1077. call tau2date( itau, idater )
  1078. else
  1079. ! read field from 3 hours before
  1080. call tau2date( itau - 3600*3, idater )
  1081. end if
  1082. do region = 1, nregions
  1083. ! entire region grid size
  1084. imr = global_lli(region)%nlon
  1085. jmr = global_lli(region)%nlat
  1086. ! for gathering
  1087. sz = shape(conv_dat(region)%dkg) ! to get 3rd dimension
  1088. if(isRoot)then
  1089. allocate(global_arr(imr,jmr,sz(3)))
  1090. else
  1091. allocate(global_arr(1,1,1))
  1092. end if
  1093. ! Create file name
  1094. call diffusion_filename(itau,region,fname)
  1095. if (isRoot) then
  1096. inquire( file=fname, exist=exist )
  1097. end if
  1098. call par_broadcast(exist, status)
  1099. if ( .not. exist ) then
  1100. write(gol,'(a,": dkg file """,a,""" not found, will create.")') rname, trim(fname)
  1101. call goPr
  1102. status = -1
  1103. return
  1104. end if
  1105. if (isRoot) then
  1106. ! Open file
  1107. ! This is a bit too verbose, leading to long output files.
  1108. ! write(*,'(a,": reading ", a)') rname, trim(fname)
  1109. call MDF_OPEN(TRIM(fname), MDF_NETCDF, MDF_READ, fid, status )
  1110. IF_NOTOK_RETURN(status=1)
  1111. ! Read dkg
  1112. CALL MDF_Inq_VarID( fid, 'dkg', varid, status )
  1113. IF_ERROR_RETURN(status=1)
  1114. CALL MDF_Get_Var( fid, varid, global_arr, status )
  1115. IF_NOTOK_RETURN(status=1)
  1116. end if
  1117. call scatter( dgrid(region), conv_dat(region)%dkg, global_arr, 0, status)
  1118. IF_NOTOK_RETURN(status=1)
  1119. if(isRoot) then
  1120. ! Read blh
  1121. CALL MDF_Inq_VarID( fid, 'blh', varid, status )
  1122. IF_ERROR_RETURN(status=1)
  1123. CALL MDF_Get_Var( fid, varid, global_arr(:,:,1), status )
  1124. IF_NOTOK_RETURN(status=1)
  1125. end if
  1126. call scatter( dgrid(region), conv_dat(region)%blh, global_arr(:,:,1), 0, status)
  1127. IF_NOTOK_RETURN(status=1)
  1128. ! Close file
  1129. if (isRoot) then
  1130. call MDF_Close( fid, status )
  1131. IF_NOTOK_RETURN(status=1)
  1132. end if
  1133. deallocate(global_arr)
  1134. end do
  1135. status = 0
  1136. end subroutine read_diffusion
  1137. !=====================================================================================================
  1138. !=====================================================================================================
  1139. subroutine write_diffusion( status )
  1140. !
  1141. ! Write vertical diffusion to files (for all regions)
  1142. !
  1143. ! Output:
  1144. ! o status 0 = ok, else not ok
  1145. !
  1146. ! 19 Dec 2012 - P. Le Sager - made a TM5-MP version , that works using MDF
  1147. ! --- modules ------------------------------
  1148. use dims, only : itau, revert, region_name, lm, okdebug
  1149. use global_data, only : conv_dat
  1150. use datetime, only : tau2date,date2tau
  1151. use MDF, only : MDF_CREATE, MDF_NETCDF4, MDF_Def_Dim, MDF_Def_Var, MDF_Close
  1152. use MDF, only : MDF_DOUBLE, MDF_NEW, MDF_Put_Var, MDF_Put_Att, MDF_EndDef, MDF_UNLIMITED
  1153. use Partools, only : isRoot
  1154. use meteodata, only : global_lli
  1155. use tm5_distgrid, only : dgrid, gather
  1156. use GO, only : pathsep
  1157. ! --- in/out ----------------------------------------------
  1158. integer, intent(out) :: status
  1159. ! --- const -----------------------------------------------
  1160. character(len=*), parameter :: rname = mname//'write_diffusion'
  1161. ! --- local -----------------------------------------------
  1162. character(len=1000) :: fname, targetdir
  1163. integer :: region, fid, dkgid, blhid, dimlon, dimlat, dimtime, dimlev, imr, jmr, sz(3)
  1164. integer :: lonid, latid, timeid
  1165. integer, dimension(6) :: idater
  1166. real, allocatable :: global_arr(:,:,:)
  1167. real :: sec_1970
  1168. integer(kind=8) :: itau_1970
  1169. integer :: ipos
  1170. ! --- begin -----------------------------------------------
  1171. if ( revert == 1 ) then
  1172. call tau2date( itau, idater )
  1173. else
  1174. ! read field from 3 hours before
  1175. call tau2date( itau - 3600*3, idater )
  1176. end if
  1177. status=0
  1178. do region = 1, nregions
  1179. ! entire region grid size
  1180. imr = global_lli(region)%nlon
  1181. jmr = global_lli(region)%nlat
  1182. ! for gathering
  1183. sz = shape(conv_dat(region)%dkg) ! to get 3rd dimension
  1184. if(isRoot)then
  1185. allocate(global_arr(imr,jmr,sz(3)))
  1186. else
  1187. allocate(global_arr(1,1,1))
  1188. end if
  1189. ! Create and define
  1190. if(isRoot) then
  1191. ! Create file name
  1192. call diffusion_filename(itau,region,fname)
  1193. ! get the target directory from fname
  1194. ipos = scan(trim(fname), pathsep, back=.true.)
  1195. targetdir = fname(1:ipos-1)
  1196. call system('mkdir -p ' // targetdir)
  1197. if(okdebug) then
  1198. write(gol,'(a,": creating ", a)') rname, trim(fname); call goPr
  1199. end if
  1200. call MDF_CREATE(TRIM(fname), MDF_NETCDF4, MDF_NEW, fid, status )
  1201. IF_NOTOK_RETURN(status=1)
  1202. call MDF_Def_Dim( fid, 'longitude', imr, dimlon, status )
  1203. IF_NOTOK_RETURN(status=1)
  1204. call MDF_Def_Dim( fid, 'latitude', jmr, dimlat, status )
  1205. IF_NOTOK_RETURN(status=1)
  1206. call MDF_Def_Dim( fid, 'time', MDF_UNLIMITED, dimtime, status )
  1207. IF_NOTOK_RETURN(status=1)
  1208. call MDF_Def_Var( fid, 'latitude', MDF_DOUBLE, (/dimlat/), latid, status=status)
  1209. IF_NOTOK_RETURN(status=1)
  1210. call MDF_Def_Var( fid, 'longitude', MDF_DOUBLE, (/dimlon/), lonid, status=status)
  1211. IF_NOTOK_RETURN(status=1)
  1212. call MDF_Put_Att(fid, lonid, "units", values = "degrees_east", status=status)
  1213. IF_NOTOK_RETURN(status=1)
  1214. call MDF_Put_Att(fid, latid, "units", values = "degrees_north", status=status)
  1215. IF_NOTOK_RETURN(status=1)
  1216. call MDF_Def_Var(fid, 'time', MDF_DOUBLE, (/dimtime/), timeid, status)
  1217. IF_NOTOK_RETURN(status=1)
  1218. call MDF_Put_Att(fid, timeid, "units", values="seconds since 1970-01-01 00:00:00", status=status)
  1219. IF_NOTOK_RETURN(status=1)
  1220. call MDF_Put_Att(fid, timeid, "comment", values="time values represent the beginning of the period", status=status)
  1221. IF_NOTOK_RETURN(status=1)
  1222. call MDF_Def_Dim( fid, 'levels', sz(3), dimlev, status )
  1223. IF_NOTOK_RETURN(status=1)
  1224. call MDF_Def_Var( fid, 'dkg', MDF_DOUBLE, (/dimlon,dimlat,dimlev,dimtime/), dkgid, &
  1225. status=status, deflate_level=3)
  1226. IF_NOTOK_RETURN(status=1)
  1227. call MDF_Put_Att(fid, dkgid, "units", values = "ks s-1", status=status)
  1228. IF_NOTOK_RETURN(status=1)
  1229. call MDF_Put_Att(fid, dkgid, "long_name", values="vertical diffusivity coefficient", status=status)
  1230. IF_NOTOK_RETURN(status=1)
  1231. call MDF_Def_Var( fid, 'blh', MDF_DOUBLE, (/dimlon,dimlat,dimtime/), blhid, &
  1232. status=status, deflate_level=3)
  1233. IF_NOTOK_RETURN(status=1)
  1234. call MDF_Put_Att(fid, blhid, "units", values = "m", status=status)
  1235. IF_NOTOK_RETURN(status=1)
  1236. call MDF_Put_Att(fid, blhid, "long_name", values="boundary layer height", status=status)
  1237. IF_NOTOK_RETURN(status=1)
  1238. call MDF_EndDef(fid, status)
  1239. end if
  1240. ! gather and write dkg/blh
  1241. call gather( dgrid(region), conv_dat(region)%dkg, global_arr, 0, status)
  1242. IF_NOTOK_RETURN(status=1)
  1243. if (isRoot) then
  1244. call MDF_Put_Var( fid, latid, global_lli(region)%lat_deg, status)
  1245. IF_NOTOK_RETURN(status=1)
  1246. call MDF_Put_Var( fid, lonid, global_lli(region)%lon_deg, status)
  1247. IF_NOTOK_RETURN(status=1)
  1248. call date2tau([1970, 1, 1, 0, 0, 0], itau_1970)
  1249. sec_1970 = real(itau) - real(itau_1970)
  1250. call MDF_Put_Var(fid, timeid, (/sec_1970/), status)
  1251. IF_NOTOK_RETURN(status=1)
  1252. call MDF_Put_Var( fid, dkgid, global_arr, status)
  1253. IF_NOTOK_RETURN(status=1)
  1254. end if
  1255. call gather( dgrid(region), conv_dat(region)%blh, global_arr(:,:,1), 0, status)
  1256. IF_NOTOK_RETURN(status=1)
  1257. if(isRoot) then
  1258. call MDF_Put_Var( fid, blhid, global_arr(:,:,1), status)
  1259. IF_NOTOK_RETURN(status=1)
  1260. ! Close file
  1261. call MDF_CLOSE( fid, status )
  1262. IF_NOTOK_RETURN(status=1)
  1263. end if
  1264. deallocate(global_arr)
  1265. end do
  1266. status = 0
  1267. end subroutine write_diffusion
  1268. !=====================================================================================================
  1269. !=====================================================================================================
  1270. END MODULE DIFFUSION