openairclim.calc_cont

Calculates the contrail response.

openairclim.calc_cont.apply_wingspan_correction(config: dict[str, Any], rf_arr: Iterable[float], ac: str) ndarray[source]

Apply wingspan correction to an array of RF values. The function is applied against a reference wingspan of 35 m. References: Bruder et al. (2025) and Megill et al. (2025).

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • rf_arr (Iterable[float]) – RF values

  • ac (str) – Aircraft identifier from config.

Raises:
  • KeyError – If wingspan “b” not defined for aircraft “ac” in config.

  • ValueError – If wingspan “b” is outside of valid range [20 m, 80 m].

Returns:

RF values with wingspan correction

Return type:

np.ndarray

openairclim.calc_cont.calc_cccov_alltau(cfdd_dict: dict[int, ndarray], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]) dict[int, ndarray][source]

Calculate contrail cirrus coverage (all tau) using the Megill et al. (2025) method.

Parameters:
  • cfdd_dict (dict[int, np.ndarray]) – Dictionary with 1D (lon) CFDD values [km/km2]. Keys are inventory years.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Returns:

Dictionary with 1D (lon) cccov (all tau) values.

Keys are inventory years

Return type:

dict[int, np.ndarray]

openairclim.calc_cont.calc_cccov_taup05(config: dict[str, Any], cccov_dict: dict[int, ndarray], ac: str) dict[int, ndarray][source]

Convert contrail cirrus coverage (all tau) to optically thick contrail cirrus coverage (tau > 0.05).

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • cccov_dict (dict[int, np.ndarray]) – Dictionary with 1D (lon) contrail cirrus coverage (all optical thicknesses). Keys are inventory years.

  • ac (str) – Aircraft identifier from config.

Raises:
  • KeyError – If “PMrel” not defined in config.

  • KeyError – If “low_soot_case” not defined in config (OpenAirClim will default to “case_mid”).

Returns:

Dictionary with 1D (lon) cccov (tau > 0.05)

values. Keys are inventory years.

Return type:

dict[int, np.ndarray]

openairclim.calc_cont.calc_cfdd(config: dict[str, Any], inv_dict: dict[int, Dataset], ds_cont: Dataset, cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]], ac: str) dict[int, ndarray][source]

Calculate the Contrail Flight Distance Density (CFDD) for each year in inv_dict. This function uses the p_pcf data calculated using ERA5 (Megill & Grewe, 2025).

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • inv_dict (dict[int, xr.Dataset]) – Dictionary of emission inventory xarray datasets. Keys are inventory years.

  • ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

  • ac (str) – Aircraft identifier from config.

Returns:

Dictionary with CFDD values [km/km2], keys are

inventory years

Return type:

dict[int, np.ndarray]

openairclim.calc_cont.calc_cont_grid_areas(lat: ndarray, lon: ndarray) ndarray[source]

Calculate the cell area of the contrail grid using a simplified method.

Parameters:
  • lat (np.ndarray) – Latitudes of the grid cells [deg].

  • lon (np.ndarray) – Longitudes of the grid cells [deg].

Returns:

Contrail grid cell areas as a function of latitude [km^2].

Return type:

np.ndarray

openairclim.calc_cont.calc_cont_rf(cccov_dict: dict[int, ndarray], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]) dict[int, ndarray][source]

Calculate contrail Radiative Forcing (RF) using the Megill et al. (2025) method.

Parameters:
  • cccov_dict (dict[int, np.ndarray]) – Dictionary with 1D (lon) contrail cirrus coverage (tau > 0.05; shape: lon), keys are inventory years

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Returns:

Dictionary with contrail RF values for all

inventory years.

Return type:

dict[int, np.ndarray]

openairclim.calc_cont.calc_contrails(ac_lst: Sequence[str], config: dict[str, Any], inv_dict: dict[int, Dataset], full_inv_dict: dict[str, dict[int, Dataset]], ds_cont: Dataset) dict[str, ndarray][source]

Contrail calculation loop for a main OpenAirClim run.

