Cantera  3.1.0a1
FlowReactor Class Reference

Adiabatic flow in a constant-area duct with homogeneous and heterogeneous reactions. More...

#include <FlowReactor.h>

Inheritance diagram for FlowReactor:
[legend]

Detailed Description

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. More...
 
bool isOde () const override
 Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs. More...
 
bool timeIsIndependent () const override
 Indicates whether the governing equations for this reactor are functions of time or a spatial variable. More...
 
void getState (double *y) override
 Not implemented; FlowReactor implements getStateDAE() instead. More...
 
void getStateDae (double *y, double *ydot) override
 Get the current state and derivative vector of the reactor for a DAE solver. More...
 
void initialize (double t0=0.0) override
 Initialize the reactor. More...
 
void syncState () override
 Set the state of the reactor to correspond to the state of the associated ThermoPhase object. More...
 
void updateState (double *y) override
 Set the state of the reactor to correspond to the state vector y. More...
 
void eval (double t, double *LHS, double *RHS) override
 Not implemented; FlowReactor implements evalDae() instead. More...
 
void evalDae (double t, double *y, double *ydot, double *residual) override
 Evaluate the reactor governing equations. More...
 
void getConstraints (double *constraints) override
 Given a vector of length neq(), mark which variables should be considered algebraic constraints. More...
 
void setMassFlowRate (double mdot)
 Set the mass flow rate through the reactor [kg/s]. More...
 
double speed () const
 The current gas speed in the reactor [m/s]. More...
 
double area () const
 The cross-sectional area of the reactor [m^2]. More...
 
void setArea (double area)
 Sets the area of the reactor [m^2]. More...
 
double surfaceAreaToVolumeRatio () const
 The ratio of the reactor's surface area to volume ratio [m^-1]. More...
 
void setSurfaceAreaToVolumeRatio (double sa_to_vol)
 Set the reactor's surface area to volume ratio [m^-1]. More...
 
double inletSurfaceAtol () const
 Get the steady state tolerances used to determine the initial state for surface coverages. More...
 
void setInletSurfaceAtol (double atol)
 Set the steady state tolerances used to determine the initial state for surface coverages. More...
 
double inletSurfaceRtol () const
 Get the steady state tolerances used to determine the initial state for surface coverages. More...
 
void setInletSurfaceRtol (double rtol)
 Set the steady state tolerances used to determine the initial state for surface coverages. More...
 
double inletSurfaceMaxSteps () const
 Get the steady state tolerances used to determine the initial state for surface coverages. More...
 
void setInletSurfaceMaxSteps (int max_steps)
 Set the steady state tolerances used to determine the initial state for surface coverages. More...
 
double inletSurfaceMaxErrorFailures () const
 Get the steady state tolerances used to determine the initial state for surface coverages. More...
 
void setInletSurfaceMaxErrorFailures (int max_fails)
 Set the steady state tolerances used to determine the initial state for surface coverages. More...
 
size_t componentIndex (const string &nm) const override
 Return the index in the solution vector for this reactor of the component named nm. More...
 
string componentName (size_t k) override
 Return the name of the solution component with index i. More...
 
void updateSurfaceState (double *y) override
 Update the state of SurfPhase objects attached to this reactor. More...
 
- Public Member Functions inherited from IdealGasReactor
string type () const override
 String indicating the reactor model implemented. More...
 
void setThermoMgr (ThermoPhase &thermo) override
 Specify the mixture contained in the reactor. More...
 
void getState (double *y) override
 Get the the current state of the reactor. More...
 
void initialize (double t0=0.0) override
 Initialize the reactor. More...
 
void eval (double t, double *LHS, double *RHS) override
 Evaluate the reactor governing equations. More...
 
void updateState (double *y) override
 Set the state of the reactor to correspond to the state vector y. More...
 
size_t componentIndex (const string &nm) const override
 Return the index in the solution vector for this reactor of the component named nm. More...
 
string componentName (size_t k) override
 Return the name of the solution component with index i. More...
 
