Chemical Kinetics

Kinetics Managers

Kinetics

class cantera.Kinetics(*args, **kwargs)

Bases: cantera._cantera._SolutionBase

Instances of class Kinetics are responsible for evaluating reaction rates of progress, species production rates, and other quantities pertaining to a reaction mechanism.

add_reaction(self, Reaction rxn)

Add a new reaction to this phase.

creation_rates

Creation rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

creation_rates_ddC

Calculate derivatives of species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

creation_rates_ddP

Calculate derivatives of species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions.

creation_rates_ddT

Calculate derivatives of species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions.

creation_rates_ddX

Calculate derivatives for species creation rates with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

delta_enthalpy

Change in enthalpy for each reaction [J/kmol].

delta_entropy

Change in entropy for each reaction [J/kmol/K].

delta_gibbs

Change in Gibbs free energy for each reaction [J/kmol].

delta_standard_enthalpy

Change in standard-state enthalpy (independent of composition) for each reaction [J/kmol].

delta_standard_entropy

Change in standard-state entropy (independent of composition) for each reaction [J/kmol/K].

delta_standard_gibbs

Change in standard-state Gibbs free energy (independent of composition) for each reaction [J/kmol].

derivative_settings

Property setting behavior of derivative evaluation.

For GasKinetics, the following keyword/value pairs are supported:

  • skip-third-bodies (boolean) … if False (default), third body concentrations are considered for the evaluation of derivatives

  • skip-falloff (boolean) … if True (default), third-body effects on reaction rates are not considered.

  • rtol-delta (double) … relative tolerance used to perturb properties when calculating numerical derivatives. The default value is 1e-8.

Derivative settings are updated using a dictionary:

>>> gas.derivative_settings = {"skip-falloff": True}

Passing an empty dictionary will reset all values to their defaults.

destruction_rates

Destruction rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

destruction_rates_ddC

Calculate derivatives of species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

destruction_rates_ddP

Calculate derivatives of species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions.

destruction_rates_ddT

Calculate derivatives of species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions.

destruction_rates_ddX

Calculate derivatives for species destruction rates with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

equilibrium_constants

Equilibrium constants in concentration units for all reactions.

forward_rate_constants

Forward rate constants for all reactions. The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Deprecated since version 2.6: Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin, and Glarborg, Chemically Reacting Flow, Wiley Interscience, 2003). To switch to new behavior, run ct.use_legacy_rate_constants(False).

forward_rate_constants_ddC

Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

forward_rate_constants_ddP

Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions.

forward_rate_constants_ddT

Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions.

forward_rates_of_progress

Forward rates of progress for the reactions. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

forward_rates_of_progress_ddC

Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

forward_rates_of_progress_ddP

Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

forward_rates_of_progress_ddT

Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

forward_rates_of_progress_ddX

Calculate derivatives for forward rates-of-progress with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

heat_production_rates

Get the volumetric heat production rates [W/m^3] on a per-reaction basis. The sum over all reactions results in the total volumetric heat release rate. Example: C. K. Law: Combustion Physics (2006), Fig. 7.8.6

>>> gas.heat_production_rates[1]  # heat production rate of the 2nd reaction
heat_release_rate

Get the total volumetric heat release rate [W/m^3].

is_reversible(self, int i_reaction)

True if reaction i_reaction is reversible.

Deprecated since version 2.6: Replaced by property Reaction.reversible. Example: gas.is_reversible(0) is replaced by gas.reaction(0).reversible

kinetics_model

Return type of kinetics.

kinetics_species_index(self, species, int phase=0)

The index of species species of phase phase within arrays returned by methods of class Kinetics. If species is a string, the phase argument is unused.

kinetics_species_name(self, k)

Name of the species with index k in the arrays returned by methods of class Kinetics.

kinetics_species_names

A list of all species names, corresponding to the arrays returned by methods of class Kinetics.

modify_reaction(self, int irxn, Reaction rxn)

Modify the Reaction with index irxn to have the same rate parameters as rxn. rxn must have the same reactants and products and be of the same type (for example, ElementaryReaction, FalloffReaction, PlogReaction, etc.) as the existing reaction. This method does not modify the third-body efficiencies, reaction orders, or reversibility of the reaction.

multiplier(self, int i_reaction)

A scaling factor applied to the rate coefficient for reaction i_reaction. Can be used to carry out sensitivity analysis or to selectively disable a particular reaction. See set_multiplier.

n_phases

Number of phases in the reaction mechanism.

n_reactions

Number of reactions in the reaction mechanism.

n_total_species

Total number of species in all phases participating in the kinetics mechanism.

net_production_rates

Net production rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

net_production_rates_ddC

Calculate derivatives of species net production rates with respect to molar density at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

net_production_rates_ddP

Calculate derivatives of species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions.

net_production_rates_ddT

Calculate derivatives of species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions.

net_production_rates_ddX

Calculate derivatives for species net production rates with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

net_rates_of_progress

Net rates of progress for the reactions. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

net_rates_of_progress_ddC

Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

net_rates_of_progress_ddP

Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

net_rates_of_progress_ddT

Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

net_rates_of_progress_ddX

Calculate derivatives for net rates-of-progress with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

product_stoich_coeff(self, k_spec, int i_reaction)

The stoichiometric coefficient of species k_spec as a product in reaction i_reaction.

product_stoich_coeffs(self)

The array of product stoichiometric coefficients. Element [k,i] of this array is the product stoichiometric coefficient of species k in reaction i.

Deprecated since version 2.6: Behavior to change after Cantera 2.6; for new behavior, see property Kinetics.reactant_stoich_coeffs3.

product_stoich_coeffs3

The array of product stoichiometric coefficients. Element [k,i] of this array is the product stoichiometric coefficient of species k in reaction i.

For sparse output, set ct.use_sparse(True).

product_stoich_coeffs_reversible

The array of product stoichiometric coefficients of reversible reactions. Element [k,i] of this array is the product stoichiometric coefficient of species k in reaction i.

For sparse output, set ct.use_sparse(True).

products(self, int i_reaction)

