Cantera 2.6.0
Public Member Functions | Protected Attributes | List of all members
IdealGasReactor Class Reference

Class IdealGasReactor is a class for stirred reactors that is specifically optimized for ideal gases. More...

#include <IdealGasReactor.h>

Inheritance diagram for IdealGasReactor:
[legend]
Collaboration diagram for IdealGasReactor:
[legend]

Public Member Functions

virtual std::string typeStr () const
 String indicating the reactor model implemented. More...
 
virtual std::string type () const
 String indicating the reactor model implemented. More...
 
virtual void setThermoMgr (ThermoPhase &thermo)
 Specify the mixture contained in the reactor. More...
 
virtual void getState (doublereal *y)
 Get the the current state of the reactor. More...
 
virtual void initialize (doublereal t0=0.0)
 Initialize the reactor. More...
 
virtual void eval (double t, double *LHS, double *RHS)
 Evaluate the reactor governing equations. More...
 
virtual void updateState (doublereal *y)
 Set the state of the reactor to correspond to the state vector y. More...
 
virtual size_t componentIndex (const std::string &nm) const
 Return the index in the solution vector for this reactor of the component named nm. More...
 
std::string componentName (size_t k)
 Return the name of the solution component with index i. More...
 
- Public Member Functions inherited from Reactor
template<class G >
void insert (G &contents)
 Insert something into the reactor. More...
 
void insert (shared_ptr< Solution > sol)
 
virtual void setKineticsMgr (Kinetics &kin)
 Specify chemical kinetics governing the reactor. More...
 
void setChemistry (bool cflag=true)
 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)
 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...
 
virtual void evalEqs (double t, double *y, double *ydot, double *params)
 
virtual void syncState ()
 Set the state of the reactor to correspond to the state of the associated ThermoPhase object. More...
 
virtual size_t nSensParams ()
 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 ()
 Check whether Reactor object uses advance limits. More...
 
bool getAdvanceLimits (double *limits)
 Retrieve absolute step size limits during advance. More...
 
void setAdvanceLimit (const std::string &nm, const double limit)
 Set individual step size limit for component name nm 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...
 
- Public Member Functions inherited from ReactorBase
 ReactorBase (const std::string &name="(none)")
 
 ReactorBase (const ReactorBase &)=delete
 
ReactorBaseoperator= (const ReactorBase &)=delete
 
std::string name () const
 Return the name of this reactor. More...
 
