Warning

This documentation is for an old version of Cantera. You can find docs for newer versions here.

Zero-Dimensional Reactor Networks

Defining Functions

class cantera.Func1

Bases: object

This class is used as a wrapper for a function of one variable, i.e. \(y = f(t)\), that is defined in Python and can be called by the Cantera C++ core. Func1 objects are constructed from callable Python objects, e.g. functions or classes which implement the __call__ method:

>>> f1 = Func1(math.sin)
>>> f1(math.pi/4)
0.7071067811865475

>>> f2 = Func1(lambda t: t**2 + 1)
>>> f2(3)
10

>>> class Multiplier(object):
...     def __init__(self, factor):
...         self.factor = factor
...     def __call__(self, t):
...         return self.factor * t
>>> f3 = Func1(Multiplier(5))
>>> f3(6)
30.0

For simplicity, constant functions can be defined by passing the constant value directly:

>>> f4 = Func1(2.5)
>>> f4(0.1)
2.5

Note that all methods which accept Func1 objects will also accept the callable object and create the wrapper on their own, so it is generally unnecessary to explicitly create a Func1 object.

Base Classes

ReactorBase

class cantera.ReactorBase(contents=None, name=None)

Bases: object

ReactorBase(ThermoPhase contents=None, name=None, volume=None, *)

Common base class for reactors and reservoirs.

T

The temperature [K] of the reactor’s contents.

Y

The mass fractions of the reactor’s contents.

density

The density [kg/m^3 or kmol/m^3] of the reactor’s contents.

inlets

List of flow devices installed as inlets to this reactor

insert(self, _SolutionBase solution)

Set solution to be the object used to compute thermodynamic properties and kinetic rates for this reactor.

mass

The mass of the reactor’s contents.

name

The name of the reactor.

outlets

List of flow devices installed as outlets to this reactor

reactor_type = 'None'
syncState(self)

Set the state of the Reactor to match that of the associated ThermoPhase object. After calling syncState(), call ReactorNet.reinitialize() before further integration.

thermo

The ThermoPhase object representing the reactor’s contents.

volume

The volume [m^3] of the reactor.

walls

List of walls installed on this reactor

FlowDevice

class cantera.FlowDevice(upstream, downstream, *, name=None)

Bases: object

FlowDevice(upstream, downstream, name=None, *)

Base class for devices that allow flow between reactors.

FlowDevice objects are assumed to be adiabatic, non-reactive, and have negligible internal volume, so that they are internally always in steady-state even if the upstream and downstream reactors are not. The fluid enthalpy, chemical composition, and mass flow rate are constant across a FlowDevice, and the pressure difference equals the difference in pressure between the upstream and downstream reactors.

mdot(self, double t)

The mass flow rate [kg/s] through this device at time t [s].

Reactor Networks

class cantera.ReactorNet(reactors=())

Bases: object

ReactorNet(reactors=())

Networks of reactors. ReactorNet objects are used to simultaneously advance the state of one or more coupled reactors.

Example:

>>> r1 = Reactor(gas1)
>>> r2 = Reactor(gas2)
>>> <... install walls, inlets, outlets, etc...>
>>> reactor_network = ReactorNet([r1, r2])
>>> reactor_network.advance(time)
add_reactor(self, Reactor r)

Add a reactor to the network.

advance(self, double t)

Advance the state of the reactor network in time from the current time to time t [s], taking as many integrator timesteps as necessary.

advance_to_steady_state(self, int max_steps=10000, double residual_threshold=0., double atol=0., bool return_residuals=False)

Advance the reactor network in time until steady state is reached.

The steady state is defined by requiring that the state of the system only changes below a certain threshold. The residual is computed using feature scaling:

\[r = \left| \frac{x(t + \Delta t) - x(t)}{\text{max}(x) + \text{atol}} \right| \cdot \frac{1}{\sqrt{n_x}} \]
Parameters:
  • max_steps – Maximum number of steps to be taken
  • residual_threshold – Threshold below which the feature-scaled residual r should drop such that the network is defines as steady state. By default, residual_threshold is 10 times the solver rtol.
  • atol – The smallest expected value of interest. Used for feature scaling. By default, this atol is identical to the solver atol.
  • return_residuals – If set to True, this function returns the residual time series as a vector with length max_steps.
atol

The absolute error tolerance used while integrating the reactor equations.

atol_sensitivity

The absolute error tolerance for sensitivity analysis.

component_name(self, int i)

Return the name of the i-th component of the global state vector. The name returned includes both the name of the reactor and the specific component, e.g. 'reactor1: CH4'.