The products portion of the reaction equation

Deprecated since version 2.6: Replaced by property Reaction.products. Example: gas.products(0) is replaced by gas.reaction(0).products

reactant_stoich_coeff(self, k_spec, int i_reaction)

The stoichiometric coefficient of species k_spec as a reactant in reaction i_reaction.

reactant_stoich_coeffs(self)

The array of reactant stoichiometric coefficients. Element [k,i] of this array is the reactant stoichiometric coefficient of species k in reaction i.

Deprecated since version 2.6: Behavior to change after Cantera 2.6; for new behavior, see property Kinetics.reactant_stoich_coeffs3.

reactant_stoich_coeffs3

The array of reactant stoichiometric coefficients. Element [k,i] of this array is the reactant stoichiometric coefficient of species k in reaction i.

For sparse output, set ct.use_sparse(True).

reactants(self, int i_reaction)

The reactants portion of the reaction equation

Deprecated since version 2.6: Replaced by property Reaction.reactants. Example: gas.reactants(0) is replaced by gas.reaction(0).reactants

reaction(self, int i_reaction)

Return a Reaction object representing the reaction with index i_reaction. Changes to this object do not affect the Kinetics or Solution object until the modify_reaction function is called.

reaction_equation(self, int i_reaction)

The equation for the specified reaction. See also reaction_equations.

Deprecated since version 2.6: Replaced by property Reaction.equation. Example: gas.reaction_equation(0) is replaced by gas.reaction(0).equation

reaction_equations(self, indices=None)

Returns a list containing the reaction equation for all reactions in the mechanism if indices is unspecified, or the equations for each reaction in the sequence indices. For example:

>>> gas.reaction_equations()
['2 O + M <=> O2 + M', 'O + H + M <=> OH + M', 'O + H2 <=> H + OH', ...]
>>> gas.reaction_equations([2,3])
['O + H + M <=> OH + M', 'O + H2 <=> H + OH']

See also reaction_equation.

reaction_phase_index

The index of the phase where the reactions occur.

reaction_type(self, int i_reaction)

Type code of reaction i_reaction.

Deprecated since version 2.6: Replaced by properties Reaction.type and Reaction.rate.type. Example: gas.reaction_type(0) is replaced by gas.reaction(0).reaction_type and gas.reaction(0).rate.type

reactions(self)

Return a list of all Reaction objects. Changes to these objects do not affect the Kinetics or Solution object until the modify_reaction function is called.

reverse_rate_constants

Reverse rate constants for all reactions. The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Deprecated since version 2.6: Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin and Glarborg, Chemically Reacting Flow, Wiley Interscience, 2003). To switch to new behavior, run ct.use_legacy_rate_constants(False).

reverse_rates_of_progress

Reverse rates of progress for the reactions. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

reverse_rates_of_progress_ddC

Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

reverse_rates_of_progress_ddP

Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

reverse_rates_of_progress_ddT

Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

reverse_rates_of_progress_ddX

Calculate derivatives for reverse rates-of-progress with respect to species concentrations at constant temperature, pressure and molar concentration. For sparse output, set ct.use_sparse(True).

Note that for derivatives with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

set_multiplier(self, double value, int i_reaction=-1)

Set the multiplier for for reaction i_reaction to value. If i_reaction is not specified, then the multiplier for all reactions is set to value. See multiplier.

third_body_concentrations

Effective third-body concentrations used by individual reactions; values are only defined for reactions involving third-bodies and are set to not-a-number otherwise.

InterfaceKinetics

class cantera.InterfaceKinetics(infile='', name='', adjacent=(), *args, **kwargs)

Bases: cantera._cantera.Kinetics

A kinetics manager for heterogeneous reaction mechanisms. The reactions are assumed to occur at an interface between bulk phases.

advance_coverages(self, double dt, double rtol=1e-7, double atol=1e-14, double max_step_size=0, int max_steps=20000, int max_error_test_failures=7)

This method carries out a time-accurate advancement of the surface coverages for a specified amount of time.

advance_coverages_to_steady_state(self)

This method advances the surface coverages to steady state.

get_creation_rates(self, phase)

Creation rates for each species in phase phase. Use the creation_rates property to get the creation rates for species in all phases.

get_destruction_rates(self, phase)

Destruction rates for each species in phase phase. Use the destruction_rates property to get the destruction rates for species in all phases.

get_net_production_rates(self, phase)

Net production rates for each species in phase phase. Use the net_production_rates property to get the net production rates for species in all phases.

phase_index(self, phase)

Get the index of the phase phase, where phase may specified using the phase object, the name, or the index itself.

write_yaml(self, filename, phases=None, units=None, precision=None, skip_user_defined=None)

See _SolutionBase.write_yaml.

Reactions

These classes contain the definition of a single reaction, independent of a specific Kinetics object. For legacy objects (CTI/XML input), each class integrates associated rate expressions, whereas for the new, YAML-based implementation, reaction rate evaluation is handled by dedicated ReactionRate objects.

Reaction

class cantera.Reaction(reactants=None, products=None, rate=None, equation=None, *, init=True, legacy=False, **kwargs)

Bases: object

A class which stores data about a reaction and its rate parameterization so that it can be added to a Kinetics object.

Parameters
  • reactants – Value used to set reactants

  • products – Value used to set products

  • rate

    The rate parameterization for the reaction, given as one of the following:

    • a ReactionRate object

    • a dict containing the parameters needed to construct a ReactionRate object, with keys corresponding to the YAML format

    • a dict containing Arrhenius parameters (A, b, and Ea)

  • equation – The reaction equation, used to set the reactants and products if values for those arguments are not provided.

Examples:

R = ct.Reaction({"O": 1, "H2": 1}, {"H": 1, "OH": 1},
                ct.ArrheniusRate(38.7, 2.7, 26191840.0))
R = ct.Reaction(equation="O + H2 <=> H + OH",
                rate={"A": 38.7, "b", 2.7, "Ea": 26191840.0})
R = ct.Reaction(equation="HO2 <=> OH + O", rate=ChebyshevRate(...))

