emission_terp.F90 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854
  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_TERP
  14. !
  15. ! !DESCRIPTION: data and methods for terpene emissions.
  16. !\\
  17. !\\
  18. ! !INTERFACE:
  19. !
  20. MODULE EMISSION_TERP
  21. !
  22. ! !USES:
  23. !
  24. use GO, only : gol, goPr, goErr, goBug
  25. use tm5_distgrid, only : dgrid, get_distgrid, scatter, gather
  26. use partools, only : isRoot, par_broadcast
  27. use dims, only : nregions, idate, okdebug
  28. use global_types, only : emis_data, d3_data, isop_data
  29. use emission_read, only : used_providers_terp, has_terp_emis
  30. IMPLICIT NONE
  31. PRIVATE
  32. !
  33. ! !PUBLIC MEMBER FUNCTIONS:
  34. !
  35. public :: Emission_terp_Init
  36. public :: Emission_terp_Done
  37. public :: emission_terp_Declare
  38. public :: Emission_terp_Apply
  39. !
  40. ! !PRIVATE DATA MEMBERS:
  41. !
  42. character(len=*), parameter :: mname = 'emission_terp'
  43. type(isop_data),dimension(nregions),target :: terp_dat
  44. !flags the daily calculation of terpene distribution func.
  45. integer,dimension(nregions) :: day_terp
  46. type( emis_data ), dimension(:,:), allocatable :: terp_emis_2d
  47. type( d3_data ), dimension(:,:), allocatable :: terp_emis_3d
  48. integer :: terp_2dsec, terp_3dsec
  49. logical, allocatable :: has_data_2d(:), has_data_3d(:)
  50. !
  51. ! !REVISION HISTORY:
  52. ! 1 Oct 2010 - Achim Strunk - overhaul for AR5
  53. ! 28 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  54. !
  55. ! !REMARKS:
  56. !
  57. !EOP
  58. !------------------------------------------------------------------------
  59. CONTAINS
  60. !--------------------------------------------------------------------------
  61. ! TM5 !
  62. !--------------------------------------------------------------------------
  63. !BOP
  64. !
  65. ! !IROUTINE: EMISSION_TERP_INIT
  66. !
  67. ! !DESCRIPTION: Allocate memory
  68. !\\
  69. !\\
  70. ! !INTERFACE:
  71. !
  72. SUBROUTINE EMISSION_TERP_INIT( status )
  73. !
  74. ! !USES:
  75. !
  76. use dims, only : ndyn_max, nsrce, tref, im, jm, lm, idate, nlon360, nlat180
  77. use dims, only : sec_month
  78. use global_data, only : rcfile
  79. use emission_read, only : providers_def, numb_providers, sectors_def
  80. !
  81. ! !OUTPUT PARAMETERS:
  82. !
  83. integer, intent(out) :: status
  84. !
  85. ! !REVISION HISTORY:
  86. ! 02 Apr 2014 - v0
  87. !
  88. ! !REMARKS:
  89. !
  90. !EOP
  91. !------------------------------------------------------------------------
  92. !BOC
  93. character(len=*), parameter :: rname = mname//'/emission_terp_init'
  94. integer :: region, i1, i2, j1, j2
  95. integer :: lmr, lsec, ntim, lprov
  96. real :: dtime
  97. ! --- begin -----------------------------------------
  98. status = 0
  99. if(.not. has_terp_emis) return
  100. terp_2dsec = 0
  101. terp_3dsec = 0
  102. do lprov = 1, numb_providers
  103. if (count(used_providers_terp == providers_def(lprov)%name)/=0) then
  104. ! limit AR5 to the 2D biomass burning sector
  105. if ( trim(providers_def(lprov)%name) == 'AR5' ) then
  106. terp_2dsec = terp_2dsec + count( sectors_def%prov == 'AR5 ' &
  107. .and. sectors_def%catname == 'biomassburning')
  108. else
  109. ! use all sectors
  110. terp_2dsec = terp_2dsec + providers_def(lprov)%nsect2d
  111. terp_3dsec = terp_3dsec + providers_def(lprov)%nsect3d
  112. end if
  113. endif
  114. enddo
  115. allocate( terp_emis_2d( nregions, terp_2dsec ) )
  116. allocate( terp_emis_3d( nregions, terp_3dsec ) )
  117. allocate( has_data_2d(terp_2dsec)) ; has_data_2d=.false.
  118. allocate( has_data_3d(terp_3dsec)) ; has_data_3d=.false.
  119. ! allocate information arrays (2d and 3d)
  120. do region=1,nregions
  121. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  122. lmr = lm(region)
  123. do lsec=1,terp_2dsec
  124. allocate( terp_emis_2d(region,lsec)%surf(i1:i2, j1:j2) )
  125. end do
  126. do lsec=1,terp_3dsec
  127. allocate( terp_emis_3d(region,lsec)%d3(i1:i2, j1:j2,lmr) )
  128. end do
  129. dtime=float(ndyn_max)/(2*tref(region)) ! timestep emissions
  130. ntim=86400/nint(dtime) ! number of timesteps in 24 hours for this region
  131. allocate(terp_dat(region)%scalef_isop(ntim,j1:j2)) ! allocate enough memory!
  132. enddo
  133. ! ok
  134. status = 0
  135. END SUBROUTINE EMISSION_TERP_INIT
  136. !EOC
  137. !--------------------------------------------------------------------------
  138. ! TM5 !
  139. !--------------------------------------------------------------------------
  140. !BOP
  141. !
  142. ! !IROUTINE: EMISSION_terp_DONE
  143. !
  144. ! !DESCRIPTION: Free memory
  145. !\\
  146. !\\
  147. ! !INTERFACE:
  148. !
  149. SUBROUTINE EMISSION_TERP_DONE( STATUS )
  150. !
  151. ! !OUTPUT PARAMETERS:
  152. !
  153. integer, intent(out) :: status
  154. !
  155. ! !REVISION HISTORY:
  156. ! 1 Oct 2010 - Achim Strunk - adapted to new structures
  157. ! 28 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  158. !
  159. !EOP
  160. !------------------------------------------------------------------------
  161. !BOC
  162. character(len=*), parameter :: rname = mname//'/Emission_terp_Done'
  163. integer :: region, lsec
  164. ! --- begin --------------------------------------
  165. status = 0
  166. if(.not. has_terp_emis) return
  167. do region = 1, nregions
  168. do lsec=1,terp_2dsec
  169. deallocate( terp_emis_2d(region,lsec)%surf )
  170. end do
  171. do lsec=1,terp_3dsec
  172. deallocate( terp_emis_3d(region,lsec)%d3 )
  173. end do
  174. deallocate(terp_dat(region)%scalef_isop )
  175. end do
  176. deallocate( terp_emis_2d, terp_emis_3d )
  177. deallocate( has_data_2d, has_data_3d )
  178. status = 0
  179. END SUBROUTINE EMISSION_TERP_DONE
  180. !EOC
  181. !--------------------------------------------------------------------------
  182. ! TM5 !
  183. !--------------------------------------------------------------------------
  184. !BOP
  185. !
  186. ! !IROUTINE: EMISSION_TERP_DECLARE
  187. !
  188. ! !DESCRIPTION: Opens, reads and evaluates input files (per month).
  189. ! Provides emissions on 2d/3d-arrays which are then added
  190. ! to mixing ratios in routine *apply.
  191. !\\
  192. !\\
  193. ! !INTERFACE:
  194. !
  195. SUBROUTINE EMISSION_TERP_DECLARE( status )
  196. !
  197. ! !USES:
  198. !
  199. use toolbox, only : coarsen_emission
  200. use dims, only : im, jm, lm, idate, sec_month, nlon360, nlat180, iglbsfc
  201. use chem_param, only : xmterp
  202. use emission_data, only : msg_emis, LCMIP6BMB, LAR5BMB, LMEGAN
  203. ! ---------------- CMIP6 - AR5 - MACC --------------------
  204. use emission_data, only : emis_input_year_nmvoc, emis_input_year_nat
  205. use emission_data, only : emis_input_dir_mac
  206. use emission_data, only : emis_input_dir_megan
  207. use emission_data, only : emis_input_dir_retro
  208. use emission_data, only : emis_input_dir_gfed
  209. use emission_read, only : emission_ar5_regrid_aircraft
  210. use emission_read, only : emission_cmip6_ReadSector
  211. use emission_read, only : emission_cmip6bmb_ReadSector
  212. use emission_read, only : emission_ar5_ReadSector
  213. use emission_read, only : emission_macc_ReadSector
  214. use emission_read, only : emission_megan_ReadSector
  215. use emission_read, only : emission_gfed_ReadSector
  216. use emission_read, only : emission_retro_ReadSector
  217. use emission_read, only : sectors_def, numb_sectors
  218. use emission_read, only : ar5_dim_3ddata
  219. !
  220. ! !OUTPUT PARAMETERS:
  221. !
  222. integer, intent(out) :: status
  223. !
  224. ! !REVISION HISTORY:
  225. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  226. !
  227. !EOP
  228. !------------------------------------------------------------------------
  229. !BOC
  230. character(len=*), parameter :: rname = mname//'/emission_terp_declare'
  231. integer :: region, hasData
  232. integer, parameter :: add_field=0
  233. integer, parameter :: amonth=2
  234. integer :: imr, jmr, lmr, lsec
  235. ! AR5
  236. real,dimension(:,:,:),allocatable :: field3d, field3d2
  237. type(d3_data) :: emis3d, work(nregions)
  238. type(emis_data) :: wrk2D(nregions)
  239. integer :: seccount2d, seccount3d
  240. ! --- begin -----------------------------------------
  241. status = 0
  242. if(.not. has_terp_emis) return
  243. write(gol,'(" EMISS-INFO ------------- read TERPENE emissions -------------")'); call goPr
  244. ! flag calculation of daily terpene emission distribution
  245. day_terp(1:nregions) = -1
  246. do region = 1, nregions
  247. do lsec=1,terp_2dsec
  248. terp_emis_2d(region,lsec)%surf = 0.0
  249. end do
  250. do lsec=1,terp_3dsec
  251. terp_emis_3d(region,lsec)%d3 = 0.0
  252. end do
  253. end do
  254. ! global arrays for coarsening
  255. do region = 1, nregions
  256. if (isRoot)then
  257. allocate(work(region)%d3(im(region),jm(region),ar5_dim_3ddata))
  258. else
  259. allocate(work(region)%d3(1,1,1))
  260. end if
  261. enddo
  262. do region = 1, nregions
  263. wrk2D(region)%surf => work(region)%d3(:,:,1)
  264. end do
  265. ! --------------------------------
  266. ! do a loop over available sectors
  267. ! --------------------------------
  268. ! count 2d and 3d sectors
  269. seccount2d = 0
  270. seccount3d = 0
  271. ! always allocate here 3d data set (for 2d sectors it will be filled in first layer only)
  272. if (isRoot) then
  273. allocate( field3d( nlon360, nlat180, ar5_dim_3ddata ) ) ; field3d = 0.0
  274. else
  275. allocate( field3d( 1, 1, 1 ) )
  276. end if
  277. sec : do lsec = 1, numb_sectors
  278. ! screen out unused providers
  279. if (count(used_providers_terp.eq.sectors_def(lsec)%prov).eq.0) cycle
  280. ! Explicitly skip AR5 non-fires and EDGAR and HYMN inventories
  281. ! Same w/ ED41 and ED42 (but already screened out when building providers list anyway)
  282. if ( sectors_def(lsec)%prov == 'AR5 ' .and. sectors_def(lsec)%catname /= 'biomassburning') cycle
  283. if ( sectors_def(lsec)%prov == 'ED41 ' .or. &
  284. sectors_def(lsec)%prov == 'ED42 ' .or. &
  285. sectors_def(lsec)%prov == 'HYMN ') cycle
  286. field3d = 0.0
  287. if( sectors_def(lsec)%f3d ) then
  288. seccount3d = seccount3d + 1
  289. else
  290. seccount2d = seccount2d + 1
  291. end if
  292. if (isRoot) then ! READ
  293. select case( trim(sectors_def(lsec)%prov) )
  294. case( 'CMIP6BMB' )
  295. call emission_cmip6bmb_ReadSector( 'NMVOC-'//'C10H16', emis_input_year_nmvoc, idate(2), lsec, field3d, status )
  296. IF_NOTOK_RETURN(status=1;deallocate(field3d))
  297. case( 'AR5' )
  298. call emission_ar5_ReadSector( 'terpenes', emis_input_year_nmvoc, idate(2), lsec, field3d, status )
  299. IF_NOTOK_RETURN(status=1)
  300. case( 'MACC' )
  301. ! skip if MEGAN is used, screen out sectors w/o terpenes else
  302. if( (trim(sectors_def(lsec)%name) .eq. 'emiss_bio') .and. (.not. LMEGAN) ) then
  303. if (trim(sectors_def(lsec)%catname) .eq. 'natural') then
  304. call emission_macc_ReadSector( emis_input_dir_mac, 'TERPENES', emis_input_year_nat, idate(2), &
  305. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  306. IF_NOTOK_RETURN(status=1)
  307. else
  308. call emission_macc_ReadSector( emis_input_dir_mac, 'TERPENES', emis_input_year_nmvoc, idate(2), &
  309. '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status )
  310. IF_NOTOK_RETURN(status=1)
  311. endif
  312. endif
  313. case('GFEDv3')
  314. call emission_gfed_ReadSector( emis_input_dir_gfed, 'terpenes', emis_input_year_nmvoc, idate(2), &
  315. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  316. IF_NOTOK_RETURN(status=1)
  317. case('RETRO')
  318. call emission_retro_ReadSector( emis_input_dir_retro, 'MONOTERPENES', emis_input_year_nmvoc, idate(2), &
  319. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  320. IF_NOTOK_RETURN(status=1)
  321. case('MEGAN')
  322. call emission_megan_ReadSector( emis_input_dir_megan, 'monoterpenes', emis_input_year_nat, idate(2), &
  323. sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status )
  324. IF_NOTOK_RETURN(status=1)
  325. case('DUMMY')
  326. case default
  327. write(gol,*) "Error in buidling list of providers USED_PROVIDERS_TERP" ; call goErr
  328. status=1; TRACEBACK; return
  329. end select
  330. ! nothing found???
  331. if( sum(field3d) < 100.*TINY(1.0) ) then
  332. if (okdebug) then
  333. write(gol,'("EMISS-INFO - no terp emissions found for ",a," ",a," for month ",i2 )') &
  334. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  335. endif
  336. hasData=0
  337. else
  338. if (okdebug) then
  339. write(gol,'("EMISS-INFO - found terp emissions for ",a," ",a," for month ",i2 )') &
  340. trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr
  341. endif
  342. ! scale from kg/s to kg/month
  343. field3d = field3d * sec_month ! kg / month
  344. hasData=1
  345. end if
  346. end if
  347. call Par_broadcast(hasData, status)
  348. IF_NOTOK_RETURN(status=1)
  349. if (hasData == 0) then
  350. cycle sec
  351. else
  352. if ( sectors_def(lsec)%f3d ) then
  353. has_data_3d(seccount3d)=.true.
  354. else
  355. has_data_2d(seccount2d)=.true.
  356. end if
  357. end if
  358. ! Distinguish 2d/3d sectors
  359. if( sectors_def(lsec)%f3d ) then
  360. write(gol,'("EMISS-ERROR - Unexpected 3D data: code needed")'); call goErr
  361. status=1; TRACEBACK; return
  362. else ! ar5_sector is 2d
  363. ! ---------------------------
  364. ! 2d data (Anthropogenic, Ships, Biomassburning)
  365. ! ---------------------------
  366. if (isRoot) then ! print total & regrid
  367. call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'TERP', xmterp, sum(field3d(:,:,1)) )
  368. call coarsen_emission( 'terp '//sectors_def(lsec)%name, nlon360, nlat180, field3d(:,:,1), wrk2D, add_field, status )
  369. IF_NOTOK_RETURN(status=1)
  370. end if
  371. do region = 1, nregions
  372. call scatter(dgrid(region), terp_emis_2d(region,seccount2d)%surf, work(region)%d3(:,:,1), 0, status)
  373. IF_NOTOK_RETURN(status=1)
  374. end do
  375. end if
  376. end do sec ! sectors
  377. deallocate( field3d )
  378. do region = 1, nregions
  379. if (associated(wrk2D(region)%surf)) nullify(wrk2D(region)%surf)
  380. deallocate( work(region)%d3 )
  381. end do
  382. ! ok
  383. status = 0
  384. END SUBROUTINE EMISSION_TERP_DECLARE
  385. !EOC
  386. !--------------------------------------------------------------------------
  387. ! TM5 !
  388. !--------------------------------------------------------------------------
  389. !BOP
  390. !
  391. ! !IROUTINE: EMISSION_TERP_APPLY
  392. !
  393. ! !DESCRIPTION: Take monthly emissions, and
  394. ! - split them vertically
  395. ! - apply time splitting factors
  396. ! - add them up (add_terp)
  397. !\\
  398. !\\
  399. ! !INTERFACE:
  400. !
  401. SUBROUTINE EMISSION_TERP_APPLY( region, status )
  402. !
  403. ! !USES:
  404. !
  405. use dims, only : itaur, nsrce, tref
  406. use dims, only : im, jm, lm, mlen, ndyn_max
  407. use dims, only : yref, dy, ybeg
  408. use datetime, only : tau2date, get_day
  409. use emission_data, only : emission_vdist_by_sector, bb_cycle
  410. use emission_data, only : emis_bb_trop_cycle
  411. use chem_param, only : iterp, xmterp
  412. use emission_data, only : do_add_3d, do_add_3d_cycle
  413. use emission_read, only : sectors_def, numb_sectors
  414. use tracer_data, only : tracer_print
  415. USE partools, ONLY : isRoot, par_reduce, par_reduce_element
  416. !
  417. ! !INPUT PARAMETERS:
  418. !
  419. integer, intent(in) :: region
  420. !
  421. ! !OUTPUT PARAMETERS:
  422. !
  423. integer, intent(out) :: status
  424. !
  425. ! !REVISION HISTORY:
  426. ! 1 Oct 2010 - Achim Strunk - adapted for AR5
  427. ! 27 Apr 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  428. ! 2 Apr 2014 - J. E. Williams - terpene emissions introduced
  429. !
  430. ! !REMARKS:
  431. !
  432. !EOP
  433. !------------------------------------------------------------------------
  434. !BOC
  435. character(len=*), parameter :: rname = mname//'/emission_terp_apply'
  436. integer, dimension(6) :: idater
  437. real :: dtime, fraction
  438. integer :: imr, jmr, lmr, lsec
  439. integer :: seccount2d, seccount3d
  440. type(d3_data) :: emis3d
  441. real :: dlatr
  442. integer :: day, ntim, i1, i2, j1, j2
  443. ! --- begin -----------------------------------------
  444. status = 0
  445. if(.not. has_terp_emis) return
  446. if( okdebug ) then
  447. write(gol,*) 'start of emission_apply_terp'; call goPr
  448. end if
  449. ! --- WHERE, WHEN
  450. CALL TAU2DATE( itaur(region), idater)
  451. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  452. day = GET_DAY(idater(2),idater(3),mlen)
  453. dtime = float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX.
  454. if(okdebug) then
  455. write(gol,*)'emission_terp_apply in region ',region,' at date: ',idater, ' with time step:', dtime
  456. call goPr
  457. end if
  458. ! --- DIURNAL CYCLE
  459. if ( day /= day_terp(region) ) then
  460. dlatr = dy/yref(region) ! should be the same as lli%dlat_deg
  461. dtime=float(ndyn_max)/(2*tref(region)) ! timestep emissions
  462. ntim=86400/nint(dtime) ! number of timesteps in 24 hours.
  463. if ( okdebug ) then
  464. write (gol,*) 'day',day,'ntim',ntim,'jm',jm(region),'ybeg',&
  465. ybeg(region),'dlatr',dlatr,'sum(mlen)',sum(mlen),nsrce,dtime
  466. call goPr
  467. end if
  468. ! NOTE : latitude_start must be an integer. This is ok as long as
  469. ! the regions ybeg are integer.
  470. call scale_terp(day, ntim, j1, terp_dat(region)%scalef_isop, ybeg(region), dlatr, sum(mlen))
  471. day_terp(region) = day !avoid recalculation during this day!
  472. endif
  473. ! --- ADD EMISSIONS
  474. ! get a working structure for 3d emissions
  475. allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  476. ! default: no additional splitting
  477. fraction = 1.0
  478. ! count 2d and 3d sectors
  479. seccount2d = 0
  480. seccount3d = 0
  481. ! cycle over sectors
  482. do lsec = 1, numb_sectors
  483. ! screen out unused providers
  484. if (count(used_providers_terp.eq.sectors_def(lsec)%prov).eq.0) cycle
  485. ! add up contribution from all proc
  486. ! screen out sectors which do not contain terpene emissions
  487. !
  488. if ( sectors_def(lsec)%prov == 'AR5 ' .and. sectors_def(lsec)%catname /= 'biomassburning') cycle
  489. if ( sectors_def(lsec)%prov == 'ED41 ' .or. &
  490. sectors_def(lsec)%prov == 'ED42 ' .or. &
  491. sectors_def(lsec)%prov == 'HYMN ') cycle
  492. ! distinguish between 2d/3d sectors
  493. if( sectors_def(lsec)%f3d ) then
  494. seccount3d = seccount3d + 1
  495. if (.not.has_data_3d(seccount3d)) cycle
  496. emis3d%d3 = terp_emis_3d(region,seccount3d)%d3
  497. else
  498. seccount2d = seccount2d + 1
  499. if (.not.has_data_2d(seccount2d)) cycle
  500. emis3d%d3 = 0.0
  501. ! vertically distribute according to sector
  502. call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'TERP', region, terp_emis_2d(region,seccount2d), &
  503. emis3d, status )
  504. IF_NOTOK_RETURN(status=1)
  505. endif
  506. ! add dataset according to sector and category
  507. if( trim(sectors_def(lsec)%catname) == "natural" ) then
  508. call do_add_terp( region, i1, j1, emis3d%d3, status) ! special routine
  509. IF_NOTOK_RETURN(status=1)
  510. elseif ( trim(sectors_def(lsec)%catname) == "biomassburning" ) then
  511. if( emis_bb_trop_cycle) then
  512. call do_add_3d_cycle( region, iterp, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmterp, xmterp, status, fraction )
  513. IF_NOTOK_RETURN(status=1)
  514. else
  515. call do_add_3d( region, iterp, i1, j1, emis3d%d3, xmterp, xmterp, status)
  516. IF_NOTOK_RETURN(status=1)
  517. end if
  518. endif
  519. end do
  520. deallocate( emis3d%d3 )
  521. if(okdebug) then
  522. write(gol,*) 'end of emission_terp_apply' ; call goPr
  523. end if
  524. ! OK
  525. status = 0
  526. END SUBROUTINE EMISSION_TERP_APPLY
  527. !EOC
  528. !--------------------------------------------------------------------------
  529. ! TM5 !
  530. !--------------------------------------------------------------------------
  531. !BOP
  532. !
  533. ! !IROUTINE: SCALE_TERP
  534. !
  535. ! !DESCRIPTION: calculates factors to be multiplied with an emission field
  536. ! (e.g. terp) in order to simulate a diurnal cycle in
  537. ! emissions (caused by solar dependency), e.g. :
  538. !
  539. ! rm(i,j,1,terp) = rm(i,j,1,terp) + e_terp(i,j)*scalef(ipos,j)*dt
  540. !
  541. ! with ipos depending on time and longitude.
  542. !
  543. ! The scalefactor is calculated for -180 longitude and the
  544. ! mean value for ntim timesteps during a day is scaled to 1.
  545. ! The routine should be called every day, since the position
  546. ! relative to the sun changes. The scaling is based on the
  547. ! cos(zenith) which is 1 for overhead sun and 0 at sunset.
  548. ! For the polar night, the values are set to 1.0 (even distribution)
  549. ! One thing is sure afterwards: we know that e_terp kg/s is added.
  550. !\\
  551. !\\
  552. ! !INTERFACE:
  553. !
  554. SUBROUTINE SCALE_TERP(day, ntim, j1, scalef, lat_start, dlat, idayy)
  555. !
  556. ! !INPUT PARAMETERS:
  557. !
  558. integer, intent(in) :: day ! day number based on 365 days a year
  559. integer, intent(in) :: ntim ! number of timesteps during one day
  560. integer, intent(in) :: j1 ! start index for latitudes
  561. integer, intent(in) :: lat_start ! southern edge of the domain (integer!)
  562. real, intent(in) :: dlat ! increment (real!)
  563. integer, intent(in) :: idayy ! # days this year
  564. !
  565. ! !OUTPUT PARAMETERS:
  566. !
  567. real, dimension(:,j1:), intent(out) :: scalef
  568. !
  569. ! !REVISION HISTORY:
  570. ! 10 May 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  571. !
  572. ! !REMARKS:
  573. !
  574. !EOP
  575. !------------------------------------------------------------------------
  576. !BOC
  577. real :: pi, piby, obliq, deday, delta, dt, lon, hr, lat, houra, smm,scale_test
  578. integer :: i, j
  579. pi = acos(-1.0)
  580. piby = pi/180.0
  581. obliq = 23.45 * piby
  582. deday = 4.88 + 2*pi/idayy * day
  583. delta = asin(sin(obliq)*sin(deday))
  584. dt = 24./ntim ! timestep in hours
  585. lon = -180.0*piby ! calculate for dateline
  586. do j=j1,ubound(scalef,2) ! over latitudes
  587. hr = - 0.5*dt ! shift half a interval back
  588. lat = (lat_start + (j-0.5)*dlat)*piby
  589. smm = 0.0
  590. do i=1,ntim
  591. hr = hr + dt
  592. houra = lon - pi + hr * (2.*pi/24.)
  593. !checking if the sun is above horizon
  594. scale_test = max(sin(delta)*sin(lat) + cos(delta)*cos(lat)*cos(houra),0.0)
  595. scalef(i,j) = cos(houra)
  596. !smm = smm+scalef(i,j)/ntim
  597. smm = smm+scale_test/ntim
  598. end do
  599. if ( smm > 1e-5 ) then
  600. scalef(1:ntim,j) = scalef(1:ntim,j)+1
  601. else
  602. scalef(1:ntim,j) = 1.0
  603. end if
  604. end do
  605. END SUBROUTINE SCALE_TERP
  606. !--------------------------------------------------------------------------
  607. ! TM5 !
  608. !--------------------------------------------------------------------------
  609. !BOP
  610. !
  611. ! !IROUTINE: DO_ADD_terp
  612. !
  613. ! !DESCRIPTION: Add local time splitted emissions to chemical array. The
  614. ! difference with "do_add_3d" is the number of levels this is
  615. ! applied to: 4 instead of all levels.
  616. !\\
  617. !\\
  618. ! !INTERFACE:
  619. !
  620. SUBROUTINE DO_ADD_TERP(region, i1, j1, em_terp, status)
  621. !
  622. ! !USES:
  623. !
  624. use dims, only : tref, nsrce, sec_month,dx,xref,xbeg, ndyn_max,itaur
  625. use global_data, only : mass_dat
  626. use chem_param, only : iterp, xmterp
  627. #ifdef with_budgets
  628. use budget_global, only : budg_dat, nzon_vg
  629. use emission_data, only : sum_emission, budemi_dat
  630. #endif
  631. use datetime, only : tau2date
  632. !
  633. ! !INPUT PARAMETERS:
  634. !
  635. integer, intent(in) :: region, i1, j1
  636. real, intent(in) :: em_terp(i1:,j1:,:)
  637. !
  638. ! !OUTPUT PARAMETERS:
  639. !
  640. integer, intent(out) :: status
  641. !
  642. ! !REVISION HISTORY:
  643. ! 03 Apr 2014 - J. E. Williams - copied from do_add_isop and adapted for terpenes
  644. !
  645. ! !REMARKS:
  646. ! - may be be wise to add a level_limit keyword to the original "do_add_3d"
  647. !
  648. !EOP
  649. !------------------------------------------------------------------------
  650. !BOC
  651. character(len=*), parameter :: rname = mname//'/do_add_terp'
  652. real :: dtime, month2dt, x, xlon, dtime2
  653. integer :: i, j, l, nzone, nzone_v, i2, j2
  654. integer :: sec_in_day,ntim,itim,ipos
  655. integer,dimension(6)::idater
  656. real, dimension(:,:,:,:), pointer :: rm
  657. ! --- begin -----------------------------------------
  658. rm => mass_dat(region)%rm
  659. CALL GET_DISTGRID( dgrid(region), I_STOP=i2, J_STOP=j2 )
  660. dtime = float(nsrce) / (2*tref(region)) ! timestep emissions
  661. dtime2 = float(ndyn_max) / (2*tref(region)) ! timestep emissions based on non-CFL interference (CMK 05/2006)
  662. month2dt = dtime/sec_month ! conversion to emission per timestep (CFl induced timestep)
  663. ntim = 86400/nint(dtime2) ! number of timesteps in 24 hours based on non-CFL interference (CMK 05/2006)
  664. call tau2date(itaur(region),idater)
  665. sec_in_day = idater(4)*3600 + idater(5)*60 + idater(6) ! elapsd seconds this day
  666. itim = sec_in_day/nint(dtime2)+1 ! # time interval CMK 05/2006
  667. !
  668. ! Terpene emissions only at the surface
  669. !
  670. do l= 1,4
  671. do j=j1,j2
  672. do i=i1,i2
  673. xlon = xbeg(region) + (i-0.5)*dx/xref(region)
  674. ! itim = 1 and lon = -180 --->position 1
  675. ! itim = ntim ant lon = -180 --->position ntim, etc.
  676. ! itim = 1 and lon = 0 ---->position ntim/2
  677. ipos = 1 + mod(int((xlon+180.)*ntim/360.) + (itim-1),ntim) !position in array depending on time and lon.
  678. x=em_terp(i,j,l)*terp_dat(region)%scalef_isop(ipos,j)*month2dt
  679. !x=em_terp(i,j,l)*month2dt
  680. !write(1111,*) terp_dat(region)%scalef_isop(ipos,j),x
  681. ! x = em_terp(i,j,l)*month2dt
  682. rm(i,j,l,iterp)=rm(i,j,l,iterp)+x
  683. #ifdef with_budgets
  684. nzone=budg_dat(region)%nzong(i,j) !global budget
  685. nzone_v=nzon_vg(l)
  686. budemi_dat(region)%emi(i,j,nzone_v,iterp) = &
  687. budemi_dat(region)%emi(i,j,nzone_v,iterp) + x/xmterp*1e3
  688. if(iterp ==1) sum_emission(region) = sum_emission(region) + x !in kg
  689. #endif
  690. enddo
  691. enddo
  692. enddo
  693. nullify(rm)
  694. status=0
  695. END SUBROUTINE DO_ADD_TERP
  696. !EOC
  697. END MODULE EMISSION_TERP