get_state(self)

Get the combined state vector of the reactor network.

The combined state vector consists of the concatenated state vectors of all entities contained.

max_err_test_fails

The maximum number of error test failures permitted by the CVODES integrator in a single time step.

n_sensitivity_params

The number of registered sensitivity parameters.

n_vars

The number of state variables in the system. This is the sum of the number of variables for each Reactor and Wall in the system. Equal to:

Reactor and IdealGasReactor: n_species + 3 (mass, volume, internal energy or temperature).

ConstPressureReactor and IdealGasConstPressureReactor: n_species + 2 (mass, enthalpy or temperature).

Wall: number of surface species

reinitialize(self)

Reinitialize the integrator after making changing to the state of the system. Changes to Reactor contents will automatically trigger reinitialization.

rtol

The relative error tolerance used while integrating the reactor equations.

rtol_sensitivity

The relative error tolerance for sensitivity analysis.

sensitivities(self)

Returns the sensitivities of all of the solution variables with respect to all of the registered parameters. The normalized sensitivity coefficient \(S_{ki}\) of the solution variable \(y_k\) with respect to sensitivity parameter \(p_i\) is defined as:

\[S_{ki} = \frac{p_i}{y_k} \frac{\partial y_k}{\partial p_i} \]

For reaction sensitivities, the parameter is a multiplier on the forward rate constant (and implicitly on the reverse rate constant for reversible reactions).

The sensitivities are returned in an array with dimensions (n_vars, n_sensitivity_params), unless no timesteps have been taken, in which case the shape is (0, n_sensitivity_params). The order of the variables (i.e., rows) is:

Reactor or IdealGasReactor:

  • 0 - mass
  • 1 - volume
  • 2 - internal energy or temperature
  • 3+ - mass fractions of the species

ConstPressureReactor or IdealGasConstPressureReactor:

  • 0 - mass
  • 1 - enthalpy or temperature
  • 2+ - mass fractions of the species
sensitivity(self, component, int p, int r=0)

Returns the sensitivity of the solution variable component in reactor r with respect to the parameter p. component can be a string or an integer. See component_index and sensitivities to determine the integer index for the variables and the definition of the resulting sensitivity coefficient. If it is not given, r defaults to the first reactor. Returns an empty array until the first time step is taken.

sensitivity_parameter_name(self, int p)

Name of the sensitivity parameter with index p.

set_initial_time(self, double t)

Set the initial time. Restarts integration from this time using the current state as the initial condition. Default: 0.0 s.

set_max_time_step(self, double t)

Set the maximum time step t [s] that the integrator is allowed to use.

step(self)

Take a single internal time step. The time after taking the step is returned.

time

The current time [s].

verbose

If True, verbose debug information will be printed during integration. The default is False.

Reactors

Reservoir

class cantera.Reservoir(contents=None, name=None)

Bases: cantera._cantera.ReactorBase

A reservoir is a reactor with a constant state. The temperature, pressure, and chemical composition in a reservoir never change from their initial values.

reactor_type = 'Reservoir'

Reactor

class cantera.Reactor(contents=None, *, name=None, energy='on')

Bases: cantera._cantera.ReactorBase

Reactor(contents=None, name=None, *, energy=’on’, **kwargs)

A homogeneous zero-dimensional reactor. By default, they are closed (no inlets or outlets), have fixed volume, and have adiabatic, chemically-inert walls. These properties may all be changed by adding appropriate components, e.g. Wall, MassFlowController and Valve.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents.
  • name – Used only to identify this reactor in output. If not specified, defaults to 'Reactor_n', where n is an integer assigned in the order Reactor objects are created.
  • energy – Set to 'on' or 'off'. If set to 'off', the energy equation is not solved, and the temperature is held at its initial value..

Some examples showing how to create Reactor objects are shown below.

>>> gas = Solution('gri30.xml')
>>> r1 = Reactor(gas)

This is equivalent to:

>>> r1 = Reactor()
>>> r1.insert(gas)

Arguments may be specified using keywords in any order:

>>> r2 = Reactor(contents=gas, energy='off',
...              name='isothermal_reactor')
>>> r3 = Reactor(name='adiabatic_reactor', contents=gas)
add_sensitivity_reaction(self, m)

Specifies that the sensitivity of the state variables with respect to reaction m should be computed. m is the 0-based reaction index. The reactor must be part of a network first. Specifying the same reaction more than one time raises an exception.

add_sensitivity_species_enthalpy(self, k)

