Cantera  3.2.0a4
Loading...
Searching...
No Matches

Class Reactor is a general-purpose class for stirred reactors. More...

#include <Reactor.h>

Inheritance diagram for Reactor:
[legend]

Detailed Description

Class Reactor is a general-purpose class for stirred reactors.

The reactor may have an arbitrary number of inlets and outlets, each of which may be connected to a "flow device" such as a mass flow controller, a pressure regulator, etc. Additional reactors may be connected to the other end of the flow device, allowing construction of arbitrary reactor networks.

The reactor class integrates the same governing equations no matter what type of reactor is simulated. The differences among reactor types are completely specified by the attached flow devices and the time-dependent user-specified boundary conditions.

If an instance of class Reactor is used directly, it will simulate an adiabatic, constant volume reactor with gas-phase chemistry but no surface chemistry. Other reactor types may be simulated by deriving a class from Reactor. This method allows specifying the following in terms of the instantaneous reactor state:

  • rate of change of the total volume (m^3/s)
  • surface heat loss rate (W)
  • species surface production rates (kmol/s)

See the Science Reference for the governing equations of class Reactor.

Definition at line 46 of file Reactor.h.

Public Member Functions

 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.
 
virtual bool timeIsIndependent () const
 Indicates whether the governing equations for this reactor are functions of time or a spatial variable.
 
void setInitialVolume (double vol) override
 Set the initial reactor volume.
 
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.
 
virtual vector< size_t > steadyConstraints () const
 Get the indices of equations that are algebraic constraints when solving the steady-state problem.
 
virtual void updateState (double *y)
 Set the state of the reactor to correspond to the state vector y.
 
size_t nSensParams () const override
 Number of sensitivity parameters associated with this reactor.
 
void addSensitivityReaction (size_t rxn) override
 Add a sensitivity parameter associated with the reaction number rxn
 
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.
 
virtual double upperBound (size_t k) const
 Get the upper bound on the k-th component of the local state vector.
 
virtual double lowerBound (size_t k) const
 Get the lower bound on the k-th component of the local state vector.
 