The static methods list_from_file, list_from_yaml, listFromCti, and listFromXml can be used to create lists of Reaction objects from existing definitions in the YAML, CTI, or XML formats. All of the following will produce a list of the 325 reactions which make up the GRI 3.0 mechanism:

R = ct.Reaction.list_from_file("gri30.yaml", gas)
R = ct.Reaction.listFromCti(open("path/to/gri30.cti").read())
R = ct.Reaction.listFromXml(open("path/to/gri30.xml").read())

where gas is a Solution object with the appropriate thermodynamic model, which is the ideal-gas model in this case.

The static method list_from_yaml can be used to create lists of Reaction objects from a YAML list:

rxns = '''
  - equation: O + H2 <=> H + OH
    rate-constant: {A: 3.87e+04, b: 2.7, Ea: 6260.0}
  - equation: O + HO2 <=> OH + O2
    rate-constant: {A: 2.0e+13, b: 0.0, Ea: 0.0}
'''
R = ct.Reaction.list_from_yaml(rxns, gas)

The methods from_yaml, fromCti, and fromXml can be used to create individual Reaction objects from definitions in these formats. In the case of using YAML or CTI definitions, it is important to verify that either the pre-exponential factor and activation energy are supplied in SI units, or that they have their units specified:

R = ct.Reaction.from_yaml('''{equation: O + H2 <=> H + OH,
        rate-constant: {A: 3.87e+04 cm^3/mol/s, b: 2.7, Ea: 6260 cal/mol}}''',
        gas)

R = ct.Reaction.fromCti('''reaction('O + H2 <=> H + OH',
        [3.87e1, 2.7, 2.619184e7])''')

R = ct.Reaction.fromCti('''reaction('O + H2 <=> H + OH',
                [(3.87e4, 'cm3/mol/s'), 2.7, (6260, 'cal/mol')])''')
ID

Get/Set the identification string for the reaction, which can be used in filtering operations.

Pmax

Maximum pressure [K] for the Chebyshev fit

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.pressure_range[1] for reactions where the rate is a a ChebyshevRate.

Pmin

Minimum pressure [Pa] for the Chebyshev fit

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.pressure_range[0] for reactions where the rate is a a ChebyshevRate.

Tmax

Maximum temperature [K] for the Chebyshev fit

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.temperature_range[1] for reactions where the rate is a ChebyshevRate.

Tmin

Minimum temperature [K] for the Chebyshev fit

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.temperature_range[0] for reactions where the rate is a ChebyshevRate.

allow_negative_orders

Get/Set a flag which is True if negative reaction orders are allowed. Default is False.

allow_negative_pre_exponential_factor

Get/Set whether the rate coefficient is allowed to have a negative pre-exponential factor.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ArrheniusRateBase.allow_negative_pre_exponential_factor.

allow_nonreactant_orders

Get/Set a flag which is True if reaction orders can be specified for non-reactant species. Default is False.

clear_user_data(self)

Clear all saved input data, so that the data given by input_data or Solution.write_yaml will only include values generated by Cantera based on the current object state.

coeffs

2D array of Chebyshev coefficients of size (n_temperature, n_pressure).

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.data for reactions where the rate is a ChebyshevRate.

coverage_deps

Get/Set a dict containing adjustments to the Arrhenius rate expression dependent on surface species coverages. The keys of the dict are species names, and the values are tuples specifying the three coverage parameters (a, m, E) which are the modifiers for the pre-exponential factor [m, kmol, s units], the temperature exponent [nondimensional], and the activation energy [J/kmol], respectively.

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated InterfaceReaction class. Replaced by Reaction.rate.coverage_dependencies for reactions where the rate is a InterfaceArrheniusRate or StickingArrheniusRate.

duplicate

Get/Set a flag which is True if this reaction is marked as a duplicate or False otherwise.

equation

A string giving the chemical equation for this reaction. Determined automatically based on reactants and products.

static fromCti(text)

Create a Reaction object from its CTI string representation.

Deprecated since version 2.5: The CTI input format is deprecated and will be removed in Cantera 3.0.

static fromXml(text)

Create a Reaction object from its XML string representation.

Deprecated since version 2.5: The XML input format is deprecated and will be removed in Cantera 3.0.

fromYaml(type cls, text, Kinetics kinetics=None)

Create a Reaction object from its YAML string representation.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by Reaction.from_yaml.

from_dict(type cls, data, Kinetics kinetics)

Create a Reaction object from a dictionary corresponding to its YAML representation.

An example for the creation of a Reaction from a dictionary is:

rxn = Reaction.from_dict(
    {"equation": "O + H2 <=> H + OH",
     "rate-constant": {"A": 38.7, "b": 2.7, "Ea": 26191840.0}},
    kinetics=gas)

In the example, gas is a Kinetics (or Solution) object.

Parameters
  • data – A dictionary corresponding to the YAML representation.

  • kinetics – A Kinetics object whose associated phase(s) contain the species involved in the reaction.

from_yaml(type cls, text, Kinetics kinetics)

Create a Reaction object from its YAML string representation.

An example for the creation of a Reaction from a YAML string is:

rxn = Reaction.from_yaml('''
    equation: O + H2 <=> H + OH
    rate-constant: {A: 38.7, b: 2.7, Ea: 6260.0 cal/mol}
    ''', kinetics=gas)

In the example, gas is a Kinetics (or Solution) object.

Parameters
  • text – The YAML reaction string.

  • kinetics – A Kinetics object whose associated phase(s) contain the species involved in the reaction.

input_data

Get input data for this reaction with its current parameter values, along with any user-specified data provided with its input (YAML) definition.

is_sticking_coefficient

Get/Set a boolean indicating if the rate coefficient for this reaction is expressed as a sticking coefficient rather than the forward rate constant.

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated InterfaceReaction class. Replaced by dedicated rate objects InterfaceArrheniusRate and StickingArrheniusRate.

static listFromCti(text)

Create a list of Reaction objects from all the reactions defined in a CTI string.

