Cantera
3.0.0
|
Adiabatic flow in a constant-area duct with homogeneous and heterogeneous reactions. More...
#include <FlowReactor.h>
Adiabatic flow in a constant-area duct with homogeneous and heterogeneous reactions.
Definition at line 16 of file FlowReactor.h.
Public Member Functions | |
string | type () const override |
String indicating the reactor model implemented. | |
bool | isOde () const override |
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs. | |
bool | timeIsIndependent () const override |
Indicates whether the governing equations for this reactor are functions of time or a spatial variable. | |
void | getState (double *y) override |
Not implemented; FlowReactor implements getStateDAE() instead. | |
void | getStateDae (double *y, double *ydot) override |
Get the current state and derivative vector of the reactor for a DAE solver. | |
void | initialize (double t0=0.0) override |
Initialize the reactor. | |
void | syncState () override |
Set the state of the reactor to correspond to the state of the associated ThermoPhase object. | |
void | updateState (double *y) override |
Set the state of the reactor to correspond to the state vector y. | |
void | eval (double t, double *LHS, double *RHS) override |
Not implemented; FlowReactor implements evalDae() instead. | |
void | evalDae (double t, double *y, double *ydot, double *residual) override |
Evaluate the reactor governing equations. | |
void | getConstraints (double *constraints) override |
Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
void | setMassFlowRate (double mdot) |
Set the mass flow rate through the reactor [kg/s]. | |
double | speed () const |
The current gas speed in the reactor [m/s]. | |
double | area () const |
The cross-sectional area of the reactor [m^2]. | |
double | distance () const |
void | setArea (double area) |
Sets the area of the reactor [m^2]. | |
double | surfaceAreaToVolumeRatio () const |
The ratio of the reactor's surface area to volume ratio [m^-1]. | |
void | setSurfaceAreaToVolumeRatio (double sa_to_vol) |
Set the reactor's surface area to volume ratio [m^-1]. | |
double | inletSurfaceAtol () const |
Get the steady state tolerances used to determine the initial state for surface coverages. | |
void | setInletSurfaceAtol (double atol) |
Set the steady state tolerances used to determine the initial state for surface coverages. | |
double | inletSurfaceRtol () const |
Get the steady state tolerances used to determine the initial state for surface coverages. | |
void | setInletSurfaceRtol (double rtol) |
Set the steady state tolerances used to determine the initial state for surface coverages. | |
double | inletSurfaceMaxSteps () const |
Get the steady state tolerances used to determine the initial state for surface coverages. | |
void | setInletSurfaceMaxSteps (int max_steps) |
Set the steady state tolerances used to determine the initial state for surface coverages. | |
double | inletSurfaceMaxErrorFailures () const |
Get the steady state tolerances used to determine the initial state for surface coverages. | |
void | setInletSurfaceMaxErrorFailures (int max_fails) |
Set the steady state tolerances used to determine the initial state for surface coverages. | |
size_t | componentIndex (const string &nm) const override |
Return the index in the solution vector for this reactor of the component named nm. | |
string | componentName (size_t k) override |
Return the name of the solution component with index i. | |
void | updateSurfaceState (double *y) override |
Update the state of SurfPhase objects attached to this reactor. | |
Public Member Functions inherited from IdealGasReactor | |
string | type () const override |
String indicating the reactor model implemented. | |
void | setThermoMgr (ThermoPhase &thermo) override |
Specify the mixture contained in the reactor. | |
void | getState (double *y) override |
Get the the current state of the reactor. | |
void | initialize (double t0=0.0) override |
Initialize the reactor. | |
void | eval (double t, double *LHS, double *RHS) override |
Evaluate the reactor governing equations. | |
void | updateState (double *y) override |
Set the state of the reactor to correspond to the state vector y. | |
size_t | componentIndex (const string &nm) const override |
Return the index in the solution vector for this reactor of the component named nm. | |
string | componentName (size_t k) override |
Return the name of the solution component with index i. | |
Public Member Functions inherited from Reactor | |
string | type () const override |
String indicating the reactor model implemented. | |
virtual bool | isOde () const |
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs. | |
virtual bool | timeIsIndependent () const |
Indicates whether the governing equations for this reactor are functions of time or a spatial variable. | |
template<class G > | |
void | insert (G &contents) |
Insert something into the reactor. | |
void | insert (shared_ptr< Solution > sol) |
void | setKineticsMgr (Kinetics &kin) override |
Specify chemical kinetics governing the reactor. | |
void | setChemistry (bool cflag=true) override |
Enable or disable changes in reactor composition due to chemical reactions. | |
bool | chemistryEnabled () const |
Returns true if changes in the reactor composition due to chemical reactions are enabled. | |
void | setEnergy (int eflag=1) override |
Set the energy equation on or off. | |
bool | energyEnabled () const |
Returns true if solution of the energy equation is enabled. | |
size_t | neq () |
Number of equations (state variables) for this reactor. | |
virtual void | getState (double *y) |
Get the the current state of the reactor. | |
virtual void | getStateDae (double *y, double *ydot) |
Get the current state and derivative vector of the reactor for a DAE solver. | |
void | initialize (double t0=0.0) override |
Initialize the reactor. | |
virtual void | eval (double t, double *LHS, double *RHS) |
Evaluate the reactor governing equations. | |
virtual void | evalDae (double t, double *y, double *ydot, double *residual) |
Evaluate the reactor governing equations. | |
virtual void | getConstraints (double *constraints) |
Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
void | syncState () override |
Set the state of the reactor to correspond to the state of the associated ThermoPhase object. | |
virtual void | updateState (double *y) |
Set the state of the reactor to correspond to the state vector y. | |
virtual size_t | nSensParams () const |
Number of sensitivity parameters associated with this reactor (including walls) | |
virtual void | addSensitivityReaction (size_t rxn) |
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase). | |
virtual void | addSensitivitySpeciesEnthalpy (size_t k) |
Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous phase) | |
virtual size_t | componentIndex (const string &nm) const |
Return the index in the solution vector for this reactor of the component named nm. | |
virtual string | componentName (size_t k) |
Return the name of the solution component with index i. | |
void | setAdvanceLimits (const double *limits) |
Set absolute step size limits during advance. | |
bool | hasAdvanceLimits () const |
Check whether Reactor object uses advance limits. | |
bool | getAdvanceLimits (double *limits) const |
Retrieve absolute step size limits during advance. | |
void | setAdvanceLimit (const string &nm, const double limit) |
Set individual step size limit for component name nm | |
virtual Eigen::SparseMatrix< double > | jacobian () |
Calculate the Jacobian of a specific Reactor specialization. | |
Eigen::SparseMatrix< double > | finiteDifferenceJacobian () |
Calculate the reactor-specific Jacobian using a finite difference method. | |
virtual void | setDerivativeSettings (AnyMap &settings) |
Use this to set the kinetics objects derivative settings. | |
virtual void | applySensitivity (double *params) |
Set reaction rate multipliers based on the sensitivity variables in params. | |
virtual void | resetSensitivity (double *params) |
Reset the reaction rate multipliers. | |
virtual bool | preconditionerSupported () const |
Return a false if preconditioning is not supported or true otherwise. | |
Public Member Functions inherited from ReactorBase | |
ReactorBase (const string &name="(none)") | |
ReactorBase (const ReactorBase &)=delete | |
ReactorBase & | operator= (const ReactorBase &)=delete |
virtual string | type () const |
String indicating the reactor model implemented. | |
string | name () const |
Return the name of this reactor. | |
void | setName (const string &name) |
Set the name of this reactor. | |
void | restoreState () |
Set the state of the Phase object associated with this reactor to the reactor's current state. | |
virtual void | syncState () |
Set the state of the reactor to correspond to the state of the associated ThermoPhase object. | |
ThermoPhase & | contents () |
return a reference to the contents. | |
const ThermoPhase & | contents () const |
double | residenceTime () |
Return the residence time (s) of the contents of this reactor, based on the outlet mass flow rates and the mass of the reactor contents. | |
ReactorNet & | network () |
The ReactorNet that this reactor belongs to. | |
void | setNetwork (ReactorNet *net) |
Set the ReactorNet that this reactor belongs to. | |
void | setInitialVolume (double vol) |
Set the initial reactor volume. By default, the volume is 1.0 m^3. | |
void | addInlet (FlowDevice &inlet) |
Connect an inlet FlowDevice to this reactor. | |
void | addOutlet (FlowDevice &outlet) |
Connect an outlet FlowDevice to this reactor. | |
FlowDevice & | inlet (size_t n=0) |
Return a reference to the n-th inlet FlowDevice connected to this reactor. | |
FlowDevice & | outlet (size_t n=0) |
Return a reference to the n-th outlet FlowDevice connected to this reactor. | |
size_t | nInlets () |
Return the number of inlet FlowDevice objects connected to this reactor. | |
size_t | nOutlets () |
Return the number of outlet FlowDevice objects connected to this reactor. | |
size_t | nWalls () |
Return the number of Wall objects connected to this reactor. | |
void | addWall (WallBase &w, int lr) |
Insert a Wall between this reactor and another reactor. | |
WallBase & | wall (size_t n) |
Return a reference to the n-th Wall connected to this reactor. | |
virtual void | addSurface (ReactorSurface *surf) |
ReactorSurface * | surface (size_t n) |
Return a reference to the n-th ReactorSurface connected to this reactor. | |
virtual size_t | nSurfs () |
Return the number of surfaces in a reactor. | |
double | volume () const |
Returns the current volume (m^3) of the reactor. | |
double | density () const |
Returns the current density (kg/m^3) of the reactor's contents. | |
double | temperature () const |
Returns the current temperature (K) of the reactor's contents. | |
double | enthalpy_mass () const |
Returns the current enthalpy (J/kg) of the reactor's contents. | |
double | intEnergy_mass () const |
Returns the current internal energy (J/kg) of the reactor's contents. | |
double | pressure () const |
Returns the current pressure (Pa) of the reactor. | |
double | mass () const |
Returns the mass (kg) of the reactor's contents. | |
const double * | massFractions () const |
Return the vector of species mass fractions. | |
double | massFraction (size_t k) const |
Return the mass fraction of the k-th species. | |
Protected Attributes | |
double | m_rho = NAN |
Density [kg/m^3]. First component of the state vector. | |
double | m_u = -1.0 |
Axial velocity [m/s]. Second component of the state vector. | |
double | m_P = NAN |
Pressure [Pa]. Third component of the state vector. | |
double | m_T = NAN |
Temperature [K]. Fourth component of the state vector. | |
const size_t | m_offset_Y = 4 |
offset to the species equations | |
double | m_area = 1.0 |
reactor area [m^2] | |
double | m_sa_to_vol = -1.0 |
reactor surface area to volume ratio [m^-1] | |
vector< double > | m_sdot_temp |
temporary storage for surface species production rates | |
vector< double > | m_hk |
temporary storage for species partial molar enthalpies | |
double | m_ss_rtol = 1e-7 |
steady-state relative tolerance, used to determine initial surface coverages | |
double | m_ss_atol = 1e-14 |
steady-state absolute tolerance, used to determine initial surface coverages | |
int | m_max_ss_steps = 20000 |
maximum number of steady-state coverage integrator-steps | |
int | m_max_ss_error_fails = 10 |
maximum number of steady-state integrator error test failures | |
Protected Attributes inherited from IdealGasReactor | |
vector< double > | m_uk |
Species molar internal energies. | |
Protected Attributes inherited from Reactor | |
Kinetics * | m_kin = nullptr |
Pointer to the homogeneous Kinetics object that handles the reactions. | |
double | m_vdot = 0.0 |
net rate of volume change from moving walls [m^3/s] | |
double | m_Qdot = 0.0 |
net heat transfer into the reactor, through walls [W] | |
double | m_mass = 0.0 |
total mass | |
vector< double > | m_work |
vector< double > | m_sdot |
Production rates of gas phase species on surfaces [kmol/s]. | |
vector< double > | m_wdot |
Species net molar production rates. | |
vector< double > | m_uk |
Species molar internal energies. | |
bool | m_chem = false |
bool | m_energy = true |
size_t | m_nv = 0 |
size_t | m_nv_surf |
vector< double > | m_advancelimits |
!< Number of variables associated with reactor surfaces | |
vector< SensitivityParameter > | m_sensParams |
vector< Eigen::Triplet< double > > | m_jac_trips |
Vector of triplets representing the jacobian. | |
Protected Attributes inherited from ReactorBase | |
size_t | m_nsp = 0 |
Number of homogeneous species in the mixture. | |
ThermoPhase * | m_thermo = nullptr |
double | m_vol = 1.0 |
Current volume of the reactor [m^3]. | |
double | m_enthalpy = 0.0 |
Current specific enthalpy of the reactor [J/kg]. | |
double | m_intEnergy = 0.0 |
Current internal energy of the reactor [J/kg]. | |
double | m_pressure = 0.0 |
Current pressure in the reactor [Pa]. | |
vector< double > | m_state |
vector< FlowDevice * > | m_inlet |
vector< FlowDevice * > | m_outlet |
vector< WallBase * > | m_wall |
vector< ReactorSurface * > | m_surfaces |
vector< int > | m_lr |
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wall. | |
string | m_name |
ReactorNet * | m_net = nullptr |
The ReactorNet that this reactor is part of. | |
Additional Inherited Members | |
Protected Member Functions inherited from Reactor | |
virtual size_t | speciesIndex (const string &nm) const |
Return the index in the solution vector for this reactor of the species named nm, in either the homogeneous phase or a surface phase, relative to the start of the species terms. | |
virtual void | evalWalls (double t) |
Evaluate terms related to Walls. | |
virtual void | evalSurfaces (double *LHS, double *RHS, double *sdot) |
Evaluate terms related to surface reactions. | |
virtual void | evalSurfaces (double *RHS, double *sdot) |
virtual void | updateSurfaceState (double *y) |
Update the state of SurfPhase objects attached to this reactor. | |
virtual void | updateConnected (bool updatePressure) |
Update the state information needed by connected reactors, flow devices, and reactor walls. | |
virtual void | getSurfaceInitialConditions (double *y) |
Get initial conditions for SurfPhase objects attached to this reactor. | |
|
inlineoverridevirtual |
String indicating the reactor model implemented.
Usually corresponds to the name of the derived class.
Reimplemented from ReactorBase.
Definition at line 21 of file FlowReactor.h.
|
inlineoverridevirtual |
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs.
In the first case, this class implements the eval() method. In the second case, this class implements the evalDae() method.
Reimplemented from Reactor.
Definition at line 25 of file FlowReactor.h.
|
inlineoverridevirtual |
Indicates whether the governing equations for this reactor are functions of time or a spatial variable.
All reactors in a network must have the same value.
Reimplemented from Reactor.
Definition at line 29 of file FlowReactor.h.
|
inlineoverridevirtual |
Not implemented; FlowReactor implements getStateDAE() instead.
Reimplemented from Reactor.
Definition at line 34 of file FlowReactor.h.
|
overridevirtual |
Get the current state and derivative vector of the reactor for a DAE solver.
[out] | y | state vector representing the initial state of the reactor |
[out] | ydot | state vector representing the initial derivatives of the reactor |
Reimplemented from Reactor.
Definition at line 20 of file FlowReactor.cpp.
|
overridevirtual |
Initialize the reactor.
Called automatically by ReactorNet::initialize.
Reimplemented from ReactorBase.
Definition at line 155 of file FlowReactor.cpp.
|
overridevirtual |
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
This is the inverse of restoreState(). Calling this will trigger integrator reinitialization.
Reimplemented from ReactorBase.
Definition at line 185 of file FlowReactor.cpp.
|
overridevirtual |
Set the state of the reactor to correspond to the state vector y.
Reimplemented from Reactor.
Definition at line 193 of file FlowReactor.cpp.
|
inlineoverridevirtual |
Not implemented; FlowReactor implements evalDae() instead.
Reimplemented from Reactor.
Definition at line 44 of file FlowReactor.h.
|
overridevirtual |
Evaluate the reactor governing equations.
Called by ReactorNet::eval.
[in] | t | time. |
[in] | y | solution vector, length neq() |
[in] | ydot | rate of change of solution vector, length neq() |
[out] | residual | residuals vector, length neq() |
Reimplemented from Reactor.
Definition at line 257 of file FlowReactor.cpp.
|
overridevirtual |
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Reimplemented from Reactor.
Definition at line 346 of file FlowReactor.cpp.
void setMassFlowRate | ( | double | mdot | ) |
Set the mass flow rate through the reactor [kg/s].
Definition at line 210 of file FlowReactor.cpp.
|
inline |
The current gas speed in the reactor [m/s].
Definition at line 56 of file FlowReactor.h.
|
inline |
The cross-sectional area of the reactor [m^2].
Definition at line 61 of file FlowReactor.h.
double distance | ( | ) | const |
Definition at line 216 of file FlowReactor.cpp.
void setArea | ( | double | area | ) |
Sets the area of the reactor [m^2].
Definition at line 227 of file FlowReactor.cpp.
double surfaceAreaToVolumeRatio | ( | ) | const |
The ratio of the reactor's surface area to volume ratio [m^-1].
Definition at line 233 of file FlowReactor.cpp.
|
inline |
Set the reactor's surface area to volume ratio [m^-1].
Definition at line 78 of file FlowReactor.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 84 of file FlowReactor.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 90 of file FlowReactor.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 96 of file FlowReactor.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 102 of file FlowReactor.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 108 of file FlowReactor.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 114 of file FlowReactor.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 120 of file FlowReactor.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 126 of file FlowReactor.h.
|
overridevirtual |
Return the index in the solution vector for this reactor of the component named nm.
Possible values for nm are "density", "speed", "pressure", "temperature", the name of a homogeneous phase species, or the name of a surface species.
Reimplemented from Reactor.
Definition at line 353 of file FlowReactor.cpp.
|
overridevirtual |
Return the name of the solution component with index i.
Reimplemented from Reactor.
Definition at line 371 of file FlowReactor.cpp.
|
overridevirtual |
Update the state of SurfPhase objects attached to this reactor.
Reimplemented from Reactor.
Definition at line 244 of file FlowReactor.cpp.
|
protected |
Density [kg/m^3]. First component of the state vector.
Definition at line 142 of file FlowReactor.h.
|
protected |
Axial velocity [m/s]. Second component of the state vector.
Definition at line 144 of file FlowReactor.h.
|
protected |
Pressure [Pa]. Third component of the state vector.
Definition at line 146 of file FlowReactor.h.
|
protected |
Temperature [K]. Fourth component of the state vector.
Definition at line 148 of file FlowReactor.h.
|
protected |
offset to the species equations
Definition at line 150 of file FlowReactor.h.
|
protected |
reactor area [m^2]
Definition at line 152 of file FlowReactor.h.
|
protected |
reactor surface area to volume ratio [m^-1]
Definition at line 154 of file FlowReactor.h.
|
protected |
temporary storage for surface species production rates
Definition at line 156 of file FlowReactor.h.
|
protected |
temporary storage for species partial molar enthalpies
Definition at line 158 of file FlowReactor.h.
|
protected |
steady-state relative tolerance, used to determine initial surface coverages
Definition at line 160 of file FlowReactor.h.
|
protected |
steady-state absolute tolerance, used to determine initial surface coverages
Definition at line 162 of file FlowReactor.h.
|
protected |
maximum number of steady-state coverage integrator-steps
Definition at line 164 of file FlowReactor.h.
|
protected |
maximum number of steady-state integrator error test failures
Definition at line 166 of file FlowReactor.h.