Cantera  3.1.0
Loading...
Searching...
No Matches
Reactor.h
Go to the documentation of this file.
1//! @file Reactor.h
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
6#ifndef CT_REACTOR_H
7#define CT_REACTOR_H
8
9#include "ReactorBase.h"
10#include "cantera/numerics/eigen_sparse.h"
11
12
13namespace Cantera
14{
15
16class Solution;
17class AnyMap;
18
19/**
20 * Class Reactor is a general-purpose class for stirred reactors. The reactor
21 * may have an arbitrary number of inlets and outlets, each of which may be
22 * connected to a "flow device" such as a mass flow controller, a pressure
23 * regulator, etc. Additional reactors may be connected to the other end of
24 * the flow device, allowing construction of arbitrary reactor networks.
25 *
26 * The reactor class integrates the same governing equations no matter what
27 * type of reactor is simulated. The differences among reactor types are
28 * completely specified by the attached flow devices and the time-dependent
29 * user-specified boundary conditions.
30 *
31 * If an instance of class Reactor is used directly, it will simulate an
32 * adiabatic, constant volume reactor with gas-phase chemistry but no surface
33 * chemistry. Other reactor types may be simulated by deriving a class from
34 * Reactor. This method allows specifying the following in terms of the
35 * instantaneous reactor state:
36 *
37 * - rate of change of the total volume (m^3/s)
38 * - surface heat loss rate (W)
39 * - species surface production rates (kmol/s)
40 *
41 * See the [Science Reference](../reference/reactors/controlreactor.html) for
42 * the governing equations of class Reactor.
43 *
44 * @ingroup reactorGroup
45 */
46class Reactor : public ReactorBase
47{
48public:
49 Reactor(shared_ptr<Solution> sol, const string& name="(none)");
50 using ReactorBase::ReactorBase; // inherit constructors
51
52 string type() const override {
53 return "Reactor";
54 }
55
56 //! Indicate whether the governing equations for this reactor type are a system of
57 //! ODEs or DAEs. In the first case, this class implements the eval() method. In the
58 //! second case, this class implements the evalDae() method.
59 virtual bool isOde() const {
60 return true;
61 }
62
63 //! Indicates whether the governing equations for this reactor are functions of time
64 //! or a spatial variable. All reactors in a network must have the same value.
65 virtual bool timeIsIndependent() const {
66 return true;
67 }
68
69 /**
70 * Insert something into the reactor. The 'something' must belong to a class
71 * that is a subclass of both ThermoPhase and Kinetics.
72 * @deprecated Unused; to be removed after %Cantera 3.1.
73 */
74 template<class G>
75 void insert(G& contents) {
76 warn_deprecated("Reactor::insert", "Unused; to be removed after Cantera 3.1.");
79 }
80
82
83 void setChemistry(bool cflag=true) override {
84 m_chem = cflag;
85 }
86
87 //! Returns `true` if changes in the reactor composition due to chemical reactions
88 //! are enabled.
89 bool chemistryEnabled() const {
90 return m_chem;
91 }
92
93 void setEnergy(int eflag=1) override {
94 if (eflag > 0) {
95 m_energy = true;
96 } else {
97 m_energy = false;
98 }
99 }
100
101 //! Returns `true` if solution of the energy equation is enabled.
102 bool energyEnabled() const {
103 return m_energy;
104 }
105
106 //! Number of equations (state variables) for this reactor
107 size_t neq() {
108 if (!m_nv) {
109 initialize();
110 }
111 return m_nv;
112 }
113
114 //! Get the the current state of the reactor.
115 /*!
116 * @param[out] y state vector representing the initial state of the reactor
117 */
118 virtual void getState(double* y);
119
120 //! Get the current state and derivative vector of the reactor for a DAE solver
121 /*!
122 * @param[out] y state vector representing the initial state of the reactor
123 * @param[out] ydot state vector representing the initial derivatives of the
124 * reactor
125 */
126 virtual void getStateDae(double* y, double* ydot) {
127 throw NotImplementedError("Reactor::getStateDae(y, ydot)");
128 }
129
130 void initialize(double t0=0.0) override;
131
132 //! Evaluate the reactor governing equations. Called by ReactorNet::eval.
133 //! @param[in] t time.
134 //! @param[out] LHS pointer to start of vector of left-hand side
135 //! coefficients for governing equations, length m_nv, default values 1
136 //! @param[out] RHS pointer to start of vector of right-hand side
137 //! coefficients for governing equations, length m_nv, default values 0
138 virtual void eval(double t, double* LHS, double* RHS);
139
140 /**
141 * Evaluate the reactor governing equations. Called by ReactorNet::eval.
142 * @param[in] t time.
143 * @param[in] y solution vector, length neq()
144 * @param[in] ydot rate of change of solution vector, length neq()
145 * @param[out] residual residuals vector, length neq()
146 */
147 virtual void evalDae(double t, double* y, double* ydot, double* residual) {
148 throw NotImplementedError("Reactor::evalDae");
149 }
150
151 //! Given a vector of length neq(), mark which variables should be
152 //! considered algebraic constraints
153 virtual void getConstraints(double* constraints) {
154 throw NotImplementedError("Reactor::getConstraints");
155 }
156
157 void syncState() override;
158
159 //! Set the state of the reactor to correspond to the state vector *y*.
160 virtual void updateState(double* y);
161
162 //! Number of sensitivity parameters associated with this reactor
163 //! (including walls)
164 virtual size_t nSensParams() const;
165
166 //! Add a sensitivity parameter associated with the reaction number *rxn*
167 //! (in the homogeneous phase).
168 virtual void addSensitivityReaction(size_t rxn);
169
170 //! Add a sensitivity parameter associated with the enthalpy formation of
171 //! species *k* (in the homogeneous phase)
172 virtual void addSensitivitySpeciesEnthalpy(size_t k);
173
174 //! Return the index in the solution vector for this reactor of the
175 //! component named *nm*. Possible values for *nm* are "mass", "volume",
176 //! "int_energy", the name of a homogeneous phase species, or the name of a
177 //! surface species.
178 virtual size_t componentIndex(const string& nm) const;
179
180 //! Return the name of the solution component with index *i*.
181 //! @see componentIndex()
182 virtual string componentName(size_t k);
183
184 //! Set absolute step size limits during advance
185 //! @param limits array of step size limits with length neq
186 void setAdvanceLimits(const double* limits);
187
188 //! Check whether Reactor object uses advance limits
189 //! @returns True if at least one limit is set, False otherwise
190 bool hasAdvanceLimits() const {
191 return !m_advancelimits.empty();
192 }
193
194 //! Retrieve absolute step size limits during advance
195 //! @param[out] limits array of step size limits with length neq
196 //! @returns True if at least one limit is set, False otherwise
197 bool getAdvanceLimits(double* limits) const;
198
199 //! Set individual step size limit for component name *nm*
200 //! @param nm component name
201 //! @param limit value for step size limit
202 void setAdvanceLimit(const string& nm, const double limit);
203
204 //! Calculate the Jacobian of a specific Reactor specialization.
205 //! @warning Depending on the particular implementation, this may return an
206 //! approximate Jacobian intended only for use in forming a preconditioner for
207 //! iterative solvers.
208 //! @ingroup derivGroup
209 //!
210 //! @warning This method is an experimental part of the %Cantera
211 //! API and may be changed or removed without notice.
212 virtual Eigen::SparseMatrix<double> jacobian() {
213 throw NotImplementedError("Reactor::jacobian");
214 }
215
216 //! Calculate the reactor-specific Jacobian using a finite difference method.
217 //!
218 //! This method is used only for informational purposes. Jacobian calculations
219 //! for the full reactor system are handled internally by CVODES.
220 //!
221 //! @warning This method is an experimental part of the %Cantera
222 //! API and may be changed or removed without notice.
223 Eigen::SparseMatrix<double> finiteDifferenceJacobian();
224
225 //! Use this to set the kinetics objects derivative settings
226 virtual void setDerivativeSettings(AnyMap& settings);
227
228 //! Set reaction rate multipliers based on the sensitivity variables in
229 //! *params*.
230 virtual void applySensitivity(double* params);
231
232 //! Reset the reaction rate multipliers
233 virtual void resetSensitivity(double* params);
234
235 //! Return a false if preconditioning is not supported or true otherwise.
236 //!
237 //! @warning This method is an experimental part of the %Cantera
238 //! API and may be changed or removed without notice.
239 //!
240 //! @since New in %Cantera 3.0
241 //!
242 virtual bool preconditionerSupported() const {return false;};
243
244protected:
245 void setKinetics(Kinetics& kin) override;
246
247 //! Return the index in the solution vector for this reactor of the species
248 //! named *nm*, in either the homogeneous phase or a surface phase, relative
249 //! to the start of the species terms. Used to implement componentIndex for
250 //! specific reactor implementations.
251 virtual size_t speciesIndex(const string& nm) const;
252
253 //! Evaluate terms related to Walls. Calculates #m_vdot and #m_Qdot based on
254 //! wall movement and heat transfer.
255 //! @param t the current time
256 virtual void evalWalls(double t);
257
258 //! Evaluate terms related to surface reactions.
259 //! @param[out] LHS Multiplicative factor on the left hand side of ODE for surface
260 //! species coverages
261 //! @param[out] RHS Right hand side of ODE for surface species coverages
262 //! @param[out] sdot array of production rates of bulk phase species on surfaces
263 //! [kmol/s]
264 virtual void evalSurfaces(double* LHS, double* RHS, double* sdot);
265
266 virtual void evalSurfaces(double* RHS, double* sdot);
267
268 //! Update the state of SurfPhase objects attached to this reactor
269 virtual void updateSurfaceState(double* y);
270
271 //! Update the state information needed by connected reactors, flow devices,
272 //! and reactor walls. Called from updateState().
273 //! @param updatePressure Indicates whether to update #m_pressure. Should
274 //! `true` for reactors where the pressure is a dependent property,
275 //! calculated from the state, and `false` when the pressure is constant
276 //! or an independent variable.
277 virtual void updateConnected(bool updatePressure);
278
279 //! Get initial conditions for SurfPhase objects attached to this reactor
280 virtual void getSurfaceInitialConditions(double* y);
281
282 //! Pointer to the homogeneous Kinetics object that handles the reactions
283 Kinetics* m_kin = nullptr;
284
285 double m_vdot = 0.0; //!< net rate of volume change from moving walls [m^3/s]
286
287 double m_Qdot = 0.0; //!< net heat transfer into the reactor, through walls [W]
288
289 double m_mass = 0.0; //!< total mass
290 vector<double> m_work;
291
292 //! Production rates of gas phase species on surfaces [kmol/s]
293 vector<double> m_sdot;
294
295 vector<double> m_wdot; //!< Species net molar production rates
296 vector<double> m_uk; //!< Species molar internal energies
297 bool m_chem = false;
298 bool m_energy = true;
299 size_t m_nv = 0;
300 size_t m_nv_surf; //!!< Number of variables associated with reactor surfaces
301
302 vector<double> m_advancelimits; //!< Advance step limit
303
304 // Data associated each sensitivity parameter
305 vector<SensitivityParameter> m_sensParams;
306
307 //! Vector of triplets representing the jacobian
308 vector<Eigen::Triplet<double>> m_jac_trips;
309};
310}
311
312#endif
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
Public interface for kinetics managers.
Definition Kinetics.h:125
An error indicating that an unimplemented function has been called.
Base class for stirred reactors.
Definition ReactorBase.h:49
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.
Definition ReactorBase.h:68
Class Reactor is a general-purpose class for stirred reactors.
Definition Reactor.h:47
virtual string componentName(size_t k)
Return the name of the solution component with index i.
Definition Reactor.cpp:474
bool chemistryEnabled() const
Returns true if changes in the reactor composition due to chemical reactions are enabled.
Definition Reactor.h:89
virtual void evalSurfaces(double *LHS, double *RHS, double *sdot)
Evaluate terms related to surface reactions.
Definition Reactor.cpp:297
virtual size_t componentIndex(const string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
Definition Reactor.cpp:458
virtual void getStateDae(double *y, double *ydot)
Get the current state and derivative vector of the reactor for a DAE solver.
Definition Reactor.h:126
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
Definition Reactor.cpp:185
void insert(G &contents)
Insert something into the reactor.
Definition Reactor.h:75
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Definition Reactor.cpp:500
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:283
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:295
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
Definition Reactor.cpp:221
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Definition Reactor.cpp:285
bool energyEnabled() const
Returns true if solution of the energy equation is enabled.
Definition Reactor.h:102
Eigen::SparseMatrix< double > finiteDifferenceJacobian()
Calculate the reactor-specific Jacobian using a finite difference method.
Definition Reactor.cpp:326
size_t neq()
Number of equations (state variables) for this reactor.
Definition Reactor.h:107
virtual void addSensitivitySpeciesEnthalpy(size_t k)
Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous p...
Definition Reactor.cpp:421
string type() const override
String indicating the reactor model implemented.
Definition Reactor.h:52
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:287
void setKinetics(Kinetics &kin) override
Specify the kinetics manager for the reactor.
Definition Reactor.cpp:49
vector< double > m_advancelimits
!< Number of variables associated with reactor surfaces
Definition Reactor.h:302
void setEnergy(int eflag=1) override
Set the energy equation on or off.
Definition Reactor.h:93
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase).
Definition Reactor.cpp:408
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor (including walls)
Definition Reactor.cpp:122
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
Definition Reactor.cpp:40
vector< double > m_uk
Species molar internal energies.
Definition Reactor.h:296
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
Definition Reactor.cpp:137
virtual void evalDae(double t, double *y, double *ydot, double *residual)
Evaluate the reactor governing equations.
Definition Reactor.h:147
void setChemistry(bool cflag=true) override
Enable or disable changes in reactor composition due to chemical reactions.
Definition Reactor.h:83
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
Definition Reactor.cpp:543
double m_mass
total mass
Definition Reactor.h:289
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
Definition Reactor.h:308
void setAdvanceLimit(const string &nm, const double limit)
Set individual step size limit for component name nm
Definition Reactor.cpp:569
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
Definition Reactor.cpp:522
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition Reactor.h:293
virtual void getState(double *y)
Get the the current state of the reactor.
Definition Reactor.cpp:59
bool hasAdvanceLimits() const
Check whether Reactor object uses advance limits.
Definition Reactor.h:190
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:285
void syncState() override
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
Definition Reactor.cpp:131
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
Definition Reactor.cpp:85
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:94
virtual void getConstraints(double *constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Definition Reactor.h:153
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
Definition Reactor.h:242
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
Definition Reactor.cpp:558
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...
Definition Reactor.cpp:436
virtual bool isOde() const
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs.
Definition Reactor.h:59
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
Definition Reactor.h:65
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
Definition Reactor.cpp:194
virtual Eigen::SparseMatrix< double > jacobian()
Calculate the Jacobian of a specific Reactor specialization.
Definition Reactor.h:212
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