- Public Member Functions inherited from Reactor
string type () const override
 String indicating the reactor model implemented. More...
 
template<class G >
void insert (G &contents)
 Insert something into the reactor. More...
 
void insert (shared_ptr< Solution > sol)
 
void setKineticsMgr (Kinetics &kin) override
 Specify chemical kinetics governing the reactor. More...
 
void setChemistry (bool cflag=true) override
 Enable or disable changes in reactor composition due to chemical reactions. More...
 
bool chemistryEnabled () const
 Returns true if changes in the reactor composition due to chemical reactions are enabled. More...
 
void setEnergy (int eflag=1) override
 Set the energy equation on or off. More...
 
bool energyEnabled () const
 Returns true if solution of the energy equation is enabled. More...
 
size_t neq ()
 Number of equations (state variables) for this reactor. More...
 
void initialize (double t0=0.0) override
 Initialize the reactor. More...
 
void syncState () override
 Set the state of the reactor to correspond to the state of the associated ThermoPhase object. More...
 
virtual size_t nSensParams () const
 Number of sensitivity parameters associated with this reactor (including walls) More...
 
virtual void addSensitivityReaction (size_t rxn)
 Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase). More...
 
virtual void addSensitivitySpeciesEnthalpy (size_t k)
 Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous phase) More...
 
void setAdvanceLimits (const double *limits)
 Set absolute step size limits during advance. More...
 
bool hasAdvanceLimits () const
 Check whether Reactor object uses advance limits. More...
 
bool getAdvanceLimits (double *limits) const
 Retrieve absolute step size limits during advance. More...
 
void setAdvanceLimit (const string &nm, const double limit)
 Set individual step size limit for component name nm More...
 
virtual Eigen::SparseMatrix< double > jacobian ()
 Calculate the Jacobian of a specific Reactor specialization. More...
 
Eigen::SparseMatrix< double > finiteDifferenceJacobian ()
 Calculate the reactor-specific Jacobian using a finite difference method. More...
 
virtual void setDerivativeSettings (AnyMap &settings)
 Use this to set the kinetics objects derivative settings. More...
 
virtual void applySensitivity (double *params)
 Set reaction rate multipliers based on the sensitivity variables in params. More...
 
virtual void resetSensitivity (double *params)
 Reset the reaction rate multipliers. More...
 
virtual bool preconditionerSupported () const
 Return a false if preconditioning is not supported or true otherwise. More...
 
- Public Member Functions inherited from ReactorBase
 ReactorBase (const string &name="(none)")
 
 ReactorBase (const ReactorBase &)=delete
 
ReactorBaseoperator= (const ReactorBase &)=delete
 
string name () const
 Return the name of this reactor. More...
 
void setName (const string &name)
 Set the name of this reactor. More...
 
void restoreState ()
 Set the state of the Phase object associated with this reactor to the reactor's current state. More...
 
ThermoPhasecontents ()
 return a reference to the contents. More...
 
const ThermoPhasecontents () 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. More...
 
ReactorNetnetwork ()
 The ReactorNet that this reactor belongs to. More...
 
void setNetwork (ReactorNet *net)
 Set the ReactorNet that this reactor belongs to. More...
 
void setInitialVolume (double vol)
 Set the initial reactor volume. By default, the volume is 1.0 m^3. More...
 
void addInlet (FlowDevice &inlet)
 Connect an inlet FlowDevice to this reactor. More...
 
void addOutlet (FlowDevice &outlet)
 Connect an outlet FlowDevice to this reactor. More...
 
FlowDeviceinlet (size_t n=0)
 Return a reference to the n-th inlet FlowDevice connected to this reactor. More...
 
FlowDeviceoutlet (size_t n=0)
 Return a reference to the n-th outlet FlowDevice connected to this reactor. More...
 
size_t nInlets ()
 Return the number of inlet FlowDevice objects connected to this reactor. More...
 
size_t nOutlets ()
 Return the number of outlet FlowDevice objects connected to this reactor. More...
 
