Zero-Dimensional Reactor Networks¶
Defining Functions¶
- class cantera.Func1(c)¶
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 aFunc1
object.
Base Classes¶
ReactorBase¶
- class cantera.ReactorBase(ThermoPhase contents=None, name=None, volume=None, *)¶
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.
- 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'¶
- surfaces¶
List of reacting surfaces installed on this reactor
- 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.
- 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(upstream, downstream, name=None, *)¶
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.
- flowdevice_type = 'None'¶
- mass_flow_rate¶
Get the mass flow rate [kg/s] through this device at the current reactor network time.
- set_pressure_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_pressure_function(lambda dP: dP**2)
where FlowDevice is either a Valve or PressureController object.
- set_time_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_time_function(lambda t: exp(-10 * (t - 0.5)**2))
where FlowDevice is either a Valve or MassFlowController object.
- type¶
The type of the flow device.
Reactor Networks¶
- class cantera.ReactorNet(reactors=())¶
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(self, Reactor r)¶
Add a reactor to the network.
- advance(self, double t, bool apply_limit=True)¶
Advance the state of the reactor network in time from the current time towards time
t
in seconds, taking as many integrator time steps as necessary. Ifapply_limit
is true and an advance limit is specified, the reactor state at the end of the timestep is estimated prior to advancing. If the difference exceed limits, the end time is reduced by half until the projected end state remains within specified limits. Returns the time 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(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 lengthmax_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, for example
'reactor1: CH4'
.
- get_derivative(self, k)¶
Get the k-th time derivative of the state vector of the reactor network.
- 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.
- global_component_index(self, name, int reactor)¶
Returns the index of a component named
name
of a reactor with indexreactor
within the global state vector. That is, this determines the absolute index of the component, wherereactor
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.
- initialize(self)¶
Force initialization of the integrator after initial setup.
- max_err_test_fails¶
The maximum number of error test failures permitted by the CVODES integrator in a single time step.
- max_steps¶
The maximum number of internal integration time-steps that CVODES is allowed to take before reaching the next output time.
- 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
andWall
in the system. Equal to:Reactor
andIdealGasReactor
:n_species
+ 3 (mass, volume, internal energy or temperature).ConstPressureReactor
andIdealGasConstPressureReactor
: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 (that is, rows) is:
0 - mass
1 - volume
2 - internal energy or temperature
3+ - mass fractions of the species
ConstPressureReactor
orIdealGasConstPressureReactor
: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 reactorr
with respect to the parameterp
.component
can be a string or an integer. Seecomponent_index
andsensitivities
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.
- step(self)¶
Take a single internal time step. The time after taking the step is returned.
- time¶
The current time [s].
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', **kwargs)¶
Bases:
cantera._cantera.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
andValve
.- 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 orderReactor
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.yaml') >>> 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 isFalse
, the reactor composition is held constant.
- component_index(self, 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(self, int i)¶
Returns the name of the component with index
i
within the array of variables returned byget_state
. This is the inverse ofcomponent_index
.
- energy_enabled¶
True
when the energy equation is being solved for this reactor. When this isFalse
, the reactor temperature is held constant.
- get_state(self)¶
Get the state vector of the reactor.
The order of the variables (that is, rows) is:
0 - mass
1 - volume
2 - internal energy or temperature
3+ - mass fractions of the species
ConstPressureReactor
orIdealGasConstPressureReactor
: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, orcomponent_name
to determine the name from the index.
- insert(self, _SolutionBase solution)¶
Set
solution
to be the object used to compute thermodynamic properties and kinetic rates for this reactor.
- n_vars¶
The number of state variables in the reactor. Equal to:
Reactor
andIdealGasReactor
:n_species
+ 3 (mass, volume, internal energy or temperature).ConstPressureReactor
andIdealGasConstPressureReactor
:n_species
+ 2 (mass, enthalpy or temperature).
- reactor_type = 'Reactor'¶
- set_advance_limit(self, name, limit)¶
Limit absolute change of component
name
duringReactorNet.advance
. (positivelimit
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.).
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 orderReactor
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.yaml') >>> 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 orderReactor
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.yaml') >>> 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 orderReactor
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.yaml') >>> 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 orderReactor
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.yaml') >>> 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
ExtensibleReactor¶
- class cantera.ExtensibleReactor(contents=None, *, name=None, energy='on')¶
Bases:
cantera._cantera.Reactor
ExtensibleReactor(*args, **kwargs)
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_
, orreplace_
, indicating whether the 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 thanNone
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 anafter
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
(lengthn_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
(lengthn_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 ofy
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 ofy
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 componenti
of the state vector, the time derivativedy[i]/dt
is calculated asLHS[i] * dy[i]/dt = RHS[i]
.LHS
andRHS
are arrays of lengthn_vars
.eval_walls(self, t : double) -> None
Responsible for calculating the net rate of volume change
vdot
and the net rate of heat transferqdot
caused by walls connected to this reactor.eval_surfaces(LHS : double[:], RHS : double[:], sdot : double[:]) -> None
Responsible for calculating the
LHS
andRHS
(length: total number of surface species in all surfaces) of the ODEs for surface species coverages, and the molar production rate of bulk phase speciessdot
(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*)')}¶
- n_vars¶
Get/Set the number of state variables in the reactor.
- qdot¶
Get/Set the net heat transfer rate (for example, through walls) [W]
- reactor_type = 'ExtensibleReactor'¶
- restore_surface_state(self, n)¶
Set the state of the thermo object for surface
n
to correspond to the state of that surface
- restore_thermo_state(self)¶
Set the state of the thermo object to correspond to the state of the reactor.
- vdot¶
Get/Set the net rate of volume change (for example, from moving walls) [m^3/s]
ExtensibleIdealGasReactor¶
- class cantera.ExtensibleIdealGasReactor(contents=None, *, name=None, energy='on')¶
Bases:
cantera._cantera.ExtensibleReactor
A variant of
ExtensibleReactor
where the base behavior corresponds to theIdealGasReactor
class.- reactor_type = 'ExtensibleIdealGasReactor'¶
ExtensibleConstPressureReactor¶
- class cantera.ExtensibleConstPressureReactor(contents=None, *, name=None, energy='on')¶
Bases:
cantera._cantera.ExtensibleReactor
A variant of
ExtensibleReactor
where the base behavior corresponds to theConstPressureReactor
class.- reactor_type = 'ExtensibleConstPressureReactor'¶
ExtensibleIdealGasConstPressureReactor¶
- class cantera.ExtensibleIdealGasConstPressureReactor(contents=None, *, name=None, energy='on')¶
Bases:
cantera._cantera.ExtensibleReactor
A variant of
ExtensibleReactor
where the base behavior corresponds to theIdealGasConstPressureReactor
class.- reactor_type = 'ExtensibleIdealGasConstPressureReactor'¶
Walls¶
Wall¶
- class cantera.Wall(left, right, *, name=None, A=None, K=None, U=None, Q=None, velocity=None)¶
Bases:
cantera._cantera.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\).
- 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]
- 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
.
- wall_type = 'Wall'¶
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
- 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)¶
Add this
ReactorSurface
to the specifiedReactor
- kinetics¶
The
InterfaceKinetics
object used for calculating reaction rates on this surface.
Flow Controllers¶
MassFlowController¶
- class cantera.MassFlowController(upstream, downstream, name=None, *, mdot=1.)¶
Bases:
cantera._cantera.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 the property
mass_flow_coeff
and the methodset_time_function
, respectively. The propertymass_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.
- 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
set_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 calls theset_time_function
method.>>> mfc.mass_flow_rate = 0.3 >>> mfc.mass_flow_rate = lambda t: 2.5 * exp(-10 * (t - 0.5)**2)
- set_pressure_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_pressure_function(lambda dP: dP**2)
where FlowDevice is either a Valve or PressureController object.
- set_time_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_time_function(lambda t: exp(-10 * (t - 0.5)**2))
where FlowDevice is either a Valve or MassFlowController object.
- type¶
The type of the flow device.
Valve¶
- class cantera.Valve(upstream, downstream, name=None, *, K=1.)¶
Bases:
cantera._cantera.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 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 = 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 method
set_time_function
, such that\[\dot m = K_v*g(t)*f(P_1 - P_2) \]See the documentation for the
valve_coeff
property as well as theset_pressure_function
andset_time_function
methods 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 aValve
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.- flowdevice_type = 'Valve'¶
- mass_flow_rate¶
Get the mass flow rate [kg/s] through this device at the current reactor network time.
- set_pressure_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_pressure_function(lambda dP: dP**2)
where FlowDevice is either a Valve or PressureController object.
- set_time_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_time_function(lambda t: exp(-10 * (t - 0.5)**2))
where FlowDevice is either a Valve or MassFlowController object.
- type¶
The type of the flow device.
- 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(upstream, downstream, name=None, *, master=None, K=1.)¶
Bases:
cantera._cantera.FlowDevice
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 correspondingPressureController
is installed on the outlet of the reactor. ThePressureController
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). \]As an alternative, an arbitrary function of pressure differential can be specified using the method
set_pressure_function
, such that\[\dot m = \dot m_{\rm master} + K_v*f(P_1 - P_2) \]where \(f\) is the arbitrary function of a single argument.
- 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.
- set_master(self, FlowDevice d)¶
Set the “master”
FlowDevice
used to compute this device’s mass flow rate.
- set_pressure_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_pressure_function(lambda dP: dP**2)
where FlowDevice is either a Valve or PressureController object.
- set_time_function(self, k)¶
Set 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. The calculation of mass flow rate depends to the flow device.
>>> F = FlowDevice(res1, reactor1) >>> F.set_time_function(lambda t: exp(-10 * (t - 0.5)**2))
where FlowDevice is either a Valve or MassFlowController object.
- type¶
The type of the flow device.