Parameters:
  • ac_lst (Sequence[str]) – List of aircraft identifiers to be used.

  • config (dict[str, Any]) – Configuration dictionary from config file.

  • inv_dict (dict[int, xr.Dataset]) – Dictionary of emission inventory xarray datasets. Keys are inventory years.

  • full_inv_dict (dict[str, dict[int, xr.Dataset]]) – Nested dictionary of emission inventories. Keys are aircraft identifier, followed by year.

  • ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

Returns:

Dictionary of RF values over time for each

aircraft identifier.

Return type:

dict[str, np.ndarray]

openairclim.calc_cont.calc_ppcf(config: dict[str, Any], ds_cont: Dataset, ac: str) DataArray[source]

Calculate Potential Persistent Contrail Formation (p_PCF) using the precalculated contrail data from the Limiting Factors study (Megill & Grewe, 2025; default).

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

  • ac (str) – Aircraft identifier from config.

Raises:
  • KeyError – If key “formation_method” is not defined.

  • ValueError – If “formation_method” key is unknown.

Returns:

Interpolated p_PCF on precalculated contrail data grid

Return type:

xr.DataArray

openairclim.calc_cont.calc_ppcf_megill(config: dict[str, Any], ds_cont: Dataset, ac: str) DataArray[source]

Calculate the Potential Persistent Contrail Formation (p_PCF) using the Megill & Grewe (2025) method and precalculated data from ERA5.

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

  • ac (str) – Aircraft identifier from config.

Returns:

Interpolated p_PCF on precalculated contrail data grid.

Return type:

xr.DataArray

openairclim.calc_cont.calc_sac_slope(p: float, sac_eq: str, q_h: float, eta: float | None = None, eta_elec: float | None = None, ei_h2o: float | None = None, r: float | None = None) float[source]

Calculates the slope of the SAC mixing line.

Parameters:
  • p (float) – Ambient pressure [Pa]

  • sac_eq (str) – SAC equation to be used, one of “CON”, “HYB”, “H2C”, “H2FC”

  • q_h (float) – Lower heating value of the fuel [J/kg] for “CON”, “HYB” and “H2C”; formation enthalpy of water vapour [J/mol] for “H2FC”

  • eta (float, optional) – Overall propulsion efficiency of the liquid fuel system [-]

  • eta_elec (float, optional) – Overall propulsion efficiency of the electric/fuel cell system [-]

  • ei_h2o (float, optional) – Emission index of water vapour [kg/kg]

  • r (float, optional) – Degree of hybridisation. r=1 is pure liquid fuel operation; r=0 pure electric operation

Raises:
Returns:

Slope of the SAC mixing line [Pa/K].

Return type:

float

openairclim.calc_cont.calc_total_over_ac(data: dict[str, dict[int, Any]], ac_lst: Iterable[str]) dict[str, dict[int, Any]][source]

Add a “TOTAL” entry to data by summing the per-year values across all aircraft identifiers.

Parameters:
  • data (dict[str, dict[int, Any]]) – Nested dictionary with keys ac, then yr

  • ac_lst (Iterable[str]) – List of aircraft identifiers to be used

Returns:

Shallow copy of data with TOTAL added or

overwritten

Return type:

dict[str, dict[int, Any]]

openairclim.calc_cont.cfdd_to_1d(cfdd_dict: dict[str, dict[int, ndarray]], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]) dict[str, dict[int, ndarray]][source]

Convert 3D CFDD to 1D (lon axis) CFDD to match contrail cirrus coverage using a vertical sum and area-weighting to remove latitude-dependence.

Parameters:
  • cfdd_dict (dict[str, dict[int, np.ndarray]]) – Dictionary with CFDD values [km/km2] in 3D (plev, lat, lon). Keys are inventory years.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Returns:

Dictionary with CFDD values in 1D (lon).

Return type:

dict[str, dict[int, np.ndarray]]

openairclim.calc_cont.check_cont_input(ds_cont: Dataset) None[source]

Checks the input data for the contrail module.

Parameters:

ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

openairclim.calc_cont.check_plev_range(inv_dict: dict[int, Dataset], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]], clamp: bool = True) dict[int, Dataset][source]