size_t nWalls ()
 Return the number of Wall objects connected to this reactor. More...
 
void addWall (WallBase &w, int lr)
 Insert a Wall between this reactor and another reactor. More...
 
WallBasewall (size_t n)
 Return a reference to the n-th Wall connected to this reactor. More...
 
virtual void addSurface (ReactorSurface *surf)
 
ReactorSurfacesurface (size_t n)
 Return a reference to the n-th ReactorSurface connected to this reactor. More...
 
virtual size_t nSurfs ()
 Return the number of surfaces in a reactor. More...
 
double volume () const
 Returns the current volume (m^3) of the reactor. More...
 
double density () const
 Returns the current density (kg/m^3) of the reactor's contents. More...
 
double temperature () const
 Returns the current temperature (K) of the reactor's contents. More...
 
double enthalpy_mass () const
 Returns the current enthalpy (J/kg) of the reactor's contents. More...
 
double intEnergy_mass () const
 Returns the current internal energy (J/kg) of the reactor's contents. More...
 
double pressure () const
 Returns the current pressure (Pa) of the reactor. More...
 
double mass () const
 Returns the mass (kg) of the reactor's contents. More...
 
const double * massFractions () const
 Return the vector of species mass fractions. More...
 
double massFraction (size_t k) const
 Return the mass fraction of the k-th species. More...
 

Protected Attributes

double m_rho = NAN
 Density [kg/m^3]. First component of the state vector. More...
 
double m_u = -1.0
 Axial velocity [m/s]. Second component of the state vector. More...
 
double m_P = NAN
 Pressure [Pa]. Third component of the state vector. More...
 
double m_T = NAN
 Temperature [K]. Fourth component of the state vector. More...
 
const size_t m_offset_Y = 4
 offset to the species equations More...
 
double m_area = 1.0
 reactor area [m^2] More...
 
double m_sa_to_vol = -1.0
 reactor surface area to volume ratio [m^-1] More...
 
vector< double > m_sdot_temp
 temporary storage for surface species production rates More...
 
vector< double > m_hk
 temporary storage for species partial molar enthalpies More...
 
double m_ss_rtol = 1e-7
 steady-state relative tolerance, used to determine initial surface coverages More...
 
double m_ss_atol = 1e-14
 steady-state absolute tolerance, used to determine initial surface coverages More...
 
int m_max_ss_steps = 20000
 maximum number of steady-state coverage integrator-steps More...
 
int m_max_ss_error_fails = 10
 maximum number of steady-state integrator error test failures More...
 
- Protected Attributes inherited from IdealGasReactor
vector< double > m_uk
 Species molar internal energies. More...
 
- Protected Attributes inherited from Reactor
Kineticsm_kin = nullptr
 Pointer to the homogeneous Kinetics object that handles the reactions. More...
 
double m_vdot = 0.0
 net rate of volume change from moving walls [m^3/s] More...
 
double m_Qdot = 0.0
 net heat transfer into the reactor, through walls [W] More...
 
double m_mass = 0.0
 total mass More...
 
vector< double > m_work
 
vector< double > m_sdot
 Production rates of gas phase species on surfaces [kmol/s]. More...
 
vector< double > m_wdot
 Species net molar production rates. More...
 
vector< double > m_uk
 Species molar internal energies. More...
 
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 More...
 
vector< SensitivityParameter > m_sensParams
 
vector< Eigen::Triplet< double > > m_jac_trips
 Vector of triplets representing the jacobian. More...
 
- Protected Attributes inherited from ReactorBase
size_t m_nsp = 0
 Number of homogeneous species in the mixture. More...
 
ThermoPhasem_thermo = nullptr
 
double m_vol = 1.0
 Current volume of the reactor [m^3]. More...
 
double m_enthalpy = 0.0
 Current specific enthalpy of the reactor [J/kg]. More...
 
double m_intEnergy = 0.0
 Current internal energy of the reactor [J/kg]. More...
 
double m_pressure = 0.0
 Current pressure in the reactor [Pa]. More...
 
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. More...
 
