emission_sox.F90 32 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910
  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_SOX
  14. !
  15. ! !DESCRIPTION: data and methods for SOx emissions.
  16. !
  17. !
  18. ! AR5 (August 2010):
  19. ! - SO2 emissions are used and split into SO2/SO4
  20. ! - in case of M7, additional partitioning is done into given modes
  21. ! E.V.: 2.5 % of SO2 emissions are supposed to be SO4 (AEROCOM) to
  22. ! take into account the sulphate formation in the plumes,
  23. ! as sub-grid effect. For SO2 the remaining 97.5 % are
  24. ! calculated in the do_add routine. For SO4 they are calculated
  25. ! in the main emission_sox routine
  26. !
  27. ! SO4 emissions:
  28. ! industrial: (From Stier et al. 2004)
  29. ! 100% accumulation mode number
  30. ! median radius = 0.075 um sigma = 1.59
  31. ! domestic - transport - biomass burning :
  32. ! 50% Aitken mode number
  33. ! median radius = 0.03 um sigma = 1.59
  34. ! 50% Accumulation mode number
  35. ! median radius = 0.075 um sigma = 1.59
  36. !>>> TvN
  37. ! The size distributions have been changed following
  38. ! the AeroCom recommendations of Dentener et al. (2006),
  39. ! adapted to the M7 modes (Stier et al., ACP, 2005).
  40. ! For emissions from industry, power plants and shipping
  41. ! a number median radius of 0.5 um was recommended,
  42. ! at the boundary between accumulation and coarse mode.
  43. ! Therefore, the emitted mass from these sources
  44. ! is equally divided between these modes,
  45. ! with median radii of 0.075 and 0.75 um, respectively.
  46. ! Emissions from road and off-road transport and the domestic sector
  47. ! are assumed to be in the Aitken mode,
  48. ! and from biomass burning in the accumulation mode.
  49. ! Since only 2.5% of SOx emissions are emitted as particles,
  50. ! these choices have only minor impacts.
  51. ! See also comments in mo_aero.F90.
  52. !<<< TvN
  53. !\\
  54. !\\
  55. ! !INTERFACE:
  56. !
  57. MODULE EMISSION_SOX
  58. !
  59. ! !USES:
  60. !
  61. use GO, only : gol, goErr, goPr
  62. use dims, only : nregions, idate, okdebug
  63. use global_types, only : emis_data, d3_data
  64. use emission_read, only : used_providers, has_emis
  65. use tm5_distgrid, only : dgrid, get_distgrid, scatter
  66. use partools, only : isRoot, par_broadcast
  67. implicit none
  68. private
  69. !
  70. ! !PUBLIC MEMBER FUNCTIONS:
  71. !
  72. public :: Emission_SOx_Init ! allocate dataset
  73. public :: Emission_SOx_Done ! deallocate dataset
  74. public :: Emission_sox_declare ! read monthly input
  75. public :: Emission_sox_apply ! distribute & add emissions to tracer array
  76. !
  77. ! !PRIVATE DATA MEMBERS:
  78. !
  79. character(len=*), parameter :: mname = 'emission_sox'
  80. type( emis_data ), dimension(:,:), allocatable :: so2_emis_2d
  81. type( d3_data ), dimension(:,:), allocatable :: so2_emis_3d
  82. integer :: so2_3dsec, so2_2dsec
  83. logical, allocatable :: has_data_3d(:), has_data_2d(:)
  84. !
  85. ! !REVISION HISTORY:
  86. ! 2005 - Elisabetta Vignati - changed for coupling with M7
  87. ! 1 Oct 2010 - Achim Strunk - AR5
  88. ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4
  89. ! 27 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  90. !
  91. ! !REMARKS:
  92. !
  93. !EOP
  94. !------------------------------------------------------------------------
  95. CONTAINS
  96. !--------------------------------------------------------------------------
  97. ! TM5 !
  98. !--------------------------------------------------------------------------
  99. !BOP
  100. !
  101. ! !IROUTINE: EMISSION_SOX_INIT
  102. !
  103. ! !DESCRIPTION: allocate memory
  104. !\\
  105. !\\
  106. ! !INTERFACE:
  107. !
  108. SUBROUTINE EMISSION_SOX_INIT( status )
  109. !
  110. ! !USES:
  111. !
  112. use dims, only : lm
  113. use emission_read, only : providers_def, numb_providers, ed42_nsect_so2
  114. use emission_data, only : LAR5BMB
  115. use emission_read, only : n_ar5_ant_sec, n_ar5_shp_sec, n_ar5_air_sec, n_ar5_bmb_sec
  116. use emission_read, only : ar5_cat_ant, ar5_cat_shp, ar5_cat_air, ar5_cat_bmb
  117. !
  118. ! !OUTPUT PARAMETERS:
  119. !
  120. integer, intent(out) :: status
  121. !
  122. ! !REVISION HISTORY:
  123. ! 1 Oct 2010 - Achim Strunk - AR5 format
  124. ! 27 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  125. !
  126. !EOP
  127. !------------------------------------------------------------------------
  128. !BOC
  129. character(len=*), parameter :: rname = mname//'/Emission_sox_Init'
  130. integer :: region, lsec
  131. integer :: lmr, lprov, i1, i2, j1, j2
  132. ! --- begin --------------------------------------
  133. status = 0
  134. if(.not. has_emis) return
  135. ! nb of sectors
  136. so2_2dsec = 0
  137. so2_3dsec = 0
  138. do lprov = 1, numb_providers
  139. if (count(used_providers.eq.providers_def(lprov)%name)/=0) then
  140. if (trim(providers_def(lprov)%name) .eq. 'AR5') then
  141. ! nb of available sectors in AR5 depends on category
  142. so2_2dsec = so2_2dsec + n_ar5_ant_sec*count('SO2'.eq.ar5_cat_ant) + &
  143. n_ar5_shp_sec*count('SO2'.eq.ar5_cat_shp)
  144. if (LAR5BMB) so2_2dsec = so2_2dsec + n_ar5_bmb_sec*count('SO2'.eq.ar5_cat_bmb)
  145. so2_3dsec = so2_3dsec + n_ar5_air_sec*count('SO2'.eq.ar5_cat_air)
  146. ! so2_2dsec = so2_2dsec + providers_def(lprov)%nsect2d
  147. ! so2_3dsec = so2_3dsec + count('SO2'.eq.ar5_cat_air)
  148. elseif (trim(providers_def(lprov)%name) .eq. 'ED42') then
  149. so2_2dsec = so2_2dsec + ed42_nsect_so2
  150. ! no 3d sectors in EDGAR 4.2
  151. else
  152. so2_2dsec = so2_2dsec + providers_def(lprov)%nsect2d
  153. so2_3dsec = so2_3dsec + providers_def(lprov)%nsect3d
  154. endif
  155. endif
  156. enddo
  157. allocate( so2_emis_2d( nregions, so2_2dsec ) )
  158. allocate( so2_emis_3d( nregions, so2_3dsec ) )
  159. allocate( has_data_2d(so2_2dsec)) ; has_data_2d=.false.
  160. allocate( has_data_3d(so2_3dsec)) ; has_data_3d=.false.
  161. ! allocate information arrays (2d and 3d)
  162. do region=1,nregions
  163. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  164. lmr = lm(region)
  165. do lsec=1,so2_2dsec
  166. allocate( so2_emis_2d(region,lsec)%surf(i1:i2,j1:j2) )
  167. end do
  168. do lsec=1,so2_3dsec
  169. allocate( so2_emis_3d(region,lsec)%d3(i1:i2,j1:j2,lmr) )
  170. end do
  171. enddo
  172. ! ok
  173. status = 0
  174. END SUBROUTINE EMISSION_SOX_INIT
  175. !EOC
  176. !--------------------------------------------------------------------------
  177. ! TM5 !
  178. !--------------------------------------------------------------------------
  179. !BOP
  180. !
  181. ! !IROUTINE: EMISSION_SOX_DONE
  182. !
  183. ! !DESCRIPTION: Free memory.
  184. !\\
  185. !\\
  186. ! !INTERFACE:
  187. !
  188. SUBROUTINE EMISSION_SOX_DONE( status )
  189. !
  190. ! !OUTPUT PARAMETERS:
  191. !
  192. integer, intent(out) :: status
  193. !
  194. ! !REVISION HISTORY:
  195. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  196. !
  197. !EOP
  198. !------------------------------------------------------------------------
  199. !BOC
  200. character(len=*), parameter :: rname = mname//'/Emission_SOx_Done'
  201. integer :: region, lsec
  202. ! --- begin -----------------------------------------
  203. status = 0
  204. if(.not. has_emis) return
  205. do region = 1, nregions
  206. do lsec=1,so2_2dsec
  207. deallocate( so2_emis_2d(region,lsec)%surf )
  208. end do
  209. do lsec=1,so2_3dsec
  210. deallocate( so2_emis_3d(region,lsec)%d3 )
  211. end do
  212. end do
  213. deallocate( so2_emis_2d, so2_emis_3d )
  214. deallocate( has_data_2d, has_data_3d )
  215. status = 0
  216. END SUBROUTINE EMISSION_SOX_DONE
  217. !EOC
  218. !--------------------------------------------------------------------------
  219. ! TM5 !
  220. !--------------------------------------------------------------------------
  221. !BOP
  222. !
  223. ! !IROUTINE: EMISSION_SOX_DECLARE
  224. !
  225. ! !DESCRIPTION: Opens, reads and evaluates input files (per month).
  226. ! Provides emissions on 2d/3d-arrays which are then added
  227. ! to tracers in routine *apply.
  228. !\\
  229. !\\
  230. ! !INTERFACE:
  231. !
  232. SUBROUTINE EMISSION_SOX_DECLARE( status )
  233. !
  234. ! !USES:
  235. !
  236. use binas, only : pi
  237. use toolbox, only : coarsen_emission
  238. use dims, only : im, jm, lm, idate, sec_month, nlon360, nlat180, iglbsfc
  239. use chem_param, only : xms, xmso2
  240. use emission_data, only : msg_emis, LAR5BMB
  241. #ifdef with_m7
  242. use chem_param, only : mode_ais, mode_acs, mode_cos
  243. use chem_param, only : iso4ais, iso4acs, iso4cos
  244. use chem_param, only : iais_n, iacs_n, icos_n
  245. use chem_param, only : xmso4, xmnumb, sigma_lognormal
  246. use chem_param, only : so4_density, h2so4_factor
  247. use emission_data, only : emis_mass, emis_number
  248. use emission_data, only : rad_so4_ait, rad_so4_acc, rad_so4_coa
  249. use emission_data, only : frac_so4
  250. use emission_data, only : emission_vdist_by_sector
  251. #endif
  252. ! ---------------- AR5 - EDGAR 4 --------------------
  253. use emission_data, only : emis_input_year_sox, emis_input_year_nat
  254. use emission_data, only : emis_input_dir_ed4, emis_input_dir_mac
  255. use emission_data, only : emis_input_dir_gfed, emis_input_dir_retro
  256. use emission_read, only : emission_ar5_regrid_aircraft
  257. use emission_read, only : emission_cmip6_ReadSector
  258. use emission_read, only : emission_cmip6bmb_ReadSector
  259. use emission_read, only : emission_ar5_ReadSector
  260. use emission_read, only : emission_macc_ReadSector
  261. use emission_read, only : emission_ed4_ReadSector
  262. use emission_read, only : emission_gfed_ReadSector
  263. use emission_read, only : emission_retro_ReadSector
  264. use emission_read, only : sectors_def, numb_sectors
  265. use emission_read, only : ar5_dim_3ddata
  266. use emission_read, only : ed42_so2_sectors
  267. !
  268. ! !OUTPUT PARAMETERS:
  269. !
  270. integer, intent(out) :: status
  271. !
  272. ! !REVISION HISTORY:
  273. ! 1 Oct 2010 - Achim Strunk - revamped for AR5
  274. ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4
  275. ! 27 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  276. !
  277. ! !REMARKS:
  278. !
  279. !EOP
  280. !------------------------------------------------------------------------
  281. !BOC
  282. character(len=*), parameter :: rname = mname//'/emission_sox_declare'
  283. integer :: region, hasData
  284. integer, parameter :: add_field=0
  285. integer, parameter :: ayear=1, amonth=2
  286. integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2
  287. type(d3_data), dimension(nregions) :: emis3d, work, work3d
  288. type(emis_data) :: wrk2D(nregions)
  289. #ifdef with_m7
  290. real :: mass2numb_acs, mass2numb_ais, mass2numb_cos
  291. real :: numbscale, sfacacs, sfacais, sfaccos
  292. #endif
  293. ! ---------------------------------------------------------------
  294. real, dimension(:,:,:), allocatable :: field3d, field3d2
  295. integer :: seccount2d, seccount3d
  296. ! --- begin ----------------------------
  297. status = 0
  298. if(.not. has_emis) return
  299. write(gol,'(" EMISS-INFO ------------- read SOx emissions -------------")'); call goPr
  300. do region = 1, nregions
  301. do lsec=1,so2_2dsec
  302. so2_emis_2d(region,lsec)%surf = 0.0
  303. end do
  304. do lsec=1,so2_3dsec
  305. so2_emis_3d(region,lsec)%d3 = 0.0
  306. end do
  307. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  308. lmr = lm(region)
  309. allocate( work3d(region)%d3 (i1:i2,j1:j2, ar5_dim_3ddata) ) ; work3d(region)%d3 = 0.0
  310. allocate( emis3d(region)%d3 (i1:i2,j1:j2, lmr ) ) ; emis3d(region)%d3 = 0.0
  311. end do
  312. ! global arrays for coarsening
  313. do region = 1, nregions
  314. if (isRoot)then
  315. allocate(work(region)%d3(im(region),jm(region),ar5_dim_3ddata))
  316. else
  317. allocate(work(region)%d3(1,1,1))
  318. end if
  319. enddo
  320. do region = 1, nregions
  321. wrk2D(region)%surf => work(region)%d3(:,:,1)
  322. end do
  323. ! --------------------------------
  324. ! do a loop over available sectors
  325. ! --------------------------------
  326. ! count 2d and 3d sectors
  327. seccount2d = 0
  328. seccount3d = 0
  329. ! always allocate here 3d data set (for 2d sectors it will be filled in first layer only)
  330. if (isRoot) then
  331. allocate( field3d( nlon360, nlat180, ar5_dim_3ddata ) ) ; field3d = 0.0
  332. else
  333. allocate( field3d( 1, 1, 1 ) )
  334. end if
  335. sec : do lsec = 1, numb_sectors
  336. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  337. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_so2_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  338. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  339. if (count('SO2'.eq.sectors_def(lsec)%species).eq.0) cycle
  340. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  341. endif
  342. field3d = 0.0
  343. if( sectors_def(lsec)%f3d ) then
  344. seccount3d = seccount3d + 1
  345. else
  346. seccount2d = seccount2d + 1
  347. end if
  348. if (isRoot) then ! READ
  349. select case( trim(sectors_def(lsec)%prov) )
  350. case( 'CMIP6' )
  351. call emission_cmip6_ReadSector( 'SO2', emis_input_year_sox, idate(2), lsec, field3d, status )
  352. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  353. case( 'CMIP6BMB' )
  354. call emission_cmip6bmb_ReadSector( 'SO2', emis_input_year_sox, idate(2), lsec, field3d, status )
  355. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  356. case( 'AR5' )
  357. ! Screen out agricultural and solvent sectors for SO2,
  358. ! because they are zero in the RCPs
  359. ! and not present in the historical files.
  360. if (trim(sectors_def(lsec)%name) .ne. 'emiss_agr' .and. &
  361. trim(sectors_def(lsec)%name) .ne. 'emiss_slv') then
  362. call emission_ar5_ReadSector( 'SO2', emis_input_year_sox, idate(2), lsec, field3d, status )
  363. IF_NOTOK_RETURN(status=1)
  364. endif
  365. case( 'MACC' )
  366. ! screen out sectors w/o SO2 (soil, bio, oc)
  367. if ( ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_soil')) .and. &
  368. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_bio') ) .and. &
  369. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_air') ) .and. &
  370. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_oc') ) ) then
  371. if (trim(sectors_def(lsec)%catname) .eq. 'natural') then
  372. call emission_macc_ReadSector( emis_input_dir_mac, 'SO2', emis_input_year_nat, idate(2), &
  373. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  374. IF_NOTOK_RETURN(status=1)
  375. else
  376. call emission_macc_ReadSector( emis_input_dir_mac, 'SO2', emis_input_year_sox, idate(2), &
  377. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  378. IF_NOTOK_RETURN(status=1)
  379. endif
  380. endif
  381. case( 'ED41' )
  382. select case(trim(sectors_def(lsec)%name))
  383. case ('1A3b_c_e','1A3d_SHIP','1A3d1')
  384. call emission_ed4_ReadSector( emis_input_dir_ed4, 'SO2', 'so2', 2005, idate(2), &
  385. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  386. IF_NOTOK_RETURN(status=1)
  387. field3d = field3d * xmso2/xms
  388. end select
  389. case( 'ED42' )
  390. ! biomass burning (GFED/RETRO/AR5BMB) and transport (ED41) are excluded through ED42_SO2_SECTORS definition
  391. call emission_ed4_ReadSector( emis_input_dir_ed4, 'SO2', 'so2', emis_input_year_sox, idate(2), &
  392. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  393. IF_NOTOK_RETURN(status=1)
  394. field3d = field3d * xmso2/xms
  395. case('GFEDv3')
  396. call emission_gfed_ReadSector( emis_input_dir_gfed, 'so2', emis_input_year_sox, idate(2), &
  397. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  398. IF_NOTOK_RETURN(status=1)
  399. case('RETRO')
  400. call emission_retro_ReadSector( emis_input_dir_retro, 'SO2', emis_input_year_sox, idate(2), &
  401. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  402. IF_NOTOK_RETURN(status=1)
  403. case('MEGAN')
  404. !
  405. ! No biogenic emission of SO2 in the MEGAN inventory
  406. !
  407. case('DUMMY')
  408. case default
  409. write(gol,*) "Error in buidling list of providers USED_PROVIDERS"; call goErr
  410. status=1; TRACEBACK; return
  411. end select
  412. ! nothing found???
  413. if( sum(field3d) < 100.*TINY(1.0) ) then
  414. if (okdebug) then
  415. write(gol,'("EMISS-INFO - no SO2 emissions found for ",a," ",a," for month ",i2 )') &
  416. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  417. endif
  418. hasData=0
  419. else
  420. if (okdebug) then
  421. write(gol,'("EMISS-INFO - found SO2 emissions for ",a," ",a," for month ",i2 )') &
  422. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  423. endif
  424. field3d = field3d * sec_month ! from kg(SO2)/s to kg(SO2)/month
  425. hasData=1
  426. end if
  427. end if
  428. call Par_broadcast(hasData, status)
  429. IF_NOTOK_RETURN(status=1)
  430. if (hasData == 0) then
  431. cycle sec
  432. else
  433. if ( sectors_def(lsec)%f3d ) then
  434. has_data_3d(seccount3d)=.true.
  435. else
  436. has_data_2d(seccount2d)=.true.
  437. end if
  438. end if
  439. ! Distinguish 2d/3d sectors
  440. if( sectors_def(lsec)%f3d ) then
  441. ! ---------------------------------------
  442. ! 3d data (AIRCRAFT), available for CMIP6
  443. ! ---------------------------------------
  444. if (isRoot) then
  445. ! write some numbers
  446. call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'SO2', xmso2, sum(field3d) )
  447. ! distribute to work arrays in regions
  448. call Coarsen_Emission( 'SO2 '//trim(sectors_def(lsec)%name), nlon360, nlat180, ar5_dim_3ddata, &
  449. field3d, work, add_field, status )
  450. IF_NOTOK_RETURN(status=1)
  451. end if
  452. ! scatter, sum up on target array
  453. do region = 1, nregions
  454. call scatter(dgrid(region), work3d(region)%d3, work(region)%d3, 0, status)
  455. IF_NOTOK_RETURN( status=1 )
  456. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, J_STRT=j1)
  457. ! aircraft data: regrid vertically to model layers
  458. call emission_ar5_regrid_aircraft( region, i1, j1, work3d(region)%d3, emis3d(region)%d3, status )
  459. IF_NOTOK_RETURN( status=1 )
  460. so2_emis_3d(region,seccount3d)%d3 = so2_emis_3d(region,seccount3d)%d3 + emis3d(region)%d3
  461. end do
  462. else ! ar5_sector is 2d
  463. ! ---------------------------
  464. ! 2d data (Anthropogenic, Ships, Biomassburning)
  465. ! ---------------------------
  466. if (isRoot) then ! print total & regrid
  467. call msg_emis( amonth, trim(sectors_def(lsec)%prov),sectors_def(lsec)%name, 'SO2', xmso2, sum(field3d(:,:,1)) )
  468. call coarsen_emission( 'SO2 '//sectors_def(lsec)%name, &
  469. nlon360, nlat180, field3d(:,:,1), wrk2D, add_field, status )
  470. IF_NOTOK_RETURN(status=1)
  471. end if
  472. do region = 1, nregions
  473. call scatter(dgrid(region), so2_emis_2d(region,seccount2d)%surf, work(region)%d3(:,:,1), 0, status)
  474. IF_NOTOK_RETURN(status=1)
  475. end do
  476. end if ! 2D/3D
  477. end do sec ! sectors
  478. ! Cleanup
  479. deallocate( field3d )
  480. do region = 1, nregions
  481. if (associated(wrk2D(region)%surf)) nullify(wrk2D(region)%surf)
  482. deallocate( work(region)%d3 )
  483. deallocate( work3d(region)%d3 )
  484. end do
  485. ! check sectors found
  486. if( seccount2d /= so2_2dsec ) then
  487. write(gol,'(80("-"))') ; call goPr
  488. write(gol,'("ERROR: 2d sectors do not equal total number:",i4," /= ",i4," !")') seccount2d, so2_2dsec ; call goErr
  489. write(gol,'(80("-"))') ; call goPr
  490. status=1; return
  491. end if
  492. if( seccount3d /= so2_3dsec ) then
  493. write(gol,'(80("-"))') ; call goPr
  494. write(gol,'("ERROR: 3d sectors do not equal total number:",i4," /= ",i4," !")') seccount3d, so2_3dsec ; call goErr
  495. write(gol,'(80("-"))') ; call goPr
  496. status=1; return
  497. end if
  498. #ifdef with_m7
  499. ! do SO4 partitioning from SO2 emissions
  500. ! --------------------------------
  501. ! do a loop over available sectors
  502. ! --------------------------------
  503. ! count 2d and 3d sectors
  504. seccount2d = 0
  505. seccount3d = 0
  506. do lsec = 1, numb_sectors
  507. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  508. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_so2_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  509. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  510. if (count('SO2'.eq.sectors_def(lsec)%species).eq.0) cycle
  511. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  512. endif
  513. if( sectors_def(lsec)%f3d ) then
  514. seccount3d = seccount3d + 1
  515. if (.not.has_data_3d(seccount3d)) cycle
  516. else
  517. seccount2d = seccount2d + 1
  518. if (.not.has_data_2d(seccount2d)) cycle
  519. end if
  520. !>>> TvN
  521. select case( trim(sectors_def(lsec)%vdisttype) )
  522. !case( 'industry' )
  523. case ( 'industry', 'combenergy' )
  524. sfacais = 0.0 * xmso4 / xmso2 * frac_so4
  525. !sfacacs = 1.0 * xmso4 / xmso2 * frac_so4
  526. sfacacs = 0.5 * xmso4 / xmso2 * frac_so4
  527. sfaccos = 0.5 * xmso4 / xmso2 * frac_so4
  528. case ( 'volcanic' )
  529. sfacais = 0.5 * xmso4 / xmso2 * frac_so4
  530. sfacacs = 0.5 * xmso4 / xmso2 * frac_so4
  531. sfaccos = 0.0 * xmso4 / xmso2 * frac_so4
  532. case default
  533. select case ( trim(sectors_def(lsec)%catname) )
  534. case ('ships')
  535. ! as industry and power plants
  536. sfacais = 0.0 * xmso4 / xmso2 * frac_so4
  537. sfacacs = 0.5 * xmso4 / xmso2 * frac_so4
  538. sfaccos = 0.5 * xmso4 / xmso2 * frac_so4
  539. case ( 'biomassburning')
  540. sfacais = 0.0 * xmso4 / xmso2 * frac_so4
  541. sfacacs = 1.0 * xmso4 / xmso2 * frac_so4
  542. sfaccos = 0.0 * xmso4 / xmso2 * frac_so4
  543. case default
  544. ! e.g. traffic, domestic
  545. sfacais = 1.0 * xmso4 / xmso2 * frac_so4
  546. sfacacs = 0.0 * xmso4 / xmso2 * frac_so4
  547. sfaccos = 0.0 * xmso4 / xmso2 * frac_so4
  548. end select
  549. end select
  550. !<<< TvN
  551. ! mass to number factors for these two modes
  552. !>>> TvN
  553. ! The mass2numb conversion factors have been
  554. ! slightly enhanced to account for the small contribution
  555. ! from the hydrogen ions to the mass,
  556. ! since so4_density is the density of sulfuric acid.
  557. numbscale = rad_so4_ait*EXP(1.5*(LOG(sigma_lognormal(mode_ais)))**2)
  558. mass2numb_ais = h2so4_factor*3./(4.*pi*numbscale**3*so4_density)
  559. numbscale = rad_so4_acc*EXP(1.5*(LOG(sigma_lognormal(mode_acs)))**2)
  560. mass2numb_acs = h2so4_factor*3./(4.*pi*numbscale**3*so4_density)
  561. numbscale = rad_so4_coa*EXP(1.5*(LOG(sigma_lognormal(mode_cos)))**2)
  562. mass2numb_cos = h2so4_factor*3./(4.*pi*numbscale**3*so4_density)
  563. !<<< TvN
  564. do region = 1, nregions
  565. ! allocate helper array
  566. call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  567. lmr = lm(region)
  568. allocate( field3d(i1:i2,j1:j2,lmr) )
  569. ! distinguish 2d/3d sectors
  570. if( sectors_def(lsec)%f3d ) then
  571. field3d = so2_emis_3d(region,seccount3d)%d3
  572. else
  573. emis3d(region)%d3 = 0.0
  574. ! do vertical distribution
  575. call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'SO2', region, so2_emis_2d(region,seccount2d), &
  576. emis3d(region), status )
  577. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  578. field3d = emis3d(region)%d3
  579. end if
  580. ! finally add it to the emis target arrays according to the splitting above
  581. ! AIS
  582. emis_mass (region,mode_ais)%d4(:,:,:,1) = &
  583. emis_mass (region,mode_ais)%d4(:,:,:,1) + field3d * sfacais
  584. ! scale mass to get number
  585. emis_number(region,mode_ais)%d4(:,:,:,1) = &
  586. emis_number(region,mode_ais)%d4(:,:,:,1) + field3d * sfacais * mass2numb_ais
  587. ! ACS
  588. emis_mass (region,mode_acs)%d4(:,:,:,1) = &
  589. emis_mass (region,mode_acs)%d4(:,:,:,1) + field3d * sfacacs
  590. ! scale mass to get number
  591. emis_number(region,mode_acs)%d4(:,:,:,1) = &
  592. emis_number(region,mode_acs)%d4(:,:,:,1) + field3d * sfacacs * mass2numb_acs
  593. !>>> TvN
  594. ! COS
  595. emis_mass (region,mode_cos)%d4(:,:,:,1) = &
  596. emis_mass (region,mode_cos)%d4(:,:,:,1) + field3d * sfaccos
  597. ! scale mass to get number
  598. emis_number(region,mode_cos)%d4(:,:,:,1) = &
  599. emis_number(region,mode_cos)%d4(:,:,:,1) + field3d * sfaccos * mass2numb_cos
  600. !<<< TvN
  601. deallocate( field3d )
  602. end do ! regions
  603. end do ! sectors
  604. do region = 1, nregions
  605. deallocate( emis3d(region)%d3 )
  606. end do
  607. #endif /* ifdef M7 */
  608. END SUBROUTINE EMISSION_SOX_DECLARE
  609. !EOC
  610. !--------------------------------------------------------------------------
  611. ! TM5 !
  612. !--------------------------------------------------------------------------
  613. !BOP
  614. !
  615. ! !IROUTINE: EMISSION_SOX_APPLY
  616. !
  617. ! !DESCRIPTION: Take monthly emissions, and
  618. ! - split them vertically
  619. ! - apply time splitting factors
  620. ! - add them up (add_3d)
  621. !\\
  622. !\\
  623. ! !INTERFACE:
  624. !
  625. SUBROUTINE EMISSION_SOX_APPLY( region, status )
  626. !
  627. ! !USES:
  628. !
  629. use dims, only : itaur, nsrce, tref
  630. use dims, only : im, jm, lm
  631. use datetime, only : tau2date
  632. use emission_data, only : frac_so4
  633. use emission_data, only : emission_vdist_by_sector, LAR5BMB
  634. use emission_data, only : do_add_3d, do_add_3d_cycle, bb_cycle
  635. use emission_data, only : emis_bb_trop_cycle
  636. use emission_read , only : ed42_so2_sectors
  637. use chem_param, only : inh3, xmnh3, xmn
  638. use emission_read, only : sectors_def, numb_sectors
  639. use chem_param, only : xmso2, xms, iso2, iso4
  640. use chem_param, only : xmh2so4
  641. !
  642. ! !INPUT PARAMETERS:
  643. !
  644. integer, intent(in) :: region
  645. !
  646. ! !OUTPUT PARAMETERS:
  647. !
  648. integer, intent(out) :: status
  649. !
  650. ! !REVISION HISTORY:
  651. ! 1 Oct 2010 - Achim Strunk - AR5 version
  652. ! 27 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  653. !
  654. ! !REMARKS:
  655. !
  656. !EOP
  657. !------------------------------------------------------------------------
  658. !BOC
  659. character(len=*), parameter :: rname = mname//'/emission_sox_apply'
  660. integer, dimension(6) :: idater
  661. real :: dtime, fraction, sfrac
  662. integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2
  663. integer :: seccount2d, seccount3d
  664. type(d3_data) :: emis3d
  665. ! --- begin -----------------------------------------
  666. status = 0
  667. if(.not. has_emis) return
  668. if( okdebug ) then
  669. write(gol,*) 'start of emission_sox_apply'; call goPr
  670. end if
  671. call tau2date(itaur(region),idater)
  672. dtime=float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX.
  673. ! get a working structure for 3d emissions
  674. call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  675. allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  676. ! count 2d and 3d sectors
  677. seccount2d = 0
  678. seccount3d = 0
  679. ! cycle over sectors
  680. do lsec = 1, numb_sectors
  681. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  682. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_so2_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  683. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  684. if (count('SO2'.eq.sectors_def(lsec)%species).eq.0) cycle
  685. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  686. endif
  687. ! default: no additional splitting
  688. fraction = 1.0
  689. ! ----------------------------------------------------------------------------------------
  690. ! distinguish here between sectors and whether they should have additional splitting, e.g.,
  691. ! if( sectors_def(lsec)%catname == 'biomassburning' ) fraction = fraction * bb_frac etc...
  692. ! ----------------------------------------------------------------------------------------
  693. ! distinguish between 2d/3d sectors
  694. if( sectors_def(lsec)%f3d ) then
  695. seccount3d = seccount3d + 1
  696. if (.not.has_data_3d(seccount3d)) cycle
  697. emis3d%d3 = so2_emis_3d(region,seccount3d)%d3
  698. else
  699. seccount2d = seccount2d + 1
  700. if (.not.has_data_2d(seccount2d)) cycle
  701. emis3d%d3 = 0.0
  702. ! vertically distribute according to sector
  703. call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'SO2', region, so2_emis_2d(region,seccount2d), emis3d, status)
  704. IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3))
  705. end if
  706. ! SO2: sfrac = fraction * (1.-frac_so4)
  707. sfrac = fraction * (1.0 - frac_so4)
  708. ! add dataset according to sector and category
  709. if( emis_bb_trop_cycle .and. trim(sectors_def(lsec)%catname) == "biomassburning" ) then
  710. call do_add_3d_cycle( region, iso2, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmso2, xmso2, status, sfrac )
  711. IF_NOTOK_RETURN(status=1)
  712. else
  713. call do_add_3d( region, iso2, i1, j1, emis3d%d3, xmso2, xmso2, status, sfrac )
  714. IF_NOTOK_RETURN(status=1)
  715. end if
  716. #ifndef with_m7
  717. ! SO4: sfrac = fraction * frac_so4
  718. sfrac = fraction * frac_so4
  719. ! add dataset according to sector and category
  720. if( emis_bb_trop_cycle .and. trim(sectors_def(lsec)%catname) == "biomassburning" ) then
  721. call do_add_3d_cycle( region, iso4, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmh2so4, xmso2, status, sfrac )
  722. IF_NOTOK_RETURN(status=1)
  723. else
  724. call do_add_3d( region, iso4, i1, j1, emis3d%d3, xmh2so4, xmso2, status, sfrac )
  725. IF_NOTOK_RETURN(status=1)
  726. end if
  727. #endif
  728. end do ! sectors_def
  729. deallocate( emis3d%d3 )
  730. if(okdebug) then
  731. write(gol,*) 'end of emission_sox_apply'; call goPr
  732. endif
  733. ! OK
  734. status = 0
  735. END SUBROUTINE EMISSION_SOX_APPLY
  736. !EOC
  737. END MODULE EMISSION_SOX