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