API reference

fair.fair

Initialiser module for FaIR.

class fair.FAIR(n_gasboxes=4, n_layers=3, iirf_max=100, br_cl_ods_potential=45, ghg_method='meinshausen2020', ch4_method='leach2021', temperature_prescribed=False)

Initialise FaIR.

Parameters:
  • n_gasboxes (int) – the number of atmospheric greenhouse gas boxes to run the model with

  • n_layers (int) – the number of ocean layers in the energy balance or impulse response model to run with

  • iirf_max (float) – limit for time-integral of greenhouse gas impulse response function.

  • br_cl_ods_potential (float) – factor describing the ratio of efficiency that each bromine atom has as an ozone depleting substance relative to each chlorine atom.

  • ghg_method (str) – method to use for calculating greenhouse gas forcing from CO2, CH4 and N2O. Valid options are {“leach2021”, “meinshausen2020”, “etminan2016”, “myhre1998”}

  • ch4_method (str) – method to use for calculating methane lifetime change. Valid options are {“leach2021”, “thornhill2021”}.

  • temperature_prescribed (bool) – Run FaIR with temperatures prescribed.

Raises:

ValueError – if ghg_method or ch4_method given are not valid options.

allocate()

Create xarray DataArrays of data input and output.

calculate_concentration_per_emission(mass_atmosphere=5.1352e+18, molecular_weight_air=28.97)

Calculate change in atmospheric concentration for each unit emission.

calculate_g(iirf_horizon=100)

Calculate lifetime scaling parameters.

Parameters:

iirf_horizon (float) – time horizon for time-integrated airborne fraction (yr).

calculate_iirf0(iirf_horizon=100)

Convert greenhouse gas lifetime to time-integrated airborne fraction.

iirf_0 is the 100-year time-integrated airborne fraction to a pulse emission. We know that the gas’s atmospheric airborne fraction \(a(t)\) for a gas with lifetime \(\tau\) after time \(t\) is therefore

\[a(t) = \exp(-t/tau)\]

and integrating this for \(H\) years after a pulse emissions gives us:

\[r_0(t) = \int_0^{H} \exp(-t/\tau) dt = \tau (1 - \exp (-H / \tau)).\]

\(H\) = 100 years is the default time horizon in FaIR but this can be set to any value.

Parameters:

iirf_horizon (float) – time horizon for time-integrated airborne fraction (yr).

property ch4_method

Return methane lifetime method.

define_configs(configs)

Define configs to analyse in FaIR.

Parameters:

configs (list) – config names to run

define_scenarios(scenarios)

Define scenarios to analyse in FaIR.

Parameters:

scenarios (list) – scenario names to run

define_species(species, properties)

Define species to run in FaIR.

Parameters:
  • species (list) – names of species to include in FaIR

  • properties (dict) –

    mapping of each specie to particular run properties. This is a nested dict, where each dict key contains a dict of 5 keys as follows:

    typestr

    the type of specie that is being provided. Valid inputs are “co2 ffi”, “co2 afolu”, “co2”, “ch4”, “n2o”, “cfc-11”, “other halogen”, “f-gas”, “sulfur”, “black carbon”, “organic carbon”, “other slcf”, “nox aviation”, “eesc”, “ozone”, “ari”, “aci”, “contrails”, “lapsi”, “h2o stratospheric”, “land use”, “volcanic”, “solar”, “unspecified”,

    input_mode{‘emissions’, ‘concentration’, ‘forcing’, ‘calculated’}

    describes how the specie is input into the model.

    greenhouse_gasbool

    is the specie a greenhouse gas?

    aerosol_chemistry_from_emissionsbool

    does the specie’s emissions affect aerosols and/or chemistry?

    aerosol_chemistry_from_concentrationbool

    does the specie’s concentration affect aerosols and/or chemistry?

Raises:
  • ValueError – if a specie in species does not have a matching key in properties.

  • ValueError – if an invalid species type is specified.

  • ValueError – if an invalid input_type (driving mode) is provided for a particular type.

  • ValueError – if duplicate species types are provided for types that must be unique.

define_time(start, end, step)

Define timebounds vector to run FaIR.

Parameters:
  • start (float) – first timebound of the model (year)

  • end (float) – last timebound of the model (year)

  • step (float) – timestep (year)

fill_from_rcmip()

Fill emissions, concentrations and/or forcing from RCMIP scenarios.

fill_species_configs(filename='/home/docs/checkouts/readthedocs.org/user_builds/fair/envs/latest/lib/python3.11/site-packages/fair/defaults/data/ar6/species_configs_properties.csv')

Fill the species_configs with values from a CSV file.

Parameters:

filename (str) – Path of the CSV file to read the species configs from. If omitted, the default configs file will be read.

property ghg_method

Return greenhouse gas forcing method.

run(progress=True, suppress_warnings=True)

Run the FaIR model.

Parameters:
  • progress (bool) – Display progress bar.

  • suppress_warnings (bool) – Hide warnings relating to covariance in energy balance matrix.

to_netcdf(filename)

Write out FaIR scenario data to a netCDF file.

Parameters:

filename (str) – file path of the file to write.

fair.energy_balance_model

n-layer energy balance representation of Earth’s climate.

