Cantera  4.0.0a1
Loading...
Searching...
No Matches
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

 FlowReactor (shared_ptr< Solution > sol, const string &name="(none)")
 
 FlowReactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 
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 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.
 
vector< size_t > steadyConstraints () const override
 Get the indices of equations that are algebraic constraints when solving the steady-state problem.
 
double massFlowRate ()
 Mass flow rate through the reactor [kg/s].
 
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 override
 The cross-sectional area of the reactor [m²].
 
void setArea (double area) override
 Sets the area of the reactor [m²].
 
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 IdealGasReactor
string type () const override
 String indicating the reactor model implemented.
 
void getState (double *y) override
 Get 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.
 
vector< size_t > steadyConstraints () const override
 Get the indices of equations that are algebraic constraints when solving the steady-state problem.
 
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.
 
double upperBound (size_t k) const override
 Get the upper bound on the k-th component of the local state vector.
 
double lowerBound (size_t k) const override
 Get the lower bound on the k-th component of the local state vector.
 
 Reactor (shared_ptr< Solution > sol, const string &name="(none)")
 
 Reactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 
- Public Member Functions inherited from Reactor
 Reactor (shared_ptr< Solution > sol, const string &name="(none)")
 
 Reactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 
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.
 
void setInitialVolume (double vol) override
 Set the initial reactor volume.
 
void setChemistryEnabled (bool cflag=true) override
 Enable or disable changes in reactor composition due to chemical reactions.
 
bool chemistryEnabled () const override
 Returns true if changes in the reactor composition due to chemical reactions are enabled.
 
void setEnergyEnabled (bool eflag=true) override
 Set the energy equation on or off.
 
bool energyEnabled () const override
 Returns true if solution of the energy equation is enabled.
 
void getState (double *y) override
 Get 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.
 
vector< size_t > steadyConstraints () const override
 Get the indices of equations that are algebraic constraints when solving the steady-state problem.
 
void updateState (double *y) override
 Set the state of the reactor to correspond to the state vector y.
 
void addSensitivityReaction (size_t rxn) override
 Add a sensitivity parameter associated with the reaction number rxn
 
void addSensitivitySpeciesEnthalpy (size_t k) override
 Add a sensitivity parameter associated with the enthalpy formation of species k.
 
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.
 
double upperBound (size_t k) const override
 Get the upper bound on the k-th component of the local state vector.
 
double lowerBound (size_t k) const override
 Get the lower bound on the k-th component of the local state vector.
 
void resetBadValues (double *y) override
 Reset physically or mathematically problematic values, such as negative species concentrations.
 
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
 
Eigen::SparseMatrix< double > finiteDifferenceJacobian ()
 Calculate the reactor-specific Jacobian using a finite difference method.
 
void setDerivativeSettings (AnyMap &settings) override
 Use this to set the kinetics objects derivative settings.
 
void applySensitivity (double *params) override
 Set reaction rate multipliers based on the sensitivity variables in params.
 
void resetSensitivity (double *params) override
 Reset the reaction rate multipliers.
 
