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