class fair.energy_balance_model.EnergyBalanceModel(ocean_heat_capacity, ocean_heat_transfer, deep_ocean_efficacy=1, forcing_4co2=8, stochastic_run=False, sigma_eta=0.5, sigma_xi=0.5, gamma_autocorrelation=2, seed=None, timestep=1, n_timesteps=1)

Energy balance model that converts forcing to temperature.

The energy balance model is converted to an impulse-response formulation (hence the IR part of FaIR) to allow efficient evaluation. The benefits of this are increased as once derived, the “layers” of the energy balance model do not communicate with each other. The model description can be found in [Leach2021], [Cummins2020], [Tsutsui2017] and [Geoffroy2013].

Parameters:
  • ocean_heat_capacity (np.ndarray) – Ocean heat capacity of each layer (top first), W m-2 yr K-1

  • ocean_heat_transfer (np.ndarray) – Heat exchange coefficient between ocean layers (top first). The first element of this array is akin to the climate feedback parameter, with the convention that stabilising feedbacks are positive (opposite to most climate sensitivity literature). W m-2 K-1

  • deep_ocean_efficacy (float) – efficacy of deepest ocean layer. See e.g. [Geoffroy2013].

  • forcing_4co2 (float) – effective radiative forcing from a quadrupling of atmospheric CO2 concentrations above pre-industrial.

  • stochastic_run (bool) – Activate the stochastic variability component from [Cummins2020].

  • sigma_eta (float) – Standard deviation of stochastic forcing component from [Cummins2020].

  • sigma_xi (float) – Standard deviation of stochastic disturbance applied to surface layer. See [Cummins2020].

  • gamma_autocorrelation (float) – Stochastic forcing continuous-time autocorrelation parameter. See [Cummins2020].

  • seed (int or None) – Random seed to use for stochastic variability.

  • timestep (float) – Time interval of the model (yr)

Raises:
  • ValueError – if the shapes of ocean_heat_capacity and ocean_heat_transfer differ.

  • ValueError – if there are not at least two layers in the energy balance model.

add_forcing(forcing, timestep)

Add a forcing time series to EnergyBalanceModel.

Parameters:
  • forcing (np.ndarray) – time series of [effective] radiative forcing

  • timestep (float) – Model timestep, in years

property eb_matrix_d

Return the discretised matrix exponential.

emergent_parameters(forcing_2co2_4co2_ratio=0.5)

Calculate emergent parameters from the energy balance parameters.

Parameters:

forcing_2co2_4co2_ratio (float) – ratio of (effective) radiative forcing converting a quadrupling of CO2 to a doubling of CO2.

property forcing_vector_d

Return the discretised forcing vector.

impulse_response()

Convert the energy balance to impulse response representation.

run()

Run the EnergyBalanceModel.

property stochastic_d

Return the stochastic matrix.

fair.energy_balance_model.calculate_toa_imbalance_postrun(state, forcing, ocean_heat_transfer, deep_ocean_efficacy)

Calculate top of atmosphere energy imbalance.

The calculation is performed after the scenario has been run to avoid looping, since no dynamic state changes affect the calculation.

Parameters:
  • state (np.ndarray) – stacked arrays of forcing and temperature of layers across the run

  • forcing (np.ndarray) – stacked arrays of [effective] radiative forcing across the run

  • ocean_heat_transfer (np.ndarray) – Heat exchange coefficient between ocean layers (top first). The first element of this array is akin to the climate feedback parameter, with the convention that stabilising feedbacks are positive (opposite to most climate sensitivity literature). W m-2 K-1

  • deep_ocean_efficacy (np.ndarray) – efficacy of deepest ocean layer.

Returns:

toa_imbalance – Top of atmsophere energy imbalance.

Return type:

np.ndarray

fair.energy_balance_model.multi_ebm(configs, ocean_heat_capacity, ocean_heat_transfer, deep_ocean_efficacy, stochastic_run, sigma_eta, sigma_xi, gamma_autocorrelation, seed, use_seed, forcing_4co2, timestep, timebounds)

Create several instances of the EnergyBalanceModel.

This allows efficient parallel implementation in FaIR. We have to use a for loop in this function as is does not look like the linear algebra functions in scipy are naturally parallel.

Parameters:
  • configs (list) – A named list of climate configurations.

  • ocean_heat_capacity (np.ndarray) – Ocean heat capacity of each layer (top first), W m-2 yr K-1

  • ocean_heat_transfer (np.ndarray) – Heat exchange coefficient between ocean layers (top first). The first element of this array is akin to the climate feedback parameter, with the convention that stabilising feedbacks are positive (opposite to most climate sensitivity literature). W m-2 K-1

  • deep_ocean_efficacy (float) – efficacy of deepest ocean layer. See e.g. [Geoffroy2013].

  • stochastic_run (bool) – Activate the stochastic variability component from [Cummins2020].

  • sigma_eta (float) – Standard deviation of stochastic forcing component from [Cummins2020].

  • sigma_xi (float) – Standard deviation of stochastic disturbance applied to surface layer. See [Cummins2020].

  • gamma_autocorrelation (float) – Stochastic forcing continuous-time autocorrelation parameter. See [Cummins2020].

  • seed (int or None) – Random seed to use for stochastic variability.

  • use_seed (bool) – Whether or not to use the random seed.

  • forcing_4co2 (float) – effective radiative forcing from a quadrupling of atmospheric CO2 concentrations above pre-industrial.

  • timestep (float) – Time interval of the model (yr)

  • timebounds (np.ndarray) – Vector representing the time snapshots to calculate temperature on.