Deprecated since version 2.5: The CTI input format is deprecated and will be removed in Cantera 3.0.

static listFromFile(filename, Kinetics kinetics=None, section=u'reactions')

Create a list of Reaction objects from all of the reactions defined in a YAML, CTI, or XML file.

For YAML input files, a Kinetics object is required as the second argument, and reactions from the section section will be returned.

Directories on Cantera’s input file path will be searched for the specified file.

In the case of an XML file, the <reactions> nodes are assumed to be children of the <reactionsData> node in a document with a <ctml> root node, as in the XML files produced by conversion from CTI files.

Deprecated since version 2.5: The CTI and XML input formats are deprecated and will be removed in Cantera 3.0.

Deprecated since version 2.6: To be removed after Cantera 2.6. Replaced by Reaction.list_from_file.

static listFromXml(text)

Create a list of Reaction objects from all the reaction defined in an XML string. The <reaction> nodes are assumed to be children of the <reactionData> node in a document with a <ctml> root node, as in the XML files produced by conversion from CTI files.

Deprecated since version 2.5: The XML input format is deprecated and will be removed in Cantera 3.0.

static listFromYaml(text, Kinetics kinetics)

Create a list of Reaction objects from all the reactions defined in a YAML string.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by Reaction.list_from_yaml.

static list_from_file(filename, Kinetics kinetics, section=u'reactions')

Create a list of Reaction objects from all of the reactions defined in a YAML file. Reactions from the section section will be returned.

Directories on Cantera’s input file path will be searched for the specified file.

static list_from_yaml(text, Kinetics kinetics)

Create a list of Reaction objects from all the reactions defined in a YAML string.

nPressure

Number of pressures over which the Chebyshev fit is computed

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.n_pressure for reactions where the rate is a a ChebyshevRate.

nTemperature

Number of temperatures over which the Chebyshev fit is computed

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by Reaction.rate.n_temperature for reactions where the rate is a ChebyshevRate.

orders

Get/Set the reaction order with respect to specific species as a dict with species names as the keys and orders as the values, or as a composition string. By default, mass-action kinetics is assumed, with the reaction order for each reactant species equal to each its stoichiometric coefficient.

product_string

A string representing the products side of the chemical equation for this reaction. Determined automatically based on products.

products

Get/Set the products in this reaction as a dict where the keys are species names and the values, are the stoichiometric coefficients, for example {'CH3':1, 'H2O':1}, or as a composition string, for example 'CH3:1, H2O:1'.

rate

Get/Set the ArrheniusRate rate coefficient for this reaction.

rate_coeff_units

Get reaction rate coefficient units

rates

For reactions with Plog rates, get/set the rate coefficients for this reaction as a list of (pressure, rate) tuples.

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated PlogReaction class. Replaced by Reaction.rate.rates for reactions where the rate is a PlogRate.

reactant_string

A string representing the reactants side of the chemical equation for this reaction. Determined automatically based on reactants.

reactants

Get/Set the reactants in this reaction as a dict where the keys are species names and the values, are the stoichiometric coefficients, for example {'CH4':1, 'OH':1}, or as a composition string, for example 'CH4:1, OH:1'.

reaction_type

Retrieve the native type name of the reaction.

reversible

Get/Set a flag which is True if this reaction is reversible or False otherwise.

set_parameters(self, Tmin, Tmax, Pmin, Pmax, coeffs)

For Chebyshev reactions, simultaneously set values for Tmin, Tmax, Pmin, Pmax, and coeffs.

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated ChebyshevReaction class. Replaced by ChebyshevRate constructor for reactions where the rate is a ChebyshevRate.

sticking_species

The name of the sticking species. Needed only for reactions with multiple non-surface reactant species, where the sticking species is ambiguous.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property StickingArrheniusRate.sticking_species.

update_user_data(self, data)

Add the contents of the provided dict as additional fields when generating YAML phase definition files with Solution.write_yaml or in the data returned by input_data. Existing keys with matching names are overwritten.

use_motz_wise_correction

Get/Set a boolean indicating whether to use the correction factor developed by Motz & Wise for reactions with high (near-unity) sticking coefficients when converting the sticking coefficient to a rate coefficient.

Deprecated since version 2.6: This property is for temporary backwards-compatibility with the deprecated InterfaceReaction class. Replaced by Reaction.rate.motz_wise_correction for reactions where the rate is a StickingArrheniusRate.

uses_legacy

Indicate whether reaction uses a legacy implementation

ElementaryReaction

class cantera.ElementaryReaction(reactants=None, products=None, rate=None, equation=None, *, Kinetics kinetics=None, init=True, **kwargs)

Bases: cantera._cantera.Reaction

A reaction which follows mass-action kinetics with a modified Arrhenius reaction rate.

An example for the definition of an ElementaryReaction object is given as:

rxn = ElementaryReaction(
    equation="O + H2 <=> H + OH",
    rate={"A": 38.7, "b": 2.7, "Ea": 2.619184e+07},
    kinetics=gas)

The YAML description corresponding to this reaction is:

equation: O + H2 <=> H + OH
rate-constant: {A: 3.87e+04 cm^3/mol/s, b: 2.7, Ea: 6260.0 cal/mol}

Deprecated since version 2.6: To be removed after Cantera 2.6. Capabilities merged directly into the base Reaction class.

allow_negative_pre_exponential_factor

Get/Set whether the rate coefficient is allowed to have a negative pre-exponential factor.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ArrheniusRateBase.allow_negative_pre_exponential_factor.

rate

Get/Set the ArrheniusRate rate coefficient for this reaction.

ThreeBodyReaction

class cantera.ThreeBodyReaction(reactants=None, products=None, rate=None, equation=None, *, efficiencies=None, Kinetics kinetics=None, legacy=False, init=True, **kwargs)

Bases: cantera._cantera.ElementaryReaction

A reaction with a non-reacting third body “M” that acts to add or remove energy from the reacting species.

An example for the definition of an ThreeBodyReaction object is given as:

rxn = ThreeBodyReaction(
    equation="2 O + M <=> O2 + M",
    rate={"A": 1.2e+17, "b": -1.0, "Ea": 0.0},
    efficiencies={"H2": 2.4, "H2O": 15.4, "AR": 0.83},
    kinetics=gas)

The YAML description corresponding to this reaction is:

equation: 2 O + M <=> O2 + M
type: three-body
rate-constant: {A: 1.2e+17 cm^6/mol^2/s, b: -1.0, Ea: 0.0 cal/mol}
efficiencies: {H2: 2.4, H2O: 15.4, AR: 0.83}
default_efficiency

Get/Set the default third-body efficiency for this reaction, used for species used for species not in efficiencies.

efficiencies

Get/Set a dict defining non-default third-body efficiencies for this reaction, where the keys are the species names and the values are the efficiencies.

efficiency(self, species)

Get the efficiency of the third body named species considering both the default efficiency and species-specific efficiencies.

FalloffReaction

class cantera.FalloffReaction(reactants=None, products=None, rate=None, equation=None, *, efficiencies=None, Kinetics kinetics=None, init=True, legacy=False, **kwargs)

Bases: cantera._cantera.Reaction

A reaction that is first-order in [M] at low pressure, like a third-body reaction, but zeroth-order in [M] as pressure increases.

An example for the definition of a FalloffReaction object is given as:

rxn = FalloffReaction(
    equation="2 OH (+ M) <=> H2O2 (+ M)",
    rate=ct.TroeRate(low=ct.Arrhenius(2.3e+12, -0.9, -7112800.0),
                     high=ct.Arrhenius(7.4e+10, -0.37, 0),
                     falloff_coeffs=[0.7346, 94.0, 1756.0, 5182.0]),
    efficiencies={"AR": 0.7, "H2": 2.0, "H2O": 6.0},
    kinetics=gas)

The YAML description corresponding to this reaction is:

equation: 2 OH (+ M) <=> H2O2 (+ M)  # Reaction 3
type: falloff
low-P-rate-constant: {A: 2.3e+12, b: -0.9, Ea: -1700.0 cal/mol}
high-P-rate-constant: {A: 7.4e+10, b: -0.37, Ea: 0.0 cal/mol}
Troe: {A: 0.7346, T3: 94.0, T1: 1756.0, T2: 5182.0}
efficiencies: {AR: 0.7, H2: 2.0, H2O: 6.0}
allow_negative_pre_exponential_factor

Get/Set whether the rate coefficient is allowed to have a negative pre-exponential factor.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property FalloffRate.allow_negative_pre_exponential_factor.

default_efficiency

Get/Set the default third-body efficiency for this reaction, used for species used for species not in efficiencies.

efficiencies

Get/Set a dict defining non-default third-body efficiencies for this reaction, where the keys are the species names and the values are the efficiencies.

efficiency(self, species)

Get the efficiency of the third body named species considering both the default efficiency and species-specific efficiencies.

falloff

Get/Set the Falloff function used to blend the high- and low-pressure rate coefficients

high_rate

Get/Set the Arrhenius rate constant in the high-pressure limit

low_rate

Get/Set the Arrhenius rate constant in the low-pressure limit

rate

Get/Set the FalloffRate rate coefficients for this reaction.

ChemicallyActivatedReaction

class cantera.ChemicallyActivatedReaction

Bases: cantera._cantera.FalloffReaction

A reaction where the rate decreases as pressure increases due to collisional stabilization of a reaction intermediate. Like a FalloffReaction, except that the forward rate constant is written as being proportional to the low- pressure rate constant.

PlogReaction

class cantera.PlogReaction(reactants=None, products=None, rate=None, equation=None, *, Kinetics kinetics=None, init=True, **kwargs)

Bases: cantera._cantera.Reaction

A pressure-dependent reaction parameterized by logarithmically interpolating between Arrhenius rate expressions at various pressures.

An example for the definition of a PlogReaction object is given as:

rxn = PlogReaction(
    equation="H2 + O2 <=> 2 OH",
    rate=[(1013.25, Arrhenius(1.2124e+16, -0.5779, 45491376.8)),
          (101325., Arrhenius(4.9108e+31, -4.8507, 103649395.2)),
          (1013250., Arrhenius(1.2866e+47, -9.0246, 166508556.0)),
          (10132500., Arrhenius(5.9632e+56, -11.529, 220076726.4))],
    kinetics=gas)

The YAML description corresponding to this reaction is:

equation: H2 + O2 <=> 2 OH
type: pressure-dependent-Arrhenius
rate-constants:
- {P: 0.01 atm, A: 1.2124e+16, b: -0.5779, Ea: 1.08727e+04 cal/mol}
- {P: 1.0 atm, A: 4.9108e+31, b: -4.8507, Ea: 2.47728e+04 cal/mol}
- {P: 10.0 atm, A: 1.2866e+47, b: -9.0246, Ea: 3.97965e+04 cal/mol}
- {P: 100.0 atm, A: 5.9632e+56, b: -11.529, Ea: 5.25996e+04 cal/mol}or.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Implemented by the Reaction class with a PlogRate reaction rate.

rates

Get/Set the rate coefficients for this reaction, which are given as a list of (pressure, Arrhenius) tuples.

ChebyshevReaction

class cantera.ChebyshevReaction(reactants=None, products=None, rate=None, equation=None, *, Kinetics kinetics=None, init=True, **kwargs)

Bases: cantera._cantera.Reaction

A pressure-dependent reaction parameterized by a bivariate Chebyshev polynomial in temperature and pressure.

An example for the definition of a ChebyshevReaction object is given as:

rxn = ChebyshevReaction(
    equation="HO2 <=> OH + O",
    rate={"temperature-range": [290.0, 3000.0],
          "pressure-range": [1e3, 1e8],
          "data": [[8.2883, -1.1397, -0.12059, 0.016034],
                   [1.9764, 1.0037, 7.2865e-03, -0.030432],
                   [0.3177, 0.26889, 0.094806, -7.6385e-03]]},
    kinetics=gas)

The YAML description corresponding to this reaction is:

