emission_co.F90 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733
  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
  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_ar5_ReadSector
  214. use emission_read, only : emission_macc_ReadSector
  215. use emission_read, only : emission_ed4_ReadSector
  216. use emission_read, only : emission_gfed_ReadSector
  217. use emission_read, only : emission_retro_ReadSector
  218. use emission_read, only : emission_megan_ReadSector
  219. use emission_read, only : sectors_def, numb_sectors
  220. use emission_read, only : ar5_dim_3ddata
  221. use emission_read, only : ed42_co_sectors
  222. !
  223. ! !OUTPUT PARAMETERS:
  224. !
  225. integer, intent(out) :: status
  226. !
  227. ! !REVISION HISTORY:
  228. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  229. ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4
  230. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  231. !
  232. ! !REMARKS:
  233. !
  234. !EOP
  235. !------------------------------------------------------------------------
  236. !BOC
  237. character(len=*), parameter :: rname = mname//'/emission_co_declare'
  238. integer :: region, hasData
  239. integer, parameter :: add_field=0
  240. integer, parameter :: amonth=2
  241. integer :: imr, jmr, lmr, lsec
  242. ! AR5
  243. real,dimension(:,:,:), allocatable :: field3d, field3d2
  244. type(d3_data) :: emis3d, work(nregions)
  245. type(emis_data) :: wrk2D(nregions)
  246. integer :: seccount2d, seccount3d
  247. ! --- begin -----------------------------------------
  248. status = 0
  249. if(.not. has_emis) return
  250. write(gol,'(" EMISS-INFO ------------- read CO emissions -------------")'); call goPr
  251. do region = 1, nregions
  252. do lsec=1,co_2dsec
  253. co_emis_2d(region,lsec)%surf = 0.0
  254. end do
  255. do lsec=1,co_3dsec
  256. co_emis_3d(region,lsec)%d3 = 0.0
  257. end do
  258. end do
  259. ! global arrays for coarsening
  260. do region = 1, nregions
  261. if (isRoot)then
  262. allocate(work(region)%d3(im(region),jm(region),lm(region)))
  263. else
  264. allocate(work(region)%d3(1,1,1))
  265. end if
  266. enddo
  267. do region = 1, nregions
  268. wrk2D(region)%surf => work(region)%d3(:,:,1)
  269. end do
  270. ! --------------------------------
  271. ! do a loop over available sectors
  272. ! --------------------------------
  273. ! count 2d and 3d sectors
  274. seccount2d = 0
  275. seccount3d = 0
  276. ! always allocate here 3d data set (for 2d sectors it will be filled in first layer only)
  277. if (isRoot) then
  278. allocate( field3d( nlon360, nlat180, ar5_dim_3ddata ) ) ; field3d = 0.0
  279. else
  280. allocate( field3d( 1, 1, 1 ) )
  281. end if
  282. sec : do lsec = 1, numb_sectors
  283. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  284. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  285. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  286. if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle
  287. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  288. endif
  289. field3d = 0.0
  290. if( sectors_def(lsec)%f3d ) then
  291. seccount3d = seccount3d + 1
  292. else
  293. seccount2d = seccount2d + 1
  294. end if
  295. if (isRoot) then ! READ
  296. select case( trim(sectors_def(lsec)%prov) )
  297. case( 'AR5' )
  298. ! Screen out agricultural sector for CO,
  299. ! because it is zero in the RCPs
  300. ! and not present in the historical files.
  301. if (trim(sectors_def(lsec)%name) .ne. 'emiss_agr') then
  302. call emission_ar5_ReadSector( 'CO', emis_input_year, idate(2), lsec, field3d, status )
  303. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  304. endif
  305. case( 'MACC' )
  306. ! Screen out 'soil', 'nat', and 'air' since they are not available for CO.
  307. ! If using MEGAN skip also biogenic source.
  308. if ( ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_soil')) .and. &
  309. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_nat') ) .and. &
  310. ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_air') ) .and. &
  311. ( .not. (LMEGAN .and. (trim(sectors_def(lsec)%name) .eq. 'emiss_bio'))) ) then
  312. call emission_macc_ReadSector( emis_input_dir_mac, 'CO', emis_input_year, idate(2), &
  313. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  314. IF_NOTOK_RETURN(status=1)
  315. end if
  316. case( 'ED41' )
  317. select case(trim(sectors_def(lsec)%name))
  318. case ('1A3b_c_e','1A3d_SHIP','1A3d1')
  319. call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', 2005, idate(2), &
  320. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  321. IF_NOTOK_RETURN(status=1)
  322. end select
  323. case( 'ED42' )
  324. ! biomass burning (GFED/RETRO/AR5BMB) and transport (ED41) are excluded through ED42_CO_SECTORS definition
  325. call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', emis_input_year, idate(2), &
  326. lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status )
  327. IF_NOTOK_RETURN(status=1)
  328. case('GFEDv3')
  329. call emission_gfed_ReadSector( emis_input_dir_gfed, 'co', emis_input_year, idate(2), &
  330. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  331. IF_NOTOK_RETURN(status=1)
  332. case('RETRO')
  333. call emission_retro_ReadSector( emis_input_dir_retro, 'CO', emis_input_year, idate(2), &
  334. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  335. IF_NOTOK_RETURN(status=1)
  336. case('MEGAN')
  337. call emission_megan_ReadSector( emis_input_dir_megan, 'CO', emis_input_year, idate(2), &
  338. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  339. IF_NOTOK_RETURN(status=1)
  340. case('DUMMY')
  341. case default
  342. write(gol,*) "Error in buidling list of providers USED_PROVIDERS"; call goErr
  343. status=1; TRACEBACK; return
  344. end select
  345. ! nothing found???
  346. if( sum(field3d) < 100.*TINY(1.0) ) then
  347. if (okdebug) then
  348. write(gol,'("EMISS-INFO - no CO emissions found for ",a," ",a," for month ",i2 )') &
  349. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  350. endif
  351. hasData=0
  352. else
  353. if (okdebug) then
  354. write(gol,'("EMISS-INFO - found CO emissions for ",a," ",a," for month ",i2 )') &
  355. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  356. endif
  357. ! scale from kg/s to kg/month
  358. field3d = field3d * sec_month ! kg / month
  359. hasData=1
  360. end if
  361. end if
  362. call Par_broadcast(hasData, status)
  363. IF_NOTOK_RETURN(status=1)
  364. if (hasData == 0) then
  365. cycle sec
  366. else
  367. if ( sectors_def(lsec)%f3d ) then
  368. has_data_3d(seccount3d)=.true.
  369. else
  370. has_data_2d(seccount2d)=.true.
  371. end if
  372. end if
  373. ! Distinguish 2d/3d sectors
  374. if( sectors_def(lsec)%f3d ) then
  375. write(gol,'("EMISS-ERROR - Unexpected 3D data: Uncomment code below ")'); call goErr
  376. status=1; TRACEBACK; return
  377. ! ---------------------------
  378. ! 3d data (AIRCRAFT)
  379. ! ---------------------------
  380. ! if (isRoot) then ! REGRID
  381. ! ! helper array for regridding
  382. ! allocate( field3d2( nlon360,nlat180,lm(1) ) ) ; field3d2 = 0.0
  383. ! ! aircraft data: regrid vertically to model layers
  384. ! call emission_ar5_regrid_aircraft( iglbsfc, field3d, nlon360, nlat180, ar5_dim_3ddata, lm(1), field3d2, status )
  385. ! IF_NOTOK_RETURN(status=1;deallocate(field3d,field3d2))
  386. ! ! write some numbers
  387. ! call msg_emis( amonth, trim(sectors_def(lsec)%prov),sectors_def(lsec)%name, 'CO', xmco, sum(field3d2) )
  388. ! ! distribute to co_emis in regions
  389. ! call Coarsen_Emission( 'CO '//trim(sectors_def(lsec)%name), nlon360, nlat180, lm(1), &
  390. ! field3d2, work, add_field, status )
  391. ! IF_NOTOK_RETURN(status=1;deallocate(field3d,field3d2))
  392. ! deallocate( field3d2 )
  393. ! end if
  394. ! do region = 1, nregions
  395. ! call scatter(dgrid(region), co_emis_3d(region,seccount3d)%d3, work(region)%d3, 0, status)
  396. ! IF_NOTOK_RETURN(status=1)
  397. ! end do
  398. else ! ar5_sector is 2d
  399. ! ---------------------------
  400. ! 2d data (Anthropogenic, Ships, Biomassburning)
  401. ! ---------------------------
  402. if (isRoot) then ! print total & regrid
  403. call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'CO', xmco, sum(field3d(:,:,1)) )
  404. call coarsen_emission( 'CO '//sectors_def(lsec)%name, &
  405. nlon360, nlat180, field3d(:,:,1), wrk2D, add_field, status )
  406. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  407. end if
  408. do region = 1, nregions
  409. call scatter(dgrid(region), co_emis_2d(region,seccount2d)%surf, work(region)%d3(:,:,1), 0, status)
  410. IF_NOTOK_RETURN(status=1)
  411. end do
  412. end if
  413. end do sec ! sectors
  414. deallocate( field3d )
  415. do region = 1, nregions
  416. if (associated(wrk2D(region)%surf)) nullify(wrk2D(region)%surf)
  417. deallocate( work(region)%d3 )
  418. end do
  419. ! check sectors found
  420. if( seccount2d /= co_2dsec ) then
  421. write(gol,'(80("-"))') ; call goPr
  422. write(gol,'("ERROR: 2d sectors do not equal total number:",i4," /= ",i4," !")') seccount2d, co_2dsec ; call goErr
  423. write(gol,'(80("-"))') ; call goPr
  424. status=1; return
  425. end if
  426. if( seccount3d /= co_3dsec ) then
  427. write(gol,'(80("-"))') ; call goPr
  428. write(gol,'("ERROR: 3d sectors do not equal total number:",i4," /= ",i4," !")') seccount3d, co_3dsec ; call goErr
  429. write(gol,'(80("-"))') ; call goPr
  430. status=1; return
  431. end if
  432. status = 0
  433. END SUBROUTINE EMISSION_CO_DECLARE
  434. !EOC
  435. !--------------------------------------------------------------------------
  436. ! TM5 !
  437. !--------------------------------------------------------------------------
  438. !BOP
  439. !
  440. ! !IROUTINE: EMISSION_CO_APPLY
  441. !
  442. ! !DESCRIPTION: Take monthly emissions, and
  443. ! - split them vertically
  444. ! - apply time splitting factors
  445. ! - add them to tracers (add_3d)
  446. !\\
  447. !\\
  448. ! !INTERFACE:
  449. !
  450. SUBROUTINE EMISSION_CO_APPLY( region, status )
  451. !
  452. ! !USES:
  453. !
  454. use dims, only : itaur, nsrce, tref
  455. use dims, only : im, jm, lm
  456. use datetime, only : tau2date
  457. use emission_data, only : emission_vdist_by_sector, LAR5BMB
  458. use emission_data, only : do_add_3d, do_add_3d_cycle, bb_cycle
  459. use emission_data, only : emis_bb_trop_cycle
  460. use chem_param, only : ico, xmco
  461. use emission_read, only : sectors_def, numb_sectors
  462. use emission_read, only : ed42_co_sectors
  463. !
  464. ! !INPUT PARAMETERS:
  465. !
  466. integer, intent(in) :: region
  467. !
  468. ! !OUTPUT PARAMETERS:
  469. !
  470. integer, intent(out) :: status
  471. !
  472. ! !REVISION HISTORY:
  473. ! 1 Oct 2010 - Achim Strunk - AR5
  474. ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  475. !
  476. !EOP
  477. !------------------------------------------------------------------------
  478. !BOC
  479. character(len=*), parameter :: rname = mname//'/emission_co_apply'
  480. integer, dimension(6) :: idater
  481. real :: dtime, fraction
  482. integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2
  483. integer :: seccount2d, seccount3d
  484. type(d3_data) :: emis3d
  485. ! --- begin -----------------------------
  486. status = 0
  487. if(.not. has_emis) return
  488. if( okdebug ) then
  489. write(gol,*) 'start of emission_co_apply'; call goPr
  490. end if
  491. call tau2date(itaur(region),idater)
  492. dtime=float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX.
  493. if(okdebug) then
  494. write(gol,*)'emission_co_apply in region ',region,' at date: ',idater, ' with time step:', dtime
  495. call goPr
  496. end if
  497. ! get a working structure for 3d emissions
  498. call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  499. allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  500. ! count 2d and 3d sectors
  501. seccount2d = 0
  502. seccount3d = 0
  503. ! cycle over sectors
  504. do lsec = 1, numb_sectors
  505. if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle
  506. if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle
  507. if (associated(sectors_def(lsec)%species)) then ! AR5 check
  508. if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle
  509. if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle
  510. endif
  511. ! default: no additional splitting
  512. fraction = 1.0
  513. ! ----------------------------------------------------------------------------------------
  514. ! distinguish here between sectors and whether they should have additional splitting
  515. ! if( sectors_def(lsec)%catname == 'biomassburning' ) fraction = fraction * bb_frac etc...
  516. ! ----------------------------------------------------------------------------------------
  517. ! distinguish between 2d/3d sectors
  518. if( sectors_def(lsec)%f3d ) then
  519. seccount3d = seccount3d + 1
  520. if (.not.has_data_3d(seccount3d)) cycle
  521. emis3d%d3 = co_emis_3d(region,seccount3d)%d3
  522. else
  523. seccount2d = seccount2d + 1
  524. if (.not.has_data_2d(seccount2d)) cycle
  525. emis3d%d3 = 0.0
  526. ! vertically distribute according to sector
  527. call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'CO', region, co_emis_2d(region,seccount2d), emis3d, status )
  528. IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3))
  529. end if
  530. ! add dataset according to sector and category
  531. if( emis_bb_trop_cycle .and. trim(sectors_def(lsec)%catname) == "biomassburning" ) then
  532. call do_add_3d_cycle( region, ico, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmco, xmco, status, fraction )
  533. IF_NOTOK_RETURN(status=1)
  534. else
  535. call do_add_3d( region, ico, i1, j1, emis3d%d3, xmco, xmco, status, fraction )
  536. IF_NOTOK_RETURN(status=1)
  537. end if
  538. end do
  539. deallocate( emis3d%d3 )
  540. if(okdebug) then
  541. write(gol,*) 'end of emission_co_apply'; call goPr
  542. end if
  543. ! OK
  544. status = 0
  545. END SUBROUTINE EMISSION_CO_APPLY
  546. !EOC
  547. !NOTYET !--------------------------------------------------------------------------
  548. !NOTYET ! TM5 !
  549. !NOTYET !--------------------------------------------------------------------------
  550. !NOTYET !BOP
  551. !NOTYET !
  552. !NOTYET ! !IROUTINE: DISTRIBUTE_HIGH_LAT_CO
  553. !NOTYET !
  554. !NOTYET ! !DESCRIPTION: Attribute seasonal cycle to the hight latitude CO industrial
  555. !NOTYET ! /fossil fuel CO emissions.
  556. !NOTYET ! This simplified distribution function should be replaced if
  557. !NOTYET ! more detailed info becomes available.
  558. !NOTYET !\\
  559. !NOTYET !\\
  560. !NOTYET ! !INTERFACE:
  561. !NOTYET !
  562. !NOTYET SUBROUTINE DISTRIBUTE_HIGH_LAT_CO( imdim, jmdim, field2d )
  563. !NOTYET !
  564. !NOTYET ! !USES:
  565. !NOTYET !
  566. !NOTYET use dims, only : idate, sec_month, sec_day, sec_year, mlen
  567. !NOTYET !
  568. !NOTYET ! !INPUT PARAMETERS:
  569. !NOTYET !
  570. !NOTYET integer, intent(in) :: jmdim, imdim ! dimension of grid
  571. !NOTYET !
  572. !NOTYET ! !INPUT/OUTPUT PARAMETERS:
  573. !NOTYET !
  574. !NOTYET real, dimension(imdim, jmdim) :: field2d
  575. !NOTYET !
  576. !NOTYET ! !REVISION HISTORY:
  577. !NOTYET ! 24 Jun 2002 - fjd
  578. !NOTYET ! 14 Feb 2013 - Ph. Le Sager - reincorporated for ED4
  579. !NOTYET !
  580. !NOTYET ! !REMARKS:
  581. !NOTYET !
  582. !NOTYET !EOP
  583. !NOTYET !------------------------------------------------------------------------
  584. !NOTYET !BOC
  585. !NOTYET
  586. !NOTYET real :: year2month
  587. !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/)
  588. !NOTYET
  589. !NOTYET ! high latitude seasonal CO distribution function
  590. !NOTYET integer :: mo,j
  591. !NOTYET
  592. !NOTYET year2month =sec_month/sec_year !scale factor for yearly total to specific month
  593. !NOTYET mo=idate(2)
  594. !NOTYET
  595. !NOTYET do j=1,jmdim
  596. !NOTYET if (j>jmdim*3/4) then
  597. !NOTYET field2d(:,j)=field2d(:,j)*coseas(mo)/12. ! convert to month
  598. !NOTYET else
  599. !NOTYET field2d(:,j)=field2d(:,j)*year2month
  600. !NOTYET endif
  601. !NOTYET enddo
  602. !NOTYET
  603. !NOTYET END SUBROUTINE DISTRIBUTE_HIGH_LAT_CO
  604. !NOTYET !EOC
  605. END MODULE EMISSION_CO