Checks whether all pressure level values in inv_dict are within the bounds of the pre-calculated contrail grid. Logs a warning if any values are found and automatically clamps values into the allowed range.

Parameters:
  • inv_dict (dict[int, xr.Dataset]) – Dictionary of emission inventory xarray datasets. Keys are inventory years.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

  • clamp (bool, optional) – Whether the values should be clamped to within the allowed range. Defaults to True.

Returns:

Dictionary of emission inventory xarray datasets

clamped to within the allowed plev range.

Return type:

dict[int, xr.Dataset]

openairclim.calc_cont.contrail_attribution(input_dict: dict[int, ndarray], ac_dict: dict[int, ndarray], total_dict: dict[int, ndarray]) dict[int, ndarray][source]

Use proportional attribution to split the input into the contribution from ac_dict (single aircraft identifier). The keys of all inputs must match.

Parameters:
  • input_dict (dict[int, np.ndarray]) – Dictionary to be split, keys are inventory years.

  • ac_dict (dict[int, np.ndarray]) – Dictionary with values for a single aircraft identifier, could be CFDD or coverage values for example. Keys are inventory years.

  • total_dict (dict[int, np.ndarray]) – Dictionary with values for all aircraft identifiers, could be CFDD or coverage values for example (must match ac_dict). Keys are inventory years

Returns:

Dictionary with proportionally attributed values.

Keys are years.

Return type:

dict[int, np.ndarray]

openairclim.calc_cont.get_cont_grid(ds_cont: Dataset) tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]][source]

Get contrail grid from ds_cont.

Parameters:

ds_cont (xr.Dataset) – Dataset of precalculated contrail data.

Returns:

Tuple (lon, lat, plev); each is 1-D float array with

shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Return type:

ContGrid

openairclim.calc_cont.interp_base_inv_dict(inv_yrs: Sequence[int], base_inv_dict: dict[int, Dataset], intrp_vars: Iterable[str], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]) dict[int, Dataset][source]

Create base emission inventories for years in inv_yrs that do not exist in base_inv_dict.

Parameters:
  • inv_yrs (Sequence[int]) – Dictionary of emission inventory xarrays, keys are inventory years.

  • base_inv_dict (dict[int, xr.Dataset]) – Dictionary of base emission inventory xarrays. Keys are inventory years.

  • intrp_vars (Iterable[str]) – List of strings of data variables in base_inv_dict that are to be included in the missing base inventories, e.g. [“distance”] (for contrail calculations).

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Returns:

Dictionary of base emission inventory xarrays

including any missing years compared to inv_dict. Keys are inventory years.

Return type:

dict[int, xr.Dataset]

Note

A custom nearest neighbour method is used for regridding and a linear interpolation method for calculating data in missing years. In future versions, the user will be able to select methods for both.

openairclim.calc_cont.interp_ppcf(inv: Dataset, p_pcf: DataArray, cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]])[source]

Interpolate p_PCF onto contrail grid.

Parameters:
  • inv (xr.Dataset) – Emission inventory for a given year.

  • p_pcf (xr.DataArray) – p_PCF on precalculated contrail data grid.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Returns:

Interpolated p_PCF; tuple of indices.

Return type:

(np.ndarray, tuple)

openairclim.calc_cont.load_base_inventories(config: dict[str, Any], inv_yrs: Sequence[int], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]]) dict[str, dict[int, Dataset]][source]

Load the base emission inventories. These must at least span the same time range as the input emission inventories, but can also be wider. The base emission inventories are linearly interpolated onto years that are defined by the input emission inventories if those years do not otherwise exist.

Parameters:
  • config (dict[str, Any]) – Configuration dictionary from config file.

  • inv_yrs (Sequence[int]) – Years for which the input emission inventories are defined.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

Raises:

ValueError – If the base emission inventories do not at least span the input emission inventories (given by inv_yrs).

Returns:

Full base emission inventory. First-

level keys are “ac”, second-level years.

Return type:

dict[str, dict[int, xr.Dataset]]

openairclim.calc_cont.logistic(x: ArrayLike, l: float, k: float, x0: float) ndarray[source]