equation: HO2 <=> OH + O
type: Chebyshev
temperature-range: [290.0, 3000.0]
pressure-range: [1.e-03 bar, 10. bar]
data:
- [8.2883, -1.1397, -0.12059, 0.016034]
- [1.9764, 1.0037, 7.2865e-03, -0.030432]
- [0.3177, 0.26889, 0.094806, -7.6385e-03]

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Implemented by the Reaction class with a ChebyshevRate reaction rate.

Pmax

Maximum pressure [Pa] for the Chebyshev fit

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.pressure_range[1].

Pmin

Minimum pressure [Pa] for the Chebyshev fit

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.pressure_range[0].

Tmax

Maximum temperature [K] for the Chebyshev fit

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.temperature_range[1].

Tmin

Minimum temperature [K] for the Chebyshev fit

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.temperature_range[0].

coeffs

2D array of Chebyshev coefficients of size (n_temperature, n_pressure).

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.data.

nPressure

Number of pressures over which the Chebyshev fit is computed

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.n_pressure.

nTemperature

Number of temperatures over which the Chebyshev fit is computed

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property ChebyshevRate.n_temperature.

set_parameters(self, Tmin, Tmax, Pmin, Pmax, coeffs)

Simultaneously set values for Tmin, Tmax, Pmin, Pmax, and coeffs.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by ChebyshevRate constructor.

InterfaceReaction

class cantera.InterfaceReaction(equation=None, rate=None, Kinetics kinetics=None, init=True, legacy=False, **kwargs)

Bases: cantera._cantera.ElementaryReaction

A reaction occurring on an Interface (that is, a surface or an edge)

rxn = InterfaceReaction(

equation=”H(S) + O(S) <=> OH(S) + PT(S)”, rate={“A”: 3.7e+20, “b”: 0, “Ea”: 1.15e7}, kinetics=surf)

The YAML description corresponding to this reaction is:

equation: H(S) + O(S) <=> OH(S) + PT(S)
rate-constant: {A: 3.7e+20, b: 0, Ea: 11500 J/mol}

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Implemented by the Reaction class with either InterfaceArrheniusRate or StickingArrheniusRate reaction rate.

coverage_deps

Get/Set a dict containing adjustments to the Arrhenius rate expression dependent on surface species coverages. The keys of the dict are species names, and the values are tuples specifying the three coverage parameters (a, m, E) which are the modifiers for the pre-exponential factor [m, kmol, s units], the temperature exponent [nondimensional], and the activation energy [J/kmol], respectively.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property InterfaceRateBase.coverage_dependencies.

is_sticking_coefficient

Get/Set a boolean indicating if the rate coefficient for this reaction is expressed as a sticking coefficient rather than the forward rate constant.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by dedicated classes StickRateBase.

sticking_species

The name of the sticking species. Needed only for reactions with multiple non-surface reactant species, where the sticking species is ambiguous.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property stickRateBase.sticking_species.

use_motz_wise_correction

Get/Set a boolean indicating whether to use the correction factor developed by Motz & Wise for reactions with high (near-unity) sticking coefficients when converting the sticking coefficient to a rate coefficient.

Deprecated since version 2.6: To be deprecated with version 2.6, and removed thereafter. Replaced by property stickRateBase.mote_wise_correction.

Reaction Rates

ReactionRate

class cantera.ReactionRate

Bases: object

Base class for ReactionRate objects.

ReactionRate objects are used to calculate reaction rates and are associated with a Reaction object.

from_dict(type cls, data)

Create a ReactionRate object from a dictionary corresponding to its YAML representation.

An example for the creation of a ReactionRate from a dictionary is:

rate = ReactionRate.from_dict(
    {"rate-constant": {"A": 38.7, "b": 2.7, "Ea": 26191840.0}})
Parameters

data – A dictionary corresponding to the YAML representation.

from_yaml(type cls, text)

Create a ReactionRate object from its YAML string representation.

An example for the creation of a ReactionRate from a YAML string is:

rate = ReactionRate.from_yaml(
    "rate-constant: {A: 38.7, b: 2.7, Ea: 6260.0 cal/mol}")

Units for A require a unit system with length in m and quantity in kmol (standard Cantera units).

Parameters

text – The YAML reaction rate string.

input_data

Get input data for this reaction rate with its current parameter values.

type

Get the C++ ReactionRate type

ArrheniusRateBase

class cantera.ArrheniusRateBase(input_data)

Bases: cantera._cantera.ReactionRate

Base class collecting commonly used features of Arrhenius-type rate objects. Objects should be instantiated by specialized classes, for example ArrheniusRate, BlowersMaselRate and TwoTempPlasmaRate.

activation_energy

The activation energy E [J/kmol].

allow_negative_pre_exponential_factor

Get/Set whether the rate coefficient is allowed to have a negative pre-exponential factor.

pre_exponential_factor

The pre-exponential factor A in units of m, kmol, and s raised to powers depending on the reaction order.

temperature_exponent

The temperature exponent b.

ArrheniusRate

class cantera.ArrheniusRate(A, b, Ea)

Bases: cantera._cantera.ArrheniusRateBase

A reaction rate coefficient which depends on temperature only and follows the modified Arrhenius form:

\[k_f = A T^b \exp(-\tfrac{E_a}{RT})\]

where A is the pre_exponential_factor, b is the temperature_exponent, and Ea is the activation_energy.

BlowersMaselRate

class cantera.BlowersMaselRate(A, b, Ea, w)

Bases: cantera._cantera.ArrheniusRateBase

A reaction rate coefficient which depends on temperature and enthalpy change of the reaction follows the Blowers-Masel approximation and modified Arrhenius form described in ArrheniusRate.

bond_energy

Average bond dissociation energy of the bond being formed and broken in the reaction E [J/kmol].

delta_enthalpy

Enthalpy change of reaction deltaH [J/kmol]

The enthalpy change of reaction is a function of temperature and thus not an independent property. Accordingly, the setter should be only used for testing purposes, as any value will be overwritten by an update of the thermodynamic state.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

TwoTempPlasmaRate

