Cantera  4.0.0a1
Loading...
Searching...
No Matches
ReactorBase.cpp
Go to the documentation of this file.
1//! @file ReactorBase.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
10#include "cantera/zeroD/Wall.h"
15
16namespace Cantera
17{
18
19ReactorBase::ReactorBase(const string& name) : m_name(name)
20{
21}
22
23ReactorBase::ReactorBase(shared_ptr<Solution> sol, const string& name)
24 : ReactorBase(sol, true, name)
25{
26}
27
28ReactorBase::ReactorBase(shared_ptr<Solution> sol, bool clone, const string& name)
29 : ReactorBase(name)
30{
31 if (!sol || !(sol->thermo())) {
32 throw CanteraError("ReactorBase::ReactorBase",
33 "Missing or incomplete Solution object.");
34 }
35 if (clone) {
36 m_solution = sol->clone({}, true, false);
37 } else {
38 m_solution = sol;
39 }
40 m_solution->thermo()->addSpeciesLock();
41 m_thermo = m_solution->thermo().get();
42 m_nsp = m_thermo->nSpecies();
43 m_enthalpy = m_thermo->enthalpy_mass(); // Needed for flow and wall interactions
44 m_pressure = m_thermo->pressure(); // Needed for flow and wall interactions
45}
46
47ReactorBase::~ReactorBase()
48{
49 if (m_solution) {
50 m_solution->thermo()->removeSpeciesLock();
51 }
52}
53
54bool ReactorBase::setDefaultName(map<string, int>& counts)
55{
56 if (m_defaultNameSet) {
57 return false;
58 }
59 m_defaultNameSet = true;
60 if (m_name == "(none)" || m_name == "") {
61 m_name = fmt::format("{}_{}", type(), counts[type()]);
62 }
63 counts[type()]++;
64 return true;
65}
66
69 if (m_net) {
70 throw CanteraError("ReactorBase::addInlet",
71 "Cannot add an inlet to reactor '{}' after it has been"
72 " added to a ReactorNet.", m_name);
73 }
74 m_inlet.push_back(&inlet);
75}
76
78{
79 if (m_net) {
80 throw CanteraError("ReactorBase::addOutlet",
81 "Cannot add an outlet to reactor '{}' after it has been"
82 " added to a ReactorNet.", m_name);
83 }
84 m_outlet.push_back(&outlet);
85}
86
88{
89 if (m_net) {
90 throw CanteraError("ReactorBase::addWall",
91 "Cannot add a wall to reactor '{}' after it has been"
92 " added to a ReactorNet.", m_name);
93 }
94 m_wall.push_back(&w);
95 if (lr == 0) {
96 m_lr.push_back(0);
97 } else {
98 m_lr.push_back(1);
99 }
100}
101
103{
104 return *m_wall[n];
105}
106
108{
109 if (m_net) {
110 throw CanteraError("ReactorBase::addSurface",
111 "Cannot add a surface to reactor '{}' after it has been"
112 " added to a ReactorNet.", m_name);
113 }
114 if (find(m_surfaces.begin(), m_surfaces.end(), surf) == m_surfaces.end()) {
115 m_surfaces.push_back(surf);
116 }
117}
118
120{
121 return m_surfaces[n];
122}
123
124Eigen::SparseMatrix<double> ReactorBase::jacobian()
125{
126 vector<Eigen::Triplet<double>> trips;
127 getJacobianElements(trips);
128 for (auto& trip : trips) {
129 trip = Eigen::Triplet<double>(trip.row() - m_offset, trip.col() - m_offset, trip.value());
130 }
131 Eigen::SparseMatrix<double> J(m_nv, m_nv);
132 J.setFromTriplets(trips.begin(), trips.end());
133 return J;
134}
135
137{
138 warn_deprecated("ReactorBase::syncState",
139 "To be removed after Cantera 4.0. Use ReactorNet::reinitialize to indicate "
140 "a change in state that requires integrator reinitialization.");
141 if (m_net) {
143 }
144}
145
146void ReactorBase::updateConnected(bool updatePressure) {
147 // save parameters needed by other connected reactors
148 m_enthalpy = m_thermo->enthalpy_mass();
149 if (updatePressure) {
150 m_pressure = m_thermo->pressure();
151 }
152
153 // Update the mass flow rate of connected flow devices
154 double time = 0.0;
155 if (m_net != nullptr) {
156 time = (timeIsIndependent()) ? m_net->time() : m_net->distance();
157 }
158 for (size_t i = 0; i < m_outlet.size(); i++) {
159 m_outlet[i]->setSimTime(time);
160 m_outlet[i]->updateMassFlowRate(time);
161 }
162 for (size_t i = 0; i < m_inlet.size(); i++) {
163 m_inlet[i]->setSimTime(time);
164 m_inlet[i]->updateMassFlowRate(time);
165 }
166 for (size_t i = 0; i < m_wall.size(); i++) {
167 m_wall[i]->setSimTime(time);
168 }
169}
170
172{
173 if (m_net) {
174 return *m_net;
175 } else {
176 throw CanteraError("ReactorBase::network",
177 "Reactor is not part of a ReactorNet");
178 }
179}
180
182{
183 if (m_net) {
184 throw CanteraError("ReactorBase::setNetwork",
185 "Reactor {} is already part of a ReactorNet", m_name);
186 }
187 m_net = net;
188}
189
191{
192 double mout = 0.0;
193 for (size_t i = 0; i < m_outlet.size(); i++) {
194 mout += m_outlet[i]->massFlowRate();
195 }
196 return mass()/mout;
197}
198
200{
201 return m_thermo->density();
202}
203
205{
206 return m_thermo->temperature();
207}
208
209const double* ReactorBase::massFractions() const
210{
211 return m_thermo->massFractions();
212}
213
214double ReactorBase::massFraction(size_t k) const
215{
216 return m_thermo->massFraction(k);
217}
218
220{
221 return *m_inlet[n];
222}
224{
225 return *m_outlet[n];
226}
227
228}
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.
Definition FlowDevice.h:25
double massFraction(size_t k) const
Return the mass fraction of a single species.
Definition Phase.cpp:480
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:245
double temperature() const
Temperature (K).
Definition Phase.h:598
const double * massFractions() const
Return a const pointer to the mass fraction array.
Definition Phase.h:478
virtual double density() const
Density (kg/m^3).
Definition Phase.h:623
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Definition Phase.h:616
Base class for reactor objects.
Definition ReactorBase.h:51
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.
Definition ReactorBase.h:76
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.
Definition ReactorNet.h:30
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.
Definition ReactorNet.h:336
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.
Definition Wall.h:23
virtual Eigen::SparseMatrix< double > jacobian()
Calculate the Jacobian of a specific reactor specialization.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Definition AnyMap.cpp:1997