Loading [MathJax]/extensions/tex2jax.js
Cantera  3.2.0a1
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 void setChemistry(bool cflag=true) override {
70 m_chem = cflag;
71 }
72
73 //! Returns `true` if changes in the reactor composition due to chemical reactions
74 //! are enabled.
75 bool chemistryEnabled() const {
76 return m_chem;
77 }
78
79 void setEnergy(int eflag=1) override {
80 if (eflag > 0) {
81 m_energy = true;
82 } else {
83 m_energy = false;
84 }
85 }
86
87 //! Returns `true` if solution of the energy equation is enabled.
88 bool energyEnabled() const {
89 return m_energy;
90 }
91
92 //! Number of equations (state variables) for this reactor
93 size_t neq() {
94 if (!m_nv) {
95 initialize();
96 }
97 return m_nv;
98 }
99
100 //! Get the the current state of the reactor.
101 /*!
102 * @param[out] y state vector representing the initial state of the reactor
103 */
104 virtual void getState(double* y);
105
106 //! Get the current state and derivative vector of the reactor for a DAE solver
107 /*!
108 * @param[out] y state vector representing the initial state of the reactor
109 * @param[out] ydot state vector representing the initial derivatives of the
110 * reactor
111 */
112 virtual void getStateDae(double* y, double* ydot) {
113 throw NotImplementedError("Reactor::getStateDae(y, ydot)");
114 }
115
116 void initialize(double t0=0.0) override;
117
118 //! Evaluate the reactor governing equations. Called by ReactorNet::eval.
119 //! @param[in] t time.
120 //! @param[out] LHS pointer to start of vector of left-hand side
121 //! coefficients for governing equations, length m_nv, default values 1
122 //! @param[out] RHS pointer to start of vector of right-hand side
123 //! coefficients for governing equations, length m_nv, default values 0
124 virtual void eval(double t, double* LHS, double* RHS);
125
126 /**
127 * Evaluate the reactor governing equations. Called by ReactorNet::eval.
128 * @param[in] t time.
129 * @param[in] y solution vector, length neq()
130 * @param[in] ydot rate of change of solution vector, length neq()
131 * @param[out] residual residuals vector, length neq()
132 */
133 virtual void evalDae(double t, double* y, double* ydot, double* residual) {
134 throw NotImplementedError("Reactor::evalDae");
135 }
136
137 //! Given a vector of length neq(), mark which variables should be
138 //! considered algebraic constraints
139 virtual void getConstraints(double* constraints) {
140 throw NotImplementedError("Reactor::getConstraints");
141 }
142
143 void syncState() override;
144
145 //! Set the state of the reactor to correspond to the state vector *y*.
146 virtual void updateState(double* y);
147
148 //! Number of sensitivity parameters associated with this reactor
149 //! (including walls)
150 virtual size_t nSensParams() const;
151
152 //! Add a sensitivity parameter associated with the reaction number *rxn*
153 //! (in the homogeneous phase).
154 virtual void addSensitivityReaction(size_t rxn);
155
156 //! Add a sensitivity parameter associated with the enthalpy formation of
157 //! species *k* (in the homogeneous phase)
158 virtual void addSensitivitySpeciesEnthalpy(size_t k);
159
160 //! Return the index in the solution vector for this reactor of the
161 //! component named *nm*. Possible values for *nm* are "mass", "volume",
162 //! "int_energy", the name of a homogeneous phase species, or the name of a
163 //! surface species.
164 virtual size_t componentIndex(const string& nm) const;
165
166 //! Return the name of the solution component with index *i*.
167 //! @see componentIndex()
168 virtual string componentName(size_t k);
169
170 //! Set absolute step size limits during advance
171 //! @param limits array of step size limits with length neq
172 void setAdvanceLimits(const double* limits);
173
174 //! Check whether Reactor object uses advance limits
175 //! @returns True if at least one limit is set, False otherwise
176 bool hasAdvanceLimits() const {
177 return !m_advancelimits.empty();
178 }
179
180 //! Retrieve absolute step size limits during advance
181 //! @param[out] limits array of step size limits with length neq
182 //! @returns True if at least one limit is set, False otherwise
183 bool getAdvanceLimits(double* limits) const;
184
185 //! Set individual step size limit for component name *nm*
186 //! @param nm component name
187 //! @param limit value for step size limit
188 void setAdvanceLimit(const string& nm, const double limit);
189
190 //! Calculate the Jacobian of a specific Reactor specialization.
191 //! @warning Depending on the particular implementation, this may return an
192 //! approximate Jacobian intended only for use in forming a preconditioner for
193 //! iterative solvers.
194 //! @ingroup derivGroup
195 //!
196 //! @warning This method is an experimental part of the %Cantera
197 //! API and may be changed or removed without notice.
198 virtual Eigen::SparseMatrix<double> jacobian() {
199 throw NotImplementedError("Reactor::jacobian");
200 }
201
202 //! Calculate the reactor-specific Jacobian using a finite difference method.
203 //!
204 //! This method is used only for informational purposes. Jacobian calculations
205 //! for the full reactor system are handled internally by CVODES.
206 //!
207 //! @warning This method is an experimental part of the %Cantera
208 //! API and may be changed or removed without notice.
209 Eigen::SparseMatrix<double> finiteDifferenceJacobian();
210
211 //! Use this to set the kinetics objects derivative settings
212 virtual void setDerivativeSettings(AnyMap& settings);
213
214 //! Set reaction rate multipliers based on the sensitivity variables in
215 //! *params*.
216 virtual void applySensitivity(double* params);
217
218 //! Reset the reaction rate multipliers
219 virtual void resetSensitivity(double* params);
220
221 //! Return a false if preconditioning is not supported or true otherwise.
222 //!
223 //! @warning This method is an experimental part of the %Cantera
224 //! API and may be changed or removed without notice.
225 //!
226 //! @since New in %Cantera 3.0
227 //!
228 virtual bool preconditionerSupported() const {return false;};
229
230protected:
231 void setKinetics(Kinetics& kin) override;
232
233 //! Return the index in the solution vector for this reactor of the species
234 //! named *nm*, in either the homogeneous phase or a surface phase, relative
235 //! to the start of the species terms. Used to implement componentIndex for
236 //! specific reactor implementations.
237 virtual size_t speciesIndex(const string& nm) const;
238
239 //! Evaluate terms related to Walls. Calculates #m_vdot and #m_Qdot based on
240 //! wall movement and heat transfer.
241 //! @param t the current time
242 virtual void evalWalls(double t);
243
244 //! Evaluate terms related to surface reactions.
245 //! @param[out] LHS Multiplicative factor on the left hand side of ODE for surface
246 //! species coverages
247 //! @param[out] RHS Right hand side of ODE for surface species coverages
248 //! @param[out] sdot array of production rates of bulk phase species on surfaces
249 //! [kmol/s]
250 virtual void evalSurfaces(double* LHS, double* RHS, double* sdot);
251
252 virtual void evalSurfaces(double* RHS, double* sdot);
253
254 //! Update the state of SurfPhase objects attached to this reactor
255 virtual void updateSurfaceState(double* y);
256
257 //! Update the state information needed by connected reactors, flow devices,
258 //! and reactor walls. Called from updateState().
259 //! @param updatePressure Indicates whether to update #m_pressure. Should
260 //! `true` for reactors where the pressure is a dependent property,
261 //! calculated from the state, and `false` when the pressure is constant
262 //! or an independent variable.
263 virtual void updateConnected(bool updatePressure);
264
265 //! Get initial conditions for SurfPhase objects attached to this reactor
266 virtual void getSurfaceInitialConditions(double* y);
267
268 //! Pointer to the homogeneous Kinetics object that handles the reactions
269 Kinetics* m_kin = nullptr;
270
271 double m_vdot = 0.0; //!< net rate of volume change from moving walls [m^3/s]
272
273 double m_Qdot = 0.0; //!< net heat transfer into the reactor, through walls [W]
274
275 double m_mass = 0.0; //!< total mass
276 vector<double> m_work;
277
278 //! Production rates of gas phase species on surfaces [kmol/s]
279 vector<double> m_sdot;
280
281 vector<double> m_wdot; //!< Species net molar production rates
282 vector<double> m_uk; //!< Species molar internal energies
283 bool m_chem = false;
284 bool m_energy = true;
285 size_t m_nv = 0;
286 size_t m_nv_surf; //!!< Number of variables associated with reactor surfaces
287
288 vector<double> m_advancelimits; //!< Advance step limit
289
290 // Data associated each sensitivity parameter
291 vector<SensitivityParameter> m_sensParams;
292
293 //! Vector of triplets representing the jacobian
294 vector<Eigen::Triplet<double>> m_jac_trips;
295};
296}
297
298#endif
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:432
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
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:465
bool chemistryEnabled() const
Returns true if changes in the reactor composition due to chemical reactions are enabled.
Definition Reactor.h:75
virtual void evalSurfaces(double *LHS, double *RHS, double *sdot)
Evaluate terms related to surface reactions.
Definition Reactor.cpp:288
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:449
virtual void getStateDae(double *y, double *ydot)
Get the current state and derivative vector of the reactor for a DAE solver.
Definition Reactor.h:112
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
Definition Reactor.cpp:176
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Definition Reactor.cpp:491
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:269
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:281
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
Definition Reactor.cpp:212
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Definition Reactor.cpp:276
bool energyEnabled() const
Returns true if solution of the energy equation is enabled.
Definition Reactor.h:88
Eigen::SparseMatrix< double > finiteDifferenceJacobian()
Calculate the reactor-specific Jacobian using a finite difference method.
Definition Reactor.cpp:317
size_t neq()
Number of equations (state variables) for this reactor.
Definition Reactor.h:93
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:412
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:273
void setKinetics(Kinetics &kin) override
Specify the kinetics manager for the reactor.
Definition Reactor.cpp:40
vector< double > m_advancelimits
!< Number of variables associated with reactor surfaces
Definition Reactor.h:288
void setEnergy(int eflag=1) override
Set the energy equation on or off.
Definition Reactor.h:79
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase).
Definition Reactor.cpp:399
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor (including walls)
Definition Reactor.cpp:113
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
Definition Reactor.cpp:31
vector< double > m_uk
Species molar internal energies.
Definition Reactor.h:282
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
Definition Reactor.cpp:128
virtual void evalDae(double t, double *y, double *ydot, double *residual)
Evaluate the reactor governing equations.
Definition Reactor.h:133
void setChemistry(bool cflag=true) override
Enable or disable changes in reactor composition due to chemical reactions.
Definition Reactor.h:69
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
Definition Reactor.cpp:534
double m_mass
total mass
Definition Reactor.h:275
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
Definition Reactor.h:294
void setAdvanceLimit(const string &nm, const double limit)
Set individual step size limit for component name nm
Definition Reactor.cpp:560
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
Definition Reactor.cpp:513
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition Reactor.h:279
virtual void getState(double *y)
Get the the current state of the reactor.
Definition Reactor.cpp:50
bool hasAdvanceLimits() const
Check whether Reactor object uses advance limits.
Definition Reactor.h:176
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:271
void syncState() override
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
Definition Reactor.cpp:122
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
Definition Reactor.cpp:76
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:85
virtual void getConstraints(double *constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Definition Reactor.h:139
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
Definition Reactor.h:228
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
Definition Reactor.cpp:549
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:427
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:185
virtual Eigen::SparseMatrix< double > jacobian()
Calculate the Jacobian of a specific Reactor specialization.
Definition Reactor.h:198
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595