Zero-Dimensional Reactor Networks#

Defining Functions#

class cantera.Func1#

Bases: object

This class is used as a wrapper for a function of one variable, \(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, for example 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:
...     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 a 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 often not necessary to explicitly create a Func1 object.

classmethod cxx_functor(functor_type, *args)#

Retrieve a C++ Func1 functor (advanced feature).

For implemented functor types, see the Cantera C++ Func1 documentation.

New in version 3.0.

type#

Return the type of the underlying C++ functor object.

New in version 3.0.

write(name='t')#

Write a \(LaTeX\) expression representing a functor.

Parameters:

name – Name of the variable to be used.

New in version 3.0.

class cantera.Tabulated1#

Bases: Func1

A Tabulated1 object representing a tabulated function is defined by sample points and corresponding function values. Inputs are specified by two iterable objects containing sample point location and function values. Between sample points, values are evaluated based on the optional argument method; options are 'linear' (linear interpolation, default) or 'previous' (nearest previous value). Outside the sample interval, the value at the closest end point is returned.

Examples for Tabulated1 objects are:

>>> t1 = Tabulated1([0, 1, 2], [2, 1, 0])
>>> [t1(v) for v in [-0.5, 0, 0.5, 1.5, 2, 2.5]]
[2.0, 2.0, 1.5, 0.5, 0.0, 0.0]

>>> t2 = Tabulated1(np.array([0, 1, 2]), np.array([2, 1, 0]))
>>> [t2(v) for v in [-0.5, 0, 0.5, 1.5, 2, 2.5]]
[2.0, 2.0, 1.5, 0.5, 0.0, 0.0]

The optional method keyword argument changes the type of interpolation from the 'linear' default to 'previous':

>>> t3 = Tabulated1([0, 1, 2], [2, 1, 0], method='previous')
>>> [t3(v) for v in [-0.5, 0, 0.5, 1.5, 2, 2.5]]
[2.0, 2.0, 2.0, 1.0, 0.0, 0.0]

New in version 3.0.

Base Classes#

ReactorBase#

class cantera.ReactorBase#

Bases: object

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.

draw(graph=None, *, graph_attr=None, node_attr=None, print_state=False, species=None, species_units='percent')#

Draw as graphviz dot node. The node is added to an existing graph if provided. Optionally include current reactor state in the node.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the reactor is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the node method invoked to draw the reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • print_state – Whether state information of the reactor is printed into the node. Defaults to False.

  • species – If print_state is True, define how species are to be printed. Options are 'X' and 'Y' for mole and mass fractions of all species, respectively, or an iterable that contains the desired species names as strings. Defaults to None.

  • species_units – Defines the units the species are displayed in as either "percent" or "ppm". Defaults to "percent".

Returns:

graphviz.graphs.BaseGraph object with reactor

New in version 3.1.

inlets#

List of flow devices installed as inlets to this reactor

insert(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.

node_attr#

A dictionary containing draw attributes for the representation of the reactor as a graphviz node. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

New in version 3.1.

outlets#

List of flow devices installed as outlets to this reactor

reactor_type = 'none'#
surfaces#

List of reacting surfaces installed on this reactor

syncState()#

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.

type#

The type of the reactor.

volume#

The volume [m^3] of the reactor.

walls#

List of walls installed on this reactor

FlowDevice#

class cantera.FlowDevice#

Bases: object

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.

downstream#

Return the Reactor or Reservoir object downstream of the flow device.

New in version 3.1.

draw(graph=None, *, graph_attr=None, node_attr=None, edge_attr=None)#

Draw as connection between upstream and downstream reactor or reservoir using graphviz.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw this flow controller connection. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

Returns:

A graphviz.graphs.BaseGraph object depicting the connection.

New in version 3.1.

edge_attr#

A dictionary containing draw attributes for the representation of the FlowDevice as a graphviz edge.See https://graphviz.org/docs/edges/ for a list of all usable attributes.

New in version 3.1.

flowdevice_type = 'none'#
mass_flow_rate#

Get the mass flow rate [kg/s] through this device at the current reactor network time.

pressure_function#

The relationship between mass flow rate and the pressure drop across a flow device. The mass flow rate [kg/s] is calculated given the pressure drop [Pa] and a coefficient set by a flow device specific function. Unless a user-defined pressure function is provided, the function returns the pressure difference across the device. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.pressure_function = lambda dP: dP**2

where FlowDevice is either a Valve or PressureController object.

New in version 3.0.

time_function#

The time dependence of a flow device. The mass flow rate [kg/s] is calculated for a Flow device, and multiplied by a function of time, which returns 1.0 unless a user-defined function is provided. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.time_function = lambda t: exp(-10 * (t - 0.5)**2)

where FlowDevice is either a Valve or MassFlowController object.

New in version 3.0.

type#

The type of the flow device.

upstream#

Return the Reactor or Reservoir object upstream of the flow device.

New in version 3.1.

Reactor Networks#

class cantera.ReactorNet#

Bases: object

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(r)#

Add a reactor to the network.

advance(t, apply_limit=True)#

Advance the state of the reactor network from the current time/distance towards the specified value t of the independent variable, which depends on the type of reactors included in the network.

The integrator will take as many steps as necessary to reach t. If apply_limit is true and an advance limit is specified, the reactor state at the end of the step is estimated prior to advancing. If the difference exceed limits, the end value is reduced by half until the projected end state remains within specified limits. Returns the time/distance reached at the end of integration.

advance_limits#

Get or set absolute limits for state changes during ReactorNet.advance (positive values are considered; negative values disable a previously set advance limit for a solution component). Note that limits are disabled by default (with individual values set to -1.).

advance_to_steady_state(max_steps=10000, residual_threshold=0.0, atol=0.0, 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(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, for example 'reactor1: CH4'.

derivative_settings#

Apply derivative settings to all reactors in the network. See also Kinetics.derivative_settings.

distance#

The current distance[ m] along the length of the reactor network, for reactors that are solved as a function of space.

draw(*, graph_attr=None, node_attr=None, edge_attr=None, heat_flow_attr=None, mass_flow_attr=None, moving_wall_edge_attr=None, surface_edge_attr=None, print_state=False, show_wall_velocity=True, **kwargs)#

Draw as graphviz.graphs.DiGraph. Connecting flow controllers and walls are depicted as arrows.

Parameters:
  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). node_attr defined in the reactor object itself have priority. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any edge (flow controllers, walls). edge_attr defined in the connection objects (subclasses of FlowDevice or walls) themselves have priority. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

  • heat_flow_attr – Same as edge_attr but only applied to edges representing walls. Default is {"color": "red", "style": "dashed"}.

  • mass_flow_attr – Same as edge_attr but only applied to edges representing FlowDevice objects.

  • moving_wall_edge_attr – Same as edge_attr but only applied to edges representing wall movement.

  • surface_edge_attr – Same as edge_attr but only applied to edges representing connections between ReactorSurface objects and reactors. Default is {"style": "dotted", "arrowhead": "none"}.

  • print_state – Whether state information of the reactors is printed into each node. Defaults to False.

  • kwargs – Additional keywords are passed on to each call of draw_reactor, draw_surface, draw_flow_controllers, and draw_walls.

Returns:

graphviz.graphs.BaseGraph object with reactor net.

New in version 3.1.

get_derivative(k)#

Get the k-th derivative of the state vector of the reactor network with respect to the independent integrator variable (time/distance).

get_state()#

Get the combined state vector of the reactor network.

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

global_component_index(name, reactor)#

Returns the index of a component named name of a reactor with index reactor within the global state vector. That is, this determines the absolute index of the component, where reactor is the index of the reactor that holds the component. name is either a species name or the name of a reactor state variable, for example, 'int_energy', 'temperature', etc. depending on the reactor’s equations.

include_algebraic_in_error_test#

Get/Set whether to include algebraic variables in the in the local error test. Applicable only to DAE systems. The default is True.

initial_time#

The initial time of the integrator. When set, integration is restarted from this time using the current state as the initial condition. Default: 0.0 s.

New in version 3.0.

initialize()#

Force initialization of the integrator after initial setup.

linear_solver_type#

The type of linear solver used in integration.

Options for this property include:

  • "DENSE"

  • "GMRES"

  • "BAND"

  • "DIAG"

max_err_test_fails#

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

max_nonlinear_convergence_failures#

Get/Set the maximum number of nonlinear solver convergence failures permitted in one step of the SUNDIALS integrator. The default value is 10.

max_nonlinear_iterations#

Get/Set the maximum number of nonlinear solver iterations permitted by the SUNDIALS solver in one solve attempt. The default value is 4.

max_order#

Get/Set the maximum order of the linear multistep method. The default value and maximum is 5.

max_steps#

The maximum number of internal integration steps that CVODES is allowed to take before reaching the next output point.

max_time_step#

Get/set the maximum time step t [s] that the integrator is allowed to use. The default value is set to zero, so that no time step maximum is used.

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

preconditioner#

Preconditioner associated with integrator

reactors#

List of all reactors that are part of the reactor network.

New in version 3.1.

reinitialize()#

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()#

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 integration steps have been taken, in which case the shape is (0, n_sensitivity_params). The order of the variables (that is, 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(component, p, 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 integration step is taken.

sensitivity_parameter_name(p)#

Name of the sensitivity parameter with index p.

solver_stats#

ODE solver stats from integrator

step()#

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

time#

The current time [s], for reactor networks that are solved in the time domain.

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: 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', node_attr=None, group_name='')#

Bases: ReactorBase

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 such as 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. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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(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(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(name)#

Returns the index of the component named name in the system. This determines the index of the component in the vector of sensitivity coefficients. name is either a species name or the name of a reactor state variable, for example 'int_energy' or 'temperature', depending on the reactor’s equations.

component_name(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.

finite_difference_jacobian#

Get the reactor-specific Jacobian, calculated using a finite difference method.

Warning

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

get_state()#

Get the state vector of the reactor.

The order of the variables (that is, 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.

group_name#

Optional name of a grouping of reactors that will be drawn as a cluster in the graphical representation using graphviz. See https://graphviz.org/Gallery/directed/cluster.html.

New in version 3.1.

insert(solution)#

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

jacobian#

Get the local, reactor-specific Jacobian or an approximation thereof

Warning

Depending on the particular implementation, this may return an approximate Jacobian intended only for use in forming a preconditioner for iterative solvers, excluding terms that would generate a fully-dense Jacobian.

Warning

This method is an experimental part of the Cantera API and may be changed or removed without notice.

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'#
set_advance_limit(name, limit)#

Limit absolute change of component name during ReactorNet.advance. (positive limit values are considered; negative values disable a previously set advance limit for a solution component). Note that limits are disabled by default (with individual values set to -1.).

MoleReactor#

class cantera.MoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: Reactor

A homogeneous zero-dimensional reactor with a mole based state vector. 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 such as Wall, MassFlowController and Valve.

New in version 3.0.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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 = 'MoleReactor'#

IdealGasReactor#

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

Bases: 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. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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'#

IdealGasMoleReactor#

class cantera.IdealGasMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: Reactor

A constant volume, zero-dimensional reactor for ideal gas mixtures with a mole based state vector

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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 = 'IdealGasMoleReactor'#

ConstPressureReactor#

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

Bases: 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. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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'#

ConstPressureMoleReactor#

class cantera.ConstPressureMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: Reactor

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

New in version 3.0.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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 = 'ConstPressureMoleReactor'#

IdealGasConstPressureReactor#

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

Bases: 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. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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'#

IdealGasConstPressureMoleReactor#

class cantera.IdealGasConstPressureMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: 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. This reactor also uses a mole based state vector.

Parameters:
  • contents – Reactor contents. If not specified, the reactor is initially empty. In this case, call insert to specify the contents. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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 = 'IdealGasConstPressureMoleReactor'#

FlowReactor#

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

Bases: Reactor

A steady-state plug flow reactor with constant cross sectional area. 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. Providing valid contents will become mandatory after Cantera 3.1.

  • 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.

  • node_attr – Attributes to be passed to the node method invoked to draw this reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • group_name – Group reactors of the same group_name when drawn using graphviz.

New in version 3.1: Added the node_attr and group_name parameters.

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

>>> gas = Solution('gri30.yaml')
>>> r1 = Reactor(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)
area#

Get/set the area of the reactor [m^2].

When the area is changed, the flow speed is scaled to keep the total mass flow rate constant.

inlet_surface_atol#

Get/Set the steady-state tolerances used to determine the initial surface species coverages.

inlet_surface_max_error_failures#

Get/Set the maximum number of integrator error failures allowed when determining the initial surface species coverages.

inlet_surface_max_steps#

Get/Set the maximum number of integrator steps used to determine the initial surface species coverages.

inlet_surface_rtol#

Get/Set the steady-state tolerances used to determine the initial surface species coverages.

mass_flow_rate#

Mass flow rate [kg/s]

reactor_type = 'FlowReactor'#
speed#

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

surface_area_to_volume_ratio#

Get/Set the surface area to volume ratio of the reactor [m^-1]

ExtensibleReactor#

class cantera.ExtensibleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: Reactor

A base class for a reactor with delegated methods where the base functionality corresponds to the Reactor class.

The following methods of the C++ Reactor class can be modified by a Python class which inherits from this class. For each method, the name below should be prefixed with before_, after_, or replace_, indicating whether this method should be called before, after, or instead of the corresponding method from the base class.

For methods that return a value and have a before method specified, if that method returns a value other than None that value will be returned without calling the base class method; otherwise, the value from the base class method will be returned. For methods that return a value and have an after method specified, the returned value wil be the sum of the values from the supplied method and the base class method.

initialize(self, t0: double) -> None

Responsible for allocating and setting the sizes of any internal variables, initializing attached walls, and setting the total number of state variables associated with this reactor, n_vars.

Called once before the start of time integration.

sync_state(self) -> None

Responsible for setting the state of the reactor to correspond to the state of the associated ThermoPhase object.

get_state(self, y : double[:]) -> None

Responsible for populating the state vector y (length n_vars) with the initial state of the reactor.

update_state(self, y : double[:]) -> None

Responsible for setting the state of the reactor object from the values in the state vector y (length n_vars)

update_surface_state(self, y : double[:]) -> None

Responsible for setting the state of surface phases in this reactor from the values in the state vector y. The length of y is the total number of surface species in all surfaces.

get_surface_initial_conditions(self, y : double[:]) -> None

Responsible for populating the state vector y with the initial state of each surface phase in this reactor. The length of y is the total number of surface species in all surfaces.

update_connected(self, update_pressure : bool) -> None

Responsible for storing properties which may be accessed by connected reactors, and for updating the mass flow rates of connected flow devices.

eval(self, t : double, LHS : double[:], RHS : double[:]) -> None

Responsible for calculating the time derivative of the state at time t based on the current state of the reactor. For each component i of the state vector, the time derivative dy[i]/dt is calculated as LHS[i] * dy[i]/dt = RHS[i]. LHS and RHS are arrays of length n_vars.

eval_walls(self, t : double) -> None

Responsible for calculating the net rate of volume change expansion_rate and the net rate of heat transfer heat_rate caused by walls connected to this reactor.

eval_surfaces(LHS : double[:], RHS : double[:], sdot : double[:]) -> None

Responsible for calculating the LHS and RHS (length: total number of surface species in all surfaces) of the ODEs for surface species coverages, and the molar production rate of bulk phase species sdot (length: number of bulk phase species).

component_name(i : int) -> string

Returns the name of the state vector component with index i

component_index(name: string) -> int

Returns the index of the state vector component named name

species_index(name : string) -> int

Returns the index of the species named name, in either the bulk phase or a surface phase, relative to the start of the species terms in the state vector.

delegatable_methods = {'component_index': ('componentIndex', 'size_t(string)'), 'component_name': ('componentName', 'string(size_t)'), 'eval': ('eval', 'void(double, double*, double*)'), 'eval_surfaces': ('evalSurfaces', 'void(double*,double*,double*)'), 'eval_walls': ('evalWalls', 'void(double)'), 'get_state': ('getState', 'void(double*)'), 'get_surface_initial_conditions': ('getSurfaceInitialConditions', 'void(double*)'), 'initialize': ('initialize', 'void(double)'), 'species_index': ('speciesIndex', 'size_t(string)'), 'sync_state': ('syncState', 'void()'), 'update_connected': ('updateConnected', 'void(bool)'), 'update_state': ('updateState', 'void(double*)'), 'update_surface_state': ('updateSurfaceState', 'void(double*)')}#
expansion_rate#

Get/Set the net rate of volume change (for example, from moving walls) [m^3/s]

New in version 3.0.

heat_rate#

Get/Set the net heat transfer rate (for example, through walls) [W]

New in version 3.0.

n_vars#

Get/Set the number of state variables in the reactor.

reactor_type = 'ExtensibleReactor'#
restore_surface_state(n)#

Set the state of the thermo object for surface n to correspond to the state of that surface

restore_thermo_state()#

Set the state of the thermo object to correspond to the state of the reactor.

ExtensibleIdealGasReactor#

class cantera.ExtensibleIdealGasReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the IdealGasReactor class.

reactor_type = 'ExtensibleIdealGasReactor'#

ExtensibleConstPressureReactor#

class cantera.ExtensibleConstPressureReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the ConstPressureReactor class.

reactor_type = 'ExtensibleConstPressureReactor'#

ExtensibleIdealGasConstPressureReactor#

class cantera.ExtensibleIdealGasConstPressureReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the IdealGasConstPressureReactor class.

reactor_type = 'ExtensibleIdealGasConstPressureReactor'#

ExtensibleMoleReactor#

class cantera.ExtensibleMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the MoleReactor class.

New in version 3.0.

reactor_type = 'ExtensibleMoleReactor'#

ExtensibleIdealGasMoleReactor#

class cantera.ExtensibleIdealGasMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the IdealGasMoleReactor class.

New in version 3.0.

reactor_type = 'ExtensibleIdealGasMoleReactor'#

ExtensibleConstPressureMoleReactor#

class cantera.ExtensibleConstPressureMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the ConstPressureMoleReactor class.

New in version 3.0.

reactor_type = 'ExtensibleConstPressureMoleReactor'#

ExtensibleIdealGasConstPressureMoleReactor#

class cantera.ExtensibleIdealGasConstPressureMoleReactor(contents=None, *, name=None, energy='on', node_attr=None, group_name='')#

Bases: ExtensibleReactor

A variant of ExtensibleReactor where the base behavior corresponds to the IdealGasConstPressureMoleReactor class.

New in version 3.0.

reactor_type = 'ExtensibleIdealGasConstPressureMoleReactor'#

Walls#

Wall#

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

Bases: WallBase

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\).

  • edge_attr – Attributes like style when drawn as a connecting edge using graphviz’s dot language. Default is {}.

New in version 3.1: Added the edge_attr parameter.

area#

The wall area [m^2].

draw(graph=None, *, graph_attr=None, node_attr=None, edge_attr=None, moving_wall_edge_attr=None, show_wall_velocity=True)#

Draw as connection between left and right reactor or reservoir using graphviz.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw this wall connection. Default is {"color": "red", "style": "dashed"}. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

  • moving_wall_edge_attr – Same as edge_attr but only applied to edges representing wall movement. Default is {"arrowtail": "icurveteecurve", "dir": "both", "style": "dotted", "arrowhead": "icurveteecurve"}.

  • show_wall_velocity – If True, wall movement will be indicated by additional arrows with the corresponding wall velocity as a label.

Returns:

A graphviz.graphs.BaseGraph object depicting the connection.

New in version 3.1.

edge_attr#

A dictionary containing draw attributes for the representation of the WallBase as a graphviz edge.See https://graphviz.org/docs/edges/ for a list of all usable attributes.

New in version 3.1.

emissivity#

The emissivity (nondimensional)

expansion_rate#

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

New in version 3.0.

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_flux#

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

New in version 3.0.

heat_rate#

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

New in version 3.0.

heat_transfer_coeff#

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

left_reactor#

Return the Reactor or Reservoir object left of the wall.

New in version 3.1.

right_reactor#

Return the Reactor or Reservoir object right of the wall.

New in version 3.1.

type#

The type of the wall.

velocity#

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

New in version 3.0.

wall_type = 'Wall'#

Surfaces#

ReactorSurface#

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

Bases: object

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]

  • node_attr – Attributes to be passed to the node method invoked to draw this surface. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

New in version 3.1: Added the node_attr parameter.

add_sensitivity_reaction(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.

draw(graph=None, *, graph_attr=None, node_attr=None, surface_edge_attr=None, print_state=False, **kwargs)#

Draw the surface as a graphviz dot node connected to its reactor. The node is added to an existing graph if provided. Optionally include current reactor state in the reactor node.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the reactor is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the node method invoked to draw the reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • surface_edge_attr – Attributes to be passed to the edge method invoked to draw the connection between the surface and its reactor. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

  • print_state – Whether state information of the reactor is printed into its node. Defaults to False

  • kwargs – Additional keywords are passed on to ~.drawnetwork.draw_reactor.

Returns:

graphviz.graphs.BaseGraph object with surface and connected reactor.

New in version 3.1.

install(r)#

Add this ReactorSurface to the specified Reactor

kinetics#

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

node_attr#

A dictionary containing draw attributes for the representation of the reactor surface as a graphviz node. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

New in version 3.1.

reactor#

Return the Reactor object the surface is connected to.

New in version 3.1.

Flow Controllers#

MassFlowController#

class cantera.MassFlowController#

Bases: FlowDevice

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*g(t), 0.),\]

where \(\dot m_0\) is a constant value and \(g(t)\) is a function of time. Both \(\dot m_0\) and \(g(t)\) can be set individually by properties mass_flow_coeff and time_function, respectively. The property mass_flow_rate combines the former into a single interface. Note that if \(\dot m_0*g(t) < 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.

downstream#

Return the Reactor or Reservoir object downstream of the flow device.

New in version 3.1.

draw(graph=None, *, graph_attr=None, node_attr=None, edge_attr=None)#

Draw as connection between upstream and downstream reactor or reservoir using graphviz.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw this flow controller connection. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

Returns:

A graphviz.graphs.BaseGraph object depicting the connection.

New in version 3.1.

edge_attr#

A dictionary containing draw attributes for the representation of the FlowDevice as a graphviz edge.See https://graphviz.org/docs/edges/ for a list of all usable attributes.

New in version 3.1.

flowdevice_type = 'MassFlowController'#
mass_flow_coeff#

Set the mass flow rate [kg/s] through the mass flow controller as a constant, which may be modified by a function of time, see time_function.

>>> mfc = MassFlowController(res1, reactor1)
>>> mfc.mass_flow_coeff = 1e-4  # Set the flow rate to a constant
>>> mfc.mass_flow_coeff  # Get the flow rate value
mass_flow_rate#

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

Note that depending on the argument type, this method either changes the property mass_flow_coeff or updates the time_function property.

>>> mfc.mass_flow_rate = 0.3
>>> mfc.mass_flow_rate = lambda t: 2.5 * exp(-10 * (t - 0.5)**2)
pressure_function#

The relationship between mass flow rate and the pressure drop across a flow device. The mass flow rate [kg/s] is calculated given the pressure drop [Pa] and a coefficient set by a flow device specific function. Unless a user-defined pressure function is provided, the function returns the pressure difference across the device. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.pressure_function = lambda dP: dP**2

where FlowDevice is either a Valve or PressureController object.

New in version 3.0.

time_function#

The time dependence of a flow device. The mass flow rate [kg/s] is calculated for a Flow device, and multiplied by a function of time, which returns 1.0 unless a user-defined function is provided. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.time_function = lambda t: exp(-10 * (t - 0.5)**2)

where FlowDevice is either a Valve or MassFlowController object.

New in version 3.0.

type#

The type of the flow device.

upstream#

Return the Reactor or Reservoir object upstream of the flow device.

New in version 3.1.

Valve#

class cantera.Valve#

Bases: FlowDevice

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 using the valve_coeff property. 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 = K_v*f(P_1 - P_2)\]

where \(f\) is the arbitrary function that multiplies \(K_v\) given a single argument, the pressure differential. Further, a valve opening function \(g\) may be specified using the time_function property, such that

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

See the documentation for the valve_coeff property as well as the pressure_function and time_function properties for examples. 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.

downstream#

Return the Reactor or Reservoir object downstream of the flow device.

New in version 3.1.

draw(graph=None, *, graph_attr=None, node_attr=None, edge_attr=None)#

Draw as connection between upstream and downstream reactor or reservoir using graphviz.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw this flow controller connection. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

Returns:

A graphviz.graphs.BaseGraph object depicting the connection.

New in version 3.1.

edge_attr#

A dictionary containing draw attributes for the representation of the FlowDevice as a graphviz edge.See https://graphviz.org/docs/edges/ for a list of all usable attributes.

New in version 3.1.

flowdevice_type = 'Valve'#
mass_flow_rate#

Get the mass flow rate [kg/s] through this device at the current reactor network time.

pressure_function#

The relationship between mass flow rate and the pressure drop across a flow device. The mass flow rate [kg/s] is calculated given the pressure drop [Pa] and a coefficient set by a flow device specific function. Unless a user-defined pressure function is provided, the function returns the pressure difference across the device. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.pressure_function = lambda dP: dP**2

where FlowDevice is either a Valve or PressureController object.

New in version 3.0.

time_function#

The time dependence of a flow device. The mass flow rate [kg/s] is calculated for a Flow device, and multiplied by a function of time, which returns 1.0 unless a user-defined function is provided. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.time_function = lambda t: exp(-10 * (t - 0.5)**2)

where FlowDevice is either a Valve or MassFlowController object.

New in version 3.0.

type#

The type of the flow device.

upstream#

Return the Reactor or Reservoir object upstream of the flow device.

New in version 3.1.

valve_coeff#

Set valve coefficient, that is, the proportionality constant between mass flow rate and pressure drop [kg/s/Pa].

>>> v = Valve(res1, reactor1)
>>> v.valve_coeff = 1e-4  # Set the value of K to a constant
>>> v.valve_coeff  # Get the value of K

PressureController#

class cantera.PressureController#

Bases: FlowDevice

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

\[\dot m = \dot m_{\rm primary} + K_v(P_1 - P_2).\]

As an alternative, an arbitrary function of pressure differential can be specified using the pressure_function property, such that

\[\dot m = \dot m_{\rm primary} + K_v*f(P_1 - P_2)\]

where \(f\) is the arbitrary function of a single argument.

downstream#

Return the Reactor or Reservoir object downstream of the flow device.

New in version 3.1.

draw(graph=None, *, graph_attr=None, node_attr=None, edge_attr=None)#

Draw as connection between upstream and downstream reactor or reservoir using graphviz.

Parameters:
  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw this flow controller connection. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

Returns:

A graphviz.graphs.BaseGraph object depicting the connection.

New in version 3.1.

edge_attr#

A dictionary containing draw attributes for the representation of the FlowDevice as a graphviz edge.See https://graphviz.org/docs/edges/ for a list of all usable attributes.

New in version 3.1.

flowdevice_type = 'PressureController'#
mass_flow_rate#

Get the mass flow rate [kg/s] through this device at the current reactor network time.

pressure_coeff#

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

pressure_function#

The relationship between mass flow rate and the pressure drop across a flow device. The mass flow rate [kg/s] is calculated given the pressure drop [Pa] and a coefficient set by a flow device specific function. Unless a user-defined pressure function is provided, the function returns the pressure difference across the device. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.pressure_function = lambda dP: dP**2

where FlowDevice is either a Valve or PressureController object.

New in version 3.0.

primary#

Primary FlowDevice used to compute this device’s mass flow rate.

New in version 3.0.

time_function#

The time dependence of a flow device. The mass flow rate [kg/s] is calculated for a Flow device, and multiplied by a function of time, which returns 1.0 unless a user-defined function is provided. The calculation of mass flow rate depends on the flow device.

>>> f = FlowDevice(res1, reactor1)
>>> f.time_function = lambda t: exp(-10 * (t - 0.5)**2)

where FlowDevice is either a Valve or MassFlowController object.

New in version 3.0.

type#

The type of the flow device.

upstream#

Return the Reactor or Reservoir object upstream of the flow device.

New in version 3.1.

Preconditioners#

AdaptivePreconditioner#

class cantera.AdaptivePreconditioner#

Bases: PreconditionerBase

ilut_drop_tol#

Property setting the linear solvers drop tolerance.

During factorization any element below the product of the drop tolerance and average magnitude is dropped.

Update the ILUT drop tolerance to a desired value as:
>>> precon.ilut_drop_tol = 1e-10

Default is 1e-10.

ilut_fill_factor#

Property setting the linear solvers fill factor.

During factorization, after row elimination, only some of the largest elements in the L and U in addition to the diagonal element are kept. The number of elements kept is computed from the fill factor (a ratio) relative to the initial number of nonzero elements.

Update the ILUT fill factor to a desired value as:
>>> precon.ilut_fill_factor = 2

Default is the state size divided by 4.

matrix#

Property to retrieve the latest internal preconditioner matrix.

precon_linear_solver_type = 'GMRES'#
precon_type = 'Adaptive'#
print_contents()#
threshold#

The threshold of the preconditioner is used to remove or prune any off diagonal elements below the given value inside of the preconditioner. In other words, after the preconditioner is formed by P = (I - gamma * Jac), the off diagonal values within P are compared with the threshold and removed if below it.

The goal of thresholding is to improve matrix sparsity while still providing a good preconditioner for the system.

Update the threshold to a desired value as:
>>> precon.threshold = 1e-8

Default is 0.0.

Drawing Reactor Networks#

These functions provide the implementation behind the draw methods of the corresponding classes.

cantera.drawnetwork.draw_reactor(r, graph=None, graph_attr=None, node_attr=None, print_state=False, species=None, species_units='percent')#

Draw ReactorBase object as graphviz dot node. The node is added to an existing graph if provided. Optionally include current reactor state in the node.

Parameters:
  • rReactorBase object or subclass.

  • graphgraphviz.graphs.BaseGraph object to which the reactor is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the node method invoked to draw the reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • print_state – Whether state information of the reactor is printed into the node. Defaults to False.

  • species – If print_state is True, define how species are to be printed. Options are 'X' and 'Y' for mole and mass fractions of all species, respectively, or an iterable that contains the desired species names as strings. Defaults to None.

  • species_units – Defines the units the species are displayed in as either "percent" or "ppm". Defaults to "percent".

Returns:

graphviz.graphs.BaseGraph object with reactor

New in version 3.1.

cantera.drawnetwork.draw_reactor_net(n, graph_attr=None, node_attr=None, edge_attr=None, heat_flow_attr=None, mass_flow_attr=None, moving_wall_edge_attr=None, surface_edge_attr=None, print_state=False, show_wall_velocity=True, **kwargs)#

Draw ReactorNet object as graphviz.graphs.DiGraph. Connecting flow controllers and walls are depicted as arrows.

Parameters:
  • nReactorNet object

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). node_attr defined in the reactor object itself have priority. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any edge (flow controllers, walls). edge_attr defined in the connection objects (subclasses of FlowDevice or walls) themselves have priority. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

  • heat_flow_attr – Same as edge_attr but only applied to edges representing walls. Default is {"color": "red", "style": "dashed"}.

  • mass_flow_attr – Same as edge_attr but only applied to edges representing FlowDevice objects.

  • moving_wall_edge_attr – Same as edge_attr but only applied to edges representing wall movement.

  • surface_edge_attr – Same as edge_attr but only applied to edges representing connections between a surface and its reactor. Default is {"style": "dotted", "arrowhead": "none"}.

  • print_state – Whether state information of each reactor is printed into the respective node. See draw_reactor for additional keywords to control this output. Defaults to False.

  • show_wall_velocity – If True, wall movement will be indicated by additional arrows with the corresponding wall velocity as a label.

  • kwargs – Additional keywords are passed on to each call of draw_reactor, draw_surface, draw_flow_controllers, and draw_walls.

Returns:

graphviz.graphs.BaseGraph object with reactor net.

New in version 3.1.

cantera.drawnetwork.draw_surface(surface, graph=None, graph_attr=None, node_attr=None, surface_edge_attr=None, print_state=False, **kwargs)#

Draw ReactorSurface object with its connected reactor.

Parameters:
  • surfaceReactorSurface object.

  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the node method invoked to draw the reactor. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • surface_edge_attr – Attributes to be passed to the edge method invoked to draw the connection between the surface and its reactor. See https://graphviz.org/docs/edges/ for a list of all usable attributes. Default is {"style": "dotted", "arrowhead": "none"}.

  • print_state – Whether state information of the reactor is printed into the node. See draw_reactor for additional keywords to control this output. Defaults to False.

  • kwargs – Additional keywords are passed on to draw_reactor.

Returns:

A graphviz.graphs.BaseGraph object depicting the surface and its reactor.

New in version 3.1.

cantera.drawnetwork.draw_flow_controllers(flow_controllers, graph=None, graph_attr=None, node_attr=None, edge_attr=None)#

Draw flow controller connections between reactors and reservoirs.

Parameters:
  • flow_controllers – Iterable of subtypes of FlowDevice.

  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw reactor connections. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

Returns:

A graphviz.graphs.BaseGraph object depicting the connections.

New in version 3.1.

cantera.drawnetwork.draw_walls(walls, graph=None, graph_attr=None, node_attr=None, edge_attr=None, moving_wall_edge_attr=None, show_wall_velocity=True)#

Draw wall connections between reactors and reservoirs.

Parameters:
  • walls – Iterable of subtypes of WallBase that connect reactors and reservoirs.

  • graphgraphviz.graphs.BaseGraph object to which the connection is added. If not provided, a new DiGraph is created. Defaults to None.

  • graph_attr – Attributes to be passed to the graphviz.Digraph function that control the general appearance of the drawn network. Has no effect if existing graph is provided. See https://graphviz.org/docs/graph/ for a list of all usable attributes.

  • node_attr – Attributes to be passed to the graphviz.Digraph function that control the default appearance of any node (reactors, reservoirs). Has no effect if existing graph is provided. See https://graphviz.org/docs/nodes/ for a list of all usable attributes.

  • edge_attr – Attributes to be passed to the edge method invoked to draw reactor connections. Default is {"color": "red", "style": "dashed"}. See https://graphviz.org/docs/edges/ for a list of all usable attributes.

  • moving_wall_edge_attr – Same as edge_attr but only applied to edges representing wall movement. Default is {"arrowtail": "icurveteecurve", "dir": "both", "style": "dotted", "arrowhead": "icurveteecurve"}.

  • show_wall_velocity – If True, wall movement will be indicated by additional arrows with the corresponding wall velocity as a label.

Returns:

A graphviz.graphs.BaseGraph object depicting the connections.

New in version 3.1.