class cantera.TwoTempPlasmaRate(A, b, Ea_gas, Ea_electron)

Bases: cantera._cantera.ArrheniusRateBase

A reaction rate coefficient which depends on both gas and electron temperature with the form similar to the modified Arrhenius form. Specifically, the temperature exponent (b) is applied to the electron temperature instead. In addition, the exponential term with activation energy for electron is included.

\[k_f = A T_e^b \exp(-\tfrac{E_{a,g}}{RT}) \exp(\tfrac{E_{a,e}(T_e - T)}{R T T_e})\]

where \(A\) is the pre_exponential_factor, \(b\) is the temperature_exponent, \(E_{a,g}\) (Ea_gas) is the activation_energy, and \(E_{a,e}\) (Ea_electron) is the activation_electron_energy.

activation_electron_energy

The activation electron energy \(E_{a,e}\) [J/kmol].

FalloffRate

class cantera.FalloffRate

Bases: cantera._cantera.ReactionRate

Base class for parameterizations used to describe the fall-off in reaction rates due to intermolecular energy transfer. These objects are used by FalloffReaction. Note that FalloffRate is a base class for specialized fall-off parameterizations and cannot be instantiated by itself.

allow_negative_pre_exponential_factor

Get/Set whether the rate coefficient is allowed to have a negative pre-exponential factor.

chemically_activated

Get whether the object is a chemically-activated reaction rate.

falloff_coeffs

The array of coefficients used to define this falloff function.

falloff_function(self, double temperature, double conc3b)

Evaluate the falloff function based on temperature and third-body concentration.

high_rate

Get/Set the Arrhenius rate constant in the high-pressure limit

low_rate

Get/Set the Arrhenius rate constant in the low-pressure limit

LindemannRate

class cantera.LindemannRate(low, high, falloff_coeffs)

Bases: cantera._cantera.FalloffRate

The Lindemann falloff parameterization.

This class implements the simple falloff function \(F(T,P_r) = 1.0\).

TroeRate

class cantera.TroeRate(low, high, falloff_coeffs)

Bases: cantera._cantera.FalloffRate

The 3- or 4-parameter Troe falloff function.

Parameters

falloff_coeffs – An array of 3 or 4 parameters: \([a, T^{***}, T^*, T^{**}]\) where the final parameter is optional (with a default value of 0).

SriRate

class cantera.SriRate(low, high, falloff_coeffs)

Bases: cantera._cantera.FalloffRate

The 3- or 5-parameter SRI falloff function.

Parameters

falloff_coeffs – An array of 3 or 5 parameters: \([a, b, c, d, e]\) where the last two parameters are optional (with default values of 1 and 0, respectively).

TsangRate

class cantera.TsangRate(low, high, falloff_coeffs)

Bases: cantera._cantera.FalloffRate

The Tsang falloff parameterization.

PlogRate

class cantera.PlogRate(rates)

Bases: cantera._cantera.ReactionRate

A pressure-dependent reaction rate parameterized by logarithmically interpolating between Arrhenius rate expressions at various pressures.

rates

Get/Set the rate coefficients for this reaction, which are given as a list of (pressure, Arrhenius) tuples.

ChebyshevRate

class cantera.ChebyshevRate(temperature_range, pressure_range, data)

Bases: cantera._cantera.ReactionRate

A pressure-dependent reaction rate parameterized by a bivariate Chebyshev polynomial in temperature and pressure.

data

2D array of Chebyshev coefficients where rows and columns correspond to temperature and pressure dimensions over which the Chebyshev fit is computed.

n_pressure

Number of pressures over which the Chebyshev fit is computed (same as number of columns of data property).

n_temperature

Number of temperatures over which the Chebyshev fit is computed. (same as number of rows of data property).

pressure_range

Valid pressure range [Pa] for the Chebyshev fit

temperature_range

Valid temperature range [K] for the Chebyshev fit

CustomRate

class cantera.CustomRate(k)

Bases: cantera._cantera.ReactionRate

A custom rate coefficient which depends on temperature only.

The simplest way to create a CustomRate object is to use a lambda function, for example:

rr = CustomRate(lambda T: 38.7 * T**2.7 * exp(-3150.15/T))

Warning: this class is an experimental part of the Cantera API and may be changed or removed without notice.

set_rate_function(self, k)

Set the function describing a custom reaction rate:

rr = CustomRate()
rr.set_rate_function(lambda T: 38.7 * T**2.7 * exp(-3150.15/T))

InterfaceRateBase

class cantera.InterfaceRateBase

Bases: cantera._cantera.ArrheniusRateBase

Base class collecting commonly used features of Arrhenius-type rate objects that include coverage dependencies.

beta

Return the charge transfer beta parameter

coverage_dependencies

Get/set a dictionary containing adjustments to the Arrhenius rate expression dependent on surface species coverages. The keys of the dictionary are species names, and the values are dictionaries specifying the three coverage parameters a, m and E which are the modifiers for the pre-exponential factor [m, kmol, s units], the temperature exponent [nondimensional], and the activation energy [J/kmol], respectively.

set_species(self, species)

Set association with an ordered list of all species associated with an InterfaceKinetics object.

site_density

Site density [kmol/m^2]

The site density is not an independent property, as it is set by an associated InterfaceKinetics object. Accordingly, the setter should be only used for testing purposes, as the value will be overwritten by an update of the thermodynamic state.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

uses_electrochemistry

Return boolean flag indicating whether rate involves a charge transfer.

InterfaceArrheniusRate

class cantera.InterfaceArrheniusRate(A, b, Ea)

Bases: cantera._cantera.InterfaceRateBase

A reaction rate coefficient which depends on temperature and interface coverage

InterfaceBlowersMaselRate

class cantera.InterfaceBlowersMaselRate(A, b, Ea0, w)

Bases: cantera._cantera.InterfaceRateBase

A reaction rate coefficient which depends on temperature and enthalpy change of the reaction follows the Blowers-Masel approximation and modified Arrhenius form described in ArrheniusRate.

bond_energy

