Cantera  4.0.0a1
Loading...
Searching...
No Matches
ReactorBase.h
Go to the documentation of this file.
1//! @file ReactorBase.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_REACTORBASE_H
7#define CT_REACTORBASE_H
8
11#include "cantera/numerics/eigen_sparse.h"
12
13namespace Cantera
14{
15
16//! @defgroup zerodGroup Zero-Dimensional Reactor Networks
17//!
18//! @details See the [Reactor Science](../reference/reactors/index.html) section of the
19//! %Cantera website for a description of the governing equations for specific reactor
20//! types and the methods used for solving networks of interconnected reactors.
21
22class FlowDevice;
23class WallBase;
24class ReactorNet;
25class ReactorSurface;
26class Kinetics;
27class ThermoPhase;
28class Solution;
29class AnyMap;
30
31enum class SensParameterType {
32 reaction,
33 enthalpy
34};
35
37{
38 size_t local; //!< local parameter index
39 size_t global; //!< global parameter index
40 double value; //!< nominal value of the parameter
41 SensParameterType type; //!< type of sensitivity parameter
42};
43
44/**
45 * Base class for reactor objects. Allows using any substance model, with arbitrary
46 * inflow, outflow, heat loss/gain, surface chemistry, and volume change, whenever
47 * defined.
48 * @ingroup reactorGroup
49 */
50class ReactorBase : public std::enable_shared_from_this<ReactorBase>
51{
52public:
53 //! Instantiate a ReactorBase object with Solution contents.
54 //! @param sol Solution object to be set.
55 //! @param name Name of the reactor.
56 //! @since New in %Cantera 3.1.
57 ReactorBase(shared_ptr<Solution> sol, const string& name="(none)");
58
59 //! Instantiate a ReactorBase object with Solution contents.
60 //! @param sol Solution object representing the contents of this reactor
61 //! @param clone Determines whether to clone `sol` so that the internal state of
62 //! this reactor is independent of the original Solution object and any Solution
63 //! objects used by other reactors in the network.
64 //! @param name Name of the reactor.
65 //! @since Added the `clone` argument in %Cantera 3.2. If not specified, the default
66 //! behavior in %Cantera 3.2 is not to clone the Solution object. This will
67 //! change after %Cantera 3.2 to default to `true`.
68 ReactorBase(shared_ptr<Solution> sol, bool clone, const string& name="(none)");
69
70 virtual ~ReactorBase();
71 ReactorBase(const ReactorBase&) = delete;
72 ReactorBase& operator=(const ReactorBase&) = delete;
73
74 //! String indicating the reactor model implemented. Usually
75 //! corresponds to the name of the derived class.
76 virtual string type() const {
77 return "ReactorBase";
78 }
79
80 //! Return the name of this reactor
81 string name() const {
82 return m_name;
83 }
84
85 //! Set the name of this reactor
86 void setName(const string& name) {
87 m_name = name;
88 }
89
90 //! Set the default name of a reactor. Returns `false` if it was previously set.
91 bool setDefaultName(map<string, int>& counts);
92
93 //! Access the Solution object used to represent the contents of this reactor.
94 //! @since New in %Cantera 3.2
95 shared_ptr<Solution> phase() { return m_solution; }
96
97 //! Access the Solution object used to represent the contents of this reactor.
98 //! @since New in %Cantera 3.2
99 shared_ptr<const Solution> phase() const { return m_solution; }
100
101 //! Indicates whether the governing equations for this reactor are functions of time
102 //! or a spatial variable. All reactors in a network must have the same value.
103 virtual bool timeIsIndependent() const {
104 return true;
105 }
106
107 //! Number of equations (state variables) for this reactor
108 size_t neq() {
109 return m_nv;
110 }
111
112 //! @name Methods to set up a simulation
113 //! @{
114
115 //! Set the initial reactor volume.
116 virtual void setInitialVolume(double vol) {
117 throw NotImplementedError("ReactorBase::setInitialVolume",
118 "Volume is undefined for reactors of type '{}'.", type());
119 }
120
121 //! Returns an area associated with a reactor [m²].
122 //! Examples: surface area of ReactorSurface or cross section area of FlowReactor.
123 virtual double area() const {
124 throw NotImplementedError("ReactorBase::area",
125 "Area is undefined for reactors of type '{}'.", type());
126 }
127
128 //! Evaluate contributions from walls connected to this reactor
129 virtual void evalWalls(double t) {
130 throw NotImplementedError("ReactorBase::evalWalls",
131 "Wall evaluation is undefined for reactors of type '{}'.", type());
132 }
133
134 //! Set an area associated with a reactor [m²].
135 //! Examples: surface area of ReactorSurface or cross section area of FlowReactor.
136 virtual void setArea(double a) {
137 throw NotImplementedError("ReactorBase::setArea",
138 "Area is undefined for reactors of type '{}'.", type());
139 }
140
141 //! Returns `true` if changes in the reactor composition due to chemical reactions
142 //! are enabled.
143 //! @since New in %Cantera 3.2.
144 virtual bool chemistryEnabled() const {
145 throw NotImplementedError("ReactorBase::chemistryEnabled",
146 "Not implemented for reactor type '{}'.", type());
147 }
148
149 //! Enable or disable changes in reactor composition due to chemical reactions.
150 //! @since New in %Cantera 3.2.
151 virtual void setChemistryEnabled(bool cflag = true) {
152 throw NotImplementedError("ReactorBase::setChemistryEnabled",
153 "Not implemented for reactor type '{}'.", type());
154 }
155
156 //! Returns `true` if solution of the energy equation is enabled.
157 //! @since New in %Cantera 3.2.
158 virtual bool energyEnabled() const {
159 throw NotImplementedError("ReactorBase::energyEnabled",
160 "Not implemented for reactor type '{}'.", type());
161 }
162
163 //! Set the energy equation on or off.
164 //! @since New in %Cantera 3.2.
165 virtual void setEnergyEnabled(bool eflag = true) {
166 throw NotImplementedError("ReactorBase::setEnergyEnabled",
167 "Not implemented for reactor type '{}'.", type());
168 }
169
170 //! Connect an inlet FlowDevice to this reactor
171 virtual void addInlet(FlowDevice& inlet);
172
173 //! Connect an outlet FlowDevice to this reactor
174 virtual void addOutlet(FlowDevice& outlet);
175
176 //! Return a reference to the *n*-th inlet FlowDevice connected to this reactor.
177 FlowDevice& inlet(size_t n = 0);
178
179 //! Return a reference to the *n*-th outlet FlowDevice connected to this reactor.
180 FlowDevice& outlet(size_t n = 0);
181
182 //! Return the number of inlet FlowDevice objects connected to this reactor.
183 size_t nInlets() {
184 return m_inlet.size();
185 }
186
187 //! Return the number of outlet FlowDevice objects connected to this reactor.
188 size_t nOutlets() {
189 return m_outlet.size();
190 }
191
192 //! Return the number of Wall objects connected to this reactor.
193 size_t nWalls() {
194 return m_wall.size();
195 }
196
197 //! Insert a Wall between this reactor and another reactor.
198 /*!
199 * `lr` = 0 if this reactor is to the left of the wall and `lr` = 1 if
200 * this reactor is to the right of the wall. This method is called
201 * automatically for both the left and right reactors by WallBase::install.
202 */
203 virtual void addWall(WallBase& w, int lr);
204
205 //! Return a reference to the *n*-th Wall connected to this reactor.
206 WallBase& wall(size_t n);
207
208 //! Add a ReactorSurface object to a Reactor object.
209 //! @attention This method should generally not be called directly by users.
210 //! Reactor and ReactorSurface objects should be connected by providing adjacent
211 //! reactors to the newReactorSurface factory function.
212 virtual void addSurface(ReactorSurface* surf);
213
214 //! Return a reference to the *n*-th ReactorSurface connected to this reactor.
215 ReactorSurface* surface(size_t n);
216
217 //! Return the number of surfaces in a reactor
218 virtual size_t nSurfs() const {
219 return m_surfaces.size();
220 }
221
222 //! Production rates on surfaces.
223 //!
224 //! For bulk reactors, this contains the total production rates [kmol/s] of bulk
225 //! phase species due to reactions on all adjacent species.
226 //!
227 //! For surfaces, this contains the production rates [kmol/m²/s] of species on the
228 //! surface and all adjacent phases, in the order defined by the InterfaceKinetics
229 //! object.
230 span<const double> surfaceProductionRates() const {
231 return m_sdot;
232 }
233
234 /**
235 * Initialize the reactor. Called automatically by ReactorNet::initialize.
236 */
237 virtual void initialize(double t0 = 0.0) {
238 throw NotImplementedError("ReactorBase::initialize");
239 }
240
241 //! Get the current state of the reactor.
242 /*!
243 * @param[out] y state vector representing the initial state of the reactor
244 */
245 virtual void getState(double* y) {
246 throw NotImplementedError("ReactorBase::getState",
247 "Not implemented for reactor type '{}'.", type());
248 }
249
250 //! Get the current state and derivative vector of the reactor for a DAE solver
251 /*!
252 * @param[out] y state vector representing the initial state of the reactor
253 * @param[out] ydot state vector representing the initial derivatives of the
254 * reactor
255 */
256 virtual void getStateDae(double* y, double* ydot) {
257 throw NotImplementedError("ReactorBase::getStateDae(y, ydot)");
258 }
259
260 //! Evaluate the reactor governing equations. Called by ReactorNet::eval.
261 //! @param[in] t time.
262 //! @param[out] LHS pointer to start of vector of left-hand side
263 //! coefficients for governing equations, length m_nv, default values 1
264 //! @param[out] RHS pointer to start of vector of right-hand side
265 //! coefficients for governing equations, length m_nv, default values 0
266 virtual void eval(double t, double* LHS, double* RHS) {
267 throw NotImplementedError("ReactorBase::eval");
268 }
269
270 /**
271 * Evaluate the reactor governing equations. Called by ReactorNet::eval.
272 * @param[in] t time.
273 * @param[in] y solution vector, length neq()
274 * @param[in] ydot rate of change of solution vector, length neq()
275 * @param[out] residual residuals vector, length neq()
276 */
277 virtual void evalDae(double t, double* y, double* ydot, double* residual) {
278 throw NotImplementedError("ReactorBase::evalDae");
279 }
280
281 //! Evaluate the governing equations with modifications for the steady-state solver.
282 //!
283 //! This method calls the standard eval() method then modifies elements of `RHS`
284 //! that correspond to algebraic constraints.
285 //!
286 //! @since New in %Cantera 4.0.
287 virtual void evalSteady(double t, double* LHS, double* RHS) {
288 throw NotImplementedError("ReactorBase::evalSteady",
289 "Not implemented for reactor type '{}'.", type());
290 }
291
292 //! Given a vector of length neq(), mark which variables should be
293 //! considered algebraic constraints
294 virtual void getConstraints(double* constraints) {
295 throw NotImplementedError("ReactorBase::getConstraints");
296 }
297
298 //! Initialize the reactor before solving a steady-state problem.
299 //!
300 //! This method is responsible for storing the initial value for any algebraic
301 //! constraints and returning the indices of those constraints.
302 //!
303 //! @return Indices of equations that are algebraic constraints when solving the
304 //! steady-state problem.
305 //!
306 //! @warning This method is an experimental part of the %Cantera API and may be
307 //! changed or removed without notice.
308 //! @since New in %Cantera 3.2. Renamed from `steadyConstraints` in %Cantera 4.0.
309 virtual vector<size_t> initializeSteady() {
310 throw NotImplementedError("ReactorBase::initializeSteady");
311 }
312
313 //! Set the state of the reactor to correspond to the state vector *y*.
314 virtual void updateState(double* y) {
315 throw NotImplementedError("ReactorBase::updateState");
316 }
317
318 //! Add a sensitivity parameter associated with the enthalpy formation of
319 //! species *k*.
320 virtual void addSensitivitySpeciesEnthalpy(size_t k) {
321 throw NotImplementedError("ReactorBase::addSensitivitySpeciesEnthalpy");
322 }
323
324 //! Return the index in the solution vector for this reactor of the
325 //! component named *nm*.
326 virtual size_t componentIndex(const string& nm) const {
327 throw NotImplementedError("ReactorBase::componentIndex");
328 }
329
330 //! Return the name of the solution component with index *i*.
331 //! @see componentIndex()
332 virtual string componentName(size_t k) {
333 throw NotImplementedError("ReactorBase::componentName");
334 }
335
336 //! Get the upper bound on the k-th component of the local state vector.
337 virtual double upperBound(size_t k) const {
338 throw NotImplementedError("ReactorBase::upperBound");
339 }
340
341 //! Get the lower bound on the k-th component of the local state vector.
342 virtual double lowerBound(size_t k) const {
343 throw NotImplementedError("ReactorBase::lowerBound");
344 }
345
346 //! Reset physically or mathematically problematic values, such as negative species
347 //! concentrations.
348 //! @param[inout] y current state vector, to be updated; length neq()
349 virtual void resetBadValues(double* y) {
350 throw NotImplementedError("ReactorBase::resetBadValues");
351 }
352
353 //! Get Jacobian elements for this reactor within the full reactor network.
354 //!
355 //! Indices within `trips` are global indices within the full reactor network. The
356 //! reactor is responsible for providing all elements of the Jacobian in the rows
357 //! corresponding to its state variables, that is, all derivatives of its state
358 //! variables with respect to all state variables in the network.
359 //!
360 //! @warning This method is an experimental part of the %Cantera API and may be
361 //! changed or removed without notice.
362 virtual void getJacobianElements(vector<Eigen::Triplet<double>>& trips) {};
363
364 //! Calculate the Jacobian of a specific reactor specialization.
365 //! @warning Depending on the particular implementation, this may return an
366 //! approximate Jacobian intended only for use in forming a preconditioner for
367 //! iterative solvers.
368 //! @ingroup derivGroup
369 //!
370 //! @warning This method is an experimental part of the %Cantera
371 //! API and may be changed or removed without notice.
372 virtual Eigen::SparseMatrix<double> jacobian();
373
374 //! Use this to set the kinetics objects derivative settings
375 virtual void setDerivativeSettings(AnyMap& settings) {
376 throw NotImplementedError("ReactorBase::setDerivativeSettings");
377 }
378
379 //! Set reaction rate multipliers based on the sensitivity variables in
380 //! *params*.
381 virtual void applySensitivity(double* params) {
382 throw NotImplementedError("ReactorBase::applySensitivity");
383 }
384
385 //! Reset the reaction rate multipliers
386 virtual void resetSensitivity(double* params) {
387 throw NotImplementedError("ReactorBase::resetSensitivity");
388 }
389
390 //! Return a false if preconditioning is not supported or true otherwise.
391 //!
392 //! @warning This method is an experimental part of the %Cantera
393 //! API and may be changed or removed without notice.
394 //!
395 //! @since New in %Cantera 3.0
396 //!
397 virtual bool preconditionerSupported() const { return false; };
398
399 //! @}
400
401 //! Set the state of the reactor to the associated ThermoPhase object.
402 //! This method will trigger integrator reinitialization.
403 //! @deprecated To be removed after %Cantera 4.0. Use ReactorNet::reinitialize to
404 //! indicate a change in state that requires integrator reinitialization.
405 virtual void syncState();
406
407 //! Update state information needed by connected reactors, flow devices, and walls.
408 //!
409 //! Called from updateState() for normal reactor types, and from
410 //! ReactorNet::updateState for Reservoir.
411 //!
412 //! @param updatePressure Indicates whether to update #m_pressure. Should
413 //! `true` for reactors where the pressure is a dependent property,
414 //! calculated from the state, and `false` when the pressure is constant
415 //! or an independent variable.
416 virtual void updateConnected(bool updatePressure);
417
418 //! Return the residence time (s) of the contents of this reactor, based
419 //! on the outlet mass flow rates and the mass of the reactor contents.
420 double residenceTime();
421
422 //! @name Solution components
423 //!
424 //! The values returned are those after the last call to ReactorNet::advance
425 //! or ReactorNet::step.
426 //! @{
427
428 //! Returns the current volume (m^3) of the reactor.
429 double volume() const {
430 return m_vol;
431 }
432
433 //! Returns the current density (kg/m^3) of the reactor's contents.
434 double density() const;
435
436 //! Returns the current temperature (K) of the reactor's contents.
437 double temperature() const;
438
439 //! Returns the current enthalpy (J/kg) of the reactor's contents.
440 double enthalpy_mass() const {
441 return m_enthalpy;
442 }
443
444 //! Returns the current pressure (Pa) of the reactor.
445 double pressure() const {
446 return m_pressure;
447 }
448
449 //! Returns the mass (kg) of the reactor's contents.
450 double mass() const {
451 return m_mass;
452 }
453
454 //! Return the vector of species mass fractions.
455 const double* massFractions() const;
456
457 //! Return the mass fraction of the *k*-th species.
458 double massFraction(size_t k) const;
459
460 //! @}
461
462 //! The ReactorNet that this reactor belongs to.
464
465 //! Set the ReactorNet that this reactor belongs to.
466 void setNetwork(ReactorNet* net);
467
468 //! Get the starting offset for this reactor's state variables within the global
469 //! state vector of the ReactorNet.
470 size_t offset() const { return m_offset; }
471
472 //! Set the starting offset for this reactor's state variables within the global
473 //! state vector of the ReactorNet.
474 void setOffset(size_t offset) { m_offset = offset; }
475
476 //! Offset of the first species in the local state vector
477 size_t speciesOffset() const {
478 return m_nv - m_nsp;
479 }
480
481 //! Get scaling factors for the Jacobian matrix terms proportional to
482 //! @f$ d\dot{n}_k/dC_j @f$.
483 //!
484 //! Used to determine contribution of surface phases to the Jacobian.
485 //!
486 //! @param f_species Scaling factor for derivatives appearing in the species
487 //! equations. Equal to $1/V$.
488 //! @param f_energy Scaling factor for each species term appearing in the energy
489 //! equation.
490 virtual void getJacobianScalingFactors(double& f_species, double* f_energy) {
491 throw NotImplementedError("ReactorBase::getJacobianScalingFactors");
492 }
493
494 //! Add a sensitivity parameter associated with the reaction number *rxn*
495 virtual void addSensitivityReaction(size_t rxn) {
496 throw NotImplementedError("ReactorBase::addSensitivityReaction");
497 }
498
499 //! Number of sensitivity parameters associated with this reactor.
500 virtual size_t nSensParams() const {
501 return m_sensParams.size();
502 }
503
504protected:
505 explicit ReactorBase(const string& name="(none)");
506
507 //! Number of homogeneous species in the mixture
508 size_t m_nsp = 0;
509
510 ThermoPhase* m_thermo = nullptr;
511 double m_vol = 0.0; //!< Current volume of the reactor [m^3]
512 double m_mass = 0.0; //!< Current mass of the reactor [kg]
513 double m_enthalpy = 0.0; //!< Current specific enthalpy of the reactor [J/kg]
514 double m_pressure = 0.0; //!< Current pressure in the reactor [Pa]
515 vector<double> m_state;
516 vector<FlowDevice*> m_inlet, m_outlet;
517
518 vector<WallBase*> m_wall;
519 vector<ReactorSurface*> m_surfaces;
520 vector<double> m_sdot; //!< species production rates on surfaces
521
522 //! Vector of length nWalls(), indicating whether this reactor is on the left (0)
523 //! or right (1) of each wall.
524 vector<int> m_lr;
525 string m_name; //!< Reactor name.
526 bool m_defaultNameSet = false; //!< `true` if default name has been previously set.
527 size_t m_nv = 0; //!< Number of state variables for this reactor
528
529 ReactorNet* m_net = nullptr; //!< The ReactorNet that this reactor is part of
530 size_t m_offset = 0; //!< Offset into global ReactorNet state vector
531
532 //! Composite thermo/kinetics/transport handler
533 shared_ptr<Solution> m_solution;
534
535 // Data associated each sensitivity parameter
536 vector<SensitivityParameter> m_sensParams;
537};
538}
539
540#endif
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
Base class for 'flow devices' (valves, pressure regulators, etc.) connecting reactors.
Definition FlowDevice.h:25
An error indicating that an unimplemented function has been called.
Base class for reactor objects.
Definition ReactorBase.h:51
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
double massFraction(size_t k) const
Return the mass fraction of the k-th species.
shared_ptr< Solution > m_solution
Composite thermo/kinetics/transport handler.
virtual void getJacobianScalingFactors(double &f_species, double *f_energy)
Get scaling factors for the Jacobian matrix terms proportional to .
virtual double upperBound(size_t k) const
Get the upper bound on the k-th component of the local state vector.
virtual void addSensitivitySpeciesEnthalpy(size_t k)
Add a sensitivity parameter associated with the enthalpy formation of species k.
virtual void evalWalls(double t)
Evaluate contributions from walls connected to this reactor.
virtual void getStateDae(double *y, double *ydot)
Get the current state and derivative vector of the reactor for a DAE solver.
bool m_defaultNameSet
true if default name has been previously set.
size_t nWalls()
Return the number of Wall objects connected to this reactor.
virtual void setArea(double a)
Set an area associated with a reactor [m²].
WallBase & wall(size_t n)
Return a reference to the n-th Wall connected to this reactor.
double density() const
Returns the current density (kg/m^3) of the reactor's contents.
double pressure() const
Returns the current pressure (Pa) of the reactor.
virtual void addOutlet(FlowDevice &outlet)
Connect an outlet FlowDevice to this reactor.
double m_pressure
Current pressure in the reactor [Pa].
virtual string type() const
String indicating the reactor model implemented.
Definition ReactorBase.h:76
ReactorNet * m_net
The ReactorNet that this reactor is part of.
virtual void addWall(WallBase &w, int lr)
Insert a Wall between this reactor and another reactor.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
size_t neq()
Number of equations (state variables) for this reactor.
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
virtual void resetBadValues(double *y)
Reset physically or mathematically problematic values, such as negative species concentrations.
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn
size_t m_nv
Number of state variables for this reactor.
virtual double area() const
Returns an area associated with a reactor [m²].
virtual size_t componentIndex(const string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
void setName(const string &name)
Set the name of this reactor.
Definition ReactorBase.h:86
virtual size_t nSurfs() const
Return the number of surfaces in a reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
double temperature() const
Returns the current temperature (K) of the reactor's contents.
virtual void evalSteady(double t, double *LHS, double *RHS)
Evaluate the governing equations with modifications for the steady-state solver.
virtual bool chemistryEnabled() const
Returns true if changes in the reactor composition due to chemical reactions are enabled.
virtual vector< size_t > initializeSteady()
Initialize the reactor before solving a steady-state problem.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
double m_vol
Current volume of the reactor [m^3].
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
size_t m_offset
Offset into global ReactorNet state vector.
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor.
void setOffset(size_t offset)
Set the starting offset for this reactor's state variables within the global state vector of the Reac...
virtual void addInlet(FlowDevice &inlet)
Connect an inlet FlowDevice to this reactor.
virtual string componentName(size_t k)
Return the name of the solution component with index i.
size_t speciesOffset() const
Offset of the first species in the local state vector.
span< const double > surfaceProductionRates() const
Production rates on surfaces.
virtual void syncState()
Set the state of the reactor to the associated ThermoPhase object.
virtual void evalDae(double t, double *y, double *ydot, double *residual)
Evaluate the reactor governing equations.
double m_mass
Current mass of the reactor [kg].
virtual void setEnergyEnabled(bool eflag=true)
Set the energy equation on or off.
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
vector< double > m_sdot
species production rates on surfaces
const double * massFractions() const
Return the vector of species mass fractions.
size_t m_nsp
Number of homogeneous species in the mixture.
string m_name
Reactor name.
double mass() const
Returns the mass (kg) of the reactor's contents.
virtual void setInitialVolume(double vol)
Set the initial reactor volume.
double volume() const
Returns the current volume (m^3) of the reactor.
ReactorNet & network()
The ReactorNet that this reactor belongs to.
double residenceTime()
Return the residence time (s) of the contents of this reactor, based on the outlet mass flow rates an...
size_t offset() const
Get the starting offset for this reactor's state variables within the global state vector of the Reac...
virtual double lowerBound(size_t k) const
Get the lower bound on the k-th component of the local state vector.
size_t nOutlets()
Return the number of outlet FlowDevice objects connected to this reactor.
virtual void setChemistryEnabled(bool cflag=true)
Enable or disable changes in reactor composition due to chemical reactions.
virtual void getConstraints(double *constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
virtual void initialize(double t0=0.0)
Initialize the reactor.
size_t nInlets()
Return the number of inlet FlowDevice objects connected to this reactor.
ReactorSurface * surface(size_t n)
Return a reference to the n-th ReactorSurface connected to this reactor.
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
double m_enthalpy
Current specific enthalpy of the reactor [J/kg].
virtual bool energyEnabled() const
Returns true if solution of the energy equation is enabled.
shared_ptr< Solution > phase()
Access the Solution object used to represent the contents of this reactor.
Definition ReactorBase.h:95
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
virtual void getState(double *y)
Get the current state of the reactor.
shared_ptr< const Solution > phase() const
Access the Solution object used to represent the contents of this reactor.
Definition ReactorBase.h:99
virtual void getJacobianElements(vector< Eigen::Triplet< double > > &trips)
Get Jacobian elements for this reactor within the full reactor network.
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
string name() const
Return the name of this reactor.
Definition ReactorBase.h:81
virtual void updateConnected(bool updatePressure)
Update state information needed by connected reactors, flow devices, and walls.
double enthalpy_mass() const
Returns the current enthalpy (J/kg) of the reactor's contents.
virtual void addSurface(ReactorSurface *surf)
Add a ReactorSurface object to a Reactor object.
A class representing a network of connected reactors.
Definition ReactorNet.h:30
A surface where reactions can occur that is in contact with the bulk fluid of a Reactor.
Base class for a phase with thermodynamic properties.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
Definition Wall.h:23
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
virtual Eigen::SparseMatrix< double > jacobian()
Calculate the Jacobian of a specific reactor specialization.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
size_t global
global parameter index
Definition ReactorBase.h:39
SensParameterType type
type of sensitivity parameter
Definition ReactorBase.h:41
size_t local
local parameter index
Definition ReactorBase.h:38
double value
nominal value of the parameter
Definition ReactorBase.h:40