Computes the logistic function, a sigmoid curve, commonly used to model growth or decay. Function from Megill & Grewe (2025): https://github.com/liammegill/contrail-limiting-factors

Parameters:
  • x (npt.ArrayLike) – The input values for which the logistic function will be computed.

  • l (float) – The maximum value or carrying capacity of the function.

  • k (float) – The steepness of the curve.

  • x0 (float) – The midpoint value of x where the function reaches half of l.

Returns:

The logistic function values for the given input x.

Return type:

np.ndarray

openairclim.calc_cont.logistic_gen(x: ArrayLike, l: float, k: float, x0: float, d: float) ndarray[source]

Computes a generalized logistic function with an additional vertical shift. Function from Megill & Grewe (2025): https://github.com/liammegill/contrail-limiting-factors

Parameters:
  • x (npt.ArrayLike) – The input values for which the logistic function will be computed.

  • l (float) – The maximum value or carrying capacity of the function.

  • k (float) – The steepness of the curve.

  • x0 (float) – The midpoint value of x where the function reaches half of l.

  • d (float) – The vertical shift applied to the function.

Returns:

The values of the shifted logistic function for the input

x.

Return type:

np.ndarray

openairclim.calc_cont.pad_inv_dict(inv_yrs: Sequence[int], inv_dict: dict[int, Dataset], pad_vars: list[str], cont_grid: tuple[ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]], ndarray[tuple[Any, ...], dtype[floating]]], ac: str) dict[source]

This function checks whether all years given in inv_yrs are present in the input emission inventory inv_dict. If a year is missing, a zero dataset is added to inv_dict for each variable in pad_vars on the pre-calculated contrail grid cont_grid. The ac variable is also added since this is necessary for other functions in the contrail module.

This functionality can be necessary if a specific aircraft identifier is not included in an emission inventory passed to OpenAirClim, for example because the aircraft newly enters service at a later time.

Parameters:
  • inv_yrs (Sequence[int]) – Years for which the inv_dict emission inventory should be defined.

  • inv_dict (dict[int, xr.Dataset]) – Dictionary of emission inventory xarrays, keys are inventory years.

  • pad_vars (list[str]) – Variables to be included in the xarrays.

  • cont_grid (tuple) – Precalculated contrail grid. Shape (lon, lat, plev); each is 1-D float array with shapes (n_lon,), (n_lat,), (n_plev,). Units: lon [deg], lat [deg], plev [hPa].

  • ac (str) – Aircraft identifier from config.

Returns:

inv_dict modified in-place with zero arrays in missing years.

Return type:

dict

openairclim.calc_cont.pm_factor(x: float, ls_case: str = 'case_mid') float[source]

Calculate the nvPM factor depending on the relative nvPM emissions. The low-soot regime (x < 0.1) is currently only available with the OpenAirClim Premium license. Please contact the OpenAirClim team if you would like access to this. Note that the factor is not validated in the low-soot regime.

Reference: Megill (2026)

Parameters:
  • x (float) – relative nvPM emissions with respect to 1.5e15 kg^-1

  • ls_case (str, optional) – Descriptor of pre-calculated case. One of: “case_low”, “case_mid” or “case_high”. Defaults to “case_mid”.

Raises:
  • ValueError – For invalid nvPM emissions.

  • ImportError – For low-soot regime calls if the openairclim_premium package cannot be found.

Returns:

nvPM factor

Return type:

float

openairclim.calc_cont.pm_factor_high(x: float, params: tuple) float[source]

Calculate nvPM factor in the high-soot regime.

Parameters:
  • x (float) – relative nvPM emissions with respect to 1.5e15 kg^-1

  • params (tuple) – fit parameters

Returns:

nvPM factor

Return type:

float

openairclim.calc_cont.pm_factor_high_prime(x: float, params: tuple) float[source]

Calculate the first-order derivative of the nvPM factor in the high-soot regime.

Parameters:
  • x (float) – relative nvPM emissions with respect to 1.5e15 kg^-1

  • params (tuple) – fit parameters

Returns:

derivative of the nvPM factor

Return type:

float