virtual void resetBadValues (double *y)
 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
 
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.
 
 ReactorBase (const string &name="(none)")
 TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.
 
 ReactorBase (shared_ptr< Solution > sol, const string &name="(none)")
 TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.
 
 ReactorBase (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.
 
 ReactorBase (const ReactorBase &)=delete
 TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.
 
- Public Member Functions inherited from ReactorBase
 ReactorBase (const string &name="(none)")
 
 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.
 
void setSolution (shared_ptr< Solution > sol)
 Set the Solution specifying the ReactorBase content.
 
shared_ptr< Solutioncontents4 ()
 Access the Solution object used to represent the contents of this reactor.
 
shared_ptr< const Solutioncontents4 () const
 Access the Solution object used to represent the contents of this reactor.
 
virtual 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 the associated ThermoPhase object.
 
ThermoPhasecontents ()
 return a reference to the contents.
 
const ThermoPhasecontents () const
 return a reference to the contents.
 
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.
 
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.
 
void addSurface (shared_ptr< ReactorBase > 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.
 
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 Member Functions

void setKinetics (Kinetics &kin) override
 
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.
 
- Protected Member Functions inherited from ReactorBase
virtual void setThermo (ThermoPhase &thermo)
 Specify the mixture contained in the reactor.
 
virtual void setKinetics (Kinetics &kin)
 Specify the kinetics manager for the reactor.
 

Protected Attributes

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_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< 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.
 
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_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
 Reactor name.
 
bool m_defaultNameSet = false
 true if default name has been previously set.
 
ReactorNetm_net = nullptr
 The ReactorNet that this reactor is part of.
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 
vector< SensitivityParameterm_sensParams
 

Constructor & Destructor Documentation

◆ Reactor() [1/2]

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

Definition at line 27 of file Reactor.cpp.

◆ Reactor() [2/2]

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

Definition at line 39 of file Reactor.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 55 of file Reactor.h.

◆ isOde()

virtual bool isOde ( ) const
inlinevirtual

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 in FlowReactor.

Definition at line 62 of file Reactor.h.

◆ timeIsIndependent()

virtual bool timeIsIndependent ( ) const
inlinevirtual

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 in FlowReactor.

Definition at line 68 of file Reactor.h.

◆ setInitialVolume()

void setInitialVolume ( double  vol)
inlineoverridevirtual

Set the initial reactor volume.

Reimplemented from ReactorBase.

Definition at line 72 of file Reactor.h.

◆ setChemistry()

void setChemistry ( bool  cflag = true)
inlineoverridevirtual

Enable or disable changes in reactor composition due to chemical reactions.

Reimplemented from ReactorBase.

Definition at line 76 of file Reactor.h.

◆ chemistryEnabled()

bool chemistryEnabled ( ) const
inline

Returns true if changes in the reactor composition due to chemical reactions are enabled.

Definition at line 82 of file Reactor.h.

◆ setEnergy()

void setEnergy ( int  eflag = 1)
inlineoverridevirtual

Set the energy equation on or off.

Reimplemented from ReactorBase.

Definition at line 86 of file Reactor.h.

◆ energyEnabled()

bool energyEnabled ( ) const
inline

Returns true if solution of the energy equation is enabled.

Definition at line 95 of file Reactor.h.

◆ neq()

size_t neq ( )
inline

Number of equations (state variables) for this reactor.

Definition at line 100 of file Reactor.h.

◆ getState()

void getState ( double *  y)
virtual

Get the the current state of the reactor.

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

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 73 of file Reactor.cpp.

◆ getStateDae()

virtual void getStateDae ( double *  y,
double *  ydot 
)
inlinevirtual

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 in FlowReactor.

Definition at line 119 of file Reactor.h.

◆ initialize()

void initialize ( double  t0 = 0.0)
overridevirtual

Initialize the reactor.

Called automatically by ReactorNet::initialize.

Reimplemented from ReactorBase.

Definition at line 108 of file Reactor.cpp.

◆ 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 in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 233 of file Reactor.cpp.

◆ evalDae()

virtual void evalDae ( double  t,
double *  y,
double *  ydot,
double *  residual 
)
inlinevirtual

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 in FlowReactor.

Definition at line 140 of file Reactor.h.

◆ getConstraints()

virtual void getConstraints ( double *  constraints)
inlinevirtual

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

Reimplemented in FlowReactor.

Definition at line 146 of file Reactor.h.

◆ steadyConstraints()

vector< size_t > steadyConstraints ( ) const
virtual

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 in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, and IdealGasReactor.

Definition at line 338 of file Reactor.cpp.

◆ updateState()

void updateState ( double *  y)
virtual

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

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 145 of file Reactor.cpp.

◆ nSensParams()

size_t nSensParams ( ) const
overridevirtual

Number of sensitivity parameters associated with this reactor.

(including walls)

Reimplemented from ReactorBase.

Definition at line 136 of file Reactor.cpp.

◆ addSensitivityReaction()

void addSensitivityReaction ( size_t  rxn)
overridevirtual

Add a sensitivity parameter associated with the reaction number rxn

Reimplemented from ReactorBase.

Definition at line 436 of file Reactor.cpp.

◆ addSensitivitySpeciesEnthalpy()

void addSensitivitySpeciesEnthalpy ( size_t  k)
virtual

Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous phase)

Definition at line 449 of file Reactor.cpp.

◆ componentIndex()

size_t componentIndex ( const 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", "int_energy", the name of a homogeneous phase species, or the name of a surface species.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 486 of file Reactor.cpp.

◆ componentName()

string componentName ( size_t  k)
virtual

Return the name of the solution component with index i.

See also
componentIndex()

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 502 of file Reactor.cpp.

◆ upperBound()

double upperBound ( size_t  k) const
virtual

Get the upper bound on the k-th component of the local state vector.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 528 of file Reactor.cpp.

◆ lowerBound()

double lowerBound ( size_t  k) const
virtual

Get the lower bound on the k-th component of the local state vector.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, and MoleReactor.

Definition at line 542 of file Reactor.cpp.

◆ resetBadValues()

void resetBadValues ( double *  y)
virtual

Reset physically or mathematically problematic values, such as negative species concentrations.

Parameters
[in,out]ycurrent state vector, to be updated; length neq()

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, and MoleReactor.

Definition at line 556 of file Reactor.cpp.

◆ setAdvanceLimits()

void setAdvanceLimits ( const double *  limits)

Set absolute step size limits during advance.

Parameters
limitsarray of step size limits with length neq

Definition at line 605 of file Reactor.cpp.

◆ hasAdvanceLimits()

bool hasAdvanceLimits ( ) const
inline

Check whether Reactor object uses advance limits.

Returns
True if at least one limit is set, False otherwise

Definition at line 198 of file Reactor.h.

◆ getAdvanceLimits()

bool getAdvanceLimits ( double *  limits) const

Retrieve absolute step size limits during advance.

Parameters
[out]limitsarray of step size limits with length neq
Returns
True if at least one limit is set, False otherwise

Definition at line 620 of file Reactor.cpp.

◆ setAdvanceLimit()

void setAdvanceLimit ( const string &  nm,
const double  limit 
)

Set individual step size limit for component name nm

Parameters
nmcomponent name
limitvalue for step size limit

Definition at line 631 of file Reactor.cpp.

◆ finiteDifferenceJacobian()

Eigen::SparseMatrix< double > finiteDifferenceJacobian ( )

Calculate the reactor-specific Jacobian using a finite difference method.

This method is used only for informational purposes. Jacobian calculations for the full reactor system are handled internally by CVODES.

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

Definition at line 354 of file Reactor.cpp.

◆ setDerivativeSettings()

void setDerivativeSettings ( AnyMap settings)
virtual

Use this to set the kinetics objects derivative settings.

Definition at line 51 of file Reactor.cpp.

◆ applySensitivity()

void applySensitivity ( double *  params)
virtual

Set reaction rate multipliers based on the sensitivity variables in params.

Definition at line 562 of file Reactor.cpp.

◆ resetSensitivity()

void resetSensitivity ( double *  params)
virtual

Reset the reaction rate multipliers.

Definition at line 584 of file Reactor.cpp.

◆ preconditionerSupported()

virtual bool preconditionerSupported ( ) const
inlinevirtual

Return a false if preconditioning is not supported or true otherwise.

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

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 250 of file Reactor.h.

◆ setKinetics()

void setKinetics ( Kinetics kin)
overrideprotectedvirtual
Deprecated:
To be removed after Cantera 3.2. Use constructor with Solution object instead.

Reimplemented from ReactorBase.

Definition at line 60 of file Reactor.cpp.

◆ speciesIndex()

size_t speciesIndex ( const string &  nm) const
protectedvirtual

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.

Used to implement componentIndex for specific reactor implementations.

Definition at line 464 of file Reactor.cpp.

◆ evalWalls()

void evalWalls ( double  t)
protectedvirtual

Evaluate terms related to Walls.

Calculates m_vdot and m_Qdot based on wall movement and heat transfer.

Parameters
tthe current time

Definition at line 297 of file Reactor.cpp.

◆ evalSurfaces() [1/2]

void evalSurfaces ( double *  LHS,
double *  RHS,
double *  sdot 
)
protectedvirtual

Evaluate terms related to surface reactions.

Parameters
[out]LHSMultiplicative factor on the left hand side of ODE for surface species coverages
[out]RHSRight hand side of ODE for surface species coverages
[out]sdotarray of production rates of bulk phase species on surfaces [kmol/s]

Reimplemented in MoleReactor.

Definition at line 309 of file Reactor.cpp.

◆ evalSurfaces() [2/2]

void evalSurfaces ( double *  RHS,
double *  sdot 
)
protectedvirtual

Definition at line 407 of file Reactor.cpp.

◆ updateSurfaceState()

void updateSurfaceState ( double *  y)
protectedvirtual

Update the state of SurfPhase objects attached to this reactor.

Reimplemented in FlowReactor, and MoleReactor.

Definition at line 193 of file Reactor.cpp.

◆ updateConnected()

void updateConnected ( bool  updatePressure)
protectedvirtual

Update the state information needed by connected reactors, flow devices, and reactor walls.

Called from updateState().

Parameters
updatePressureIndicates whether to update m_pressure. Should true for reactors where the pressure is a dependent property, calculated from the state, and false when the pressure is constant or an independent variable.

Definition at line 202 of file Reactor.cpp.

◆ getSurfaceInitialConditions()

void getSurfaceInitialConditions ( double *  y)
protectedvirtual

Get initial conditions for SurfPhase objects attached to this reactor.

Reimplemented in MoleReactor.

Definition at line 99 of file Reactor.cpp.

◆ ReactorBase() [1/4]

ReactorBase ( const string &  name = "(none)")
explicit

TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.

Definition at line 52 of file ReactorBase.cpp.

◆ ReactorBase() [2/4]

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

TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.

Definition at line 58 of file ReactorBase.cpp.

◆ ReactorBase() [3/4]

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

TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.

Definition at line 69 of file ReactorBase.cpp.

◆ ReactorBase() [4/4]

ReactorBase ( const ReactorBase )
delete

TODO: Remove after Cantera 3.2 – all derived classes should use Solution-based constructors.

Member Data Documentation

◆ m_kin

Kinetics* m_kin = nullptr
protected

Pointer to the homogeneous Kinetics object that handles the reactions.

Definition at line 293 of file Reactor.h.

◆ m_vdot

double m_vdot = 0.0
protected

net rate of volume change from moving walls [m^3/s]

Definition at line 295 of file Reactor.h.

◆ m_Qdot

double m_Qdot = 0.0
protected

net heat transfer into the reactor, through walls [W]

Definition at line 297 of file Reactor.h.

◆ m_work

vector<double> m_work
protected

Definition at line 299 of file Reactor.h.

◆ m_sdot

vector<double> m_sdot
protected

Production rates of gas phase species on surfaces [kmol/s].

Definition at line 302 of file Reactor.h.

◆ m_wdot

vector<double> m_wdot
protected

Species net molar production rates.

Definition at line 304 of file Reactor.h.

◆ m_uk

vector<double> m_uk
protected

Species molar internal energies.

Definition at line 305 of file Reactor.h.

◆ m_chem

bool m_chem = false
protected

Definition at line 306 of file Reactor.h.

◆ m_energy

bool m_energy = true
protected

Definition at line 307 of file Reactor.h.

◆ m_nv

size_t m_nv = 0
protected

Definition at line 308 of file Reactor.h.

◆ m_nv_surf

size_t m_nv_surf
protected

Definition at line 309 of file Reactor.h.

◆ m_advancelimits

vector<double> m_advancelimits
protected

!< Number of variables associated with reactor surfaces

Advance step limit

Definition at line 311 of file Reactor.h.

◆ m_jac_trips

vector<Eigen::Triplet<double> > m_jac_trips
protected

Vector of triplets representing the jacobian.

Definition at line 314 of file Reactor.h.


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