string m_name
 
ReactorNetm_net = nullptr
 The ReactorNet that this reactor is part of. More...
 

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. More...
 
virtual void evalWalls (double t)
 Evaluate terms related to Walls. More...
 
virtual void evalSurfaces (double *LHS, double *RHS, double *sdot)
 Evaluate terms related to surface reactions. More...
 
virtual void evalSurfaces (double *RHS, double *sdot)
 
virtual void updateConnected (bool updatePressure)
 Update the state information needed by connected reactors, flow devices, and reactor walls. More...
 
virtual void getSurfaceInitialConditions (double *y)
 Get initial conditions for SurfPhase objects attached to this reactor. More...
 

Member Function Documentation

◆ type()

string type ( ) const
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.

◆ isOde()

bool isOde ( ) const
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.

◆ timeIsIndependent()

bool timeIsIndependent ( ) const
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.

◆ getState()

void getState ( double *  y)
inlineoverridevirtual

Not implemented; FlowReactor implements getStateDAE() instead.

Reimplemented from Reactor.

Definition at line 34 of file FlowReactor.h.

◆ getStateDae()

void getStateDae ( double *  y,
double *  ydot 
)
overridevirtual

Get the current state and derivative vector of the reactor for a DAE solver.

Parameters
[out]ystate vector representing the initial state of the reactor
[out]ydotstate vector representing the initial derivatives of the reactor

Reimplemented from Reactor.

Definition at line 20 of file FlowReactor.cpp.

◆ initialize()

void initialize ( double  t0 = 0.0)
overridevirtual

Initialize the reactor.

Called automatically by ReactorNet::initialize.

Reimplemented from ReactorBase.

Definition at line 155 of file FlowReactor.cpp.

◆ syncState()

void syncState ( )
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.

◆ updateState()

void updateState ( double *  y)
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.

◆ eval()

void eval ( double  t,
double *  LHS,
double *  RHS 
)
inlineoverridevirtual

Not implemented; FlowReactor implements evalDae() instead.

Reimplemented from Reactor.

Definition at line 44 of file FlowReactor.h.

◆ evalDae()

void evalDae ( double  t,
double *  y,
double *  ydot,
double *  residual 
)
overridevirtual

Evaluate the reactor governing equations.

Called by ReactorNet::eval.

Parameters
[in]ttime.
[in]ysolution vector, length neq()
[in]ydotrate of change of solution vector, length neq()
[out]residualresiduals vector, length neq()

Reimplemented from Reactor.

Definition at line 246 of file FlowReactor.cpp.

◆ getConstraints()

void getConstraints ( double *  constraints)
overridevirtual

Given a vector of length neq(), mark which variables should be considered algebraic constraints.

Reimplemented from Reactor.

Definition at line 333 of file FlowReactor.cpp.

◆ setMassFlowRate()

void setMassFlowRate ( double  mdot)

Set the mass flow rate through the reactor [kg/s].

Definition at line 210 of file FlowReactor.cpp.

◆ speed()

double speed ( ) const
inline

The current gas speed in the reactor [m/s].

Definition at line 56 of file FlowReactor.h.

◆ area()

double area ( ) const
inline

The cross-sectional area of the reactor [m^2].

Definition at line 61 of file FlowReactor.h.

◆ setArea()

void setArea ( double  area)

Sets the area of the reactor [m^2].

Definition at line 216 of file FlowReactor.cpp.

◆ surfaceAreaToVolumeRatio()

double surfaceAreaToVolumeRatio ( ) const

The ratio of the reactor's surface area to volume ratio [m^-1].

Note
If the surface area to volume ratio is unspecified by the user, this will be calculated assuming the reactor is a cylinder.

Definition at line 222 of file FlowReactor.cpp.

◆ setSurfaceAreaToVolumeRatio()

void setSurfaceAreaToVolumeRatio ( double  sa_to_vol)
inline

Set the reactor's surface area to volume ratio [m^-1].

Definition at line 74 of file FlowReactor.h.

