Emissions Module
AEIC.emissions.emission.Emission is the module that uses
AEIC.performance_model.PerformanceModel and a flown
AEIC.trajectories.trajectory.Trajectory to compute emissions for the
entire mission. It layers multiple methods for emission calculations
from user choices in the configuration file.
Overview
Computes trajectory, LTO, APU, GSE, and life-cycle \(\mathrm{CO_2}\).
Emits structured arrays (grams by species) plus convenience containers (
EmissionSliceandEmissionsOutput) for downstream analysis.Has helper methods such as
Emission.emit_trajectory()orEmission.emit_lto()when only a subset is needed.
Configuration Inputs
The [Emissions] section of the configuration TOML file is validated through
EmissionsConfig. Keys and meanings are summarised below.
Key |
Allowed values |
Description |
|---|---|---|
|
any fuel name matching |
Selects the fuel file used for EIs and life-cycle data. |
|
|
When |
|
|
Toggles calculation of fuel-dependent, constant EI species. |
|
|
Selects the method for \(\mathrm{NO_x}\) calculation (None disables calculation). |
|
|
Selects the method for HC/CO calculation (None disables calculation). |
|
|
Chooses the PMvol method. |
|
|
Chooses the PMnvol method. |
|
|
Pulls LTO EI/fuel-flow data from the performance tables or EDB file. |
|
path |
Required when |
|
|
Enables non-trajectory emission sources and life-cycle \(\mathrm{CO_2}\) adjustments. |
Usage Example
from AEIC.performance_model import PerformanceModel
from AEIC.trajectories.trajectory import Trajectory
from AEIC.emissions.emission import Emission
perf = PerformanceModel("IO/default_config.toml")
mission = perf.missions[0]
traj = Trajectory(perf, mission, optimize_traj=True, iterate_mass=False)
traj.fly_flight()
emitter = Emission(perf)
output = emitter.emit(traj)
print("Total CO2 (g)", output.total['CO2'])
print("Taxi NOx (g)", output.lto.emissions_g['NOx'][0])
# Need only the trajectory segment
segments = emitter.emit_trajectory(traj)
print("Per-segment PM number", segments.emissions_g['PMnvol'])
Inner Containers
The module defines dataclasses that document both inputs and outputs of the computation:
EmissionsConfig: user-facing configuration parsed from the TOML file. It validates enums (LTOInputMode,EINOxMethod,PMvolMethod,PMnvolMethod), resolves defaults, and ensures databank paths are present when required.EmissionSettings: flattened, runtime-only view of the above. It keeps booleans for metric flags, file paths, and LTO/auxiliary toggles so subsequent runs avoid re-validating the original mapping.AtmosphericState: carries temperature, pressure, and Mach arrays that emission-index models reuse when HC/CO/\(\text{NO}_x\)/PM need ambient conditions.EmissionSlice: describes any source (trajectory, LTO, APU, GSE). It storesindices(emission indices in g/kg) and the realizedemissions_g.TrajectoryEmissionSlice: extendsEmissionSlicewithfuel_burn_per_segment(kg) andtotal_fuel_burn(kg) so users can derive intensity metrics.EmissionsOutput: top-level container returned byEmission.emit(). It exposestrajectory,lto,apu,gse,total(summed structured array), and optionallifecycle_co2_g.
Computation Workflow
The Emission object is instanced once per performance model:
EmissionsConfigis materialized fromPerformanceModel.config.emissionsand converted toEmissionSettings.Fuel properties are read from
fuels/<Fuel>.toml. These provide \(\mathrm{CO_2}\)/\(\mathrm{H_2O}\)/\(\mathrm{SO_x}\) emission indices, and life-cycle factors.emit(traj)resets internal arrays sized to the trajectory stepsEmission.get_trajectory_emissions()computes EI values for each mission point:Constant EI species (\(\mathrm{CO_2}\), \(\mathrm{H_2O}\), \(\mathrm{SO}_x`\)).
Methods for HC/CO/\(\mathrm{NO_x}\)/PMvol/PMnvol applied according to user specification.
Emission.get_LTO_emissions()builds the ICAO style landing and take off emissions using either databank values (LTO_input_mode = "edb") or the per-mode inputs embedded in the performance file.AEIC.emissions.APU_emissions.get_APU_emissions()andEmission.get_GSE_emissions()contributions are added if enabled.Emission.sum_total_emissions()aggregates each pollutant intoself.summed_emission_gand, when requested, life-cycle \(\mathrm{CO_2}\) is appended viaEmission.get_lifecycle_emissions().
Structured Arrays
All emission indices and gram totals share the dtype emitted by the private
__emission_dtype helper. Each field is float64:
CO2, H2O, HC, CO, NOx, NO, NO2, HONO,
PMnvol, PMnvolGMD, PMvol, OCic, SO2, SO4.
If EI_PMnvol_method is SCOPE11 or MEEM, the additional PMnvolN
field is emitted. Metric-specific flags (see Emission.metric_flags) determine
which fields are populated; disabled species stay as 0, making it easy to filter downstream.
API Reference
- class AEIC.emissions.config.EmissionsConfig(*, fuel: str, climb_descent_usage: bool = True, co2_enabled: bool = True, h2o_enabled: bool = True, sox_enabled: bool = True, nox_method: EINOxMethod = EINOxMethod.BFFM2, hc_method: EINOxMethod = EINOxMethod.BFFM2, co_method: EINOxMethod = EINOxMethod.BFFM2, pmvol_method: PMvolMethod = PMvolMethod.FUEL_FLOW, pmnvol_method: PMnvolMethod = PMnvolMethod.MEEM, apu_enabled: bool = True, gse_enabled: bool = True, lifecycle_enabled: bool = True)
Configuration data for emissions module.
- DEFAULT_METHOD: ClassVar[EINOxMethod] = 'bffm2'
Default method for NOx, HC, and CO emissions calculations.
- apu_enabled: bool
APU emission calculation flag.
- climb_descent_usage: bool
Flag controlling flight phases for which emissions are calculated. If true, emissions are calculated for the entire trajectory (takeoff, climb, cruise, descent); if false, emissions are only calculated for cruise and LTO data is used for takeoff, climb and approach.
- co2_enabled: bool
CO2 emission calculation flag.
- property co_enabled: bool
CO emission calculation flag.
- co_method: EINOxMethod
CO emission calculation method.
- fuel: str
Fuel used (conventional Jet-A, SAF, etc.).
- property fuel_file: str
Fuel file path.
- gse_enabled: bool
GSE emission calculation flag.
- h2o_enabled: bool
H2O emission calculation flag.
- property hc_enabled: bool
HC emission calculation flag.
- hc_method: EINOxMethod
HC emission calculation method.
- lifecycle_enabled: bool
Lifecycle emission calculation flag.
- model_config: ClassVar[ConfigDict] = {'frozen': True}
Configuration is frozen after creation.
- property nox_enabled: bool
NOx emission calculation flag.
- nox_method: EINOxMethod
NOx emission calculation method.
- property pmnvol_enabled: bool
PMnvol emission calculation flag.
- pmnvol_method: PMnvolMethod
PMnvol emission calculation method.
- property pmvol_enabled: bool
PMvol emission calculation flag.
- pmvol_method: PMvolMethod
PMvol emission calculation method.
- sox_enabled: bool
SOx emission calculation flag.
- class AEIC.emissions.emission.Emission(ac_performance: PerformanceModel)
Bases:
objectStateless emissions calculator that can be reused across trajectories. Configure it with a
PerformanceModelonce, then callemitwith each trajectory to obtain anEmissionsOutputdata container.- __init__(ac_performance: PerformanceModel)
- Parameters:
ac_performance (PerformanceModel) – Aircraft performance object containing climb/cruise/descent performance and LTO data.
- emit(trajectory: Trajectory) EmissionsOutput
Compute emissions for a single trajectory and return structured results.
- emit_apu(trajectory: Trajectory) EmissionSlice | None
Convenience helper returning only the APU emissions (if enabled).
- emit_gse(trajectory: Trajectory) EmissionSlice | None
Convenience helper returning only the GSE emissions (if enabled).
- emit_lto(trajectory: Trajectory) EmissionSlice
Convenience helper returning only the LTO emissions.
- emit_trajectory(trajectory: Trajectory) TrajectoryEmissionSlice
Convenience helper returning only the trajectory portion of emissions.
- class AEIC.emissions.emission.EmissionSlice(indices: ndarray | None, emissions_g: ndarray)
Emission data for a single source such as trajectory, LTO, APU, or GSE.
indices identifies the contributing mission segments (or None for aggregate-only sources) and emissions_g is the structured per-species array of integrated masses in grams aligned with those indices.
- class AEIC.emissions.emission.TrajectoryEmissionSlice(indices: ndarray | None, emissions_g: ndarray, fuel_burn_per_segment: ndarray, total_fuel_burn: float)
EmissionSlice specialized for total trajectory contribution.
Adds per-segment fuel burn (kg) and the total fuel burned so downstream analysis can compute intensity metrics or lifecycle adjustments alongside the species masses.
- class AEIC.emissions.emission.EmissionsOutput(trajectory: TrajectoryEmissionSlice, lto: EmissionSlice, apu: EmissionSlice | None, gse: EmissionSlice | None, total: ndarray, lifecycle_co2_g: float | None = None)
Data container that stores aggregated emissions for a mission broken down by flight phase/components.
Helper Functions
- AEIC.emissions.APU_emissions.get_APU_emissions(APU_emission_indices, APU_emissions_g, LTO_emission_indices, APU_data, LTO_noProp, LTO_no2Prop, LTO_honoProp, EI_H2O, nvpm_method: PMnvolMethod = PMnvolMethod.MEEM, apu_tim=900)
Calculate APU emissions using time in modes and given APU data.
- Parameters:
APU_emission_indices (ndarray) – self.APU_emission_indices from Emissions class
APU_emissions_g (ndarray) – self.APU_emissions_g from Emissions class
LTO_emission_indices (ndarray) – self.LTO_emission_indices from Emissions class
APU_data (dict) – dictionary containing fuel flows and EIs of chosen APU
LTO_noProp (float) – NOx speciation elements from LTO analysis
LTO_no2Prop (float) – NOx speciation elements from LTO analysis
LTO_honoProp (float) – NOx speciation elements from LTO analysis
apu_time (float) – Time in mode for APU; default value = 900 seconds (Stettler et al. 2011)
- Returns:
APU_emission_indices (ndarray) – Emissions indicies for APU
APU_emissions_g (ndarray) – Emissions in g for APU
apu_fuel_burn (float) – kg of fuel burnt by APU
- class AEIC.emissions.EI_CO2.CO2EmissionResult(EI_CO2: float, nvolCarbCont: float)
Named return values for CO₂ EI calculations.
- AEIC.emissions.EI_CO2.EI_CO2(fuel: Mapping[str, Any]) CO2EmissionResult
Calculate carbon-balanced CO2 emissions index (EI).
- Parameters:
fuel (Mapping[str, Any]) – Fuel information (input from toml file)
- Returns:
Structured CO2 emissions metadata.
- Return type:
- AEIC.emissions.EI_H2O.EI_H2O(fuel)
Calculate H2O emissions index (EI).
- Parameters:
fuel (dictionary) – Fuel information (input from toml file)
- Returns:
H2O_EI – H2O emissions index [g/kg fuel]
- Return type:
float
- AEIC.emissions.EI_SOx.EI_SOx(fuel: Mapping[str, Any]) SOxEmissionResult
Calculate universal SOx emissions indices (SO2EI and SO4EI).
- Parameters:
fuel (Mapping[str, Any]) – Fuel information (input from toml file)
- Returns:
Structured SO2/SO4 emissions indices [g/kg fuel]
- Return type:
- class AEIC.emissions.EI_SOx.SOxEmissionResult(EI_SO2: float, EI_SO4: float)
Structured SOx emission indices.
- AEIC.emissions.EI_HCCO.EI_HCCO(fuelflow_evaluate: ndarray, x_EI_matrix: ndarray, fuelflow_calibrate: ndarray, Tamb: ndarray = array(1.5e-323), Pamb: ndarray = array(4.65066429e-310), cruiseCalc: bool = False) ndarray
BFFM2 bilinear HC/CO fit to SLS data
- Parameters:
fuelflow_evaluate (ndarray, shape (n_points,)) – Fuel flows [kg/s] at which to compute xEI. Must be 1D.
x_EI_matrix (ndarray, shape (4,)) – Baseline emission indices [g x / kg fuel] at four calibration fuel‐flow points.
fuelflow_calibrate (ndarray, shape (4,)) – Calibration fuel flows [kg/s] corresponding to x_EI_matrix
cruiseCalc (bool) – If True, apply cruise correction (ambient T and P) to the final xEI.
Tamb (ndarray, shape (n_points,)) – Ambient temperature [K] for cruise correction (if cruiseCalc is True).
Pamb (ndarray, shape (n_points,)) – Ambient pressure [Pa] for cruise correction (if cruiseCalc is True).
- Returns:
xEI – The HC+CO emission index [g x / kg fuel] at each fuelflow_evaluate.
- Return type:
ndarray, shape (n_points,)
- class AEIC.emissions.EI_NOx.BFFM2EINOxResult(NOxEI: ndarray, NOEI: ndarray, NO2EI: ndarray, HONOEI: ndarray, noProp: ndarray, no2Prop: ndarray, honoProp: ndarray)
Bundled NOx emissions indices and speciation data.
- AEIC.emissions.EI_NOx.BFFM2_EINOx(sls_equiv_fuel_flow: ndarray, NOX_EI_matrix: ndarray, fuelflow_performance: ndarray, Tamb: ndarray, Pamb: ndarray, cruiseCalc: bool = True) BFFM2EINOxResult
Calculate NOx, NO, NO2, and HONO emission indices All inputs are 1-dimensional arrays of equal length for calibration (fuelflow_KGperS vs. NOX_EI_matrix) and 1-dimensional for SLS_equivalent_fuel_flow (multiple evaluation points).
- Parameters:
fuelflow_trajectory (ndarray, shape (n_times,)) – Fuel flow at which to compute EI.
NOX_EI_matrix (ndarray, shape (n_cal,)) – Baseline NOx EI values [g NOx / kg fuel] corresponding to calibration fuel flows.
fuelflow_performance (ndarray, shape (n_cal,)) – Calibration fuel flow values [kg/s] for which NOX_EI_matrix is defined.
cruiseCalc (bool) – If True, apply cruise ambient corrections (temperature and pressure) to NOx EI.
Tamb (float) – Ambient temperature [K].
Pamb (float) – Ambient pressure [Pa].
- Returns:
Structured NOx EI arrays and speciation fractions.
- Return type:
- AEIC.emissions.EI_PMnvol.PMnvol_MEEM(EDB_data: dict, altitudes: ndarray, Tamb_cruise: ndarray, Pamb_cruise: ndarray, machFlight: ndarray)
Estimate non-volatile particulate matter (nvPM) emissions at cruise using the Mission Emissions Estimation Methodology (MEEM) based on Ahrens et al. (2022), SCOPE11, and the methodology of Peck et al. (2013).
This function computes: - Geometric mean diameter (GMD) of emitted particles - Mass-based emissions index (EI) in g/kg of fuel - Number-based emissions index (EI) in #/kg of fuel
- Parameters:
EDB_data (dict) – EDB data containing engine type, bypass ratio, pressure ratio, smoke number (SN) matrix, and emission indices (mass and number).
altitudes (ndarray) – Array of flight altitudes [m] over the mission trajectory.
Tamb_cruise (ndarray) – Ambient temperature [K] at each point in the trajectory.
Pamb_cruise (ndarray) – Ambient pressure [Pa] at each point in the trajectory.
machFlight (ndarray) – Mach number at each point in the trajectory.
- Returns:
EI_PMnvol_GMD (ndarray) – Geometric mean diameter of particles [nm], constant along the trajectory.
EI_PMnvol (ndarray) – Emissions index of non-volatile PM mass [g/kg fuel] along the trajectory.
EI_PMnvolN (ndarray) – Emissions index of non-volatile PM number [#/kg fuel] along the trajectory.
Notes
If nvPM_mass_matrix or nvPM_num_matrix is undefined or negative in the EDB_data, this function reconstructs the values using the SN matrix and correlations from the literature.
Adjustments for altitude and in-flight thermodynamic conditions are made using combustor inlet temperature and pressure estimates derived from ambient conditions and engine pressure ratio.
Interpolated values account for max thrust EI values where provided.
Results with invalid SN or negative EI are set to zero with a warning.
- AEIC.emissions.EI_PMnvol.calculate_PMnvolEI_scope11(SN_matrix, PR, ENGINE_TYPE, BP_Ratio)
Calculate PM non-volatile Emission Index (EI) using SCOPE11 methodology.
- Parameters:
SN_matrix (ndarray (n x 4)) – Smoke number matrix for each engine and ICAO mode.
PR (ndarray (n x 4)) – Pressure ratio matrix.
ENGINE_TYPE (list of str) – Engine types (‘TF’, ‘MTF’, etc.) for each engine.
BP_Ratio (ndarray (n,)) – Bypass ratio for each engine.
- Returns:
PMnvolEI_best_ICAOthrust – Emission index of non-volatile PM mass [g/kg_fuel] for each engine, including 0% thrust extrapolation as first column.
- Return type:
ndarray (n x 5)
- AEIC.emissions.EI_PMvol.EI_PMvol_FOA3(thrusts: ndarray, HCEI: ndarray)
Calculate volatile organic PM emissions index (PMvoloEI) and OC internal EI (OCicEI) using the FOA3.0 method (Wayson et al., 2009).
- Parameters:
thrusts (ndarray, shape (n_types, n_times)) – ICAO thrust settings (%) for each mode and time.
HCEI (ndarray, shape (n_types, n_times)) – Hydrocarbon emissions index [g/kg fuel] for each mode and time.
- Returns:
PMvoloEI (ndarray, shape (n_types, n_times)) – Emissions index for volatile organic PM [g/kg fuel].
OCicEI (ndarray, shape (n_types, n_times)) – Same as PMvoloEI (internal organic carbon component).
- AEIC.emissions.EI_PMvol.EI_PMvol_FuelFlow(fuelflow: ndarray, thrustCat: ndarray)
Calculate EI(PMvolo) and OCicEI based on fuel flow
- Parameters:
fuelflow (ndarray, shape (n_types, 4)) – Fuel flow factor per type and 4 thrust modes.
- Returns:
PMvoloEI (ndarray, shape (n_types, 4)) – Emissions index for volatile organic PM [g/kg fuel].
OCicEI (ndarray, shape (n_types, 4)) – Emissions index for organic carbon internal [g/kg fuel].
- AEIC.emissions.lifecycle_CO2.lifecycle_CO2(fuel, fuel_burn)
Calculate lifecycle CO2 emissions.
- Parameters:
fuel (dictionary) – Fuel information (input from toml file)
- Returns:
lifecycle_CO2
- Return type:
float