10#include "cantera/numerics/eigen_sparse.h"
49 Reactor(shared_ptr<Solution> sol,
const string&
name=
"(none)");
50 using ReactorBase::ReactorBase;
52 string type()
const override {
76 warn_deprecated(
"Reactor::insert",
"Unused; to be removed after Cantera 3.1.");
138 virtual void eval(
double t,
double* LHS,
double* RHS);
147 virtual void evalDae(
double t,
double* y,
double* ydot,
double* residual) {
264 virtual void evalSurfaces(
double* LHS,
double* RHS,
double* sdot);
290 vector<double> m_work;
298 bool m_energy =
true;
305 vector<SensitivityParameter> m_sensParams;
A map of string keys to values whose type can vary at runtime.
Public interface for kinetics managers.
An error indicating that an unimplemented function has been called.
Base class for stirred reactors.
virtual void setThermo(ThermoPhase &thermo)
Specify the mixture contained in the reactor.
void insert(shared_ptr< Solution > sol)
ThermoPhase & contents()
return a reference to the contents.
string name() const
Return the name of this reactor.
Class Reactor is a general-purpose class for stirred reactors.
virtual string componentName(size_t k)
Return the name of the solution component with index i.
bool chemistryEnabled() const
Returns true if changes in the reactor composition due to chemical reactions are enabled.
virtual void evalSurfaces(double *LHS, double *RHS, double *sdot)
Evaluate terms related to surface reactions.
virtual size_t componentIndex(const string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
virtual void getStateDae(double *y, double *ydot)
Get the current state and derivative vector of the reactor for a DAE solver.
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
void insert(G &contents)
Insert something into the reactor.
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
vector< double > m_wdot
Species net molar production rates.
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
virtual void evalWalls(double t)
Evaluate terms related to Walls.
bool energyEnabled() const
Returns true if solution of the energy equation is enabled.
Eigen::SparseMatrix< double > finiteDifferenceJacobian()
Calculate the reactor-specific Jacobian using a finite difference method.
size_t neq()
Number of equations (state variables) for this reactor.
virtual void addSensitivitySpeciesEnthalpy(size_t k)
Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous p...
string type() const override
String indicating the reactor model implemented.
double m_Qdot
net heat transfer into the reactor, through walls [W]
void setKinetics(Kinetics &kin) override
Specify the kinetics manager for the reactor.
vector< double > m_advancelimits
!< Number of variables associated with reactor surfaces
void setEnergy(int eflag=1) override
Set the energy equation on or off.
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase).
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor (including walls)
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
vector< double > m_uk
Species molar internal energies.
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
virtual void evalDae(double t, double *y, double *ydot, double *residual)
Evaluate the reactor governing equations.
void setChemistry(bool cflag=true) override
Enable or disable changes in reactor composition due to chemical reactions.
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
void setAdvanceLimit(const string &nm, const double limit)
Set individual step size limit for component name nm
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
virtual void getState(double *y)
Get the the current state of the reactor.
bool hasAdvanceLimits() const
Check whether Reactor object uses advance limits.
double m_vdot
net rate of volume change from moving walls [m^3/s]
void syncState() override
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
void initialize(double t0=0.0) override
Initialize the reactor.
virtual void getConstraints(double *constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
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 homog...
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 variabl...
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
virtual Eigen::SparseMatrix< double > jacobian()
Calculate the Jacobian of a specific Reactor specialization.
Namespace for the Cantera kernel.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.