Specifies that the sensitivity of the state variables with respect to species k should be computed. The reactor must be part of a network first.

chemistry_enabled

True when the reactor composition is allowed to change due to chemical reactions in this reactor. When this is False, the reactor composition is held constant.

component_index(self, name)

Returns the index of the component named name in the system. This determines the (relative) index of the component in the vector of sensitivity coefficients. name is either a species name or the name of a reactor state variable, e.g. ‘int_energy’, ‘temperature’, depending on the reactor’s equations.

component_name(self, int i)

Returns the name of the component with index i within the array of variables returned by get_state. This is the inverse of component_index.

energy_enabled

True when the energy equation is being solved for this reactor. When this is False, the reactor temperature is held constant.

get_state(self)

Get the state vector of the reactor.

The order of the variables (i.e. rows) is:

Reactor or IdealGasReactor:

  • 0 - mass
  • 1 - volume
  • 2 - internal energy or temperature
  • 3+ - mass fractions of the species

ConstPressureReactor or IdealGasConstPressureReactor:

  • 0 - mass
  • 1 - enthalpy or temperature
  • 2+ - mass fractions of the species

You can use the function component_index to determine the location of a specific component from its name, or component_name to determine the name from the index.

insert(self, _SolutionBase solution)
kinetics

The Kinetics object used for calculating kinetic rates in this reactor.

n_vars

The number of state variables in the reactor. Equal to:

Reactor and IdealGasReactor: n_species + 3 (mass, volume, internal energy or temperature).

ConstPressureReactor and IdealGasConstPressureReactor: n_species + 2 (mass, enthalpy or temperature).

reactor_type = 'Reactor'

IdealGasReactor

class cantera.IdealGasReactor(contents=None, *, name=None, energy='on')

Bases: cantera._cantera.Reactor

A constant volume, zero-dimensional reactor for ideal gas mixtures.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents.
  • name – Used only to identify this reactor in output. If not specified, defaults to 'Reactor_n', where n is an integer assigned in the order Reactor objects are created.
  • energy – Set to 'on' or 'off'. If set to 'off', the energy equation is not solved, and the temperature is held at its initial value..

Some examples showing how to create Reactor objects are shown below.

>>> gas = Solution('gri30.xml')
>>> r1 = Reactor(gas)

This is equivalent to:

>>> r1 = Reactor()
>>> r1.insert(gas)

Arguments may be specified using keywords in any order:

>>> r2 = Reactor(contents=gas, energy='off',
...              name='isothermal_reactor')
>>> r3 = Reactor(name='adiabatic_reactor', contents=gas)
reactor_type = 'IdealGasReactor'

ConstPressureReactor

class cantera.ConstPressureReactor(contents=None, *, name=None, energy='on')

Bases: cantera._cantera.Reactor

A homogeneous, constant pressure, zero-dimensional reactor. The volume of the reactor changes as a function of time in order to keep the pressure constant.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents.
  • name – Used only to identify this reactor in output. If not specified, defaults to 'Reactor_n', where n is an integer assigned in the order Reactor objects are created.
  • energy – Set to 'on' or 'off'. If set to 'off', the energy equation is not solved, and the temperature is held at its initial value..

Some examples showing how to create Reactor objects are shown below.

>>> gas = Solution('gri30.xml')
>>> r1 = Reactor(gas)

This is equivalent to:

>>> r1 = Reactor()
>>> r1.insert(gas)

Arguments may be specified using keywords in any order:

>>> r2 = Reactor(contents=gas, energy='off',
...              name='isothermal_reactor')
>>> r3 = Reactor(name='adiabatic_reactor', contents=gas)
reactor_type = 'ConstPressureReactor'

IdealGasConstPressureReactor

class cantera.IdealGasConstPressureReactor(contents=None, *, name=None, energy='on')

Bases: cantera._cantera.Reactor

A homogeneous, constant pressure, zero-dimensional reactor for ideal gas mixtures. The volume of the reactor changes as a function of time in order to keep the pressure constant.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents.
  • name – Used only to identify this reactor in output. If not specified, defaults to 'Reactor_n', where n is an integer assigned in the order Reactor objects are created.
  • energy – Set to 'on' or 'off'. If set to 'off', the energy equation is not solved, and the temperature is held at its initial value..

Some examples showing how to create Reactor objects are shown below.

>>> gas = Solution('gri30.xml')
>>> r1 = Reactor(gas)

This is equivalent to:

>>> r1 = Reactor()
>>> r1.insert(gas)

Arguments may be specified using keywords in any order:

>>> r2 = Reactor(contents=gas, energy='off',
...              name='isothermal_reactor')
>>> r3 = Reactor(name='adiabatic_reactor', contents=gas)
reactor_type = 'IdealGasConstPressureReactor'

FlowReactor

class cantera.FlowReactor(contents=None, *, name=None, energy='on')

Bases: cantera._cantera.Reactor

A steady-state plug flow reactor with constant cross sectional area. Time integration follows a fluid element along the length of the reactor. The reactor is assumed to be frictionless and adiabatic.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents.
  • name – Used only to identify this reactor in output. If not specified, defaults to 'Reactor_n', where n is an integer assigned in the order Reactor objects are created.
  • energy – Set to 'on' or 'off'. If set to 'off', the energy equation is not solved, and the temperature is held at its initial value..

Some examples showing how to create Reactor objects are shown below.

>>> gas = Solution('gri30.xml')
>>> r1 = Reactor(gas)

This is equivalent to:

>>> r1 = Reactor()
>>> r1.insert(gas)

Arguments may be specified using keywords in any order:

>>> r2 = Reactor(contents=gas, energy='off',
...              name='isothermal_reactor')
>>> r3 = Reactor(name='adiabatic_reactor', contents=gas)
distance

The distance of the fluid element from the inlet of the reactor.

mass_flow_rate

Mass flow rate per unit area [kg/m^2*s]

reactor_type = 'FlowReactor'
speed

Speed [m/s] of the flow in the reactor at the current position

Walls

Wall

class cantera.Wall(left, right, *, name=None, A=None, K=None, U=None, Q=None, velocity=None, kinetics=(None, None))

Bases: object

Wall(left, right, name=None, *, A=None, K=None, U=None, Q=None, velocity=None)

A Wall separates two reactors, or a reactor and a reservoir. A wall has a finite area, may conduct or radiate heat between the two reactors on either side, and may move like a piston.

Walls are stateless objects in Cantera, meaning that no differential equation is integrated to determine any wall property. Since it is the wall (piston) velocity that enters the energy equation, this means that it is the velocity, not the acceleration or displacement, that is specified. The wall velocity is computed from

\[v = K(P_{\rm left} - P_{\rm right}) + v_0(t), \]

where \(K\) is a non-negative constant, and \(v_0(t)\) is a specified function of time. The velocity is positive if the wall is moving to the right.

The heat flux through the wall is computed from

\[q = U(T_{\rm left} - T_{\rm right}) + \epsilon\sigma (T_{\rm left}^4 - T_{\rm right}^4) + q_0(t), \]

where \(U\) is the overall heat transfer coefficient for conduction/convection, and \(\epsilon\) is the emissivity. The function \(q_0(t)\) is a specified function of time. The heat flux is positive when heat flows from the reactor on the left to the reactor on the right.

Parameters:
  • left – Reactor or reservoir on the left. Required.
  • right – Reactor or reservoir on the right. Required.
  • name – Name string. If omitted, the name is 'Wall_n', where 'n' is an integer assigned in the order walls are created.
  • A – Wall area [m^2]. Defaults to 1.0 m^2.
  • K – Wall expansion rate parameter [m/s/Pa]. Defaults to 0.0.
  • U – Overall heat transfer coefficient [W/m^2]. Defaults to 0.0 (adiabatic wall).
  • Q – Heat flux function \(q_0(t)\) [W/m^2]. Optional. Default: \(q_0(t) = 0.0\).
  • velocity – Wall velocity function \(v_0(t)\) [m/s]. Default: \(v_0(t) = 0.0\).
area

The wall area [m^2].

emissivity

The emissivity (nondimensional)

expansion_rate_coeff

The coefficient K [m/s/Pa] that determines the velocity of the wall as a function of the pressure difference between the adjacent reactors.

heat_transfer_coeff

the overall heat transfer coefficient [W/m^2/K]

left

The left surface of this wall.

qdot(self, double t)

Total heat flux [W] through the wall at time t. A positive value corresponds to heat flowing from the left-hand reactor to the right-hand one.

right

The right surface of this wall.

set_heat_flux(self, q)

Heat flux [W/m^2] across the wall. May be either a constant or an arbitrary function of time. See Func1.

set_velocity(self, v)

The wall velocity [m/s]. May be either a constant or an arbitrary function of time. See Func1.

vdot(self, double t)

The rate of volumetric change [m^3/s] associated with the wall at time t. A positive value corresponds to the left-hand reactor volume increasing, and the right-hand reactor volume decreasing.

WallSurface

class cantera.WallSurface(wall, side)

Bases: object