- Public Member Functions inherited from ReactorBase
 ReactorBase (shared_ptr< Solution > sol, const string &name="(none)")
 Instantiate a ReactorBase object with Solution contents.
 
 ReactorBase (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 Instantiate a ReactorBase object with Solution contents.
 
 ReactorBase (const ReactorBase &)=delete
 
ReactorBaseoperator= (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.
 
bool setDefaultName (map< string, int > &counts)
 Set the default name of a reactor. Returns false if it was previously set.
 
shared_ptr< Solutionphase ()
 Access the Solution object used to represent the contents of this reactor.
 
shared_ptr< const Solutionphase () const
 Access the Solution object used to represent the contents of this reactor.
 
virtual bool timeIsIndependent () const
 Indicates whether the governing equations for this reactor are functions of time or a spatial variable.
 
size_t neq ()
 Number of equations (state variables) for this reactor.
 
virtual void syncState ()
 Set the state of the reactor to the associated ThermoPhase object.
 
virtual void updateConnected (bool updatePressure)
 Update state information needed by connected reactors, flow devices, and walls.
 
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.
 
ReactorNetnetwork ()
 The ReactorNet that this reactor belongs to.
 
void setNetwork (ReactorNet *net)
 Set the ReactorNet that this reactor belongs to.
 
size_t offset () const
 Get the starting offset for this reactor's state variables within the global state vector of the ReactorNet.
 
void setOffset (size_t offset)
 Set the starting offset for this reactor's state variables within the global state vector of the ReactorNet.
 
size_t speciesOffset () const
 Offset of the first species in the local state vector.
 
virtual void getJacobianScalingFactors (double &f_species, double *f_energy)
 Get scaling factors for the Jacobian matrix terms proportional to \( d\dot{n}_k/dC_j \).
 
virtual void addSensitivityReaction (size_t rxn)
 Add a sensitivity parameter associated with the reaction number rxn
 
virtual size_t nSensParams () const
 Number of sensitivity parameters associated with this reactor.
 
virtual void addInlet (FlowDevice &inlet)
 Connect an inlet FlowDevice to this reactor.
 
virtual void addOutlet (FlowDevice &outlet)
 Connect an outlet FlowDevice to this reactor.
 
FlowDeviceinlet (size_t n=0)
 Return a reference to the n-th inlet FlowDevice connected to this reactor.
 
FlowDeviceoutlet (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.
 
virtual void addWall (WallBase &w, int lr)
 Insert a Wall between this reactor and another reactor.
 
WallBasewall (size_t n)
 Return a reference to the n-th Wall connected to this reactor.
 
virtual void addSurface (ReactorSurface *surf)
 Add a ReactorSurface object to a Reactor object.
 
ReactorSurfacesurface (size_t n)
 Return a reference to the n-th ReactorSurface connected to this reactor.
 
virtual size_t nSurfs () const
 Return the number of surfaces in a reactor.
 
span< const double > surfaceProductionRates () const
 Production rates on surfaces.
 
virtual void getJacobianElements (vector< Eigen::Triplet< double > > &trips)
 Get Jacobian elements for this reactor within the full reactor network.
 
virtual Eigen::SparseMatrix< double > jacobian ()
 Calculate the Jacobian of a specific reactor specialization.
 
virtual bool preconditionerSupported () const
 Return a false if preconditioning is not supported or true otherwise.
 
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 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.
 
const size_t m_offset_Y = 4
 offset to the species equations
 
double m_area = 1.0
 reactor area [m^2]
 
vector< double > m_hk
 temporary storage for species partial molar enthalpies
 
- Protected Attributes inherited from IdealGasReactor
vector< double > m_uk
 Species molar internal energies.
 
- Protected Attributes inherited from Reactor
Kineticsm_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]
 
vector< double > m_wdot
 Species net molar production rates.
 
vector< double > m_uk
 Species molar internal energies.
 
vector< double > m_sdot
 Total production rate of bulk phase species on surfaces [kmol/s].
 
bool m_chem = false
 
bool m_energy = true
 
vector< double > m_advancelimits
 Advance step limit.
 
- Protected Attributes inherited from ReactorBase
size_t m_nsp = 0
 Number of homogeneous species in the mixture.
 
ThermoPhasem_thermo = nullptr
 
double m_vol = 0.0
 Current volume of the reactor [m^3].
 
double m_mass = 0.0
 Current mass of the reactor [kg].
 
double m_enthalpy = 0.0
 Current specific enthalpy 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< double > m_sdot
 species production rates on 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
 Reactor name.
 
bool m_defaultNameSet = false
 true if default name has been previously set.
 
size_t m_nv = 0
 Number of state variables for this reactor.
 
ReactorNetm_net = nullptr
 The ReactorNet that this reactor is part of.
 
size_t m_offset = 0
 Offset into global ReactorNet state vector.
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 
vector< SensitivityParameterm_sensParams
 

Additional Inherited Members

- Protected Member Functions inherited from Reactor
void updateSurfaceProductionRates ()
 Update m_sdot to reflect current production rates of bulk phase species due to reactions on adjacent surfaces.
 
void evalWalls (double t) override
 Evaluate terms related to Walls.
 
- Protected Member Functions inherited from ReactorBase
 ReactorBase (const string &name="(none)")
 

Constructor & Destructor Documentation

◆ FlowReactor() [1/2]

FlowReactor ( shared_ptr< Solution sol,
const string &  name = "(none)" 
)

Definition at line 20 of file FlowReactor.cpp.

◆ FlowReactor() [2/2]

FlowReactor ( shared_ptr< Solution sol,
bool  clone,
const string &  name = "(none)" 
)

Definition at line 25 of file FlowReactor.cpp.

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 22 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 26 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 ReactorBase.

Definition at line 30 of file FlowReactor.h.

◆ getState()

void getState ( double *  y)
inlineoverridevirtual

Not implemented; FlowReactor implements getStateDae() instead.

Reimplemented from ReactorBase.

Definition at line 35 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 ReactorBase.

Definition at line 35 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 ReactorBase.

Definition at line 150 of file FlowReactor.cpp.

◆ eval()

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

Not implemented; FlowReactor implements evalDae() instead.

Reimplemented from ReactorBase.

Definition at line 43 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 ReactorBase.

Definition at line 171 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 ReactorBase.

Definition at line 239 of file FlowReactor.cpp.

◆ steadyConstraints()

vector< size_t > steadyConstraints ( ) const
inlineoverridevirtual

Get the indices of equations that are algebraic constraints when solving the steady-state problem.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.2.

Reimplemented from ReactorBase.

Definition at line 50 of file FlowReactor.h.

◆ massFlowRate()

double massFlowRate ( )
inline

Mass flow rate through the reactor [kg/s].

Definition at line 56 of file FlowReactor.h.

◆ setMassFlowRate()

void setMassFlowRate ( double  mdot)

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

Definition at line 159 of file FlowReactor.cpp.

◆ speed()

double speed ( ) const
inline

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

Definition at line 64 of file FlowReactor.h.

◆ area()

double area ( ) const
inlineoverridevirtual

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

Reimplemented from ReactorBase.

Definition at line 69 of file FlowReactor.h.

◆ setArea()

void setArea ( double  area)
overridevirtual

Sets the area of the reactor [m²].

Reimplemented from ReactorBase.

Definition at line 165 of file FlowReactor.cpp.

◆ 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" or the name of a homogeneous phase species.

Reimplemented from ReactorBase.

Definition at line 244 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 ReactorBase.

Definition at line 263 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 85 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 87 of file FlowReactor.h.

◆ m_offset_Y

const size_t m_offset_Y = 4
protected

offset to the species equations

Definition at line 89 of file FlowReactor.h.

◆ m_area

double m_area = 1.0
protected

reactor area [m^2]

Definition at line 91 of file FlowReactor.h.

◆ m_hk

vector<double> m_hk
protected

temporary storage for species partial molar enthalpies

Definition at line 93 of file FlowReactor.h.


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