fair.energy_balance_model.step_temperature(state_old, eb_matrix_d, forcing_vector_d, stochastic_d, forcing)

Advance parallel energy balance models forward one timestep.

Parameters:
  • state_old (np.ndarray) – stacked arrays of forcing and temperature of layers in previous timestep

  • eb_matrix_d (np.ndarray) – stacked discretised energy balance matrices

  • forcing_vector_d (np.ndarray) – stacked discretised forcing vectors

  • _stochastic_d (np.ndarray) – stacked matrices of stochastic internal variability

  • forcing (np.ndarray) – stacked vectors of [effective] radiative forcing

Returns:

state_new – stacked arrays of forcing and temperature of layers in this timestep

Return type:

np.ndarray

fair.interface

Convenience functions for filling FaIR xarray instances.

fair.interface.fill(var, data, **kwargs)

Fill a FAIR variable instance.

One could work directly with the xarray DataArrays, but this function includes additional error checking and validation, and is slightly less cumbersome than the xarray method of allocation.

Parameters:
  • var (attr) – FAIR variable attribute to fill, for example fair.climate_configs[“ocean_heat_capacity”]

  • data (np.ndarray) – data to fill the variable with

  • **kwargs – the dimensions represeted in data

Raises:

ValueError – if a kwargs element provided doesn’t correspond to a dimension name in var

fair.interface.initialise(var, value, **kwargs)

Fill a fair variable instance with value in first timebound.

Otherwise identical to fill.

Parameters:
  • var (attr) – FAIR variable attribute to fill for example fair.climate_configs[“ocean_heat_capacity”]

  • value (np.ndarray) – value to fill the first timebound with

  • **kwargs – the dimensions represeted in data

fair.io

Tools for getting data into and out of FaIR.

fair.io.read_properties(filename='/home/docs/checkouts/readthedocs.org/user_builds/fair/envs/latest/lib/python3.11/site-packages/fair/defaults/data/ar6/species_configs_properties.csv', species=None)

Get a properties file.

Parameters:
  • filename (str) – path to a csv file. Default is an AR6 WG1-like config for FaIR covering all of the species considered in CMIP6.

  • species (list of str or None) – the species that are to be included in the FaIR run. All of these species should be present in the index (first column) of the csv. If None (default), return all of the species in the defaults.

Returns:

  • species (list) – a list of species names that are included in the FaIR run.

  • properties (dict) – species properties that control the FaIR run

fair.constants

General constants.

fair.constants.DOUBLING_TIME_1PCT = 69.66071689357483

Time in years for a doubling to occur at a rate of 1% per year.

fair.constants.TIME_AXIS = 0

Axis of the data arrays referencing time.

fair.constants.SCENARIO_AXIS = 1

Axis of the data arrays referencing scenario.

fair.constants.CONFIG_AXIS = 2

Axis of the data arrays referencing climate or specie configuration.

fair.constants.SPECIES_AXIS = 3

Axis of emissions, concentration and forcing data arrays representing species.

fair.constants.GASBOX_AXIS = 4

Axis of atmopsheric gas box for the gas partition data array in full emissions to concentration mode.

fair.earth_params

Parameters (not constants) that define Earth properties.

fair.earth_params.molecular_weight_air = 28.97

Molecular weight of air, g/mol.

fair.earth_params.earth_radius = 6371000

Radius of Earth, m

fair.earth_params.mass_atmosphere = 5.1352e+18

Mass of Earth’s atmosphere, kg.

fair.earth_params.seconds_per_year = 31556925.216

Length of a tropical year, s.

fair.forcing

Module for effective radiative forcing.

fair.forcing.ghg

Module for greenhouse gas forcing.

fair.forcing.ghg.etminan2016(concentration, baseline_concentration, forcing_scaling, radiative_efficiency, co2_indices, ch4_indices, n2o_indices, minor_greenhouse_gas_indices, a1=-2.4e-07, b1=0.00072, c1=-0.00021, d1=5.36, a2=-8e-06, b2=4.2e-06, c2=-4.9e-06, d2=0.117, a3=-1.3e-06, b3=-8.2e-06, d3=0.043)

Greenhouse gas forcing from concentrations.

This uses the [Etminan2016] relationships for CO2, CH4 and N2O including band overlaps between these three gases. Forcing from minor greenhouse gases are a linear function of their concentration based on their radiative efficiency.

Note that the relationship can be poorly behaved outside of the range of valid values in [Etminan2016], (180-2000 ppm CO2, 340-3500 ppb CH4, 200-525 ppb N2O) particularly for CO2 concentrations above 2000 ppm. It is recommended to use “meinshausen2020”, which is a re-fit of the coefficients and extension to be better behaved outside of the valid concentration range.

