! #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if ! #include "tm5.inc" ! !----------------------------------------------------------------------------- ! TM5 ! !----------------------------------------------------------------------------- !BOP ! ! !MODULE: EMISSION_CO ! ! !DESCRIPTION: Hold data and methods for CO emissions. !\\ !\\ ! !INTERFACE: ! MODULE EMISSION_CO ! ! !USES: ! use GO, only : gol, goErr, goPr use dims, only : nregions, idate, okdebug use global_types, only : emis_data, d3_data use emission_read, only : used_providers, has_emis use tm5_distgrid, only : dgrid, get_distgrid, scatter use partools, only : isRoot, par_broadcast implicit none private ! ! !PUBLIC MEMBER FUNCTIONS: ! public :: Emission_CO_Init ! allocate dataset public :: Emission_CO_Done ! deallocate dataset public :: Emission_CO_Declare ! read monthly input public :: Emission_CO_Apply ! distribute & add emissions to tracer array ! ! !PRIVATE DATA MEMBERS: ! character(len=*), parameter :: mname = 'emission_co' type( emis_data ), dimension(:,:), allocatable :: co_emis_2d type( d3_data ), dimension(:,:), allocatable :: co_emis_3d logical, allocatable :: has_data_3d(:), has_data_2d(:) integer :: co_2dsec, co_3dsec ! ! !REVISION HISTORY: ! 1 Oct 2010 - Achim Strunk - revamped for AR5 ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4 ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! 29 Jan 2014 - Jason Williams - added yearly specific biogenic ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ CONTAINS !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: EMISSION_CO_INIT ! ! !DESCRIPTION: Allocate memory to handle the emissions. !\\ !\\ ! !INTERFACE: ! SUBROUTINE EMISSION_CO_INIT( status ) ! ! !USES: ! use dims, only : lm use emission_read, only : providers_def, numb_providers, ed42_nsect_co use emission_data, only : LAR5BMB use emission_read, only : n_ar5_ant_sec, n_ar5_shp_sec, n_ar5_air_sec, n_ar5_bmb_sec use emission_read, only : ar5_cat_ant, ar5_cat_shp, ar5_cat_air, ar5_cat_bmb ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 1 Oct 2010 - Achim Strunk - v0 ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/Emission_CO_Init' integer :: region integer :: imr, jmr, lmr, lsec, lprov, i1, i2, j1, j2 ! --- begin -------------------------------------- status = 0 if(.not. has_emis) return ! nb of sectors co_2dsec = 0 co_3dsec = 0 do lprov = 1, numb_providers if (count(used_providers.eq.providers_def(lprov)%name)/=0) then if (trim(providers_def(lprov)%name) .eq. 'AR5') then ! nb of available sectors in AR5 depends on category co_2dsec = co_2dsec + n_ar5_ant_sec*count('CO'.eq.ar5_cat_ant) + & n_ar5_shp_sec*count('CO'.eq.ar5_cat_shp) if (LAR5BMB) co_2dsec = co_2dsec + n_ar5_bmb_sec*count('CO'.eq.ar5_cat_bmb) co_3dsec = co_3dsec + n_ar5_air_sec*count('CO'.eq.ar5_cat_air) ! co_2dsec = co_2dsec + providers_def(lprov)%nsect2d ! co_3dsec = co_3dsec + count('CO'.eq.ar5_cat_air) elseif (trim(providers_def(lprov)%name) .eq. 'ED42') then co_2dsec = co_2dsec + ed42_nsect_co ! no 3d sectors in EDGAR 4.2 else co_2dsec = co_2dsec + providers_def(lprov)%nsect2d co_3dsec = co_3dsec + providers_def(lprov)%nsect3d endif endif enddo allocate( co_emis_2d( nregions, co_2dsec ) ) allocate( co_emis_3d( nregions, co_3dsec ) ) allocate( has_data_2d(co_2dsec)) ; has_data_2d=.false. allocate( has_data_3d(co_3dsec)) ; has_data_3d=.false. ! allocate information arrays (2d and 3d) do region=1,nregions CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 ) lmr = lm(region) do lsec=1,co_2dsec allocate( co_emis_2d(region,lsec)%surf(i1:i2,j1:j2) ) end do do lsec=1,co_3dsec allocate( co_emis_3d(region,lsec)%d3(i1:i2,j1:j2,lmr) ) end do enddo ! ok status = 0 END SUBROUTINE EMISSION_CO_INIT !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: EMISSION_CO_DONE ! ! !DESCRIPTION: Free memory after handling of the emissions. !\\ !\\ ! !INTERFACE: ! SUBROUTINE EMISSION_CO_DONE( status ) ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 1 Oct 2010 - Achim Strunk - v0 ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/Emission_CO_Done' integer :: region, lsec ! --- begin -------------------------------------- status = 0 if(.not. has_emis) return do region = 1, nregions do lsec=1,co_2dsec deallocate( co_emis_2d(region,lsec)%surf ) end do do lsec=1,co_3dsec deallocate( co_emis_3d(region,lsec)%d3 ) end do end do deallocate( co_emis_2d ) deallocate( co_emis_3d ) deallocate( has_data_2d, has_data_3d) ! ok status = 0 END SUBROUTINE EMISSION_CO_DONE !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: EMISSION_CO_DECLARE ! ! !DESCRIPTION: Opens, reads and evaluates input files (per month). ! Provides emissions on 2d/3d-arrays which are then added ! to mixing ratios in routine *apply. !\\ !\\ ! !INTERFACE: ! SUBROUTINE EMISSION_CO_DECLARE( status ) ! ! !USES: ! use toolbox, only : coarsen_emission use dims, only : im, jm, lm, idate, sec_month, nlon360, nlat180, iglbsfc use chem_param, only : xmco use emission_data, only : msg_emis, LAR5BMB, LMEGAN ! ---------------- AR5 - EDGAR 4 - ETC. -------------------- use emission_data, only : emis_input_year use emission_data, only : emis_input_dir_mac use emission_data, only : emis_input_dir_megan use emission_data, only : emis_input_dir_retro use emission_data, only : emis_input_dir_gfed use emission_data, only : emis_input_dir_ed4 use emission_read, only : emission_ar5_regrid_aircraft use emission_read, only : emission_ar5_ReadSector use emission_read, only : emission_macc_ReadSector use emission_read, only : emission_ed4_ReadSector use emission_read, only : emission_gfed_ReadSector use emission_read, only : emission_retro_ReadSector use emission_read, only : emission_megan_ReadSector use emission_read, only : sectors_def, numb_sectors use emission_read, only : ar5_dim_3ddata use emission_read, only : ed42_co_sectors ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 1 Oct 2010 - Achim Strunk - adapted for AR5 ! 1 Dec 2011 - Narcisa Banda - added EDGAR 4 ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! ! !REMARKS: ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/emission_co_declare' integer :: region, hasData integer, parameter :: add_field=0 integer, parameter :: amonth=2 integer :: imr, jmr, lmr, lsec ! AR5 real,dimension(:,:,:), allocatable :: field3d, field3d2 type(d3_data) :: emis3d, work(nregions) type(emis_data) :: wrk2D(nregions) integer :: seccount2d, seccount3d ! --- begin ----------------------------------------- status = 0 if(.not. has_emis) return write(gol,'(" EMISS-INFO ------------- read CO emissions -------------")'); call goPr do region = 1, nregions do lsec=1,co_2dsec co_emis_2d(region,lsec)%surf = 0.0 end do do lsec=1,co_3dsec co_emis_3d(region,lsec)%d3 = 0.0 end do end do ! global arrays for coarsening do region = 1, nregions if (isRoot)then allocate(work(region)%d3(im(region),jm(region),lm(region))) else allocate(work(region)%d3(1,1,1)) end if enddo do region = 1, nregions wrk2D(region)%surf => work(region)%d3(:,:,1) end do ! -------------------------------- ! do a loop over available sectors ! -------------------------------- ! count 2d and 3d sectors seccount2d = 0 seccount3d = 0 ! always allocate here 3d data set (for 2d sectors it will be filled in first layer only) if (isRoot) then allocate( field3d( nlon360, nlat180, ar5_dim_3ddata ) ) ; field3d = 0.0 else allocate( field3d( 1, 1, 1 ) ) end if sec : do lsec = 1, numb_sectors if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle if (associated(sectors_def(lsec)%species)) then ! AR5 check if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle endif field3d = 0.0 if( sectors_def(lsec)%f3d ) then seccount3d = seccount3d + 1 else seccount2d = seccount2d + 1 end if if (isRoot) then ! READ select case( trim(sectors_def(lsec)%prov) ) case( 'AR5' ) ! Screen out agricultural sector for CO, ! because it is zero in the RCPs ! and not present in the historical files. if (trim(sectors_def(lsec)%name) .ne. 'emiss_agr') then call emission_ar5_ReadSector( 'CO', emis_input_year, idate(2), lsec, field3d, status ) IF_NOTOK_RETURN(status=1;deallocate(field3d)) endif case( 'MACC' ) ! Screen out 'soil', 'nat', and 'air' since they are not available for CO. ! If using MEGAN skip also biogenic source. if ( ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_soil')) .and. & ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_nat') ) .and. & ( .not. (trim(sectors_def(lsec)%name) .eq. 'emiss_air') ) .and. & ( .not. (LMEGAN .and. (trim(sectors_def(lsec)%name) .eq. 'emiss_bio'))) ) then call emission_macc_ReadSector( emis_input_dir_mac, 'CO', emis_input_year, idate(2), & '0.5x0.5_kg.nc', sectors_def(lsec)%name, 'kg / s', field3d, status ) IF_NOTOK_RETURN(status=1) end if case( 'ED41' ) select case(trim(sectors_def(lsec)%name)) case ('1A3b_c_e','1A3d_SHIP','1A3d1') call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', 2005, idate(2), & lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status ) IF_NOTOK_RETURN(status=1) end select case( 'ED42' ) ! biomass burning (GFED/RETRO/AR5BMB) and transport (ED41) are excluded through ED42_CO_SECTORS definition call emission_ed4_ReadSector( emis_input_dir_ed4, 'CO', 'co', emis_input_year, idate(2), & lsec, trim(sectors_def(lsec)%prov), 'kg / s', field3d, status ) IF_NOTOK_RETURN(status=1) case('GFEDv3') call emission_gfed_ReadSector( emis_input_dir_gfed, 'co', emis_input_year, idate(2), & sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status ) IF_NOTOK_RETURN(status=1) case('RETRO') call emission_retro_ReadSector( emis_input_dir_retro, 'CO', emis_input_year, idate(2), & sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status ) IF_NOTOK_RETURN(status=1) case('MEGAN') call emission_megan_ReadSector( emis_input_dir_megan, 'CO', emis_input_year, idate(2), & sectors_def(lsec)%name, 'kg / s', field3d(:,:,1), status ) IF_NOTOK_RETURN(status=1) case('DUMMY') case default write(gol,*) "Error in buidling list of providers USED_PROVIDERS"; call goErr status=1; TRACEBACK; return end select ! nothing found??? if( sum(field3d) < 100.*TINY(1.0) ) then if (okdebug) then write(gol,'("EMISS-INFO - no CO emissions found for ",a," ",a," for month ",i2 )') & trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr endif hasData=0 else if (okdebug) then write(gol,'("EMISS-INFO - found CO emissions for ",a," ",a," for month ",i2 )') & trim(sectors_def(lsec)%prov), trim(sectors_def(lsec)%name), idate(2) ; call goPr endif ! scale from kg/s to kg/month field3d = field3d * sec_month ! kg / month hasData=1 end if end if call Par_broadcast(hasData, status) IF_NOTOK_RETURN(status=1) if (hasData == 0) then cycle sec else if ( sectors_def(lsec)%f3d ) then has_data_3d(seccount3d)=.true. else has_data_2d(seccount2d)=.true. end if end if ! Distinguish 2d/3d sectors if( sectors_def(lsec)%f3d ) then write(gol,'("EMISS-ERROR - Unexpected 3D data: Uncomment code below ")'); call goErr status=1; TRACEBACK; return ! --------------------------- ! 3d data (AIRCRAFT) ! --------------------------- ! if (isRoot) then ! REGRID ! ! helper array for regridding ! allocate( field3d2( nlon360,nlat180,lm(1) ) ) ; field3d2 = 0.0 ! ! aircraft data: regrid vertically to model layers ! call emission_ar5_regrid_aircraft( iglbsfc, field3d, nlon360, nlat180, ar5_dim_3ddata, lm(1), field3d2, status ) ! IF_NOTOK_RETURN(status=1;deallocate(field3d,field3d2)) ! ! write some numbers ! call msg_emis( amonth, trim(sectors_def(lsec)%prov),sectors_def(lsec)%name, 'CO', xmco, sum(field3d2) ) ! ! distribute to co_emis in regions ! call Coarsen_Emission( 'CO '//trim(sectors_def(lsec)%name), nlon360, nlat180, lm(1), & ! field3d2, work, add_field, status ) ! IF_NOTOK_RETURN(status=1;deallocate(field3d,field3d2)) ! deallocate( field3d2 ) ! end if ! do region = 1, nregions ! call scatter(dgrid(region), co_emis_3d(region,seccount3d)%d3, work(region)%d3, 0, status) ! IF_NOTOK_RETURN(status=1) ! end do else ! ar5_sector is 2d ! --------------------------- ! 2d data (Anthropogenic, Ships, Biomassburning) ! --------------------------- if (isRoot) then ! print total & regrid call msg_emis( amonth, trim(sectors_def(lsec)%prov), sectors_def(lsec)%name, 'CO', xmco, sum(field3d(:,:,1)) ) call coarsen_emission( 'CO '//sectors_def(lsec)%name, & nlon360, nlat180, field3d(:,:,1), wrk2D, add_field, status ) IF_NOTOK_RETURN(status=1;deallocate(field3d)) end if do region = 1, nregions call scatter(dgrid(region), co_emis_2d(region,seccount2d)%surf, work(region)%d3(:,:,1), 0, status) IF_NOTOK_RETURN(status=1) end do end if end do sec ! sectors deallocate( field3d ) do region = 1, nregions if (associated(wrk2D(region)%surf)) nullify(wrk2D(region)%surf) deallocate( work(region)%d3 ) end do ! check sectors found if( seccount2d /= co_2dsec ) then write(gol,'(80("-"))') ; call goPr write(gol,'("ERROR: 2d sectors do not equal total number:",i4," /= ",i4," !")') seccount2d, co_2dsec ; call goErr write(gol,'(80("-"))') ; call goPr status=1; return end if if( seccount3d /= co_3dsec ) then write(gol,'(80("-"))') ; call goPr write(gol,'("ERROR: 3d sectors do not equal total number:",i4," /= ",i4," !")') seccount3d, co_3dsec ; call goErr write(gol,'(80("-"))') ; call goPr status=1; return end if status = 0 END SUBROUTINE EMISSION_CO_DECLARE !EOC !-------------------------------------------------------------------------- ! TM5 ! !-------------------------------------------------------------------------- !BOP ! ! !IROUTINE: EMISSION_CO_APPLY ! ! !DESCRIPTION: Take monthly emissions, and ! - split them vertically ! - apply time splitting factors ! - add them to tracers (add_3d) !\\ !\\ ! !INTERFACE: ! SUBROUTINE EMISSION_CO_APPLY( region, status ) ! ! !USES: ! use dims, only : itaur, nsrce, tref use dims, only : im, jm, lm use datetime, only : tau2date use emission_data, only : emission_vdist_by_sector, LAR5BMB use emission_data, only : do_add_3d, do_add_3d_cycle, bb_cycle use emission_data, only : emis_bb_trop_cycle use chem_param, only : ico, xmco use emission_read, only : sectors_def, numb_sectors use emission_read, only : ed42_co_sectors ! ! !INPUT PARAMETERS: ! integer, intent(in) :: region ! ! !OUTPUT PARAMETERS: ! integer, intent(out) :: status ! ! !REVISION HISTORY: ! 1 Oct 2010 - Achim Strunk - AR5 ! 26 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition ! !EOP !------------------------------------------------------------------------ !BOC character(len=*), parameter :: rname = mname//'/emission_co_apply' integer, dimension(6) :: idater real :: dtime, fraction integer :: imr, jmr, lmr, lsec, i1, i2, j1, j2 integer :: seccount2d, seccount3d type(d3_data) :: emis3d ! --- begin ----------------------------- status = 0 if(.not. has_emis) return if( okdebug ) then write(gol,*) 'start of emission_co_apply'; call goPr end if call tau2date(itaur(region),idater) dtime=float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX. if(okdebug) then write(gol,*)'emission_co_apply in region ',region,' at date: ',idater, ' with time step:', dtime call goPr end if ! get a working structure for 3d emissions call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 ) allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0 ! count 2d and 3d sectors seccount2d = 0 seccount3d = 0 ! cycle over sectors do lsec = 1, numb_sectors if (count(used_providers.eq.sectors_def(lsec)%prov).eq.0) cycle if ((trim(sectors_def(lsec)%prov).eq.'ED42') .and. (count(ed42_co_sectors.eq.sectors_def(lsec)%name) .eq. 0)) cycle if (associated(sectors_def(lsec)%species)) then ! AR5 check if (count('CO'.eq.sectors_def(lsec)%species).eq.0) cycle if ((trim(sectors_def(lsec)%catname) .eq. 'biomassburning').and.(.not.LAR5BMB)) cycle endif ! default: no additional splitting fraction = 1.0 ! ---------------------------------------------------------------------------------------- ! distinguish here between sectors and whether they should have additional splitting ! if( sectors_def(lsec)%catname == 'biomassburning' ) fraction = fraction * bb_frac etc... ! ---------------------------------------------------------------------------------------- ! distinguish between 2d/3d sectors if( sectors_def(lsec)%f3d ) then seccount3d = seccount3d + 1 if (.not.has_data_3d(seccount3d)) cycle emis3d%d3 = co_emis_3d(region,seccount3d)%d3 else seccount2d = seccount2d + 1 if (.not.has_data_2d(seccount2d)) cycle emis3d%d3 = 0.0 ! vertically distribute according to sector call emission_vdist_by_sector( sectors_def(lsec)%vdisttype, 'CO', region, co_emis_2d(region,seccount2d), emis3d, status ) IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3)) end if ! add dataset according to sector and category if( emis_bb_trop_cycle .and. trim(sectors_def(lsec)%catname) == "biomassburning" ) then call do_add_3d_cycle( region, ico, i1, j1, emis3d%d3, bb_cycle(region)%scalef, xmco, xmco, status, fraction ) IF_NOTOK_RETURN(status=1) else call do_add_3d( region, ico, i1, j1, emis3d%d3, xmco, xmco, status, fraction ) IF_NOTOK_RETURN(status=1) end if end do deallocate( emis3d%d3 ) if(okdebug) then write(gol,*) 'end of emission_co_apply'; call goPr end if ! OK status = 0 END SUBROUTINE EMISSION_CO_APPLY !EOC !NOTYET !-------------------------------------------------------------------------- !NOTYET ! TM5 ! !NOTYET !-------------------------------------------------------------------------- !NOTYET !BOP !NOTYET ! !NOTYET ! !IROUTINE: DISTRIBUTE_HIGH_LAT_CO !NOTYET ! !NOTYET ! !DESCRIPTION: Attribute seasonal cycle to the hight latitude CO industrial !NOTYET ! /fossil fuel CO emissions. !NOTYET ! This simplified distribution function should be replaced if !NOTYET ! more detailed info becomes available. !NOTYET !\\ !NOTYET !\\ !NOTYET ! !INTERFACE: !NOTYET ! !NOTYET SUBROUTINE DISTRIBUTE_HIGH_LAT_CO( imdim, jmdim, field2d ) !NOTYET ! !NOTYET ! !USES: !NOTYET ! !NOTYET use dims, only : idate, sec_month, sec_day, sec_year, mlen !NOTYET ! !NOTYET ! !INPUT PARAMETERS: !NOTYET ! !NOTYET integer, intent(in) :: jmdim, imdim ! dimension of grid !NOTYET ! !NOTYET ! !INPUT/OUTPUT PARAMETERS: !NOTYET ! !NOTYET real, dimension(imdim, jmdim) :: field2d !NOTYET ! !NOTYET ! !REVISION HISTORY: !NOTYET ! 24 Jun 2002 - fjd !NOTYET ! 14 Feb 2013 - Ph. Le Sager - reincorporated for ED4 !NOTYET ! !NOTYET ! !REMARKS: !NOTYET ! !NOTYET !EOP !NOTYET !------------------------------------------------------------------------ !NOTYET !BOC !NOTYET !NOTYET real :: year2month !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/) !NOTYET !NOTYET ! high latitude seasonal CO distribution function !NOTYET integer :: mo,j !NOTYET !NOTYET year2month =sec_month/sec_year !scale factor for yearly total to specific month !NOTYET mo=idate(2) !NOTYET !NOTYET do j=1,jmdim !NOTYET if (j>jmdim*3/4) then !NOTYET field2d(:,j)=field2d(:,j)*coseas(mo)/12. ! convert to month !NOTYET else !NOTYET field2d(:,j)=field2d(:,j)*year2month !NOTYET endif !NOTYET enddo !NOTYET !NOTYET END SUBROUTINE DISTRIBUTE_HIGH_LAT_CO !NOTYET !EOC END MODULE EMISSION_CO