void setName (const std::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
 
doublereal 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 (doublereal 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...
 
void addSurface (ReactorSurface *surf)
 
ReactorSurfacesurface (size_t n)
 Return a reference to the n-th ReactorSurface connected to this reactor. More...
 
doublereal volume () const
 Returns the current volume (m^3) of the reactor. More...
 
doublereal density () const
 Returns the current density (kg/m^3) of the reactor's contents. More...
 
doublereal temperature () const
 Returns the current temperature (K) of the reactor's contents. More...
 
doublereal enthalpy_mass () const
 Returns the current enthalpy (J/kg) of the reactor's contents. More...
 
doublereal intEnergy_mass () const
 Returns the current internal energy (J/kg) of the reactor's contents. More...
 
doublereal pressure () const
 Returns the current pressure (Pa) of the reactor. More...
 
doublereal mass () const
 Returns the mass (kg) of the reactor's contents. More...
 
const doublereal * massFractions () const
 Return the vector of species mass fractions. More...
 
doublereal massFraction (size_t k) const
 Return the mass fraction of the k-th species. More...
 

Protected Attributes

vector_fp m_uk
 Species molar internal energies. More...
 
- Protected Attributes inherited from Reactor
Kineticsm_kin
 Pointer to the homogeneous Kinetics object that handles the reactions. More...
 
doublereal m_vdot
 net rate of volume change from moving walls [m^3/s] More...
 
double m_Q
 net heat transfer out of the reactor, through walls [W] More...
 
double m_Qdot
 net heat transfer into the reactor, through walls [W] More...
 
doublereal m_mass
 total mass More...
 
vector_fp m_work
 
vector_fp m_sdot
 Production rates of gas phase species on surfaces [kmol/s]. More...
 
vector_fp m_wdot
 Species net molar production rates. More...
 
vector_fp m_uk
 Species molar internal energies. More...
 
bool m_chem
 
bool m_energy
 
size_t m_nv
 
size_t m_nv_surf
 
vector_fp m_advancelimits
 !< Number of variables associated with reactor surfaces More...
 
std::vector< SensitivityParameter > m_sensParams
 
- Protected Attributes inherited from ReactorBase
size_t m_nsp
 Number of homogeneous species in the mixture. More...
 
ThermoPhasem_thermo
 
double m_vol
 Current volume of the reactor [m^3]. More...
 
double m_enthalpy
 Current specific enthalpy of the reactor [J/kg]. More...
 
double m_intEnergy
 Current internal energy of the reactor [J/kg]. More...
 
double m_pressure
 Current pressure in the reactor [Pa]. More...
 
vector_fp m_state
 
std::vector< FlowDevice * > m_inlet
 
std::vector< FlowDevice * > m_outlet
 
std::vector< WallBase * > m_wall
 
std::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...
 
std::string m_name
 
ReactorNetm_net
 The ReactorNet that this reactor is part of. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Reactor
virtual size_t speciesIndex (const std::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 updateSurfaceState (double *y)
 Update the state of SurfPhase objects attached to this reactor. More...
 
virtual void updateConnected (bool updatePressure)
 Update the state information needed by connected reactors and flow devices. More...
 
virtual void getSurfaceInitialConditions (double *y)
 Get initial conditions for SurfPhase objects attached to this reactor. More...
 

Detailed Description

Class IdealGasReactor is a class for stirred reactors that is specifically optimized for ideal gases.

In this formulation, temperature replaces the total internal energy as a state variable.

Definition at line 19 of file IdealGasReactor.h.

Constructor & Destructor Documentation

◆ IdealGasReactor()

IdealGasReactor ( )
inline

Definition at line 22 of file IdealGasReactor.h.

Member Function Documentation

◆ typeStr()

virtual std::string typeStr ( ) const
inlinevirtual

String indicating the reactor model implemented.

Usually corresponds to the name of the derived class.

Deprecated:
To be removed after Cantera 2.6. Use type() instead.

Reimplemented from Reactor.

Definition at line 24 of file IdealGasReactor.h.

References Cantera::warn_deprecated().

◆ type()

virtual std::string type ( ) const
inlinevirtual

String indicating the reactor model implemented.

Usually corresponds to the name of the derived class.

Reimplemented from Reactor.

Definition at line 30 of file IdealGasReactor.h.

◆ setThermoMgr()

void setThermoMgr ( ThermoPhase thermo)
virtual

Specify the mixture contained in the reactor.

Note that a pointer to this substance is stored, and as the integration proceeds, the state of the substance is modified.

Reimplemented from ReactorBase.

Definition at line 18 of file IdealGasReactor.cpp.

References ReactorBase::setThermoMgr(), and ThermoPhase::type().

◆ getState()

void getState ( doublereal *  y)
virtual

Get the the current state of the reactor.

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

Reimplemented from Reactor.

Definition at line 29 of file IdealGasReactor.cpp.

References Phase::density(), Phase::getMassFractions(), Reactor::getSurfaceInitialConditions(), Reactor::m_mass, ReactorBase::m_nsp, ReactorBase::m_vol, Phase::restoreState(), and Phase::temperature().

◆ initialize()

void initialize ( doublereal  t0 = 0.0)
virtual

Initialize the reactor.

Called automatically by ReactorNet::initialize.

Reimplemented from Reactor.

Definition at line 55 of file IdealGasReactor.cpp.

References Reactor::initialize(), ReactorBase::m_nsp, and IdealGasReactor::m_uk.

◆ eval()

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

Evaluate the reactor governing equations.

Called by ReactorNet::eval.

Parameters
[in]ttime.
[out]LHSpointer to start of vector of left-hand side coefficients for governing equations, length m_nv, default values 1
[out]RHSpointer to start of vector of right-hand side coefficients for governing equations, length m_nv, default values 0

Reimplemented from Reactor.

Definition at line 74 of file IdealGasReactor.cpp.

References ThermoPhase::cv_mass(), Cantera::dot(), FlowDevice::enthalpy_mass(), Reactor::evalSurfaces(), Reactor::evalWalls(), Kinetics::getNetProductionRates(), ThermoPhase::getPartialMolarIntEnergies(), ReactorBase::inlet(), Reactor::m_kin, Reactor::m_mass, ReactorBase::m_nsp, ReactorBase::m_pressure, Reactor::m_Qdot, Reactor::m_sdot, IdealGasReactor::m_uk, Reactor::m_vdot, ReactorBase::m_vol, Reactor::m_wdot, FlowDevice::massFlowRate(), Phase::massFractions(), Phase::molecularWeights(), ReactorBase::outlet(), FlowDevice::outletSpeciesMassFlowRate(), and Phase::restoreState().

◆ updateState()

void updateState ( doublereal *  y)
virtual

Set the state of the reactor to correspond to the state vector y.

Reimplemented from Reactor.

Definition at line 61 of file IdealGasReactor.cpp.

References Reactor::m_mass, ReactorBase::m_nsp, ReactorBase::m_vol, Phase::setMassFractions_NoNorm(), Phase::setState_TR(), Reactor::updateConnected(), and Reactor::updateSurfaceState().

◆ componentIndex()

size_t componentIndex ( const std::string &  nm) const
virtual

Return the index in the solution vector for this reactor of the component named nm.

Possible values for nm are "mass", "volume", "temperature", the name of a homogeneous phase species, or the name of a surface species.

Reimplemented from Reactor.

Definition at line 140 of file IdealGasReactor.cpp.

References Cantera::npos, and Reactor::speciesIndex().

◆ componentName()

std::string componentName ( size_t  k)
virtual

Return the name of the solution component with index i.

See also
componentIndex()

Reimplemented from Reactor.

Definition at line 156 of file IdealGasReactor.cpp.

References Reactor::componentName().

Member Data Documentation

◆ m_uk

vector_fp m_uk
protected

Species molar internal energies.

Definition at line 52 of file IdealGasReactor.h.

Referenced by IdealGasReactor::eval(), and IdealGasReactor::initialize().


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