Parameters:
  • concentration (np.ndarray) – concentration of greenhouse gases. “CO2”, “CH4” and “N2O” must be included in units of [ppm, ppb, ppb]. Other GHGs are units of ppt.

  • baseline_concentration (np.ndarray) – pre-industrial concentration of the gases (see above).

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency to use for linear-forcing gases, in W m-2 ppb-1

  • co2_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CO2.

  • ch4_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CH4.

  • n2o_indices (np.ndarray of bool) – index along SPECIES AXIS relating to N2O.

  • minor_greenhouse_gas_indices (np.ndarray of bool) – indices of other GHGs that are not CO2, CH4 or N2O.

  • a1 (float) – fitting parameter (see [Etminan2016])

  • b1 (float) – fitting parameter (see [Etminan2016])

  • c1 (float) – fitting parameter (see [Etminan2016])

  • d1 (float) – fitting parameter (see [Etminan2016])

  • a2 (float) – fitting parameter (see [Etminan2016])

  • b2 (float) – fitting parameter (see [Etminan2016])

  • c2 (float) – fitting parameter (see [Etminan2016])

  • d2 (float) – fitting parameter (see [Etminan2016])

  • a3 (float) – fitting parameter (see [Etminan2016])

  • b3 (float) – fitting parameter (see [Etminan2016])

  • d3 (float) – fitting parameter (see [Etminan2016])

Returns:

effective_radiative_forcing – effective radiative forcing (W/m2) from greenhouse gases

Return type:

np.ndarray

fair.forcing.ghg.leach2021ghg(concentration, baseline_concentration, forcing_scaling, radiative_efficiency, co2_indices, ch4_indices, n2o_indices, minor_greenhouse_gas_indices, f1_co2=4.57, f3_co2=0.086, f3_ch4=0.038, f3_n2o=0.106)

Greenhouse gas forcing from concentrations.

This uses the [Leach2021] relationships for CO2, CH4 and N2O that do not include band overlaps between these gases, allowing single-forcing runs. This is the default treatment in FaIR2.0. This is a re-fit of the [Etminan2016] formulation with improved coefficient fits. Minor greenhouse gases are a linear function of their concentration based on their radiative efficiency.

Parameters:
  • concentration (np.ndarray) – concentration of greenhouse gases. “CO2”, “CH4” and “N2O” must be included in units of [ppm, ppb, ppb]. Other GHGs are units of ppt.

  • baseline_concentration (np.ndarray) – pre-industrial concentration of the gases (see above).

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency to use for linear-forcing gases, in W m-2 ppb-1

  • co2_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CO2.

  • ch4_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CH4.

  • n2o_indices (np.ndarray of bool) – index along SPECIES AXIS relating to N2O.

  • minor_greenhouse_gas_indices (np.ndarray of bool) – indices of other GHGs that are not CO2, CH4 or N2O.

  • f1_co2 (float) – factor relating logarithm of CO2 conentration to radiative forcing.

  • f3_co2 (float) – factor relating square root of CO2 conentration to radiative forcing.

  • f3_ch4 (float) – factor relating square root of CH4 conentration to radiative forcing.

  • f3_n2o (float) – factor relating square root of N2O conentration to radiative forcing.

Returns:

effective_radiative_forcing – effective radiative forcing (W/m2) from greenhouse gases

Return type:

np.ndarray

fair.forcing.ghg.meinshausen2020(concentration, reference_concentration, forcing_scaling, radiative_efficiency, co2_indices, ch4_indices, n2o_indices, minor_greenhouse_gas_indices, a1=-2.4785e-07, b1=0.00075906, c1=-0.0021492, d1=5.2488, a2=-0.00034197, b2=0.00025455, c2=-0.00024357, d2=0.12173, a3=-8.9603e-05, b3=-0.00012462, d3=0.045194)

Greenhouse gas forcing from concentrations.

This uses the [Meinshausen2020] relationships for CO2, CH4 and N2O including band overlaps between these three gases. This is a rescaled [Etminan2016] function with improved stability outside the range of validity in [Etminan2016]. Minor greenhouse gases are a linear function of their concentration based on their radiative efficiency.

Parameters:
  • concentration (np.ndarray) – concentration of greenhouse gases. “CO2”, “CH4” and “N2O” must be included in units of [ppm, ppb, ppb]. Other GHGs are units of ppt.

  • reference_concentration (np.ndarray) – pre-industrial concentration of the gases (see above) used as the reference to calculate the forcing.

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency to use for linear-forcing gases, in W m-2 ppb-1

  • co2_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CO2.

  • ch4_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CH4.

  • n2o_indices (np.ndarray of bool) – index along SPECIES AXIS relating to N2O.

  • minor_greenhouse_gas_indices (np.ndarray of bool) – indices of other GHGs that are not CO2, CH4 or N2O.

  • a1 (float) – fitting parameter (see [Meinshausen2020])

  • b1 (float) – fitting parameter (see [Meinshausen2020])

  • c1 (float) – fitting parameter (see [Meinshausen2020])

  • d1 (float) – fitting parameter (see [Meinshausen2020])

  • a2 (float) – fitting parameter (see [Meinshausen2020])

  • b2 (float) – fitting parameter (see [Meinshausen2020])

  • c2 (float) – fitting parameter (see [Meinshausen2020])

  • d2 (float) – fitting parameter (see [Meinshausen2020])

  • a3 (float) – fitting parameter (see [Meinshausen2020])

  • b3 (float) – fitting parameter (see [Meinshausen2020])

  • d3 (float) – fitting parameter (see [Meinshausen2020])

Returns:

effective_radiative_forcing – effective radiative forcing (W/m2) from greenhouse gases

Return type:

np.ndarray

fair.forcing.ghg.myhre1998(concentration, baseline_concentration, forcing_scaling, radiative_efficiency, co2_indices, ch4_indices, n2o_indices, minor_greenhouse_gas_indices, alpha_co2=5.35, alpha_ch4=0.036, alpha_n2o=0.12, alpha_ch4_n2o=0.47, a1=2.01e-05, exp1=0.75, a2=5.32e-15, exp2=1.52)

