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