emission_co.F90 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756
  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_CO
  14. !
  15. ! !DESCRIPTION: Hold data and methods for CO emissions.
  16. !\\
  17. !\\
  18. ! !INTERFACE:
  19. !
  20. MODULE EMISSION_CO
  21. !
  22. ! !USES:
  23. !
  24. use GO, only : gol, goErr, goPr
  25. use dims, only : nregions, idate, okdebug
  26. use global_types, only : emis_data, d3_data
  27. use emission_read, only : used_providers, has_emis
  28. use tm5_distgrid, only : dgrid, get_distgrid, scatter
  29. use partools, only : isRoot, par_broadcast
  30. implicit none
  31. private
  32. !
  33. ! !PUBLIC MEMBER FUNCTIONS:
  34. !
  35. public :: Emission_CO_Init ! allocate dataset
  36. public :: Emission_CO_Done ! deallocate dataset
  37. public :: Emission_CO_Declare ! read monthly input
  38. public :: Emission_CO_Apply ! distribute & add emissions to tracer array
  39. !
  40. ! !PRIVATE DATA MEMBERS:
  41. !
  42. character(len=*), parameter :: mname = 'emission_co'
  43. type( emis_data ), dimension(:,:), allocatable :: co_emis_2d
  44. type( d3_data ), dimension(:,:), allocatable :: co_emis_3d
  45. logical, allocatable :: has_data_3d(:), has_data_2d(:)
  46. integer :: co_2dsec, co_3dsec
  47. !
  48. ! !REVISION HISTORY:
  49. ! 1 Oct 2010 - Achim Strunk - revamped for AR5
  50. ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4
  51. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  52. ! 29 Jan 2014 - Jason Williams - added yearly specific biogenic
  53. !
  54. ! !REMARKS:
  55. !
  56. !EOP
  57. !------------------------------------------------------------------------
  58. CONTAINS
  59. !--------------------------------------------------------------------------
  60. ! TM5 !
  61. !--------------------------------------------------------------------------
  62. !BOP
  63. !
  64. ! !IROUTINE: EMISSION_CO_INIT
  65. !
  66. ! !DESCRIPTION: Allocate memory to handle the emissions.
  67. !\\
  68. !\\
  69. ! !INTERFACE:
  70. !
  71. SUBROUTINE EMISSION_CO_INIT( status )
  72. !
  73. ! !USES:
  74. !
  75. use dims, only : lm
  76. use emission_read, only : providers_def, numb_providers, ed42_nsect_co
  77. use emission_data, only : LAR5BMB
  78. use emission_read, only : n_ar5_ant_sec, n_ar5_shp_sec, n_ar5_air_sec, n_ar5_bmb_sec
  79. use emission_read, only : ar5_cat_ant, ar5_cat_shp, ar5_cat_air, ar5_cat_bmb
  80. !
  81. ! !OUTPUT PARAMETERS:
  82. !
  83. integer, intent(out) :: status
  84. !
  85. ! !REVISION HISTORY:
  86. ! 1 Oct 2010 - Achim Strunk - v0
  87. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  88. !
  89. !EOP
  90. !------------------------------------------------------------------------
  91. !BOC
  92. character(len=*), parameter :: rname = mname//'/Emission_CO_Init'
  93. integer :: region
  94. integer :: imr, jmr, lmr, lsec, lprov, i1, i2, j1, j2
  95. ! --- begin --------------------------------------
  96. status = 0
  97. if(.not. has_emis) return
  98. ! nb of sectors
  99. co_2dsec = 0
  100. co_3dsec = 0
  101. do lprov = 1, numb_providers
  102. if (count(used_providers.eq.providers_def(lprov)%name)/=0) then
  103. if (trim(providers_def(lprov)%name) .eq. 'AR5') then
  104. ! nb of available sectors in AR5 depends on category
  105. co_2dsec = co_2dsec + n_ar5_ant_sec*count('CO'.eq.ar5_cat_ant) + &
  106. n_ar5_shp_sec*count('CO'.eq.ar5_cat_shp)
  107. if (LAR5BMB) co_2dsec = co_2dsec + n_ar5_bmb_sec*count('CO'.eq.ar5_cat_bmb)
  108. co_3dsec = co_3dsec + n_ar5_air_sec*count('CO'.eq.ar5_cat_air)
  109. ! co_2dsec = co_2dsec + providers_def(lprov)%nsect2d
  110. ! co_3dsec = co_3dsec + count('CO'.eq.ar5_cat_air)
  111. elseif (trim(providers_def(lprov)%name) .eq. 'ED42') then
  112. co_2dsec = co_2dsec + ed42_nsect_co
  113. ! no 3d sectors in EDGAR 4.2
  114. else
  115. co_2dsec = co_2dsec + providers_def(lprov)%nsect2d
  116. co_3dsec = co_3dsec + providers_def(lprov)%nsect3d
  117. endif
  118. endif
  119. enddo
  120. allocate( co_emis_2d( nregions, co_2dsec ) )
  121. allocate( co_emis_3d( nregions, co_3dsec ) )
  122. allocate( has_data_2d(co_2dsec)) ; has_data_2d=.false.
  123. allocate( has_data_3d(co_3dsec)) ; has_data_3d=.false.
  124. ! allocate information arrays (2d and 3d)
  125. do region=1,nregions
  126. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  127. lmr = lm(region)
  128. do lsec=1,co_2dsec
  129. allocate( co_emis_2d(region,lsec)%surf(i1:i2,j1:j2) )
  130. end do
  131. do lsec=1,co_3dsec
  132. allocate( co_emis_3d(region,lsec)%d3(i1:i2,j1:j2,lmr) )
  133. end do
  134. enddo
  135. ! ok
  136. status = 0
  137. END SUBROUTINE EMISSION_CO_INIT
  138. !EOC
  139. !--------------------------------------------------------------------------
  140. ! TM5 !
  141. !--------------------------------------------------------------------------
  142. !BOP
  143. !
  144. ! !IROUTINE: EMISSION_CO_DONE
  145. !
  146. ! !DESCRIPTION: Free memory after handling of the emissions.
  147. !\\
  148. !\\
  149. ! !INTERFACE:
  150. !
  151. SUBROUTINE EMISSION_CO_DONE( status )
  152. !
  153. ! !OUTPUT PARAMETERS:
  154. !
  155. integer, intent(out) :: status
  156. !
  157. ! !REVISION HISTORY:
  158. ! 1 Oct 2010 - Achim Strunk - v0
  159. !
  160. !EOP
  161. !------------------------------------------------------------------------
  162. !BOC
  163. character(len=*), parameter :: rname = mname//'/Emission_CO_Done'
  164. integer :: region, lsec
  165. ! --- begin --------------------------------------
  166. status = 0
  167. if(.not. has_emis) return
  168. do region = 1, nregions
  169. do lsec=1,co_2dsec
  170. deallocate( co_emis_2d(region,lsec)%surf )
  171. end do
  172. do lsec=1,co_3dsec
  173. deallocate( co_emis_3d(region,lsec)%d3 )
  174. end do
  175. end do
  176. deallocate( co_emis_2d )
  177. deallocate( co_emis_3d )
  178. deallocate( has_data_2d, has_data_3d)
  179. ! ok
  180. status = 0
  181. END SUBROUTINE EMISSION_CO_DONE
  182. !EOC
  183. !--------------------------------------------------------------------------
  184. ! TM5 !
  185. !--------------------------------------------------------------------------
  186. !BOP
  187. !
  188. ! !IROUTINE: EMISSION_CO_DECLARE
  189. !
  190. ! !DESCRIPTION: Opens, reads and evaluates input files (per month).
  191. ! Provides emissions on 2d/3d-arrays which are then added
  192. ! to mixing ratios in routine *apply.
  193. !\\
  194. !\\
  195. ! !INTERFACE:
  196. !
  197. SUBROUTINE EMISSION_CO_DECLARE( status )
  198. !
  199. ! !USES:
  200. !
  201. use toolbox, only : coarsen_emission
  202. use dims, only : im, jm, lm, idate, sec_month, nlon360, nlat180, iglbsfc
  203. use chem_param, only : xmco
  204. use emission_data, only : msg_emis, LAR5BMB, LMEGAN
  205. ! ---------------- AR5 - EDGAR 4 - ETC. --------------------
  206. use emission_data, only : emis_input_year_co, emis_input_year_nat
  207. use emission_data, only : emis_input_dir_mac
  208. use emission_data, only : emis_input_dir_megan
  209. use emission_data, only : emis_input_dir_retro
  210. use emission_data, only : emis_input_dir_gfed
  211. use emission_data, only : emis_input_dir_ed4
  212. use emission_read, only : emission_ar5_regrid_aircraft
  213. use emission_read, only : emission_cmip6_ReadSector
  214. use emission_read, only : emission_cmip6bmb_ReadSector
  215. use emission_read, only : emission_ar5_ReadSector
  216. use emission_read, only : emission_macc_ReadSector
  217. use emission_read, only : emission_ed4_ReadSector
  218. use emission_read, only : emission_gfed_ReadSector
  219. use emission_read, only : emission_retro_ReadSector
  220. use emission_read, only : emission_megan_ReadSector
  221. use emission_read, only : sectors_def, numb_sectors
  222. use emission_read, only : ar5_dim_3ddata
  223. use emission_read, only : ed42_co_sectors
  224. !
  225. ! !OUTPUT PARAMETERS:
  226. !
  227. integer, intent(out) :: status
  228. !
  229. ! !REVISION HISTORY:
  230. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  231. ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4
  232. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  233. !
  234. ! !REMARKS:
  235. !
  236. !EOP
  237. !------------------------------------------------------------------------
  238. !BOC
  239. character(len=*), parameter :: rname = mname//'/emission_co_declare'
  240. integer :: region, hasData
  241. integer, parameter :: add_field=0
  242. integer, parameter :: amonth=2
  243. integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2
  244. ! AR5
  245. real,dimension(:,:,:), allocatable :: field3d, field3d2
  246. type(d3_data), dimension(nregions) :: emis3d, work, work3d
  247. type(emis_data) :: wrk2D(nregions)
  248. integer :: seccount2d, seccount3d
  249. ! --- begin -----------------------------------------
  250. status = 0
  251. if(.not. has_emis) return
  252. write(gol,'(" EMISS-INFO ------------- read CO emissions -------------")'); call goPr
  253. do region = 1, nregions
  254. do lsec=1,co_2dsec
  255. co_emis_2d(region,lsec)%surf = 0.0
  256. end do
  257. do lsec=1,co_3dsec
  258. co_emis_3d(region,lsec)%d3 = 0.0
  259. end do
  260. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  261. lmr = lm(region)
  262. allocate( work3d(region)%d3 (i1:i2,j1:j2, ar5_dim_3ddata) ) ; work3d(region)%d3 = 0.0
  263. allocate( emis3d(region)%d3 (i1:i2,j1:j2, lmr ) ) ; emis3d(region)%d3 = 0.0
  264. end do
  265. ! global arrays for coarsening
  266. do region = 1, nregions
  267. if (isRoot)then
  268. allocate(work(region)%d3(im(region),jm(region),ar5_dim_3ddata))
  269. else
  270. allocate(work(region)%d3(1,1,1))
  271. end if
  272. enddo
  273. do region = 1, nregions
  274. wrk2D(region)%surf => work(region)%d3(:,:,1)
  275. end do
  276. ! --------------------------------
  277. ! do a loop over available sectors
  278. ! --------------------------------
  279. ! count 2d and 3d sectors
  280. seccount2d = 0
  281. seccount3d = 0
  282. ! always allocate here 3d data set (for 2d sectors it will be filled in first layer only)
  283. if (isRoot) then
  284. allocate( field3d( nlon360, nlat180, ar5_dim_3ddata ) ) ; field3d = 0.0
  285. else
  286. allocate( field3d( 1, 1, 1 ) )
  287. end if
  288. sec : do lsec = 1, numb_sectors
  289. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  290. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  291. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  292. if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle
  293. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  294. endif
  295. field3d = 0.0
  296. if( sectors_def(lsec)%f3d ) then
  297. seccount3d = seccount3d + 1
  298. else
  299. seccount2d = seccount2d + 1
  300. end if
  301. if (isRoot) then ! READ
  302. select case( trim(sectors_def(lsec)%prov) )
  303. case( 'CMIP6' )
  304. call emission_cmip6_ReadSector( 'CO', emis_input_year_co, idate(2), lsec, field3d, status )
  305. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  306. case( 'CMIP6BMB' )
  307. call emission_cmip6bmb_ReadSector( 'CO', emis_input_year_co, idate(2), lsec, field3d, status )
  308. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  309. case( 'AR5' )
  310. ! Screen out agricultural sector for CO,
  311. ! because it is zero in the RCPs
  312. ! and not present in the historical files.
  313. if (trim(sectors_def(lsec)%name) .ne. 'emiss_agr') then
  314. call emission_ar5_ReadSector( 'CO', emis_input_year_co, idate(2), lsec, field3d, status )
  315. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  316. endif
  317. case( 'MACC' )
  318. ! Screen out 'soil', 'nat', and 'air' since they are not available for CO.
  319. ! If using MEGAN skip also biogenic source.
  320. if ( ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_soil')) .and. &
  321. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_nat') ) .and. &
  322. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_air') ) .and. &
  323. ( .not. (LMEGAN .and. (trim(sectors_def(lsec)%name) .eq. 'emiss_bio'))) ) then
  324. if (trim(sectors_def(lsec)%catname) .eq. 'natural') then
  325. call emission_macc_ReadSector( emis_input_dir_mac, 'CO', emis_input_year_nat, idate(2), &
  326. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  327. IF_NOTOK_RETURN(status=1)
  328. else
  329. call emission_macc_ReadSector( emis_input_dir_mac, 'CO', emis_input_year_co, idate(2), &
  330. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  331. IF_NOTOK_RETURN(status=1)
  332. endif
  333. end if
  334. case( 'ED41' )
  335. select case(trim(sectors_def(lsec)%name))
  336. case ('1A3b_c_e','1A3d_SHIP','1A3d1')
  337. call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', 2005, idate(2), &
  338. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  339. IF_NOTOK_RETURN(status=1)
  340. end select
  341. case( 'ED42' )
  342. ! biomass burning (GFED/RETRO/AR5BMB) and transport (ED41) are excluded through ED42_CO_SECTORS definition
  343. call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', emis_input_year_co, idate(2), &
  344. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  345. IF_NOTOK_RETURN(status=1)
  346. case('GFEDv3')
  347. call emission_gfed_ReadSector( emis_input_dir_gfed, 'co', emis_input_year_co, idate(2), &
  348. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  349. IF_NOTOK_RETURN(status=1)
  350. case('RETRO')
  351. call emission_retro_ReadSector( emis_input_dir_retro, 'CO', emis_input_year_co, idate(2), &
  352. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  353. IF_NOTOK_RETURN(status=1)
  354. case('MEGAN')
  355. call emission_megan_ReadSector( emis_input_dir_megan, 'CO', emis_input_year_nat, idate(2), &
  356. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  357. IF_NOTOK_RETURN(status=1)
  358. case('DUMMY')
  359. case default
  360. write(gol,*) "Error in buidling list of providers USED_PROVIDERS"; call goErr
  361. status=1; TRACEBACK; return
  362. end select
  363. ! nothing found???
  364. if( sum(field3d) < 100.*TINY(1.0) ) then
  365. if (okdebug) then
  366. write(gol,'("EMISS-INFO - no CO emissions found for ",a," ",a," for month ",i2 )') &
  367. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  368. endif
  369. hasData=0
  370. else
  371. if (okdebug) then
  372. write(gol,'("EMISS-INFO - found CO emissions for ",a," ",a," for month ",i2 )') &
  373. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  374. endif
  375. ! scale from kg/s to kg/month
  376. field3d = field3d * sec_month ! kg / month
  377. hasData=1
  378. end if
  379. end if
  380. call Par_broadcast(hasData, status)
  381. IF_NOTOK_RETURN(status=1)
  382. if (hasData == 0) then
  383. cycle sec
  384. else
  385. if ( sectors_def(lsec)%f3d ) then
  386. has_data_3d(seccount3d)=.true.
  387. else
  388. has_data_2d(seccount2d)=.true.
  389. end if
  390. end if
  391. ! Distinguish 2d/3d sectors
  392. if( sectors_def(lsec)%f3d ) then
  393. ! ---------------------------------------
  394. ! 3d data (AIRCRAFT), available for CMIP6
  395. ! ---------------------------------------
  396. if (isRoot) then
  397. ! write some numbers
  398. call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'CO', xmco, sum(field3d) )
  399. ! distribute to work arrays in regions
  400. call Coarsen_Emission( 'CO '//trim(sectors_def(lsec)%name), nlon360, nlat180, ar5_dim_3ddata, &
  401. field3d, work, add_field, status )
  402. IF_NOTOK_RETURN(status=1)
  403. end if
  404. ! scatter, sum up on target array
  405. do region = 1, nregions
  406. call scatter(dgrid(region), work3d(region)%d3, work(region)%d3, 0, status)
  407. IF_NOTOK_RETURN( status=1 )
  408. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, J_STRT=j1)
  409. ! aircraft data: regrid vertically to model layers
  410. call emission_ar5_regrid_aircraft( region, i1, j1, work3d(region)%d3, emis3d(region)%d3, status )
  411. IF_NOTOK_RETURN( status=1 )
  412. co_emis_3d(region,seccount3d)%d3 = co_emis_3d(region,seccount3d)%d3 + emis3d(region)%d3
  413. end do
  414. else ! ar5_sector is 2d
  415. ! ---------------------------
  416. ! 2d data (Anthropogenic, Ships, Biomassburning)
  417. ! ---------------------------
  418. if (isRoot) then ! print total & regrid
  419. call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'CO', xmco, sum(field3d(:,:,1)) )
  420. call coarsen_emission( 'CO '//sectors_def(lsec)%name, &
  421. nlon360, nlat180, field3d(:,:,1), wrk2D, add_field, status )
  422. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  423. end if
  424. do region = 1, nregions
  425. call scatter(dgrid(region), co_emis_2d(region,seccount2d)%surf, work(region)%d3(:,:,1), 0, status)
  426. IF_NOTOK_RETURN(status=1)
  427. end do
  428. end if
  429. end do sec ! sectors
  430. deallocate( field3d )
  431. do region = 1, nregions
  432. if (associated(wrk2D(region)%surf)) nullify(wrk2D(region)%surf)
  433. deallocate( work(region)%d3 )
  434. deallocate( work3d(region)%d3 )
  435. deallocate( emis3d(region)%d3 )
  436. end do
  437. ! check sectors found
  438. if( seccount2d /= co_2dsec ) then
  439. write(gol,'(80("-"))') ; call goPr
  440. write(gol,'("ERROR: 2d sectors do not equal total number:",i4," /= ",i4," !")') seccount2d, co_2dsec ; call goErr
  441. write(gol,'(80("-"))') ; call goPr
  442. status=1; return
  443. end if
  444. if( seccount3d /= co_3dsec ) then
  445. write(gol,'(80("-"))') ; call goPr
  446. write(gol,'("ERROR: 3d sectors do not equal total number:",i4," /= ",i4," !")') seccount3d, co_3dsec ; call goErr
  447. write(gol,'(80("-"))') ; call goPr
  448. status=1; return
  449. end if
  450. status = 0
  451. END SUBROUTINE EMISSION_CO_DECLARE
  452. !EOC
  453. !--------------------------------------------------------------------------
  454. ! TM5 !
  455. !--------------------------------------------------------------------------
  456. !BOP
  457. !
  458. ! !IROUTINE: EMISSION_CO_APPLY
  459. !
  460. ! !DESCRIPTION: Take monthly emissions, and
  461. ! - split them vertically
  462. ! - apply time splitting factors
  463. ! - add them to tracers (add_3d)
  464. !\\
  465. !\\
  466. ! !INTERFACE:
  467. !
  468. SUBROUTINE EMISSION_CO_APPLY( region, status )
  469. !
  470. ! !USES:
  471. !
  472. use dims, only : itaur, nsrce, tref
  473. use dims, only : im, jm, lm
  474. use datetime, only : tau2date
  475. use emission_data, only : emission_vdist_by_sector, LAR5BMB
  476. use emission_data, only : do_add_3d, do_add_3d_cycle, bb_cycle
  477. use emission_data, only : emis_bb_trop_cycle
  478. use chem_param, only : ico, xmco
  479. use emission_read, only : sectors_def, numb_sectors
  480. use emission_read, only : ed42_co_sectors
  481. !
  482. ! !INPUT PARAMETERS:
  483. !
  484. integer, intent(in) :: region
  485. !
  486. ! !OUTPUT PARAMETERS:
  487. !
  488. integer, intent(out) :: status
  489. !
  490. ! !REVISION HISTORY:
  491. ! 1 Oct 2010 - Achim Strunk - AR5
  492. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  493. !
  494. !EOP
  495. !------------------------------------------------------------------------
  496. !BOC
  497. character(len=*), parameter :: rname = mname//'/emission_co_apply'
  498. integer, dimension(6) :: idater
  499. real :: dtime, fraction
  500. integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2
  501. integer :: seccount2d, seccount3d
  502. type(d3_data) :: emis3d
  503. ! --- begin -----------------------------
  504. status = 0
  505. if(.not. has_emis) return
  506. if( okdebug ) then
  507. write(gol,*) 'start of emission_co_apply'; call goPr
  508. end if
  509. call tau2date(itaur(region),idater)
  510. dtime=float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX.
  511. if(okdebug) then
  512. write(gol,*)'emission_co_apply in region ',region,' at date: ',idater, ' with time step:', dtime
  513. call goPr
  514. end if
  515. ! get a working structure for 3d emissions
  516. call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  517. allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  518. ! count 2d and 3d sectors
  519. seccount2d = 0
  520. seccount3d = 0
  521. ! cycle over sectors
  522. do lsec = 1, numb_sectors
  523. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  524. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  525. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  526. if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle
  527. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  528. endif
  529. ! default: no additional splitting
  530. fraction = 1.0
  531. ! ----------------------------------------------------------------------------------------
  532. ! distinguish here between sectors and whether they should have additional splitting
  533. ! if( sectors_def(lsec)%catname == 'biomassburning' ) fraction = fraction * bb_frac etc...
  534. ! ----------------------------------------------------------------------------------------
  535. ! distinguish between 2d/3d sectors
  536. if( sectors_def(lsec)%f3d ) then
  537. seccount3d = seccount3d + 1
  538. if (.not.has_data_3d(seccount3d)) cycle
  539. emis3d%d3 = co_emis_3d(region,seccount3d)%d3
  540. else
  541. seccount2d = seccount2d + 1
  542. if (.not.has_data_2d(seccount2d)) cycle
  543. emis3d%d3 = 0.0
  544. ! vertically distribute according to sector
  545. call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'CO', region, co_emis_2d(region,seccount2d), emis3d, status )
  546. IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3))
  547. end if
  548. ! add dataset according to sector and category
  549. if( emis_bb_trop_cycle .and. trim(sectors_def(lsec)%catname) == "biomassburning" ) then
  550. call do_add_3d_cycle( region, ico, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmco, xmco, status, fraction )
  551. IF_NOTOK_RETURN(status=1)
  552. else
  553. call do_add_3d( region, ico, i1, j1, emis3d%d3, xmco, xmco, status, fraction )
  554. IF_NOTOK_RETURN(status=1)
  555. end if
  556. end do
  557. deallocate( emis3d%d3 )
  558. if(okdebug) then
  559. write(gol,*) 'end of emission_co_apply'; call goPr
  560. end if
  561. ! OK
  562. status = 0
  563. END SUBROUTINE EMISSION_CO_APPLY
  564. !EOC
  565. !NOTYET !--------------------------------------------------------------------------
  566. !NOTYET ! TM5 !
  567. !NOTYET !--------------------------------------------------------------------------
  568. !NOTYET !BOP
  569. !NOTYET !
  570. !NOTYET ! !IROUTINE: DISTRIBUTE_HIGH_LAT_CO
  571. !NOTYET !
  572. !NOTYET ! !DESCRIPTION: Attribute seasonal cycle to the hight latitude CO industrial
  573. !NOTYET ! /fossil fuel CO emissions.
  574. !NOTYET ! This simplified distribution function should be replaced if
  575. !NOTYET ! more detailed info becomes available.
  576. !NOTYET !\\
  577. !NOTYET !\\
  578. !NOTYET ! !INTERFACE:
  579. !NOTYET !
  580. !NOTYET SUBROUTINE DISTRIBUTE_HIGH_LAT_CO( imdim, jmdim, field2d )
  581. !NOTYET !
  582. !NOTYET ! !USES:
  583. !NOTYET !
  584. !NOTYET use dims, only : idate, sec_month, sec_day, sec_year, mlen
  585. !NOTYET !
  586. !NOTYET ! !INPUT PARAMETERS:
  587. !NOTYET !
  588. !NOTYET integer, intent(in) :: jmdim, imdim ! dimension of grid
  589. !NOTYET !
  590. !NOTYET ! !INPUT/OUTPUT PARAMETERS:
  591. !NOTYET !
  592. !NOTYET real, dimension(imdim, jmdim) :: field2d
  593. !NOTYET !
  594. !NOTYET ! !REVISION HISTORY:
  595. !NOTYET ! 24 Jun 2002 - fjd
  596. !NOTYET ! 14 Feb 2013 - Ph. Le Sager - reincorporated for ED4
  597. !NOTYET !
  598. !NOTYET ! !REMARKS:
  599. !NOTYET !
  600. !NOTYET !EOP
  601. !NOTYET !------------------------------------------------------------------------
  602. !NOTYET !BOC
  603. !NOTYET
  604. !NOTYET real :: year2month
  605. !NOTYET real, parameter, dimension(12) :: coseas=(/1.146,1.134,1.081,0.995,0.916,0.920,0.910,0.907,0.934,0.962,1.019,1.076/)
  606. !NOTYET
  607. !NOTYET ! high latitude seasonal CO distribution function
  608. !NOTYET integer :: mo,j
  609. !NOTYET
  610. !NOTYET year2month =sec_month/sec_year !scale factor for yearly total to specific month
  611. !NOTYET mo=idate(2)
  612. !NOTYET
  613. !NOTYET do j=1,jmdim
  614. !NOTYET if (j>jmdim*3/4) then
  615. !NOTYET field2d(:,j)=field2d(:,j)*coseas(mo)/12. ! convert to month
  616. !NOTYET else
  617. !NOTYET field2d(:,j)=field2d(:,j)*year2month
  618. !NOTYET endif
  619. !NOTYET enddo
  620. !NOTYET
  621. !NOTYET END SUBROUTINE DISTRIBUTE_HIGH_LAT_CO
  622. !NOTYET !EOC
  623. END MODULE EMISSION_CO