Greenhouse gas forcing from concentrations.

Band overlaps are included between CH4 and N2O. This relationship comes from [Myhre1998], and was used up until the IPCC’s Fifth Assessment Report. Minor greenhouse gases are a linear function of their concentration based on their radiative efficiency.

Parameters:
  • concentration (np.ndarray) – concentration of greenhouse gases. “CO2”, “CH4” and “N2O” must be included in units of [ppm, ppb, ppb]. Other GHGs are units of ppt.

  • baseline_concentration (np.ndarray) – pre-industrial concentration of the gases (see above).

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency to use for linear-forcing gases, in W m-2 ppb-1

  • co2_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CO2.

  • ch4_indices (np.ndarray of bool) – index along SPECIES_AXIS relating to CH4.

  • n2o_indices (np.ndarray of bool) – index along SPECIES AXIS relating to N2O.

  • minor_greenhouse_gas_indices (np.ndarray of bool) – indices of other GHGs that are not CO2, CH4 or N2O.

  • alpha_co2 (float) – factor relating logarithm of CO2 conentration to radiative forcing.

  • alpha_ch4 (float) – factor relating square root of CH4 conentration to radiative forcing.

  • alpha_n2o (float) – factor relating square root of N2O conentration to radiative forcing.

Returns:

effective_radiative_forcing – effective radiative forcing (W/m2) from greenhouse gases

Return type:

np.ndarray

fair.forcing.minor

Module for a generic linear emissions or concentration to forcing calculation.

fair.forcing.minor.calculate_linear_forcing(driver, baseline_driver, forcing_scaling, radiative_efficiency)

Calculate effective radiative forcing from a linear relationship.

Parameters:
  • driver (np.ndarray) – input emissions, concentration or forcing

  • baseline_driver (np.ndarray) – baseline, possibly pre-industrial, emissions, concentration or forcing.

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency (W m-2 (<driver unit>)-1) of each species.

Returns:

erf_out – effective radiative forcing (W m-2)

Return type:

np.ndarray

fair.forcing.ozone

Module for ozone forcing.

fair.forcing.ozone.thornhill2021(emissions, concentration, baseline_emissions, baseline_concentration, forcing_scaling, ozone_radiative_efficiency, slcf_indices, ghg_indices)

Determine ozone effective radiative forcing.

Calculates total ozone forcing from precursor emissions and concentrations based on AerChemMIP and CMIP6 Historical behaviour in [Skeie2020] and [Thornhill2021a].

The treatment is identical to FaIR2.0 [Leach2021].

Parameters:
  • emissions (np.ndarray) – emissions in timestep

  • concentration (np.ndarray) – concentrations in timestep

  • baseline_emissions (np.ndarray) – reference, possibly pre-industrial, emissions

  • baseline_concentration (np.ndarray) – reference, possibly pre-industrial, concentrations

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • ozone_radiative_efficiency (np.ndarray) – the radiative efficiency at which ozone precursor emissions or concentrations are converted to ozone radiative forcing. The unit is W m-2 (<native emissions or concentration unit>)-1, where the emissions unit is usually Mt/yr for short-lived forcers, ppb for CH4 and N2O concentrations, and ppt for halogenated species. Note this is not the same radiative efficiency that is used in converting GHG concentrations to forcing.

  • slcf_indices (list of int) – provides a mapping of which aerosol species corresponds to which emitted species index along the SPECIES_AXIS.

  • ghg_indices (list of int) – provides a mapping of which aerosol species corresponds to which atmospheric GHG concentration along the SPECIES_AXIS.

Returns:

_erf_out – ozone forcing

Return type:

np.ndarray

fair.forcing.aerosol

Module for aerosol effective radiative forcing.

fair.forcing.aerosol.erfari

Module for calculating forcing from aerosol-radiation interactions.

fair.forcing.aerosol.erfari.calculate_erfari_forcing(emissions, concentration, baseline_emissions, baseline_concentration, forcing_scaling, radiative_efficiency, emissions_indices, concentration_indices)

Calculate effective radiative forcing from aerosol-radiation interactions.

Parameters:
  • emissions (np.ndarray) – emissions in timestep

  • concentration (np.ndarray) – concentrations in timestep

  • baseline_emissions (np.ndarray) – baseline emissions to evaluate forcing against

  • baseline_concentration (np.ndarray) – baseline concentrations to evaluate forcing against

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • radiative_efficiency (np.ndarray) – radiative efficiency (W m-2 (emission_unit yr-1)-1) of each species.

  • emissions_indices (list of int) – provides a mapping of which aerosol species corresponds to which emitted species index along the SPECIES_AXIS.

  • concentration_indices (list of int) – provides a mapping of which aerosol species corresponds to which atmospheric GHG concentration along the SPECIES_AXIS.

Returns:

effective_radiative_forcing – effective radiative forcing (W m-2) from aerosol-radiation interactions

Return type:

np.ndarray

fair.forcing.aerosol.erfaci

Module for forcing from aerosol-cloud interactions.

fair.forcing.aerosol.erfaci.logsum(emissions, concentration, baseline_emissions, baseline_concentration, forcing_scaling, scale, sensitivity, slcf_indices, ghg_indices)

