emission__co2.F90 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618
  1. #define TRACEBACK write (gol,'("in ",a," (",a,i6,")")') rname, __FILE__, __LINE__ ; call goErr
  2. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  3. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  4. !
  5. #include "tm5.inc"
  6. #include "output.inc"
  7. !
  8. !------------------------------------------------------------------------------
  9. ! TM5 !
  10. !------------------------------------------------------------------------------
  11. !BOP
  12. !
  13. ! !MODULE: EMISSION
  14. !
  15. ! !DESCRIPTION: wrappers around various emissions (init/declare/apply/done)
  16. ! routines, needed for TM5 CBM4 version.
  17. ! Also hold emissions budget variables.
  18. !
  19. !\\
  20. !\\
  21. ! !INTERFACE:
  22. !
  23. MODULE EMISSION
  24. !
  25. ! !USES:
  26. !
  27. USE GO, ONLY : gol, goErr, goPr
  28. USE GO, ONLY : GO_Timer_Def, GO_Timer_End, GO_Timer_Start
  29. USE TM5_DISTGRID, ONLY : dgrid, Get_DistGrid, scatter, gather
  30. USE dims, ONLY : nregions, okdebug
  31. USE emission_data, ONLY : plandr, emis2D
  32. #ifdef with_budgets
  33. USE emission_data, ONLY : budemi_dat, budemi_data, sum_emission
  34. USE budget_global, ONLY : nbud_vg,nbudg
  35. USE chem_param, ONLY : ntracet
  36. #endif
  37. use emission_co2 , only : Emission_CO2_Init , Emission_CO2_Done , emission_co2_declare , emission_co2_apply
  38. IMPLICIT NONE
  39. PRIVATE
  40. !
  41. ! !PUBLIC MEMBER FUNCTIONS:
  42. !
  43. PUBLIC :: Emission_Init ! allocate/init budget var; call other emiss-related init
  44. PUBLIC :: Emission_Done ! gather/write final budget
  45. PUBLIC :: Declare_Emission ! allocate emiss var (new run), read emiss data (new month)
  46. PUBLIC :: Emission_Apply !
  47. !
  48. ! !PRIVATE DATA MEMBERS:
  49. !
  50. ! These budgets should be communicated to chemistry
  51. #ifdef with_budgets
  52. REAL, DIMENSION(nbudg, nbud_vg, ntracet) :: budemig
  53. REAL, DIMENSION(nbudg, nbud_vg, ntracet), PUBLIC :: budemig_all ! for buggy MPI (see budget_global.F90)
  54. #endif
  55. integer :: itim_appl, itim_co2
  56. CHARACTER(len=*), PARAMETER :: mname = 'emission'
  57. !
  58. ! !REVISION HISTORY:
  59. ! 20 Aug 2010 - A. Strunk - Adapted to AR5 emissions + various other changes.
  60. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  61. ! 14 May 2014 - T. van Noije- created from chemistry version
  62. ! production of CO2 from other C emissions is
  63. ! currently not included
  64. !
  65. ! - NO MORE REVISION HISTORY. This is gathered from repository from now on.
  66. !EOP
  67. !----------------------------------------------------------------------
  68. CONTAINS
  69. !----------------------------------------------------------------------
  70. ! TM5 !
  71. !----------------------------------------------------------------------
  72. !BOP
  73. !
  74. ! !IROUTINE: EMISSION_INIT
  75. !
  76. ! !DESCRIPTION: Initialise emission fields and parameters by reading
  77. ! the rc-file. Allocate and initialize budget variables.
  78. !\\
  79. !\\
  80. ! !INTERFACE:
  81. !
  82. SUBROUTINE EMISSION_INIT( status )
  83. !
  84. ! !USES:
  85. !
  86. use GO, only : TrcFile, Init, Done, ReadRc
  87. use meteodata, only : Set, oro_dat
  88. use dims, only : iglbsfc, nregions, lm
  89. use dims, only : idate, ndyn_max, tref
  90. use global_data, only : inputdir, rcfile
  91. use emission_data, only : emis_input_dir
  92. use emission_data, only : emis_input_year_co2
  93. use emission_data, only : LCMIP6
  94. use emission_data, only : LAR5
  95. use emission_data, only : emis_input_dir_cmip6
  96. use emission_data, only : emis_input_dir_ar5
  97. use emission_read, only : emission_read_init
  98. !
  99. ! !OUTPUT PARAMETERS:
  100. !
  101. INTEGER, INTENT(out) :: status
  102. !
  103. ! !REVISION HISTORY:
  104. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  105. !
  106. !EOP
  107. !------------------------------------------------------------------------------
  108. !BOC
  109. CHARACTER(len=*), PARAMETER :: rname = mname//'/Emission_Init'
  110. INTEGER :: region, i1, i2, j1, j2, ntim, lmr, imode
  111. TYPE(TrcFile) :: rcF
  112. REAL :: dtime
  113. integer :: emis_input_year
  114. ! -----------------------------------
  115. ! read settings from rcfile
  116. ! -----------------------------------
  117. call Init( rcF, rcfile, status )
  118. IF_NOTOK_RETURN(status=1)
  119. if (okdebug) then
  120. write(gol,*) "EMISS-INFO - running year : ", idate(1) ; call goPr
  121. end if
  122. ! emission base year (assumption: no run overlapping more than one year)
  123. call ReadRc( rcF, 'input.emis.year', emis_input_year, status, default=idate(1) )
  124. IF_ERROR_RETURN(status=1)
  125. write(gol,*) 'EMISS-INFO - Anthropogenic emissions'; call goPr
  126. write(gol,*) 'EMISS-INFO - requested for the following years:'; call goPr
  127. call ReadRc( rcF, 'input.emis.year.co2', emis_input_year_co2, status, default=emis_input_year )
  128. IF_ERROR_RETURN(status=1)
  129. write(gol,*) 'EMISS-INFO - CO2: ', emis_input_year_co2; call goPr
  130. ! default directory for emissions data is "standard input files" dir
  131. emis_input_dir=trim(inputdir)
  132. ! directory of each data provider
  133. call ReadRc( rcF, 'input.emis.dir.CMIP6', emis_input_dir_cmip6, status, default=emis_input_dir )
  134. IF_ERROR_RETURN(status=1)
  135. call ReadRc( rcF, 'input.emis.dir.AR5', emis_input_dir_ar5, status, default=emis_input_dir )
  136. IF_ERROR_RETURN(status=1)
  137. ! Flags
  138. call ReadRc( rcF, 'use_cmip6', LCMIP6, status, default=.false. )
  139. IF_ERROR_RETURN(status=1)
  140. call ReadRc( rcF, 'use_ar5', LAR5, status, default=.false. )
  141. IF_ERROR_RETURN(status=1)
  142. ! very basic checks
  143. if (count((/ LCMIP6, LAR5 /)) > 1) then
  144. write(gol,*) 'ERROR: You use more than one ANTHROPOGENIC inventory'; call goErr
  145. status=1; TRACEBACK; return
  146. end if
  147. ! init providers info
  148. call emission_read_init( rcF, status )
  149. ! used by vertical distribution:
  150. CALL Set( oro_dat(iglbsfc), status, used=.TRUE. )
  151. ! Allocate data
  152. ! -------------
  153. DO region=1,nregions
  154. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  155. lmr = lm(region)
  156. ALLOCATE( plandr(region)%surf(i1:i2, j1:j2))
  157. ALLOCATE( emis2D(region)%surf(i1:i2, j1:j2))
  158. #ifdef with_budgets
  159. ALLOCATE( budemi_dat(region)%emi(i1:i2, j1:j2, nbud_vg, ntracet) )
  160. budemi_dat(region)%emi = 0.0
  161. sum_emission(region) = 0.0
  162. #endif
  163. ENDDO
  164. ! Done
  165. ! -------------
  166. call Done( rcF, status )
  167. IF_NOTOK_RETURN(status=1)
  168. ! define timers:
  169. call GO_Timer_Def( itim_appl, 'emission appl', status )
  170. IF_NOTOK_RETURN(status=1)
  171. call GO_Timer_Def( itim_co2, 'emission co2', status )
  172. IF_NOTOK_RETURN(status=1)
  173. status = 0
  174. END SUBROUTINE EMISSION_INIT
  175. !EOC
  176. !------------------------------------------------------------------------------
  177. ! TM5 !
  178. !------------------------------------------------------------------------------
  179. !BOP
  180. !
  181. ! !IROUTINE: EMISSION_DONE
  182. !
  183. ! !DESCRIPTION: calculate and write final budgets
  184. !\\
  185. !\\
  186. ! !INTERFACE:
  187. !
  188. SUBROUTINE EMISSION_DONE( status )
  189. !
  190. ! !USES:
  191. !
  192. USE dims, ONLY : nregions, im, jm
  193. #ifdef with_budgets
  194. USE chem_param, ONLY : ntracet, names
  195. USE budget_global, ONLY : budget_file_global, nbud_vg, budg_dat, nbudg, NHAB
  196. #ifdef with_hdf4
  197. USE file_hdf, ONLY : THdfFile, TSds
  198. USE file_hdf, ONLY : Init, Done, WriteAttribute, WriteData, SetDim
  199. #endif
  200. USE Dims, ONLY : region_name
  201. USE partools, ONLY : isRoot, par_reduce, par_reduce_element
  202. #endif
  203. !
  204. ! !OUTPUT PARAMETERS:
  205. !
  206. INTEGER, INTENT(out) :: status
  207. !
  208. ! !REVISION HISTORY:
  209. ! 16 Jul 2010 - A. Strunk -
  210. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  211. !
  212. !EOP
  213. !------------------------------------------------------------
  214. !BOC
  215. CHARACTER(len=*), PARAMETER :: rname = mname//'/Emission_Done'
  216. INTEGER :: region, i1, i2, j1, j2
  217. #ifdef with_budgets
  218. #ifdef with_hdf4
  219. TYPE(THdfFile) :: io
  220. TYPE(TSds) :: sds
  221. #endif
  222. REAL, DIMENSION(:,:,:,:), ALLOCATABLE :: collect_emissions
  223. INTEGER :: nsend,j,i,n,nzone,nzone_v
  224. real, dimension(nregions) :: sum_emission_all
  225. #endif
  226. ! --- begin ---------------------------------
  227. #ifdef with_budgets
  228. ! add up contribution from all proc
  229. DO region = 1, nregions
  230. CALL PAR_REDUCE(sum_emission(region), 'SUM', sum_emission_all(region), status)
  231. IF_NOTOK_RETURN(status=1)
  232. END DO
  233. ! Write global budget of tracer #1
  234. IF ( isRoot ) THEN
  235. write (gol,'("EMISS-INFO - ----------------------------------------------")'); call goPr
  236. write (gol,'("EMISS-INFO - Budget of tracer ",a," (kg) ")') trim(names(1)) ; call goPr
  237. write (gol,'("EMISS-INFO - ----------------------------------------------")'); call goPr
  238. do region = 1, nregions
  239. write (gol,'(A,E13.6)') 'EMISS-INFO - mass emitted : ',sum_emission_all(region); call goPr
  240. enddo
  241. #ifdef with_hdf4
  242. CALL Init(io, budget_file_global, 'write', status)
  243. IF_NOTOK_RETURN(status=1)
  244. CALL WriteAttribute(io, 'sum_emission', sum_emission_all, status)
  245. IF_NOTOK_RETURN(status=1)
  246. #endif
  247. budemig = 0.0
  248. END IF
  249. ! Gather budgets
  250. REG: DO region = 1, nregions
  251. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  252. if (isRoot) then
  253. ALLOCATE(collect_emissions(im(region), jm(region), nbud_vg, ntracet))
  254. else
  255. ALLOCATE(collect_emissions(1,1,1,1) )
  256. end if
  257. CALL GATHER( dgrid(region), budemi_dat(region)%emi, collect_emissions, 0, status)
  258. IF_NOTOK_RETURN(status=1)
  259. #ifdef with_hdf4
  260. ! write Not-Horizontally-Aggregated-Budgets
  261. IF (isRoot.and.NHAB) THEN
  262. CALL Init(Sds,io, 'budemi_dat_'//region_name(region),(/im(region),jm(region),nbud_vg,ntracet/), 'real(4)', status)
  263. CALL SetDim(Sds, 0, 'im_'//region_name(region),'longitude', (/(j, j=1,im(region))/), status)
  264. CALL SetDim(Sds, 1, 'jm_'//region_name(region),'latitude', (/(j, j=1,jm(region))/), status)
  265. CALL SetDim(Sds, 2, 'nbud_vg','vertical layer', (/(j, j=1,nbud_vg)/), status)
  266. CALL SetDim(Sds, 3, 'ntracet','tracer number', (/(j, j=1,ntracet)/), status)
  267. IF_NOTOK_RETURN(status=1)
  268. CALL WriteData(Sds,collect_emissions,status)
  269. IF_NOTOK_RETURN(status=1)
  270. CALL Done(Sds, status)
  271. IF_NOTOK_RETURN(status=1)
  272. ENDIF
  273. #endif
  274. ! horizontally aggregates budgets
  275. DO n=1,ntracet
  276. DO nzone_v=1,nbud_vg
  277. DO j=j1,j2
  278. DO i=i1,i2
  279. nzone = budg_dat(region)%nzong(i,j)
  280. budemig(nzone,nzone_v,n) = budemig(nzone,nzone_v,n) + budemi_dat(region)%emi(i,j,nzone_v,n)
  281. END DO
  282. END DO !j
  283. END DO !nzone_v
  284. END DO !nt
  285. DEALLOCATE( collect_emissions )
  286. DEALLOCATE( budemi_dat(region)%emi )
  287. ENDDO REG
  288. CALL PAR_REDUCE_ELEMENT( budemig, 'SUM', budemig_all, status)
  289. IF_NOTOK_RETURN(status=1)
  290. #ifdef with_hdf4
  291. ! Write horizontally aggregated budget
  292. IF ( isRoot ) THEN
  293. CALL Init(Sds,io, 'budemi',(/nbudg,nbud_vg,ntracet/), 'real(8)', status)
  294. IF_NOTOK_RETURN(status=1)
  295. CALL SetDim(Sds, 0, 'nbudg','horizontal region', (/(j, j=1,nbudg)/), status)
  296. CALL SetDim(Sds, 1, 'nbud_vg','vertical layer', (/(j, j=1,nbud_vg)/), status)
  297. CALL SetDim(Sds, 2, 'ntracet','tracer number', (/(j, j=1,ntracet)/), status)
  298. IF_NOTOK_RETURN(status=1)
  299. CALL WriteData(Sds,budemig_all,status)
  300. IF_NOTOK_RETURN(status=1)
  301. CALL Done(Sds, status)
  302. IF_NOTOK_RETURN(status=1)
  303. CALL Done(io, status)
  304. IF_NOTOK_RETURN(status=1)
  305. ENDIF
  306. #endif
  307. #endif /* BUDGETS */
  308. ! call other emission_*_done routines
  309. CALL FREE_EMISSION(status)
  310. IF_NOTOK_RETURN(status=1)
  311. ! ok
  312. status = 0
  313. END SUBROUTINE EMISSION_DONE
  314. !EOC
  315. !---------------------------------------------------------------------------
  316. ! TM5 !
  317. !---------------------------------------------------------------------------
  318. !BOP
  319. !
  320. ! !IROUTINE: DECLARE_EMISSION
  321. !
  322. ! !DESCRIPTION: Called at run start (init/allocate emiss data) and then at
  323. ! beginning of every month to just read in data.
  324. ! Called from SS_MONTHLY_UPDATE.
  325. !\\
  326. !\\
  327. ! !INTERFACE:
  328. !
  329. SUBROUTINE DECLARE_EMISSION( status )
  330. !
  331. ! !USES:
  332. !
  333. USE Grid, ONLY : FillGrid
  334. USE MDF, ONLY : MDF_Open, MDF_NETCDF, MDF_READ, MDF_Inq_VarID, MDF_Get_Var, MDF_Close
  335. USE dims, ONLY : im, jm, lm, newsrun
  336. USE dims, ONLY : nregions, iglbsfc, nlat180, nlon360
  337. USE chem_param
  338. USE partools, ONLY : isRoot
  339. USE global_data, ONLY : emis_data
  340. USE meteodata, ONLY : global_lli
  341. ! AR5/EDGAR4
  342. use emission_data, only : emis_input_dir
  343. !
  344. ! !OUTPUT PARAMETERS:
  345. !
  346. INTEGER, INTENT(out) :: status
  347. !
  348. ! !REVISION HISTORY:
  349. ! 16 Jul 2010 - A. Strunk - Adapted to revised emission_*.F90 routines
  350. ! 27 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  351. !
  352. ! !REMARKS:
  353. ! - anything that is done only if it's a newrun, and that does not require meteo data, should go in INIT
  354. !
  355. !EOP
  356. !------------------------------------------------------------------------------
  357. !BOC
  358. CHARACTER(len=*), PARAMETER :: rname = mname//'/declare_emission'
  359. INTEGER :: region, imode, hid, varid
  360. REAL, DIMENSION(:,:), ALLOCATABLE :: pland
  361. type(emis_data), dimension(nregions) :: wrk
  362. ! ---------------------------------------------------------------
  363. ! ** land fraction
  364. ! ---------------------------------------------------------------
  365. IF( newsrun ) THEN
  366. if(isRoot)then
  367. ALLOCATE( pland(nlon360,nlat180) )
  368. DO region=1,nregions
  369. allocate( wrk(region)%surf(im(region),jm(region)) )
  370. end DO
  371. else
  372. ALLOCATE( pland(1,1) )
  373. DO region=1,nregions
  374. allocate( wrk(region)%surf(1,1))
  375. end DO
  376. end if
  377. if (isRoot)then
  378. CALL MDF_Open( TRIM(emis_input_dir)//'/land/landfraction.nc4', MDF_NETCDF, MDF_READ, hid, status )
  379. IF_NOTOK_RETURN(status=1)
  380. CALL MDF_Inq_VarID( hid, 'LANDFRACTION', varid, status )
  381. IF_NOTOK_RETURN(status=1)
  382. CALL MDF_Get_Var( hid, varid, pland, status )
  383. IF_NOTOK_RETURN(status=1)
  384. CALL MDF_Close( hid, status )
  385. IF_NOTOK_RETURN(status=1)
  386. ! coarsen or distribute to zoom regions:
  387. DO region = 1, nregions
  388. ! convert grid:
  389. CALL FillGrid( global_lli(region), 'n', wrk(region)%surf, &
  390. global_lli(iglbsfc), 'n', pland, 'area-aver', status )
  391. IF_NOTOK_RETURN(status=1)
  392. END DO
  393. end if
  394. DO region = 1, nregions
  395. call scatter( dgrid(region), plandr(region)%surf, wrk(region)%surf, 0, status)
  396. IF_NOTOK_RETURN(status=1)
  397. DEALLOCATE( wrk(region)%surf )
  398. END DO
  399. DEALLOCATE( pland )
  400. ENDIF
  401. ! ---------------------------------------------------------------
  402. ! ** init each constituent
  403. ! ---------------------------------------------------------------
  404. ! ** 1st time, initialise emissions
  405. IF ( newsrun ) THEN
  406. CALL Emission_CO2_Init( status )
  407. IF_NOTOK_RETURN(status=1)
  408. END IF
  409. ! ** every month, read and re-grid
  410. CALL emission_co2_declare( status )
  411. IF_NOTOK_RETURN(status=1)
  412. ! ok
  413. status = 0
  414. END SUBROUTINE DECLARE_EMISSION
  415. !EOC
  416. !-------------------------------------------------------------------
  417. ! TM5 !
  418. !-------------------------------------------------------------------
  419. !BOP
  420. !
  421. ! !IROUTINE: EMISSION_APPLY
  422. !
  423. ! !DESCRIPTION: Call emission_apply methods of constituent modules.
  424. ! --> add current emissions to tracers array.
  425. !\\
  426. !\\
  427. ! !INTERFACE:
  428. !
  429. SUBROUTINE EMISSION_APPLY( region, status )
  430. !
  431. ! !USES:
  432. !
  433. USE chem_param
  434. !
  435. ! !INPUT PARAMETERS:
  436. !
  437. INTEGER, INTENT(in) :: region
  438. !
  439. ! !OUTPUT PARAMETERS:
  440. !
  441. INTEGER, INTENT(out) :: status
  442. !
  443. ! !REVISION HISTORY:
  444. ! 16 Jul 2010 - A. Strunk - Adapted to revised emission_*.F90 routines
  445. ! 27 Mar 2012 - P. Le Sager - cleanup for lat-lon mpi decomposition
  446. !
  447. !EOP
  448. !-----------------------------------------------------------------
  449. !BOC
  450. CHARACTER(len=*), PARAMETER :: rname = mname//'/emission_apply'
  451. ! --- begin --------------------------------------
  452. ! start timing:
  453. call GO_Timer_Start( itim_appl, status )
  454. IF_NOTOK_RETURN(status=1)
  455. IF (okdebug) then
  456. WRITE(gol,*) 'start of emission_apply for region:',region ; call goPr
  457. END IF
  458. ! CO2 emissions
  459. call GO_Timer_Start( itim_co2, status )
  460. IF_NOTOK_RETURN(status=1)
  461. CALL emission_co2_apply( region, status )
  462. IF_NOTOK_RETURN(status=1)
  463. call GO_Timer_End( itim_co2, status )
  464. IF_NOTOK_RETURN(status=1)
  465. IF(okdebug) then
  466. WRITE(gol,*) 'End of adding emission '; call goPr
  467. END IF
  468. ! end timing:
  469. call GO_Timer_End( itim_appl, status )
  470. IF_NOTOK_RETURN(status=1)
  471. ! ok
  472. status = 0
  473. END SUBROUTINE EMISSION_APPLY
  474. !EOC
  475. !----------------------------------------------------------------------------
  476. ! TM5 !
  477. !----------------------------------------------------------------------------
  478. !BOP
  479. !
  480. ! !IROUTINE: FREE_EMISSION
  481. !
  482. ! !DESCRIPTION: Deallocate space needed to handle the emissions by calling
  483. ! *done methods of constituents' modules.
  484. !\\
  485. !\\
  486. ! !INTERFACE:
  487. !
  488. SUBROUTINE FREE_EMISSION( status )
  489. !
  490. ! !USES:
  491. !
  492. USE dims, ONLY : nregions
  493. !
  494. ! !OUTPUT PARAMETERS:
  495. !
  496. INTEGER, INTENT(out) :: status
  497. !
  498. ! !REVISION HISTORY:
  499. ! 16 Jul 2010 - A. Strunk - Adapted to revised emission_*.F90 routines
  500. ! 27 Mar 2012 - P. Le Sager - Adapted for lon-lat MPI domain decomposition
  501. !
  502. !EOP
  503. !---------------------------------------------------------------------------
  504. !BOC
  505. CHARACTER(len=*), PARAMETER :: rname = mname//'/free_emission'
  506. INTEGER :: region, imode
  507. ! --- begin -----------------------------------
  508. DO region = 1, nregions
  509. DEALLOCATE(plandr(region)%surf)
  510. DEALLOCATE(emis2D(region)%surf)
  511. ENDDO
  512. CALL Emission_CO2_Done( status )
  513. IF_NOTOK_RETURN(status=1)
  514. ! done
  515. status = 0
  516. END SUBROUTINE FREE_EMISSION
  517. !EOC
  518. END MODULE EMISSION