23ReactorBase::ReactorBase(shared_ptr<Solution> sol,
const string& name)
31 if (!sol || !(sol->thermo())) {
33 "Missing or incomplete Solution object.");
47ReactorBase::~ReactorBase()
71 "Cannot add an inlet to reactor '{}' after it has been"
72 " added to a ReactorNet.",
m_name);
74 m_inlet.push_back(&
inlet);
81 "Cannot add an outlet to reactor '{}' after it has been"
82 " added to a ReactorNet.",
m_name);
84 m_outlet.push_back(&
outlet);
91 "Cannot add a wall to reactor '{}' after it has been"
92 " added to a ReactorNet.",
m_name);
111 "Cannot add a surface to reactor '{}' after it has been"
112 " added to a ReactorNet.",
m_name);
114 if (find(m_surfaces.begin(), m_surfaces.end(), surf) == m_surfaces.end()) {
115 m_surfaces.push_back(surf);
121 return m_surfaces[n];
126 vector<Eigen::Triplet<double>> trips;
128 for (
auto& trip : trips) {
129 trip = Eigen::Triplet<double>(trip.row() -
m_offset, trip.col() -
m_offset, trip.value());
131 Eigen::SparseMatrix<double> J(
m_nv,
m_nv);
132 J.setFromTriplets(trips.begin(), trips.end());
139 "To be removed after Cantera 4.0. Use ReactorNet::reinitialize to indicate "
140 "a change in state that requires integrator reinitialization.");
149 if (updatePressure) {
155 if (
m_net !=
nullptr) {
158 for (
size_t i = 0; i < m_outlet.size(); i++) {
159 m_outlet[i]->setSimTime(time);
160 m_outlet[i]->updateMassFlowRate(time);
162 for (
size_t i = 0; i < m_inlet.size(); i++) {
163 m_inlet[i]->setSimTime(time);
164 m_inlet[i]->updateMassFlowRate(time);
166 for (
size_t i = 0; i < m_wall.size(); i++) {
167 m_wall[i]->setSimTime(time);
177 "Reactor is not part of a ReactorNet");
185 "Reactor {} is already part of a ReactorNet",
m_name);
193 for (
size_t i = 0; i < m_outlet.size(); i++) {
194 mout += m_outlet[i]->massFlowRate();
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
Header file for class ReactorSurface.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
Header file for base class WallBase.
Base class for exceptions thrown by Cantera classes.
Base class for 'flow devices' (valves, pressure regulators, etc.) connecting reactors.
double massFraction(size_t k) const
Return the mass fraction of a single species.
size_t nSpecies() const
Returns the number of species in the phase.
double temperature() const
Temperature (K).
const double * massFractions() const
Return a const pointer to the mass fraction array.
virtual double density() const
Density (kg/m^3).
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Base class for reactor objects.
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
double massFraction(size_t k) const
Return the mass fraction of the k-th species.
shared_ptr< Solution > m_solution
Composite thermo/kinetics/transport handler.
bool m_defaultNameSet
true if default name has been previously set.
WallBase & wall(size_t n)
Return a reference to the n-th Wall connected to this reactor.
double density() const
Returns the current density (kg/m^3) of the reactor's contents.
virtual void addOutlet(FlowDevice &outlet)
Connect an outlet FlowDevice to this reactor.
double m_pressure
Current pressure in the reactor [Pa].
virtual string type() const
String indicating the reactor model implemented.
ReactorNet * m_net
The ReactorNet that this reactor is part of.
virtual void addWall(WallBase &w, int lr)
Insert a Wall between this reactor and another reactor.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
size_t m_nv
Number of state variables for this reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
double temperature() const
Returns the current temperature (K) of the reactor's contents.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
size_t m_offset
Offset into global ReactorNet state vector.
virtual void addInlet(FlowDevice &inlet)
Connect an inlet FlowDevice to this reactor.
virtual void syncState()
Set the state of the reactor to the associated ThermoPhase object.
const double * massFractions() const
Return the vector of species mass fractions.
size_t m_nsp
Number of homogeneous species in the mixture.
string m_name
Reactor name.
double mass() const
Returns the mass (kg) of the reactor's contents.
ReactorBase(shared_ptr< Solution > sol, const string &name="(none)")
Instantiate a ReactorBase object with Solution contents.
ReactorNet & network()
The ReactorNet that this reactor belongs to.
double residenceTime()
Return the residence time (s) of the contents of this reactor, based on the outlet mass flow rates an...
ReactorSurface * surface(size_t n)
Return a reference to the n-th ReactorSurface connected to this reactor.
double m_enthalpy
Current specific enthalpy of the reactor [J/kg].
virtual void getJacobianElements(vector< Eigen::Triplet< double > > &trips)
Get Jacobian elements for this reactor within the full reactor network.
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 state information needed by connected reactors, flow devices, and walls.
virtual void addSurface(ReactorSurface *surf)
Add a ReactorSurface object to a Reactor object.
A class representing a network of connected reactors.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
A surface where reactions can occur that is in contact with the bulk fluid of a Reactor.
double enthalpy_mass() const
Specific enthalpy. Units: J/kg.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
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.