Calculate effective radiative forcing from aerosol-cloud interactions.

This uses the relationship to calculate ERFaci as follows

\[F = \beta \log \left( 1 + \sum_{i} s_i A_i \right)\]

where \(A_i\) is the emissions or concentration of a specie, \(\beta\) is the scale factor and \(s_i\) describes how sensitive a specie is in contributing to ERFaci.

The calculation is performed for the emissions/concentration of interest and then for the baseline. The difference between the two values is the forcing.

This relationship is a generalisation of [Stevens2015]. To recover [Stevens2015], set \(s_i\) to zero for all except SO2, and \(s_{\text{SO2}} = 1/Q_n\) where \(Q_n\) is the natural emissions source in [Stevens2015].

Parameters:
  • emissions (np.ndarray) – input emissions

  • concentration (np.ndarray) – input emissions

  • baseline_emissions (np.ndarray) – baseline emissions

  • baseline_concentration (np.ndarray) – baseline concentration

  • forcing_scaling (np.ndarray) – scaling of the calculated radiative forcing (e.g. for conversion to effective radiative forcing and forcing uncertainty).

  • scale (np.ndarray) – per-species scaling factor to apply to the logarithm

  • sensitivity (np.ndarray) – per-species sensitivity factor for the logarithm

  • slcf_indices (list of int) – provides a mapping of which aerosol species corresponds to which emitted species index along the SPECIES_AXIS.

  • ghg_indices (list of int) – provides a mapping of which aerosol species corresponds to which atmospheric GHG concentration along the SPECIES_AXIS.

Returns:

effective_radiative_forcing – effective radiative forcing (W m-2) from aerosol-cloud interactions

Return type:

np.ndarray

fair.gas_cycle

Module containing gas cycle functions.

fair.gas_cycle.ch4_lifetime

Methane lifetime definition that is based on multiple species.

fair.gas_cycle.ch4_lifetime.calculate_alpha_ch4(emissions, concentration, temperature, baseline_emissions, baseline_concentration, ch4_lifetime_chemical_sensitivity, ch4_lifetime_temperature_sensitivity, emissions_indices, concentration_indices)

Calculate methane lifetime scaling factor.

Parameters:
  • emissions (np.ndarray) – Emitted species

  • concentration (np.ndarray) – Species derived from concentration (including EESC)

  • temperature (np.ndarray) – Global mean surface level temperature anomaly

  • baseline_emissions (np.ndarray) – reference, possibly pre-industrial, emissions

  • baseline_concentration (np.ndarray) – reference, possibly pre-industrial, concentration

  • ch4_lifetime_chemical_sensitivity (np.ndarray) – per-species sensitivity of change in CH4 lifetime per unit emission.

  • ch4_lifetime_temperature_sensitivity (np.ndarray) – per-species sensitivity of change in CH4 lifetime per unit concentration increase.

  • emissions_indices (np.ndarray of bool) – Which species indices to use emissions from.

  • concentration_indices (np.ndarray of bool) – Which species indices to use concentration from.

Returns:

alpha_ch4 – Scaling factor to baseline methane lifetime.

Return type:

np.ndarray

fair.gas_cycle.eesc

Module for calculaing equivalent effective stratospheric chlorine.

fair.gas_cycle.eesc.calculate_eesc(concentration, fractional_release, cl_atoms, br_atoms, cfc_11_index, halogen_indices, br_cl_ratio)

Calculate equivalent effective stratospheric chlorine.

Parameters:
  • concentration (np.ndarray) – concentrations in timestep

  • fractional_release (np.ndarray) – fractional release describing the proportion of available ODS that actually contributes to ozone depletion.

  • cl_atoms (np.ndarray) – Chlorine atoms in each species

  • br_atoms (np.ndarray) – Bromine atoms in each species

  • cfc_11_index (int or None) – array index along SPECIES_AXIS corresponding to CFC-11.

  • halogen_indices (list of int) – provides a mapping of which halogen species corresponds to which index along the SPECIES_AXIS.

  • br_cl_ratio (float) – how much more effective bromine is as an ozone depletor than chlorine.

Returns:

eesc_out – equivalent effective stratospheric chlorine

Return type:

np.ndarray

Notes

At present, CFC-11 must be provided in the scenario, with species type cfc-11.

fair.gas_cycle.forward

Module for the forward (emissions to concentration) model.

fair.gas_cycle.forward.step_concentration(emissions, gasboxes_old, airborne_emissions_old, alpha_lifetime, baseline_concentration, baseline_emissions, concentration_per_emission, lifetime, partition_fraction, timestep)

Calculate concentrations from emissions of any greenhouse gas.

Parameters:
  • emissions (np.ndarray) – emissions rate (emissions unit per year) in timestep.

  • gas_boxes_old (np.ndarray) – the greenhouse gas atmospheric burden in each lifetime box at the end of the previous timestep.

  • airborne_emissions_old (np.ndarray) – The total airborne emissions at the beginning of the timestep. This is the concentrations above the pre-industrial control. It is also the sum of gas_boxes_old if this is an array.

  • alpha_lifetime (np.ndarray) – scaling factor for lifetime. Necessary where there is a state- dependent feedback.

  • baseline_concentration (np.ndarray) – original (possibly pre-industrial) concentration of gas(es) in question.

  • baseline_emissions (np.ndarray or float) – original (possibly pre-industrial) emissions of gas(es) in question.

  • concentration_per_emission (np.ndarray) – how much atmospheric concentrations grow (e.g. in ppm) per unit (e.g. GtCO2) emission.

  • lifetime (np.ndarray) – atmospheric burden lifetime of greenhouse gas (yr). For multiple lifetimes gases, it is the lifetime of each box.

  • partition_fraction (np.ndarray) – the partition fraction of emissions into each gas box. If array, the entries should be individually non-negative and sum to one.

  • timestep (float) – emissions timestep in years.