◆ inletSurfaceAtol()

double inletSurfaceAtol ( ) const
inline

Get the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 80 of file FlowReactor.h.

◆ setInletSurfaceAtol()

void setInletSurfaceAtol ( double  atol)
inline

Set the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 86 of file FlowReactor.h.

◆ inletSurfaceRtol()

double inletSurfaceRtol ( ) const
inline

Get the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 92 of file FlowReactor.h.

◆ setInletSurfaceRtol()

void setInletSurfaceRtol ( double  rtol)
inline

Set the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 98 of file FlowReactor.h.

◆ inletSurfaceMaxSteps()

double inletSurfaceMaxSteps ( ) const
inline

Get the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 104 of file FlowReactor.h.

◆ setInletSurfaceMaxSteps()

void setInletSurfaceMaxSteps ( int  max_steps)
inline

Set the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 110 of file FlowReactor.h.

◆ inletSurfaceMaxErrorFailures()

double inletSurfaceMaxErrorFailures ( ) const
inline

Get the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 116 of file FlowReactor.h.

◆ setInletSurfaceMaxErrorFailures()

void setInletSurfaceMaxErrorFailures ( int  max_fails)
inline

Set the steady state tolerances used to determine the initial state for surface coverages.

Definition at line 122 of file FlowReactor.h.

◆ componentIndex()

size_t componentIndex ( const string &  nm) const
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 340 of file FlowReactor.cpp.

◆ componentName()

string componentName ( size_t  k)
overridevirtual

Return the name of the solution component with index i.

See also
componentIndex()

Reimplemented from Reactor.

Definition at line 358 of file FlowReactor.cpp.

◆ updateSurfaceState()

void updateSurfaceState ( double *  y)
overridevirtual

Update the state of SurfPhase objects attached to this reactor.

Reimplemented from Reactor.

Definition at line 233 of file FlowReactor.cpp.

Member Data Documentation

◆ m_rho

double m_rho = NAN
protected

Density [kg/m^3]. First component of the state vector.

Definition at line 138 of file FlowReactor.h.

◆ m_u

double m_u = -1.0
protected

Axial velocity [m/s]. Second component of the state vector.

Definition at line 140 of file FlowReactor.h.

◆ m_P

double m_P = NAN
protected

Pressure [Pa]. Third component of the state vector.

Definition at line 142 of file FlowReactor.h.

◆ m_T

double m_T = NAN
protected

Temperature [K]. Fourth component of the state vector.

Definition at line 144 of file FlowReactor.h.

◆ m_offset_Y

const size_t m_offset_Y = 4
protected

offset to the species equations

Definition at line 146 of file FlowReactor.h.

◆ m_area

double m_area = 1.0
protected

reactor area [m^2]

Definition at line 148 of file FlowReactor.h.

◆ m_sa_to_vol

double m_sa_to_vol = -1.0
protected

reactor surface area to volume ratio [m^-1]

Definition at line 150 of file FlowReactor.h.

◆ m_sdot_temp

vector<double> m_sdot_temp
protected

temporary storage for surface species production rates

Definition at line 152 of file FlowReactor.h.

◆ m_hk

vector<double> m_hk
protected

temporary storage for species partial molar enthalpies

Definition at line 154 of file FlowReactor.h.

◆ m_ss_rtol

double m_ss_rtol = 1e-7
protected

steady-state relative tolerance, used to determine initial surface coverages

Definition at line 156 of file FlowReactor.h.

◆ m_ss_atol

double m_ss_atol = 1e-14
protected

steady-state absolute tolerance, used to determine initial surface coverages

Definition at line 158 of file FlowReactor.h.

◆ m_max_ss_steps

int m_max_ss_steps = 20000
protected

maximum number of steady-state coverage integrator-steps

Definition at line 160 of file FlowReactor.h.

◆ m_max_ss_error_fails

int m_max_ss_error_fails = 10
protected

maximum number of steady-state integrator error test failures

Definition at line 162 of file FlowReactor.h.


The documentation for this class was generated from the following files: