emission_data.F90 53 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328
  1. !
  2. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  3. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  4. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  5. !
  6. #include "tm5.inc"
  7. !
  8. !-----------------------------------------------------------------------------
  9. ! TM5 !
  10. !-----------------------------------------------------------------------------
  11. !BOP
  12. !
  13. ! !MODULE: EMISSION_DATA
  14. !
  15. ! !DESCRIPTION: Provides general methods and fields needed for the emission
  16. ! routines.
  17. !
  18. ! emission_vdist_by_sector : distribute vertically emissions according to source sector
  19. ! msg_emis : print formated emission total
  20. ! do_add_2d : add 2D emission to tracer mass array
  21. ! do_add_3d : add 3D emission to tracer mass array
  22. ! do_add_3d_cycle : add 3D emission to tracer mass array with a daily cycle
  23. ! scale_cycle : evaluate a daily cycle
  24. ! distribute_l44 : scale input field with NH (SH) fires distributed over JAS (JFM)
  25. !\\
  26. !\\
  27. ! !INTERFACE:
  28. !
  29. MODULE EMISSION_DATA
  30. !
  31. ! !USES:
  32. !
  33. use GO, only : gol, goPr, goErr
  34. use Dims, only : nregions, lm
  35. use global_types, only : emis_data, d2_data, d23_data, d3_data
  36. #ifdef with_m7
  37. use mo_aero_m7, only : nmod
  38. use mo_aero, only : cmr_ff, cmr_bb, cmr_sk, cmr_sa, cmr_sc, facso2
  39. use mo_aero, only : zbb_wsoc_perc, zbge_wsoc_perc
  40. !use mo_aero, only : zom2oc
  41. use global_types, only : modal_emissions
  42. #endif
  43. #ifdef with_budgets
  44. use budget_global,only : budg_dat, nzon_vg
  45. #endif
  46. IMPLICIT NONE
  47. !
  48. ! !PUBLIC DATA MEMBERS & TYPES:
  49. !
  50. public :: emis_data, d2_data, d23_data
  51. !
  52. type(emis_data), target :: plandr(nregions) ! Land fraction
  53. type(emis_data), target :: emis2D(nregions) ! 2D emiss prior vertical redistribution
  54. type(emis_data), target :: emis_temp(nregions) ! temp array to store emission
  55. #ifdef with_m7
  56. type(modal_emissions), dimension(nregions,nmod), target :: emis_mass
  57. type(modal_emissions), dimension(nregions,nmod), target :: emis_number
  58. #endif
  59. #ifdef with_m7
  60. ! Count median / geometric mean radii of
  61. ! primary carbonaceous and sulfate emissions
  62. ! based on AeroCom-I recommendations (Dentener et al., ACP, 2006).
  63. ! The corresponding values for the M7 modes are given by Stier et al. (ACP, 2005).
  64. ! The same values have also been adopted in GLOMAP (Mann et al., 2010).
  65. !
  66. ! Count median radii for carbonaceous aerosol emissions from Dentener et al.,
  67. ! corresponding to sigma = 1.8:
  68. !real, parameter :: rad_emi_ff = 0.015e-6
  69. !real, parameter :: rad_emi_bb = 0.04e-6
  70. ! The corresponding values for sigma = 1.59 used in M7
  71. ! are cmr_ff and cmr_bb, specified in mo_aero.F90.
  72. ! For comparison, Bond et al. (JGR, 2013) give number median radii
  73. ! between 25 and 40 nm for fresh BC in the urban areas
  74. ! of Tokyo, Nagoya, and Seoul,
  75. ! of 60 nm in plumes associated with wildfires,
  76. ! and about 15 nm from aircraft jet engines.
  77. ! These values are volume-equivalent radii (see their Fig. 4).
  78. !
  79. ! According to the original paper by Schwarz et al. (GRL, 2008)
  80. ! the corresponding geometric standard deviation
  81. ! is sigma = 1.71 for the urban BC
  82. ! and 1.43 for the biomass burning aerosol.
  83. !
  84. ! For BC in biomass burning plumes,
  85. ! Kondo et al (JGR, 2011) estimated
  86. ! number median radii in the range 68-70.5 nm (+- 6-8 nm)
  87. ! and geometric standard deviation between 1.32 and 1.36 (+- 0.01-0.04),
  88. ! for particles thickly coated by organics.
  89. !
  90. ! Janhaell et al. (ACP, 2010) have compiled measurements of
  91. ! particle size in fresh biomass burning smoke from vegetation fires.
  92. ! They mention that particles from biomass burning are dominated
  93. ! by an accumulation mode.
  94. ! They also present a relation between the geometric mean diameter Dg
  95. ! and geometric standard deviation sigma for fresh smoke:
  96. ! Dg (um) = (584 +- 5) - (269 +-1) sigma
  97. ! This gives a geometric mean radius of 78 um for sigma = 1.59,
  98. ! in close agreement with the value used by Stier et al.
  99. !
  100. ! Another way to account for differences in sigma
  101. ! is to modify the number median radius rg
  102. ! such that the number of particles (N)
  103. ! emitted for a certain mass (M) is the same.
  104. ! N is proportional to M/rv^3,
  105. ! where rv is the volume mean radius,
  106. ! which is related to rg by
  107. ! rv = rg * exp(1.5*(ln(sigma))^2).
  108. ! According to Janhall, fresh smoke has an average
  109. ! Dg of 117 +- 13 nm and sigma of 1.7 +- 0.1.
  110. ! At sigma=1.59, this would translate into Dg of about 129 nm,
  111. ! or rg of about 65 nm.
  112. ! For the estimates from Kondo et al. and Schwarz et al.,
  113. ! this would give somewhat smaller rg values
  114. ! of about 58 and 53 nm, resp.
  115. !
  116. ! Particles emitted by grass and savannah fires are generally
  117. ! somewhat smaller than those from wood burning.
  118. ! Janhaell et al. estimate that the mean emission radii
  119. ! for grass and savannah fires, resp., are 12.5 and 10 nm smaller.
  120. ! These differences is not accounted for in the model.
  121. !
  122. ! In a later version of ECHAM-HAM particles the emission radius
  123. ! for biomass burning was reduced to the value for fossil fuel
  124. ! (Zhang et al., ACP, 2012).
  125. ! However, such as a small value seems inconsistent with measurements.
  126. !
  127. ! In the CMIP6 emission data set,
  128. ! the contributions from solid biofuel combustion
  129. ! are included in the 2-D anthropogenic sectors,
  130. ! and provided separately in a supplementary data set.
  131. ! In CMIP6, biofuel is only non-zero
  132. ! for the energy, industry and residential and
  133. ! transportation sectors.
  134. ! Most of the residential emissions are due
  135. ! to solid biofuel combustion.
  136. ! The contribution to the industrial sector
  137. ! is only substantial in developing countries,
  138. ! while the contributions to the energy
  139. ! and transportation sectors are generally very small.
  140. !
  141. ! For inventories other than CMIP6,
  142. ! no distinction between fossil and biofuel emissions
  143. ! is made in the model.
  144. !
  145. ! Winijkul et al. (Atm. Env., 2015) have measured
  146. ! size distributions from energy-related combustion
  147. ! sources for the residential, industrial, power and
  148. ! transportation sectors.
  149. ! They give regional and global estimates of
  150. ! mass median diameters.
  151. ! In their supplementary material they give
  152. ! an overview of results from other studies.
  153. ! These results indicate that the size distribution
  154. ! for both biofuel and fossil fuel combustion sources
  155. ! are strongly dependent on the technique.
  156. ! Count median radii for the residential sector,
  157. ! presented in their Table S1, vary between about
  158. ! 15 nm for modern (improved) wood-fueled stoves,
  159. ! and 25-30 for fireplaces,
  160. ! to 160 nm for (regular) wood-fueled heating stoves,
  161. ! to 270 nm for traditional cookstoves.
  162. ! General fossil fuel (2-d sectors)
  163. ! emitted in the Aitken mode
  164. real, parameter :: rad_emi_ff_sol = cmr_ff
  165. real, parameter :: rad_emi_ff_insol = cmr_ff
  166. ! Energy sector
  167. ! emitted in the Aitken mode.
  168. real, parameter :: rad_emi_ene_sol = cmr_ff
  169. real, parameter :: rad_emi_ene_insol = cmr_ff
  170. ! Industry sector
  171. ! emitted in the Aitken mode.
  172. real, parameter :: rad_emi_ind_sol = cmr_ff
  173. real, parameter :: rad_emi_ind_insol = cmr_ff
  174. ! Transportation sector
  175. ! emitted in the Aitken mode.
  176. real, parameter :: rad_emi_tra_sol = cmr_ff
  177. real, parameter :: rad_emi_tra_insol = cmr_ff
  178. ! Shipping sector
  179. ! emitted in the Aitken mode.
  180. ! Currently set to cmr_ff,
  181. ! but could be changed.
  182. real, parameter :: rad_emi_shp_sol = cmr_ff
  183. real, parameter :: rad_emi_shp_insol = cmr_ff
  184. ! Aircraft sector
  185. ! emitted in the Aitken mode.
  186. ! Currently set to cmr_ff,
  187. ! but a smaller value could be used,
  188. ! e.g. 10 or 15 nm.
  189. real, parameter :: rad_emi_air_sol = cmr_ff
  190. real, parameter :: rad_emi_air_insol = cmr_ff
  191. ! Open biomass burning
  192. ! Soluble part emitted in the accumulation mode,
  193. ! insoluble part in the Aitken mode.
  194. ! Stier et al. apply cmr_ff to the insoluble part,
  195. ! but a slightly larger value seems more realistic,
  196. ! e.g. 40 nm.
  197. real, parameter :: rad_emi_bb_sol = cmr_bb
  198. real, parameter :: rad_emi_bb_insol = cmr_ff
  199. ! Solid biofuel combustion
  200. ! Currently treated as biomass burning,
  201. ! but could be changed (see comments above).
  202. ! The value cmr_bb = 75 nm corresponds
  203. ! with the value of 50 nm at sigma=2
  204. ! assumed by Kondros et al. (ACP, 2015)
  205. ! for emissions from biofuel combusion
  206. ! in their baseline run.
  207. ! In one of their sensitivity runs,
  208. ! they increase it by a factor of 2.
  209. real, parameter :: rad_emi_bf_sol = cmr_bb
  210. real, parameter :: rad_emi_bf_insol = cmr_ff
  211. ! not used anymore:
  212. !real, parameter :: rad_soa = 0.01e-6 ! soa average radius
  213. ! assuming 3nm particle formation and growth to
  214. ! that size in half an hour
  215. ! Count median radii for sulfate aerosol emissions adapted to the M7 modes:
  216. real, parameter :: rad_so4_ait = cmr_sk ! aitken mode radius
  217. real, parameter :: rad_so4_acc = cmr_sa ! accumulation mode radius
  218. real, parameter :: rad_so4_coa = cmr_sc ! coarse mode radius
  219. ! Count median dry radii for sea salt emissions.
  220. ! These values have been updated
  221. ! following Vignati et al. (Atmos. Environ., 2010)
  222. ! For further explanations, see emission_ss.F90.
  223. !real, parameter :: radius_ssa = 0.0794e-6
  224. !real, parameter :: radius_ssc = 0.63e-6
  225. real, parameter :: radius_ssa = 0.09e-6 ! accumulation mode
  226. real, parameter :: radius_ssc = 0.794e-6 ! coarse mode
  227. ! Soluble fraction of POM mass
  228. ! According to Janhall et al. (ACP, 2010)
  229. ! 40 to 80% of the organic matter from
  230. ! vegetation fires is water soluble.
  231. ! Kondros et al. (ACP, 2015) use 80%
  232. ! for POM from biofuel combustion
  233. ! in their base run.
  234. ! For the moment, emisions from
  235. ! biofuel combustion and open biomass burning
  236. ! are treated in the same way.
  237. !
  238. ! open biomass burning
  239. real, parameter :: frac_pom_sol_bb = zbb_wsoc_perc
  240. ! solid biofuel combustion
  241. real, parameter :: frac_pom_sol_bf = zbb_wsoc_perc
  242. !real, parameter :: frac_pom_sol_bf = 0.8 ! alternative value
  243. ! fossil fuel combustion
  244. !real, parameter :: frac_pom_sol_ff = 0.65 ! original value
  245. real, parameter :: frac_pom_sol_ff = 0.0 ! value since 2015 revision
  246. ! Soluble fraction of BC mass
  247. ! In the original code,
  248. ! fresh BC was assumed 100% insoluble,
  249. ! and therefore emitted into the Aitken mode.
  250. ! The new code allows to use a non-zero fraction,
  251. ! to account for non-resolved ageing close to the source.
  252. ! See e.g. Kondros et al.,
  253. ! who use 50% for biofuel combustion
  254. ! in their base run.
  255. ! In the standard EMEP MSC-W model,
  256. ! 20% of the elemental carbon (EC) from
  257. ! anthropogenic sources
  258. ! (representative of fossil fuel combustion)
  259. ! and all EC from open biomass fires
  260. ! is assumed hygroscopic.
  261. !
  262. ! open biomass burning
  263. !real, parameter :: frac_bc_sol_bb = 0.0 ! original value
  264. !real, parameter :: frac_bc_sol_bb = 0.5
  265. real, parameter :: frac_bc_sol_bb = 0.95 ! To reduce the AOD over china and outflow region of
  266. ! Africa the water soluble fraction was increasde to 95%
  267. ! in preparation for CMIP6.
  268. ! solid biofuel combustion
  269. !real, parameter :: frac_bc_sol_bf = 0.0 ! original value
  270. !real, parameter :: frac_bc_sol_bf = 0.5 ! Aerocom
  271. real, parameter :: frac_bc_sol_bf = 0.95 ! TB:
  272. ! To reduce the AOD over china and outflow region of
  273. ! Africa the water soluble fraction was increasde to 95%
  274. ! in preparation for CMIP6.
  275. !
  276. ! Some basiss for the choice can be found here:
  277. ! (e.g. Janhall et al., 2010;
  278. ! https://doi.org/10.5194/acp-10-1427-2010 ; Winijkul et al., 2015;
  279. ! https://doi.org/10.1016/j.atmosenv.2015.02.037; Li et al., 2009;
  280. ! https://pubs.acs.org/doi/abs/10.1021/es803330j).
  281. ! fossil fuel combustion
  282. real, parameter :: frac_bc_sol_ff = 0.0
  283. ! Soluble fraction of surrogate SOA emissions
  284. ! POM from SOA is considered 65% soluble,
  285. ! as recommended by AeroCom.
  286. ! The paper by Kanakidou et al. (ACP, 2004)
  287. ! (mentioned in emission_pom.F90)
  288. ! does not seem to support 100% solubility,
  289. ! as was previously assumed.
  290. !real, parameter :: frac_soa_sol = 1.0 ! original value
  291. real, parameter :: frac_soa_sol = zbge_wsoc_perc ! value since 2015 revision
  292. ! Fraction of SOx mass emitted directly as sulfate
  293. real, parameter :: frac_so4=1.-facso2
  294. #else
  295. real, parameter :: frac_so4=1.-0.975
  296. #endif
  297. #ifdef with_m7
  298. ! The value of 1.4 for the POM to OC mass ratio,
  299. ! set in mo_aero.F90, is an outdated estimate,
  300. ! see e.g. Turpin and Lim (Aerosol Sci. Technol., 2001) and
  301. ! Aiken et al. (Environ. Sci. Technol., 2008).
  302. ! Turpin and Lim estimate a ratio of 1.6 +- 0.2 for urban aerosol,
  303. ! and 2.1 +- 0.2 for aged (nonurban) aerosol;
  304. ! They also note that aerosols heavily impacted by woodsmoke can
  305. ! have an even higher ratio (2.2 to 2.6).
  306. ! According to Reid et al. (ACP, 2005),
  307. ! the POM to OC ratio in fresh biomass burning smoke is very uncertain,
  308. ! somewhere in between 1.4 and ~2.
  309. ! Aiken et al. measure ambient aerosol values
  310. ! between 1.25 for O/C = 0 to 2.44 for O/C = 1.0.
  311. ! For OA from biomass burning, they measure 1.56-1.70,
  312. ! lower than the estimates from Turpin and Lim (~2.0).
  313. ! They find the highest ratios for
  314. ! aged and freshly formed SOA (~2.4 and ~1.9, respectively)
  315. ! and lowest values for primary OA from urban combustion.
  316. ! Based on these studies, we apply different values
  317. ! for emissions from biomass burning versus other emissions,
  318. ! as is also done in some other models
  319. ! (see Table 1 in Tsigaridis et al., ACP, 2014),
  320. ! For SOA (seem emission_pom.F90)
  321. ! we apply a relatively high value valid for aged SOA.
  322. ! This compensates for the lack of SOA formation from isoprene,
  323. ! and improves the agreement with aerosol optical depth (AOD)
  324. ! derived from satellite observations (MODIS).
  325. !
  326. ! We could go even further and apply different values for the
  327. ! water soluble and insoluble fractions (Turpin and Lim).
  328. !
  329. ! It should be acknowledged that the representation of OA
  330. ! with a single tracer is very simplistic.
  331. ! In particular, increase in OA mass due to ageing
  332. ! is not properly accounted for.
  333. !
  334. !real, parameter :: oc2pom = zom2oc !factor for conversion of OC mass to POM
  335. real, parameter :: oc2pom_ff = 1.6 ! fossil fuel
  336. real, parameter :: oc2pom_bf = 1.6 ! solid biofuel combustion
  337. real, parameter :: oc2pom_bb = 1.6 ! open biomass burning
  338. real, parameter :: oc2pom_soa = 2.4 ! SOA
  339. #endif
  340. ! CMIP6 - AR5 - EDGAR4 - RETRO - GFED3 - MACC - LPJ - HYMN - MEGAN
  341. logical :: LCMIP6, LCMIP6BMB, LCMIP6_CH4
  342. logical :: LAR5, LAR5BMB, LEDGAR4, LRETROF, LGFED3, LMACC, LLPJ, LHYMN, LMACCITY, LMEGAN
  343. character(len=256) :: cmip6_ch4_dirname
  344. character(len=256) :: emis_input_dir
  345. character(len=256) :: emis_input_dir_cmip6
  346. character(len=256) :: emis_input_dir_ar5
  347. character(len=256) :: emis_input_dir_megan
  348. character(len=256) :: emis_input_dir_ed4
  349. character(len=256) :: emis_input_dir_mac
  350. integer :: emis_input_year_nat
  351. integer :: emis_input_year_bc
  352. integer :: emis_input_year_oc
  353. integer :: emis_input_year_sox
  354. integer :: emis_input_year_nh3
  355. integer :: emis_input_year_nox
  356. integer :: emis_input_year_co
  357. integer :: emis_input_year_nmvoc
  358. integer :: emis_input_year_ch4
  359. character(len=256) :: emis_input_dir_gfed
  360. character(len=256) :: emis_input_dir_retro
  361. character(len=256) :: emis_input_dir_aerocom
  362. character(len=256) :: emis_input_dir_dms
  363. character(len=256) :: emis_input_dir_rn222
  364. character(len=256) :: emis_input_dir_dust, emis_input_dust
  365. #ifdef with_ch4_emis
  366. character(len=256) :: emis_input_dir_natch4
  367. #endif
  368. logical :: emis_ch4_single
  369. logical :: emis_ch4_fix3d
  370. real :: emis_ch4_fixed_ppb
  371. character(len=256) :: emis_zch4_fname
  372. logical :: ch4_fixyear
  373. integer :: ch4_year
  374. !------------------------------
  375. ! SSP scenario name
  376. !------------------------------
  377. character(len=14) :: ssp_name
  378. ! biomass burning levels:
  379. integer, parameter :: bb_nlev = 5
  380. integer, parameter :: bb_lm = lm(1) ! to be reduced to strip stratosphere!
  381. integer, parameter :: bl_lm = lm(1) ! to be reduced to the BL (around 8)
  382. real, parameter :: bb_hlow (0:bb_nlev) = (/ 0., 100., 500., 1000., 2000., 3000./)
  383. real, parameter :: bb_hhigh(0:bb_nlev) = (/ 100., 500., 1000., 2000., 3000., 6000./)
  384. ! biomass burning levels in tropics (-20 - 20)::
  385. integer, parameter :: bb_nlev_trop = 3
  386. real, parameter :: bb_hlow_trop (0:bb_nlev_trop) = (/ 0., 100., 500., 1000./)
  387. real, parameter :: bb_hhigh_trop(0:bb_nlev_trop) = (/ 100., 500., 1000., 2000./)
  388. ! -----------------------------------------------------------------------------------
  389. ! Biomassburning time splitting factors (now globally, instead of by constituent)
  390. logical :: emis_bb_trop_cycle
  391. type bb_cycle_data
  392. real, pointer :: scalef(:)
  393. end type bb_cycle_data
  394. type(bb_cycle_data), dimension(nregions), target :: bb_cycle
  395. ! -----------------------------------------------------------------------------------
  396. ! Budgets
  397. #ifdef with_budgets
  398. real, dimension(nregions) :: sum_emission
  399. type budemi_data
  400. real,dimension(:,:,:,:),pointer :: emi ! [i, j, nbud_vg, ntracet]
  401. end type budemi_data
  402. type(budemi_data),dimension(nregions),target :: budemi_dat
  403. #endif
  404. logical :: use_tiedkte ! read from rc file, used to scale LiNOx
  405. ! ----------------------------------------------------
  406. ! Vertical splitting : used now for all sectors, for biomassburning and other categories
  407. ! (tables in Dentener, 2006, ACP)
  408. !
  409. ! ATTENTION: changes here have impact on the routine emission_vdist_by_sector (contained)
  410. ! Note: vdist_* variable could be hold private
  411. !
  412. integer, parameter :: vd_class_name_len = 12
  413. ! biomassburning
  414. integer, parameter :: vdist_bb_nlev = 6
  415. real,dimension(vdist_bb_nlev),parameter :: vdist_bb_hlow = (/ 0., 100., 500., 1000., 2000., 3000./)
  416. real,dimension(vdist_bb_nlev),parameter :: vdist_bb_hhigh = (/ 100., 500., 1000., 2000., 3000., 6000./)
  417. ! other sources (related to surface)
  418. integer, parameter :: vdist_nn_nlev = 6
  419. real,dimension(vdist_nn_nlev),parameter :: vdist_nn_hlow = (/ 0., 30. , 100., 300., 600., 1000. /)
  420. real,dimension(vdist_nn_nlev),parameter :: vdist_nn_hhigh = (/ 30., 100., 300., 600., 1000., 2000. /)
  421. !
  422. ! !PRIVATE TYPES:
  423. !
  424. character(len=*), parameter, private :: mname = 'emission_data'
  425. !
  426. ! !REVISION HISTORY:
  427. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  428. ! 28 Jun 2012 - Ph. Le Sager - made everything public by default
  429. ! - adapted for lon-lat MPI domain decomposition
  430. ! !REMARKS:
  431. !
  432. !EOP
  433. !------------------------------------------------------------------------
  434. CONTAINS
  435. !--------------------------------------------------------------------------
  436. ! TM5 !
  437. !--------------------------------------------------------------------------
  438. !BOP
  439. !
  440. ! !IROUTINE: EMISSION_VDIST_BY_SECTOR
  441. !
  442. ! !DESCRIPTION: Return vertically distributed emissions depending on given
  443. ! sector and constituent.
  444. ! Distribution is done as 'compromise' between
  445. ! 1) Dentener et al., 2006, ACP, slightly modified for
  446. ! emissions supposed to enter the surface layer.
  447. ! 2) De Meij et al., ACP, 2006
  448. ! 3) EMEP model (docu, model code, publications)
  449. ! 4) Bieser et al., GMDD, 2010
  450. ! 5) Small updates to biomassburning following
  451. ! Huijnen et al., GMDD, 2010
  452. !\\
  453. !\\
  454. ! !INTERFACE:
  455. !
  456. SUBROUTINE EMISSION_VDIST_BY_SECTOR( splitname, constituent, region, emis2D, emis3D, status )
  457. !
  458. ! !USES:
  459. !
  460. use toolbox, only : distribute_emis2D
  461. !
  462. ! !INPUT PARAMETERS:
  463. !
  464. character(len=vd_class_name_len), intent(in) :: splitname
  465. character(len=*), intent(in) :: constituent
  466. integer, intent(in) :: region
  467. type(emis_data), intent(in) :: emis2D
  468. !
  469. ! !INPUT/OUTPUT PARAMETERS:
  470. !
  471. type(d3_data), intent(inout) :: emis3D
  472. !
  473. ! !OUTPUT PARAMETERS:
  474. !
  475. integer, intent(out) :: status
  476. !
  477. ! !REVISION HISTORY:
  478. ! 1 Oct 2010 - Achim Strunk - v0
  479. !
  480. ! !REMARKS: Splitting depending on constituent is still to code.
  481. !
  482. !EOP
  483. !------------------------------------------------------------------------
  484. !BOC
  485. character(len=*), parameter :: rname = mname//'/emission_vdist_by_sector'
  486. integer :: ilev
  487. character(len=2) :: splittype
  488. integer, parameter :: maxvdist_nlev = max(vdist_bb_nlev,vdist_nn_nlev)
  489. real, dimension(maxvdist_nlev) :: sfglobal, sfintrop, sfintemp
  490. ! --- begin --------------------------------------
  491. ! control implicit settings here in this routine
  492. if( vdist_bb_nlev /= 6 .or. vdist_nn_nlev /= 6 ) then
  493. write(gol, '("ERROR: wrong number of layers: bb_nlev /= 6 or nn_lev /= 6")'); call goErr
  494. status=1; return
  495. end if
  496. ! initialise
  497. sfglobal = 0.0
  498. sfintrop = 0.0
  499. sfintemp = 0.0
  500. ! -------------------------------------------------
  501. ! This here is selection by emission source sectors
  502. ! -------------------------------------------------
  503. select case( trim(splitname) )
  504. case( 'combenergy' )
  505. ! Combustion from energy production and transformation, power-plants
  506. ! no zonal differences, no dependency on species
  507. ! --
  508. ! assumed to be related to EMEP SNAP category [1]
  509. ! AEROCOM equivalent: power-plants
  510. ! --
  511. splittype = 'nn'
  512. sfglobal = (/ 0.0, 0.1, 0.7, 0.2, 0.0, 0.0 /)
  513. case( 'combrescom' )
  514. ! Residential and commercial combustion
  515. ! no zonal differences, no dependency on species
  516. ! --
  517. ! assumed to be related to EMEP SNAP category [2]
  518. ! AEROCOM equivalent: domestic/industry
  519. ! --
  520. splittype = 'nn'
  521. sfglobal = (/ 0.4, 0.4, 0.2, 0.0, 0.0, 0.0 /)
  522. case( 'industry' )
  523. ! Industrial processes and combustion
  524. ! no zonal differences, no dependency on species
  525. ! --
  526. ! assumed to be an aggregate of EMEP SNAP categories [3,4,5]
  527. ! AEROCOM equivalent: industry
  528. ! --
  529. splittype = 'nn'
  530. sfglobal = (/ 0.1, 0.2, 0.6, 0.1, 0.0, 0.0 /)
  531. case( 'waste' )
  532. ! Waste treatment and disposal
  533. ! no zonal differences, no dependency on species
  534. ! --
  535. ! assumed to be related to EMEP SNAP category [9]
  536. ! AEROCOM equivalent: -
  537. ! --
  538. splittype = 'nn'
  539. sfglobal = (/ 0.1, 0.2, 0.4, 0.3, 0.0, 0.0 /)
  540. case( 'nearsurface' )
  541. ! Near surface emissions
  542. ! no zonal differences, no dependency on species
  543. ! --
  544. ! AR5: solvent production and use, agricultural waste burning, &
  545. ! ships, grassland fire, mineral dust (AEROCOM or Tegen-Vignati)
  546. ! EMEP equivalents: SNAP categories [6,8]
  547. ! AEROCOM equivalents: ships, agricultural waste
  548. ! --
  549. splittype = 'nn'
  550. sfglobal = (/ 0.8, 0.2, 0.0, 0.0, 0.0, 0.0 /)
  551. case( 'surface' )
  552. ! Surface emissions
  553. ! no zonal differences, no dependency on species
  554. ! --
  555. ! AR5: agriculture, transport, soil, oceanic, biogenic
  556. ! EMEP equivalents: SNAP categories [7,8,10]
  557. ! AEROCOM equivalents: surface emissions from road/off-road, ships, agricultural waste
  558. ! --
  559. splittype = 'nn'
  560. sfglobal = (/ 1.0, 0.0, 0.0, 0.0, 0.0, 0.0 /)
  561. !!$ case ( 'mindust' )
  562. !!$ ! Mineral dust emissions
  563. !!$ ! no zonal differences
  564. !!$ ! --
  565. !!$ ! AR5: surface dust emissions from either from AEROCOM or Tegen-Vignati
  566. !!$ ! assumed to be slightly higher than near-surface
  567. !!$ splittype = 'nn'
  568. !!$ sfglobal = (/ 0.8, 0.2, 0.0, 0.0, 0.0, 0.0 /)
  569. case( 'volcanic' )
  570. ! Volcanic emissions
  571. ! --
  572. ! AR5: natural SOx emissions (from MACC repository)
  573. ! EMEP equivalents: SNAP categories [11]
  574. ! AEROCOM equivalents: continuous: 2/3 to 1 of height of volcano top
  575. ! explosive: 0.5 to 1.5 km above volcano top
  576. ! Since no data is available about volcano top levels, we assume here, that
  577. ! a volcano top is usually 200-1000 m higher than the grid's surrounding terrain
  578. ! height, thus emissions are distributed from 100 - 2000 m almost homogeneously
  579. ! --
  580. splittype = 'nn'
  581. sfglobal = (/ 0.0, 0.0, 0.1, 0.3, 0.4, 0.2 /)
  582. case( 'forestfire' )
  583. ! forest fires: 6 layer model, different splitting for different regions
  584. ! dont distinguish between boreal Europe and boreal Canada, use Europe globally
  585. ! Dentener et al., modified in tropics (up to 2 km) with respect to Huijnen et al., GMD 2010, \
  586. ! who are citing Labonne et al., 2007 ofr new information based on satellite data.
  587. splittype = 'bb'
  588. sfglobal = (/ 0.1, 0.1, 0.2, 0.2, 0.4, 0.0 /)
  589. sfintemp = (/ 0.2, 0.2, 0.2, 0.4, 0.0, 0.0 /)
  590. sfintrop = (/ 0.2, 0.2, 0.2, 0.4, 0.0, 0.0 /)
  591. case default
  592. write(gol, '("ERROR: wrong specification of emission sector ",a,"in vdist_by_sector called for constituent ",a)') &
  593. trim(splitname), trim(constituent); call goErr
  594. status=1; return
  595. end select
  596. ! --------------------
  597. ! Do the splitting -
  598. ! no need to skip cases of fraction (sfglobal, sfintemp, ..) being zero:
  599. ! distribute_emis2D is optimized for those and will return right away
  600. ! --------------------
  601. select case( splittype )
  602. case( 'nn' )
  603. do ilev = 1, vdist_nn_nlev
  604. call distribute_emis2D( region, emis2D, emis3D, vdist_nn_hlow(ilev), vdist_nn_hhigh(ilev), status, sfglobal(ilev) )
  605. IF_NOTOK_RETURN(status=1)
  606. end do
  607. case( 'bb' )
  608. do ilev = 1, vdist_bb_nlev
  609. ! boreal
  610. call distribute_emis2D( region, emis2D, emis3D, vdist_bb_hlow(ilev), vdist_bb_hhigh(ilev), status, &
  611. sfglobal(ilev), -90.,-61.)
  612. IF_NOTOK_RETURN(status=1)
  613. call distribute_emis2D( region, emis2D, emis3D, vdist_bb_hlow(ilev), vdist_bb_hhigh(ilev), status, &
  614. sfglobal(ilev), 60., 90.)
  615. IF_NOTOK_RETURN(status=1)
  616. ! temperate
  617. call distribute_emis2D( region, emis2D, emis3D, vdist_bb_hlow(ilev), vdist_bb_hhigh(ilev), status, &
  618. sfintemp(ilev), -60.,-31.)
  619. IF_NOTOK_RETURN(status=1)
  620. call distribute_emis2D( region, emis2D, emis3D, vdist_bb_hlow(ilev), vdist_bb_hhigh(ilev), status, &
  621. sfintemp(ilev), 30., 59.)
  622. IF_NOTOK_RETURN(status=1)
  623. ! tropics
  624. call distribute_emis2D( region, emis2D, emis3D, vdist_bb_hlow(ilev), vdist_bb_hhigh(ilev), status, &
  625. sfintrop(ilev), -30., 29.)
  626. IF_NOTOK_RETURN(status=1)
  627. end do
  628. case default
  629. write(gol, '("ERROR: wrong specification of splitting type in vdist_by_sector.")'); call goErr
  630. status=1; return
  631. end select
  632. ! OK
  633. status=0
  634. END SUBROUTINE EMISSION_VDIST_BY_SECTOR
  635. !EOC
  636. !--------------------------------------------------------------------------
  637. ! TM5 !
  638. !--------------------------------------------------------------------------
  639. !BOP
  640. !
  641. ! !IROUTINE: MSG_EMIS
  642. !
  643. ! !DESCRIPTION: Provide output to verify the amount of emissions
  644. ! entering the calculations
  645. !\\
  646. !\\
  647. ! !INTERFACE:
  648. !
  649. SUBROUTINE MSG_EMIS(year_month, provider, sector, msg_comp, msg_mass, summed_emissions)
  650. !
  651. ! !USES:
  652. !
  653. use dims, only: idate
  654. !
  655. ! !INPUT PARAMETERS:
  656. !
  657. integer, intent(in) :: year_month ! 1: year, 2: monthly amount
  658. character(len=*), intent(in) :: provider, sector, msg_comp
  659. real, intent(in) :: msg_mass, summed_emissions
  660. !
  661. ! !REVISION HISTORY:
  662. ! 1 Oct 2010 - Achim Strunk - protex; increase space for output
  663. !
  664. !EOP
  665. !------------------------------------------------------------------------
  666. !BOC
  667. character(len=7),dimension(2), parameter:: ym=(/'Yearly ','Monthly'/)
  668. 1111 format(' EMISS-INFO - ',a7,' Emission from ',a8,' -> ',a18,1x,a6,' : mw=',f5.1,', month=',i2.2,', Sum [kg]=',1x,1pe11.4)
  669. 1112 format(' EMISS-INFO - ',a7,' Emission from ',a8,' -> ',a18,1x,a6,' : mw=',f5.1,', --------, Sum [kg]=',1x,1pe11.4)
  670. if (year_month==1) then
  671. write (gol,1112) ym(year_month), provider, sector, msg_comp, msg_mass, summed_emissions; call goPr
  672. else
  673. write (gol,1111) ym(year_month), provider, sector, msg_comp, msg_mass, idate(2), summed_emissions; call goPr
  674. endif
  675. END SUBROUTINE MSG_EMIS
  676. !EOC
  677. !--------------------------------------------------------------------------
  678. ! TM5 !
  679. !--------------------------------------------------------------------------
  680. !BOP
  681. !
  682. ! !IROUTINE: SCALE_CYCLE
  683. !
  684. ! !DESCRIPTION: Evaluates a daily cycle to be used for the biomass
  685. ! burning emissions. It is based on the isoprene daily
  686. ! cycle, except that the max peak is at 14:00h LT.
  687. ! Calculates factors to be multiplied with an
  688. ! emission field (e.g. co) in order to simulate a
  689. ! diurnal cycle in emissions (caused by solar dependency/
  690. ! human activity), e.g.
  691. ! rm(i,j,1,isop) = rm(i,j,1,isop) + e_isop(i,j)*scalef(ipos,j)*dt
  692. ! with ipos depending on time and longitude.
  693. ! The scalefactor is calculated for -180 longitude and
  694. ! the mean value for ntim timesteps during a day is scaled
  695. ! to 1.
  696. ! The routine should be called only once, since
  697. ! the position relative to the sun is not relevant here.
  698. !\\
  699. !\\
  700. ! !INTERFACE:
  701. !
  702. SUBROUTINE SCALE_CYCLE(ntim, scalef)
  703. !
  704. ! !INPUT PARAMETERS:
  705. !
  706. integer, intent(in) :: ntim
  707. !
  708. ! !OUTPUT PARAMETERS:
  709. !
  710. real, dimension(ntim), intent(out) :: scalef
  711. !
  712. ! !REVISION HISTORY:
  713. ! 1 Oct 2010 - Achim Strunk - added without_diurnal_emission_cycle
  714. ! (to be removed again :-( )
  715. ! 19 Oct 2010 - Narcisa Banda - changed without_diurnal_emission_cycle to with_....
  716. ! 31 Jan 2013 - Ph. Le Sager - get rid of with_diurnal_emission_cycle!
  717. !
  718. ! !REMARKS:
  719. ! - This routine is called only once (see emission.F90), and only if
  720. ! emis_bb_trop_cycle is .true.. The later is set with the
  721. ! input.emis.bb.dailycycle in the rc file (chem.input.rc)
  722. !
  723. !EOP
  724. !------------------------------------------------------------------------
  725. !BOC
  726. real :: pi, piby, obliq, deday, delta, dt, lon, hr, lat, houra, smm, w, sig
  727. integer :: i
  728. ! The calling routine (emission_declare) defines NTIM as follows:
  729. ! dtime=float(ndyn_max)/(2*tref(region)) ! timestep emissions
  730. ! ntim=86400/nint(dtime) ! number of timesteps in 24 hours.
  731. pi = acos(-1.0)
  732. piby = pi/180.0
  733. ! obliq = 23.45 * piby
  734. w = 0.2
  735. sig = 2.0
  736. dt = 24./ntim !timestep in hours
  737. lon = -180.0*piby !calculate for dateline
  738. hr = - 0.5*dt !shift half an interval back
  739. if (hr .gt. 12) hr = hr - 24
  740. ! hr = hr - 2 !shift two hours. This means that Local-day-time maximum will be at 14:00 h LT.
  741. smm = 0.0
  742. do i=1,ntim
  743. hr = hr + dt
  744. houra = lon - pi + hr * (2.*pi/24.)
  745. scalef(i) = w + (1-w)*24.0/(sig*sqrt(2*pi))*exp(-0.5*(((hr-1.5)/sig)**2))
  746. ! scalef(i) = max(1+(cos(houra)),0.0)
  747. smm = smm+scalef(i)/ntim
  748. end do
  749. if ( smm > 1e-5 ) then
  750. scalef(1:ntim) = scalef(1:ntim)/smm
  751. else
  752. scalef(1:ntim) = 1.0
  753. end if
  754. END SUBROUTINE SCALE_CYCLE
  755. !EOC
  756. !--------------------------------------------------------------------------
  757. ! TM5 !
  758. !--------------------------------------------------------------------------
  759. !BOP
  760. !
  761. ! !IROUTINE: DO_ADD_2D
  762. !
  763. ! !DESCRIPTION: General routine to add a 2D emission field (given in kg
  764. ! /month) into tracer mass array. The mass increase is done
  765. ! at level L_LEVEL, and for tracer I_TRACER.
  766. !
  767. ! Update accordingly the budget.
  768. !\\
  769. !\\
  770. ! !INTERFACE:
  771. !
  772. SUBROUTINE DO_ADD_2D( region, i_tracer, l_level, is, js, emis_field, &
  773. mol_mass, mol_mass_e, status, xfrac )
  774. !
  775. ! !USES:
  776. !
  777. use dims, only : CheckShape, tref, nsrce
  778. use dims, only : sec_month, im, jm, okdebug
  779. use global_data, only : mass_dat, region_dat
  780. use chem_param, only : ntracet, names
  781. use tm5_distgrid, only : get_distgrid, dgrid
  782. #ifdef with_budgets
  783. use budget_global, only : budg_dat, nzon_vg
  784. #endif
  785. !
  786. ! !INPUT PARAMETERS:
  787. !
  788. integer, intent(in) :: region
  789. integer, intent(in) :: i_tracer
  790. integer, intent(in) :: l_level
  791. integer, intent(in) :: is, js
  792. real, intent(in) :: emis_field(is:,js:)
  793. real, intent(in) :: mol_mass ! mol mass of tracer field
  794. real, intent(in) :: mol_mass_e ! mol mass of emission field (e.g. NOx emission ar in kgN/month)
  795. real, intent(in), optional :: xfrac ! partial addition
  796. !
  797. ! !OUTPUT PARAMETERS:
  798. !
  799. integer, intent(out) :: status
  800. !
  801. ! !REVISION HISTORY:
  802. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  803. !
  804. !EOP
  805. !------------------------------------------------------------------------
  806. !BOC
  807. character(len=*), parameter :: rname = mname//'/do_add_2d'
  808. integer :: i, j, nzone, nzone_v, i1, i2, j1, j2
  809. real :: x, efrac, dtime, month2dt
  810. real :: mbef, maft, sume ! mass BEFore, AFTer, SUM Emiss (debug)
  811. real, dimension(:,:,:,:), pointer :: rm, rzm
  812. ! --- begin --------------------------------------
  813. call GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  814. call CheckShape( (/i2-i1+1,j2-j1+1/), shape(emis_field), status )
  815. IF_NOTOK_RETURN(status=1)
  816. dtime=float(nsrce)/(2*tref(region)) ! timestep emissions
  817. month2dt=dtime/sec_month ! conversion to emission per timestep
  818. rm => mass_dat(region)%rm
  819. rzm => mass_dat(region)%rzm
  820. if(present(xfrac)) then
  821. efrac = xfrac
  822. else
  823. efrac = 1.0
  824. endif
  825. if (okdebug) then
  826. write (gol, '("--------------------------------------------------------------")') ; call goPr
  827. write (gol, '(" DO ADD 2D debug in region:",i3) ' ) region ; call goPr
  828. write (gol, '(" itracer :",i3,a8)' ) i_tracer, ' :'//names(i_tracer) ; call goPr
  829. write (gol, '(" xfrac :",f6.1)' ) efrac ; call goPr
  830. write (gol, '(" mol_mass :",f6.1 ) ' ) mol_mass ; call goPr
  831. write (gol, '(" mol_mass_e :",f6.1 )' ) mol_mass_e ; call goPr
  832. write (gol, '(" emisfield :",2e12.4 )') minval(emis_field),maxval(emis_field) ; call goPr
  833. sume = sum(emis_field)*month2dt*(mol_mass/mol_mass_e)*efrac
  834. write (gol, '(" sum emmis :",e12.4 )' ) sume ; call goPr
  835. mbef = sum(rm(i1:i2,j1:j2,:,i_tracer))
  836. endif
  837. do j=j1,j2
  838. do i=i1,i2
  839. x = emis_field(i,j)*month2dt*(mol_mass/mol_mass_e)*efrac
  840. rm(i,j,l_level,i_tracer) = rm(i,j,l_level,i_tracer) + x
  841. #ifdef slopes
  842. ! injection of emissions at surface, thus vertical slope more negative
  843. ! (keep concentration at top of layer the same as before emission)
  844. if ( i_tracer <= ntracet ) then
  845. rzm(i,j,l_level,i_tracer) = rzm(i,j,l_level,i_tracer) - x
  846. end if
  847. #endif
  848. #ifdef with_budgets
  849. nzone = budg_dat(region)%nzong(i,j) !global budget
  850. nzone_v = nzon_vg(l_level)
  851. #ifndef without_emission
  852. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) = &
  853. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) + x/mol_mass*1e3
  854. if(i_tracer == 1) sum_emission(region) = sum_emission(region) + x !in kg
  855. #endif
  856. #endif
  857. enddo
  858. enddo
  859. if (okdebug) then
  860. maft = sum(rm(i1:i2,j1:j2,:,i_tracer))
  861. write (gol, '(" Added amount :",e12.4 )' ) maft-mbef ; call goPr
  862. write (gol, '(" Added amount can be different when inner zoom is present!" )' ) ; call goPr
  863. if (maft-mbef-sume > 1e-8*sume) then
  864. write (gol, '(" Warning: error in emission ")' ) ; call goErr
  865. endif
  866. endif
  867. nullify(rm)
  868. nullify(rzm)
  869. status=0
  870. END SUBROUTINE DO_ADD_2D
  871. !EOC
  872. !--------------------------------------------------------------------------
  873. ! TM5 !
  874. !--------------------------------------------------------------------------
  875. !BOP
  876. !
  877. ! !IROUTINE: DO_ADD_3D
  878. !
  879. ! !DESCRIPTION: General routine to add a 3D emission field (given in kg
  880. ! /month) into tracer mass array. The mass increase is done
  881. ! for tracer I_TRACER.
  882. !
  883. ! Update accordingly the budget.
  884. !\\
  885. !\\
  886. ! !INTERFACE:
  887. !
  888. SUBROUTINE DO_ADD_3D( region, i_tracer, is, js, emis_field, mol_mass, mol_mass_e, status, xfrac)
  889. !
  890. ! !USES:
  891. !
  892. use dims, only : CheckShape, nsrce, tref
  893. use dims, only : sec_month, im, jm, lm, okdebug
  894. use global_data, only : mass_dat, region_dat
  895. use chem_param, only : ntracet, names
  896. use tm5_distgrid, only : get_distgrid, dgrid
  897. #ifdef with_budgets
  898. use budget_global, only : budg_dat, nzon_vg
  899. #endif
  900. !
  901. ! !INPUT PARAMETERS:
  902. !
  903. integer, intent(in) :: region ! region #
  904. integer, intent(in) :: i_tracer ! id of tracer to increase
  905. integer, intent(in) :: is, js
  906. real, intent(in) :: emis_field(is:,js:,:)
  907. real, intent(in) :: mol_mass, mol_mass_e
  908. real, intent(in), optional :: xfrac !partial addition
  909. !
  910. ! !OUTPUT PARAMETERS:
  911. !
  912. integer, intent(out) :: status
  913. !
  914. ! !REVISION HISTORY:
  915. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  916. !
  917. !EOP
  918. !------------------------------------------------------------------------
  919. !BOC
  920. character(len=*), parameter :: rname = mname//'/do_add_3d'
  921. integer :: i, j, l, nzone, nzone_v, i1,j1,i2,j2
  922. integer, dimension(3) :: ubound_e
  923. real :: x, efrac, dtime, month2dt
  924. real, dimension(:,:,:,:), pointer :: rm
  925. real :: mbef, maft, sume ! debug
  926. ! --- begin --------------------------------------
  927. call GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  928. status = 0
  929. ! check arguments
  930. ubound_e = ubound(emis_field)
  931. if(ubound_e(1) /= i2) then
  932. write(gol,'(A)') 'do_add_3d: im emission not consistent' ; call goPr
  933. status = 1
  934. endif
  935. if(ubound_e(2) /= j2) then
  936. write(gol,'(A)') 'do_add_3d: jm emission not consistent' ; call goPr
  937. status = 1
  938. endif
  939. if(ubound_e(3) > lm(region)) then
  940. write(gol,'(A)') 'do_add_3d: lm emission not consistent' ; call goPr
  941. status = 1
  942. endif
  943. IF_NOTOK_RETURN(status=1)
  944. rm => mass_dat(region)%rm
  945. dtime=float(nsrce)/(2*tref(region)) ! timestep emissions
  946. month2dt=dtime/sec_month ! conversion to emission per timestep
  947. if(present(xfrac)) then
  948. efrac = xfrac
  949. else
  950. efrac = 1.0
  951. endif
  952. ! if (okdebug) then
  953. ! write (gol, '("--------------------------------------------------------------")') ; call goPr
  954. ! write (gol, '(" DO ADD 3D debug in region:",i3) ' ) region ; call goPr
  955. ! write (gol, '(" itracer :",i3,a8)' ) i_tracer, ' :'//names(i_tracer) ; call goPr
  956. ! write (gol, '(" xfrac :",f6.1)' ) efrac ; call goPr
  957. ! write (gol, '(" mol_mass :",f6.1 ) ' ) mol_mass ; call goPr
  958. ! write (gol, '(" mol_mass_e :",f6.1 )' ) mol_mass_e ; call goPr
  959. ! write (gol, '(" emisfield :",2e12.4 )') minval(emis_field), maxval(emis_field) ; call goPr
  960. ! write (gol, '(" ubound(3) :",i3 )' ) ubound_e(3) ; call goPr
  961. ! sume = sum(emis_field)*month2dt*(mol_mass/mol_mass_e)*efrac
  962. ! write (gol, '(" sume :", e12.4 )') sume ; call goPr
  963. ! mbef = sum(rm(i1:i2,j1:j2,:,i_tracer))
  964. ! endif
  965. do l=1,ubound_e(3)
  966. do j=j1,j2
  967. do i=i1,i2
  968. x = emis_field(i,j,l)*month2dt*(mol_mass/mol_mass_e)*efrac
  969. rm(i,j,l,i_tracer) = rm(i,j,l,i_tracer) + x
  970. #ifdef with_budgets
  971. ! budget___
  972. nzone=budg_dat(region)%nzong(i,j) !global budget
  973. nzone_v=nzon_vg(l)
  974. #ifndef without_emission
  975. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) = &
  976. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) + x/mol_mass*1e3
  977. if(i_tracer ==1) sum_emission(region) = sum_emission(region) + x !in kg
  978. #endif
  979. #endif
  980. enddo
  981. enddo
  982. enddo !levels
  983. ! if (okdebug) then
  984. ! maft = sum(rm(i1:i2,j1:j2,:,i_tracer))
  985. ! write (gol, '(" after-before :", e12.4 )' ) maft-mbef ; call goPr
  986. ! write (gol, '(" END debug output ")' ) ; call goPr
  987. ! endif
  988. nullify(rm)
  989. END SUBROUTINE DO_ADD_3D
  990. !EOC
  991. !--------------------------------------------------------------------------
  992. ! TM5 !
  993. !--------------------------------------------------------------------------
  994. !BOP
  995. !
  996. ! !IROUTINE: DO_ADD_3D_CYCLE
  997. !
  998. ! !DESCRIPTION: Routine to add a 3D emission field (given in kg/month),
  999. ! scaled with a daily cycle. To be used for biomass burning
  1000. ! emissions.
  1001. !\\
  1002. !\\
  1003. ! !INTERFACE:
  1004. !
  1005. SUBROUTINE DO_ADD_3D_CYCLE( region, i_tracer, is, js, emis_field, curr_cycle, mol_mass, mol_mass_e, status, xfrac)
  1006. !
  1007. ! !USES:
  1008. !
  1009. use dims, only : CheckShape, tref, nsrce, itaur
  1010. use dims, only : sec_month, im,jm,lm, okdebug
  1011. use dims, only : ndyn_max,dx, xref, xbeg,dy,yref,ybeg
  1012. use global_data, only : mass_dat, region_dat
  1013. use chem_param, only : ntracet, names
  1014. use tm5_distgrid, only : get_distgrid, dgrid
  1015. #ifdef with_budgets
  1016. use budget_global, only : budg_dat, nzon_vg
  1017. #endif
  1018. use datetime, only : tau2date
  1019. use toolbox, only : dumpfield
  1020. !
  1021. ! !INPUT PARAMETERS:
  1022. !
  1023. integer, intent(in) :: region
  1024. integer, intent(in) :: i_tracer
  1025. integer, intent(in) :: is, js
  1026. real, intent(in) :: emis_field(is:,js:,:)
  1027. real, dimension(:), intent(in) :: curr_cycle
  1028. real, intent(in) :: mol_mass ! first is the mass of tracer field
  1029. real, intent(in) :: mol_mass_e ! mass of emission field (e.g. NOx emission ar in kgN/month)
  1030. real, intent(in), optional :: xfrac ! partial addition
  1031. !
  1032. ! !OUTPUT PARAMETERS:
  1033. !
  1034. integer, intent(out) :: status
  1035. !
  1036. ! !REVISION HISTORY:
  1037. ! 1 Oct 2010 - Achim Strunk - v0 based on do_add_3d
  1038. ! 20 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  1039. !
  1040. ! !REMARKS:
  1041. !
  1042. !EOP
  1043. !------------------------------------------------------------------------
  1044. !BOC
  1045. character(len=*), parameter :: rname = mname//'/do_add_3d_cycle'
  1046. integer :: i, j, l, ipos, nzone, nzone_v
  1047. integer,dimension(6) :: idater
  1048. integer :: sec_in_day, ntim, itim, i1, j1, i2, j2
  1049. integer, dimension(3) :: ubound_e
  1050. real :: x, xlon, xlat, efrac, dtime, month2dt, dtime2
  1051. real,dimension(:,:,:,:),pointer :: rm
  1052. real :: mbef, maft, sume
  1053. ! --- begin --------------------------------------
  1054. call GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  1055. status = 0
  1056. ! check arguments
  1057. ubound_e = ubound(emis_field)
  1058. if(ubound_e(1) /= i2) then
  1059. write(gol,'(A)') 'do_add_3d: im emission not consistent' ; call goPr
  1060. status = 1
  1061. endif
  1062. if(ubound_e(2) /= j2) then
  1063. write(gol,'(A)') 'do_add_3d: jm emission not consistent' ; call goPr
  1064. status = 1
  1065. endif
  1066. if(ubound_e(3) > lm(region)) then
  1067. write(gol,'(A)') 'do_add_3d: lm emission not consistent' ; call goPr
  1068. status = 1
  1069. endif
  1070. IF_NOTOK_RETURN(status=1)
  1071. rm => mass_dat(region)%rm
  1072. dtime = float(nsrce)/(2*tref(region)) !timestep emissions
  1073. month2dt = dtime/sec_month ! conversion to emission per timestep
  1074. dtime2 = float(ndyn_max)/(2*tref(region)) !timestep emissions based on non-CFL interference (CMK 05/2006)
  1075. ntim = 86400/nint(dtime2) !number of timesteps in 24 hours based on non-CFL interference (CMK 05/2006)
  1076. call tau2date(itaur(region),idater)
  1077. sec_in_day = idater(4)*3600 + idater(5)*60 + idater(6) !elapsed seconds this day
  1078. itim = sec_in_day/nint(dtime2)+1 !time interval
  1079. if(present(xfrac)) then
  1080. efrac = xfrac
  1081. else
  1082. efrac = 1.0
  1083. endif
  1084. ! if (okdebug) then
  1085. ! write (gol, '("--------------------------------------------------------------")') ; call goPr
  1086. ! write (gol, '(" DO ADD 3D CYCLE debug in region:",i3) ' ) region ; call goPr
  1087. ! write (gol, '(" itracer :",i3,a8)' ) i_tracer, ' :'//names(i_tracer) ; call goPr
  1088. ! write (gol, '(" xfrac :",f6.1)' ) efrac ; call goPr
  1089. ! write (gol, '(" mol_mass :",f6.1 ) ' ) mol_mass ; call goPr
  1090. ! write (gol, '(" mol_mass_e :", f6.1 )' ) mol_mass_e ; call goPr
  1091. ! write (gol, '(" emisfield :", 2e12.4 )') minval(emis_field),maxval(emis_field) ; call goPr
  1092. ! write (gol, '(" ubound(3) :", i3 )' ) ubound_e(3) ; call goPr
  1093. ! sume = sum(emis_field)*month2dt*(mol_mass/mol_mass_e)*efrac
  1094. ! write (gol, '(" sume :", e12.4 )' ) sume ; call goPr
  1095. ! mbef = sum(rm(i1:i2,j1:j2,:,i_tracer))
  1096. ! write (gol,*) 'lbounds:', is,i1,js,j1 ; call goPr
  1097. ! endif
  1098. do l=1,ubound_e(3)
  1099. do j=j1,j2
  1100. do i=i1,i2
  1101. xlat = ybeg(region) + (j-0.5)*dy/yref(region)
  1102. if (xlat .gt. -20 .and. xlat .lt. 20) then
  1103. ! apply emission cycle in tropics only
  1104. ! itim = 1 and lon = -180 --->position 1
  1105. ! itim = ntim ant lon = -180 --->position ntim, etc.
  1106. ! itim = 1 and lon = 0 ---->position ntim/2
  1107. xlon = xbeg(region) + (i-0.5)*dx/xref(region)
  1108. ipos = 1 + mod(int((xlon+180.)*ntim/360.) + (itim-1),ntim) !position in array depending on time and lon.
  1109. x = emis_field(i,j,l)*month2dt*(mol_mass/mol_mass_e)*efrac*curr_cycle(ipos)
  1110. else
  1111. x = emis_field(i,j,l)*month2dt*(mol_mass/mol_mass_e)*efrac
  1112. endif
  1113. rm(i,j,l,i_tracer) = rm(i,j,l,i_tracer) + x
  1114. #ifdef with_budgets
  1115. ! budget___
  1116. ! if(region==nregions) then !sub budget finest region
  1117. ! nzone=nzon(i,j)
  1118. ! nzone_v=nzon_v(l)
  1119. ! budemi(nzone,nzone_v,i_tracer)=budemi(nzone,nzone_v,i_tracer)+x/mol_mass*1e3 !mole
  1120. ! endif
  1121. nzone=budg_dat(region)%nzong(i,j) !global budget
  1122. nzone_v=nzon_vg(l)
  1123. #ifndef without_emission
  1124. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) = &
  1125. budemi_dat(region)%emi(i,j,nzone_v,i_tracer) + x/mol_mass*1e3
  1126. if(i_tracer ==1) sum_emission(region) = sum_emission(region) + x !in kg
  1127. #endif
  1128. #endif
  1129. enddo
  1130. enddo
  1131. enddo !levels
  1132. ! if (okdebug) then
  1133. ! maft = sum(rm(i1:i2,j1:j2,:,i_tracer))
  1134. ! write (gol, '(" after-before :", e12.4 )' ) maft-mbef ; call goPr
  1135. ! !if ( sume > tiny(sume) ) then
  1136. ! ! This is wrong...
  1137. ! ! if ( (maft-mbef)/sume < 0.5) then
  1138. ! ! write (gol, '(" spurious do_add_3d_cycle.... ")' ) ; call goPr
  1139. ! ! call dumpfield(1, 'spurious.hdf', im(region), jm(region), ubound_e(3), 1, emis_field, names(i_tracer))
  1140. ! ! endif
  1141. ! !endif
  1142. ! write (gol, '(" END debug output ")' ) ; call goPr
  1143. ! endif
  1144. nullify(rm)
  1145. end subroutine do_add_3d_cycle
  1146. !EOC
  1147. !--------------------------------------------------------------------------
  1148. ! TM5 !
  1149. !--------------------------------------------------------------------------
  1150. !BOP
  1151. !
  1152. ! !IROUTINE: DISTRIBUTE_L44
  1153. !
  1154. ! !DESCRIPTION: Scales input field (field2d) such that
  1155. ! 1) NH temperate vegetation fires are distributed over the
  1156. ! months: july-august-september
  1157. ! 2) SH temperate vegetation fires are distributed over the
  1158. ! months: January-February-March
  1159. !\\
  1160. !\\
  1161. ! !INTERFACE:
  1162. !
  1163. SUBROUTINE DISTRIBUTE_L44(month, imdim, jmdim, field2d)
  1164. !
  1165. ! !USES:
  1166. !
  1167. use dims, only : okdebug
  1168. !
  1169. ! !INPUT PARAMETERS:
  1170. !
  1171. integer, intent(in) :: month ! month
  1172. integer, intent(in) :: jmdim, imdim ! dimension of grid
  1173. !
  1174. ! !REVISION HISTORY:
  1175. ! 1 Oct 2010 - Achim Strunk - protex
  1176. !
  1177. !EOP
  1178. !------------------------------------------------------------------------
  1179. !BOC
  1180. real,dimension(imdim,jmdim) :: field2d
  1181. real :: bb_nh,bb_sh ! multiplication factor for yearly temperate fires
  1182. integer :: j,i
  1183. bb_nh=0.
  1184. bb_sh=0.
  1185. if (month==7.or.month==8.or.month==9) bb_nh=3./12.
  1186. if (month==1.or.month==2.or.month==3) bb_sh=3./12.
  1187. do i=1, imdim
  1188. do j=1, jmdim/2
  1189. field2d(i,j)=field2d(i,j)*bb_sh
  1190. enddo
  1191. do j=jmdim/2+1,jmdim
  1192. field2d(i,j)=field2d(i,j)*bb_nh
  1193. enddo
  1194. enddo
  1195. if(okdebug) then
  1196. write(gol,*) 'l44 is distributed' ; call goPr
  1197. endif
  1198. END SUBROUTINE DISTRIBUTE_L44
  1199. !EOC
  1200. END MODULE EMISSION_DATA