Notes

Emissions are given in time intervals and concentrations are also reported on the same time intervals: the airborne_emissions values are on time boundaries and these are averaged before being returned.

Where array input is taken, the arrays always have the dimensions of (scenario, species, time, gas_box). Dimensionality can be 1, but we retain the singleton dimension in order to preserve clarity of calculation and speed.

Returns:

  • concentration_out (np.ndarray) – greenhouse gas concentrations at the centre of the timestep.

  • gas_boxes_new (np.ndarray) – the greenhouse gas atmospheric burden in each lifetime box at the end of the timestep.

  • airborne_emissions_new (np.ndarray) – airborne emissions (concentrations above pre-industrial control level) at the end of the timestep.

fair.gas_cycle.inverse

Module for deriving emissions from concentration.

fair.gas_cycle.inverse.unstep_concentration(concentration, gasboxes_old, airborne_emissions_old, alpha_lifetime, baseline_concentration, baseline_emissions, concentration_per_emission, lifetime, partition_fraction, timestep)

Calculate emissions from concentrations of any greenhouse gas.

Parameters:
  • concentration (np.ndarray) – greenhouse gas concentration at the centre of the timestep.

  • gas_boxes_old (np.ndarray) – the greenhouse gas atmospheric burden in each lifetime box at the end of the previous timestep.

  • airborne_emissions_old (np.ndarray) – The total airborne emissions at the beginning of the timestep. This is the concentrations above the pre-industrial control. It is also the sum of gas_boxes_old if this is an array.

  • alpha_lifetime (np.ndarray) – scaling factor for lifetime. Necessary where there is a state- dependent feedback.

  • baseline_concentration (np.ndarray) – baseline (possibly pre-industrial) concentration of gas in question.

  • baseline_emissions (np.ndarray) – baseline (possibly pre-industrial) emissions of gas in question.

  • concentration_per_emission (np.ndarray) – how much atmospheric concentrations grow (e.g. in ppm) per unit (e.g. GtCO2) emission.

  • lifetime (np.ndarray) – atmospheric burden lifetime of greenhouse gas (yr). For multiple lifetimes gases, it is the lifetime of each box.

  • partition_fraction (np.ndarray) – the partition fraction of emissions into each gas box. If array, the entries should be individually non-negative and sum to one.

  • timestep (float) – emissions timestep in years.

Notes

Emissions are given in time intervals and concentrations are also reported on the same time intervals: the airborne_emissions values are on time boundaries. Therefore it is not actually possible to provide the exact emissions that would reproduce the concentrations without using a slower root-finding mechanism (that was present in v1.4) and will always be half a time step out.

Where array input is taken, the arrays always have the dimensions of (scenario, species, time, gas_box). Dimensionality can be 1, but we retain the singleton dimension in order to preserve clarity of calculation and speed.

Returns:

  • emissions_out (np.ndarray) – emissions in timestep.

  • gas_boxes_new (np.ndarray) – the greenhouse gas atmospheric burden in each lifetime box at the end of the timestep.

  • airborne_emissions_new (np.ndarray) – airborne emissions (concentrations above pre-industrial control level) at the end of the timestep.

fair.structure

Module for organising FaIR’s structure.

fair.structure.species

Define the types of beasts in FaIR, how they roar, and whether they are friendly.

fair.structure.species.multiple_allowed = {'aci': False, 'ari': False, 'black carbon': False, 'cfc-11': False, 'ch4': False, 'co2': False, 'co2 afolu': False, 'co2 ffi': False, 'contrails': False, 'eesc': False, 'f-gas': True, 'h2o stratospheric': False, 'land use': False, 'lapsi': False, 'n2o': False, 'nox aviation': False, 'organic carbon': False, 'other halogen': True, 'other slcf': True, 'ozone': False, 'solar': False, 'sulfur': False, 'unspecified': True, 'volcanic': False}

Whether multiple instances are allowed for each specie type.

fair.structure.species.species_types = ['co2 ffi', 'co2 afolu', 'co2', 'ch4', 'n2o', 'cfc-11', 'other halogen', 'f-gas', 'sulfur', 'black carbon', 'organic carbon', 'other slcf', 'nox aviation', 'eesc', 'ozone', 'ari', 'aci', 'contrails', 'lapsi', 'h2o stratospheric', 'land use', 'volcanic', 'solar', 'unspecified']

Species types recognised by FaIR.

fair.structure.species.valid_input_modes = {'aci': ['calculated', 'forcing'], 'ari': ['calculated', 'forcing'], 'black carbon': ['emissions'], 'cfc-11': ['emissions', 'concentration', 'forcing'], 'ch4': ['emissions', 'concentration', 'forcing'], 'co2': ['emissions', 'calculated', 'concentration', 'forcing'], 'co2 afolu': ['emissions'], 'co2 ffi': ['emissions'], 'contrails': ['calculated', 'forcing'], 'eesc': ['concentration', 'calculated'], 'f-gas': ['emissions', 'concentration', 'forcing'], 'h2o stratospheric': ['calculated', 'forcing'], 'land use': ['calculated', 'forcing'], 'lapsi': ['calculated', 'forcing'], 'n2o': ['emissions', 'concentration', 'forcing'], 'nox aviation': ['emissions'], 'organic carbon': ['emissions'], 'other halogen': ['emissions', 'concentration', 'forcing'], 'other slcf': ['emissions'], 'ozone': ['calculated', 'forcing'], 'solar': ['forcing'], 'sulfur': ['emissions'], 'unspecified': ['forcing'], 'volcanic': ['forcing']}

Valid run modes for each species.

fair.structure.units

Define species units and conversions.

fair.structure.units.desired_concentration_units = {'C2F6': 'ppt', 'C3F8': 'ppt', 'C4F10': 'ppt', 'C5F12': 'ppt', 'C6F14': 'ppt', 'C7F16': 'ppt', 'C8F18': 'ppt', 'CCl4': 'ppt', 'CF4': 'ppt', 'CFC-11': 'ppt', 'CFC-113': 'ppt', 'CFC-114': 'ppt', 'CFC-115': 'ppt', 'CFC-12': 'ppt', 'CH2Cl2': 'ppt', 'CH3Br': 'ppt', 'CH3CCl3': 'ppt', 'CH3Cl': 'ppt', 'CH4': 'ppb', 'CHCl3': 'ppt', 'CO2': 'ppm', 'Equivalent effective stratospheric chlorine': 'ppt', 'HCFC-141b': 'ppt', 'HCFC-142b': 'ppt', 'HCFC-22': 'ppt', 'HFC-125': 'ppt', 'HFC-134a': 'ppt', 'HFC-143a': 'ppt', 'HFC-152a': 'ppt', 'HFC-227ea': 'ppt', 'HFC-23': 'ppt', 'HFC-236fa': 'ppt', 'HFC-245fa': 'ppt', 'HFC-32': 'ppt', 'HFC-365mfc': 'ppt', 'HFC-4310mee': 'ppt', 'Halon-1202': 'ppt', 'Halon-1211': 'ppt', 'Halon-1301': 'ppt', 'Halon-2402': 'ppt', 'N2O': 'ppb', 'NF3': 'ppt', 'SF6': 'ppt', 'SO2F2': 'ppt', 'c-C4F8': 'ppt'}

Desired concentration units for each default specie.

fair.structure.units.desired_emissions_units = {'BC': 'Mt BC/yr', 'C2F6': 'kt C2F6/yr', 'C3F8': 'kt C3F8/yr', 'C4F10': 'kt C4F10/yr', 'C5F12': 'kt C5F12/yr', 'C6F14': 'kt C6F14/yr', 'C7F16': 'kt C7F16/yr', 'C8F18': 'kt C8F18/yr', 'CCl4': 'kt CCl4/yr', 'CF4': 'kt CF4/yr', 'CFC-11': 'kt CFC11/yr', 'CFC-113': 'kt CFC113/yr', 'CFC-114': 'kt CFC114/yr', 'CFC-115': 'kt CFC115/yr', 'CFC-12': 'kt CFC12/yr', 'CH2Cl2': 'kt CH2Cl2/yr', 'CH3Br': 'kt CH3Br/yr', 'CH3CCl3': 'kt CH3CCl3/yr', 'CH3Cl': 'kt CH3Cl/yr', 'CH4': 'Mt CH4/yr', 'CHCl3': 'kt CHCl3/yr', 'CO': 'Mt CO/yr', 'CO2': 'Gt CO2/yr', 'CO2 AFOLU': 'Gt CO2/yr', 'CO2 FFI': 'Gt CO2/yr', 'HCFC-141b': 'kt HCFC141b/yr', 'HCFC-142b': 'kt HCFC142b/yr', 'HCFC-22': 'kt HCFC22/yr', 'HFC-125': 'kt HFC125/yr', 'HFC-134a': 'kt HFC134a/yr', 'HFC-143a': 'kt HFC143a/yr', 'HFC-152a': 'kt HFC152a/yr', 'HFC-227ea': 'kt HFC227ea/yr', 'HFC-23': 'kt HFC23/yr', 'HFC-236fa': 'kt HFC236fa/yr', 'HFC-245fa': 'kt HFC245fa/yr', 'HFC-32': 'kt HFC32/yr', 'HFC-365mfc': 'kt HFC365mfc/yr', 'HFC-4310mee': 'kt HFC4310mee/yr', 'Halon-1202': 'kt Halon1202/yr', 'Halon-1211': 'kt Halon1211/yr', 'Halon-1301': 'kt Halon1301/yr', 'Halon-2402': 'kt Halon2402/yr', 'N2O': 'Mt N2O/yr', 'NF3': 'kt NF3/yr', 'NH3': 'Mt NH3/yr', 'NOx': 'Mt NO2/yr', 'NOx aviation': 'Mt NO2/yr', 'OC': 'Mt OC/yr', 'SF6': 'kt SF6/yr', 'SO2F2': 'kt SO2F2/yr', 'Sulfur': 'Mt SO2/yr', 'VOC': 'Mt VOC/yr', 'c-C4F8': 'kt cC4F8/yr'}

Desired emissions units for each specie.