Average bond dissociation energy of the bond being formed and broken in the reaction E [J/kmol].

delta_enthalpy

Enthalpy change of reaction deltaH [J/kmol]

The enthalpy change of reaction is a function of temperature and thus not an independent property. Accordingly, the setter should be only used for testing purposes, as any value will be overwritten by an update of the thermodynamic state.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

StickRateBase

class cantera.StickRateBase

Bases: cantera._cantera.InterfaceRateBase

Base class collecting commonly used features of Arrhenius-type sticking rate objects that include coverage dependencies.

motz_wise_correction

Get/Set a boolean indicating whether to use the correction factor developed by Motz & Wise for reactions with high (near-unity) sticking coefficients when converting the sticking coefficient to a rate coefficient.

sticking_order

The exponent applied to site density (sticking order).

The sticking order is not an independent property and is detected automatically by Cantera. Accordingly, the setter should be only used for testing purposes.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

sticking_species

The name of the sticking species. Needed only for reactions with multiple non-surface reactant species, where the sticking species is ambiguous.

sticking_weight

The molecular weight of the sticking species.

The sticking weight is not an independent property and is detected automatically by Cantera. Accordingly, the setter should be only used for testing purposes.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

StickingArrheniusRate

class cantera.StickingArrheniusRate(A, b, Ea)

Bases: cantera._cantera.StickRateBase

A surface sticking rate expression based on the Arrhenius parameterization

StickingBlowersMaselRate

class cantera.StickingBlowersMaselRate(A, b, Ea0, w)

Bases: cantera._cantera.StickRateBase

A surface sticking rate expression based on the Blowers-Masel parameterization

bond_energy

Average bond dissociation energy of the bond being formed and broken in the reaction E [J/kmol].

delta_enthalpy

Enthalpy change of reaction deltaH [J/kmol]

The enthalpy change of reaction is a function of temperature and thus not an independent property. Accordingly, the setter should be only used for testing purposes, as any value will be overwritten by an update of the thermodynamic state.

Warning: this property is an experimental part of the Cantera API and may be changed or removed without notice.

Auxiliary Reaction Data (legacy only)

Arrhenius

class cantera.Arrhenius(A, b, E)

Bases: object

A reaction rate coefficient which depends on temperature only and follows the modified Arrhenius form:

\[k_f = A T^b \exp(-\tfrac{E}{RT})\]

where A is the pre_exponential_factor, b is the temperature_exponent, and E is the activation_energy.

activation_energy

The activation energy E [J/kmol].

pre_exponential_factor

The pre-exponential factor A in units of m, kmol, and s raised to powers depending on the reaction order.

temperature_exponent

The temperature exponent b.

Falloff

class cantera.Falloff(params=(), init=True)

Bases: object

A parameterization used to describe the fall-off in reaction rate constants due to intermolecular energy transfer. These functions are used by reactions defined using the FalloffReaction and ChemicallyActivatedReaction classes.

This base class implements the simple falloff function \(F(T,P_r) = 1.0\).

Parameters
  • params – Not used for the “simple” falloff parameterization.

  • init – Used internally when wrapping FalloffRate objects returned from C++.

Deprecated since version 2.6: To be removed after Cantera 2.6. Capabilities merged into the LindemannRate class.

parameters

The array of parameters used to define this falloff function.

type

A string defining the type of the falloff parameterization

TroeFalloff

class cantera.TroeFalloff(params=(), init=True)

Bases: cantera._cantera.Falloff

The 3- or 4-parameter Troe falloff function.

Parameters

params – An array of 3 or 4 parameters: \([a, T^{***}, T^*, T^{**}]\) where the final parameter is optional (with a default value of 0).

Deprecated since version 2.6: To be removed after Cantera 2.6. Capabilities merged into the TroeRate class.

SriFalloff

class cantera.SriFalloff(params=(), init=True)

Bases: cantera._cantera.Falloff

The 3- or 5-parameter SRI falloff function.

Parameters

params – An array of 3 or 5 parameters: \([a, b, c, d, e]\) where the last two parameters are optional (with default values of 1 and 0, respectively).

Deprecated since version 2.6: To be removed after Cantera 2.6. Capabilities merged into the SriRate class.

Reaction Path Analysis

ReactionPathDiagram

class cantera.ReactionPathDiagram(Kinetics kin, str element)

Bases: object

ReactionPathDiagram(Kinetics kin, unicode element)

Create a reaction path diagram for the fluxes of the element element according the the net reaction rates determined by the Kinetics object kin.

add(self, ReactionPathDiagram other)

Add fluxes from other to this diagram

arrow_width

Get/Set the arrow width. If < 0, then scale with flux value.

bold_color

Get/Set the color for bold lines

bold_threshold

Get/Set the minimum relative flux for bold lines

build(self, verbose=False)

Build the reaction path diagram. Called automatically by methods which return representations of the diagram, for example write_dot().

dashed_color

Get/Set the color for dashed lines

display_only(self, int k)

Include only species and fluxes that are directly connected to the species with index k. Set to -1 to include all species.

dot_options

Get/Set options for the ‘dot’ program

flow_type

Get/Set the way flows are drawn. Either ‘NetFlow’ or ‘OneWayFlow’

font

Get/Set the name of the font used

get_data(self)

Get a (roughly) human-readable representation of the reaction path diagram.

get_dot(self)

Return a string containing the reaction path diagram formatted for use by Graphviz’s ‘dot’ program.

label_threshold

Get/Set the minimum relative flux for labels

log

Logging messages generated while building the reaction path diagram

normal_color

Get/Set the color for normal-weight lines

normal_threshold

Get/Set the maximum relative flux for dashed lines

scale

Get/Set the scaling factor for the fluxes. Set to -1 to normalize by the maximum net flux.

show_details

Get/Set whether to show the details of which reactions contribute to the flux.

threshold

Get/Set the threshold for the minimum flux relative value that will be plotted.

title

Get/Set the diagram title

write_dot(self, filename)

Write the reaction path diagram formatted for use by Graphviz’s ‘dot’ program to the file named filename.