Surfaces

ReactorSurface

class cantera.ReactorSurface(kin=None, r=None, *, A=None)

Bases: object

ReactorSurface(kin=None, Reactor r=None, A=None, *)

Represents a surface in contact with the contents of a reactor.

Parameters:
  • kin – The Kinetics or Interface object representing reactions on this surface.
  • r – The Reactor into which this surface should be installed.
  • A – The area of the reacting surface [m^2]
add_sensitivity_reaction(self, int m)

Specifies that the sensitivity of the state variables with respect to reaction m should be computed. m is the 0-based reaction index. The Surface must be installed on a reactor and part of a network first.

area

Area on which reactions can occur [m^2]

coverages

The fraction of sites covered by each surface species.

install(self, Reactor r)
kinetics

The InterfaceKinetics object used for calculating reaction rates on this surface.

Flow Controllers

MassFlowController

class cantera.MassFlowController(upstream, downstream, *, name=None, mdot=None)

Bases: cantera._cantera.FlowDevice

MassFlowController(upstream, downstream, name=None, *, mdot=None)

A mass flow controller maintains a specified mass flow rate independent of upstream and downstream conditions. The equation used to compute the mass flow rate is

\[\dot m = \max(\dot m_0, 0.0),\]

where \(\dot m_0\) is either a constant value or a function of time. Note that if \(\dot m_0 < 0\), the mass flow rate will be set to zero, since reversal of the flow direction is not allowed.

Unlike a real mass flow controller, a MassFlowController object will maintain the flow even if the downstream pressure is greater than the upstream pressure. This allows simple implementation of loops, in which exhaust gas from a reactor is fed back into it through an inlet. But note that this capability should be used with caution, since no account is taken of the work required to do this.

set_mass_flow_rate(self, m)

Set the mass flow rate [kg/s] through this controller to be either a constant or an arbitrary function of time. See Func1.

>>> mfc.set_mass_flow_rate(0.3)
>>> mfc.set_mass_flow_rate(lambda t: 2.5 * exp(-10 * (t - 0.5)**2))

Valve

class cantera.Valve(upstream, downstream, *, name=None, K=None)

Bases: cantera._cantera.FlowDevice

Valve(upstream, downstream, name=None, *, K=None)

In Cantera, a Valve is a flow device with mass flow rate that is a function of the pressure drop across it. The default behavior is linear:

\[\dot m = K_v*(P_1 - P_2) \]

where \(K_v\) is a constant set by the set_valve_coeff method. Note that \(P_1\) must be greater than \(P_2\); otherwise, \(\dot m = 0\). However, an arbitrary function can also be specified, such that

\[\dot m = f(P_1 - P_2) \]

where \(f\) is the arbitrary function that returns the mass flow rate given a single argument, the pressure differential. See the documentation for the set_valve_coeff method for an example. Note that it is never possible for the flow to reverse and go from the downstream to the upstream reactor/reservoir through a line containing a Valve object.

Valve objects are often used between an upstream reactor and a downstream reactor or reservoir to maintain them both at nearly the same pressure. By setting the constant \(K_v\) to a sufficiently large value, very small pressure differences will result in flow between the reactors that counteracts the pressure difference.

set_valve_coeff(self, k)

Set the relationship between mass flow rate and the pressure drop across the valve. If a number is given, it is the proportionality constant [kg/s/Pa]. If a function is given, it should compute the mass flow rate [kg/s] given the pressure drop [Pa].

>>> V = Valve(res1, reactor1)
>>> V.set_valve_coeff(1e-4)  # Set the value of K to a constant
>>> V.set_valve_coeff(lambda dP: (1e-5 * dP)**2)  # Set the value of K to a function

PressureController

class cantera.PressureController(upstream, downstream, *, name=None, master=None, K=None)

Bases: cantera._cantera.FlowDevice

PressureController(upstream, downstream, name=None, *, master=None, K=None)

A PressureController is designed to be used in conjunction with another ‘master’ flow controller, typically a MassFlowController. The master flow controller is installed on the inlet of the reactor, and the corresponding PressureController is installed on on outlet of the reactor. The PressureController mass flow rate is equal to the master mass flow rate, plus a small correction dependent on the pressure difference:

\[\dot m = \dot m_{\rm master} + K_v(P_1 - P_2). \]
set_master(self, FlowDevice d)

Set the “master” FlowDevice used to compute this device’s mass flow rate.

set_pressure_coeff(self, double k)

Set the proportionality constant \(K_v\) [kg/s/Pa] between the pressure drop and the mass flow rate.