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(span<double> y) {
246 throw NotImplementedError("ReactorBase::getState",
247 "Not implemented for reactor type '{}'.", type());
248 }
249
250 //! Set absolute tolerances for this reactor's state variables.
251 //!
252 //! @param atol User-specified absolute tolerances. The length must be equal
253 //! to neq(), and entries are ordered according to componentName().
254 //!
255 //! ReactorNet uses these values for this reactor's portion of the state vector
256 //! instead of the network-level scalar absolute tolerance or reactor-specific
257 //! default tolerances.
258 //! @since New in %Cantera 4.0.
259 void setAbsoluteTolerances(span<const double> atol);
260
261 //! Get the absolute tolerances set for this reactor's state variables.
262 //!
263 //! @param[out] atol User-specified absolute tolerances. The length must be
264 //! equal to neq().
265 //! @return `true` if user-specified absolute tolerances are set for this
266 //! reactor; `false` otherwise.
267 //! @since New in %Cantera 4.0.
268 bool getAbsoluteTolerances(span<double> atol) const;
269
270 //! Clear user-specified absolute tolerances for this reactor.
271 //! @since New in %Cantera 4.0.
273
274 //! Returns `true` if user-specified absolute tolerances are set for this reactor.
275 //! @since New in %Cantera 4.0.
276 bool hasUserTolerances() const {
277 return !m_userAtol.empty();
278 }
279
280 //! Update default absolute tolerances based on this reactor's state variables.
281 //!
282 //! @param[in,out] atol Mutable slice of ReactorNet's absolute tolerance
283 //! vector for this reactor. The length is neq(), and entries are ordered
284 //! according to componentName(). On entry, all values are set to baseAtol.
285 //! @param baseAtol The network-level scalar default absolute tolerance.
286 //!
287 //! ReactorNet calls this hook only when neither this reactor nor the network
288 //! has user-specified absolute tolerances.
289 //! @since New in %Cantera 4.0.
290 virtual void updateDefaultTolerances(span<double> atol, double baseAtol) {}
291
292 //! Get the current state and derivative vector of the reactor for a DAE solver
293 /*!
294 * @param[out] y state vector representing the initial state of the reactor
295 * @param[out] ydot state vector representing the initial derivatives of the
296 * reactor
297 */
298 virtual void getStateDae(span<double> y, span<double> ydot) {
299 throw NotImplementedError("ReactorBase::getStateDae(y, ydot)");
300 }
301
302 //! Evaluate the reactor governing equations. Called by ReactorNet::eval.
303 //! @param[in] t time.
304 //! @param[out] LHS pointer to start of vector of left-hand side
305 //! coefficients for governing equations, length m_nv, default values 1
306 //! @param[out] RHS pointer to start of vector of right-hand side
307 //! coefficients for governing equations, length m_nv, default values 0
308 virtual void eval(double t, span<double> LHS, span<double> RHS) {
309 throw NotImplementedError("ReactorBase::eval");
310 }
311
312 /**
313 * Evaluate the reactor governing equations. Called by ReactorNet::eval.
314 * @param[in] t time.
315 * @param[in] y solution vector, length neq()
316 * @param[in] ydot rate of change of solution vector, length neq()
317 * @param[out] residual residuals vector, length neq()
318 */
319 virtual void evalDae(double t, span<const double> y, span<const double> ydot,
320 span<double> residual) {
321 throw NotImplementedError("ReactorBase::evalDae");
322 }
323
324 //! Evaluate the governing equations with modifications for the steady-state solver.
325 //!
326 //! This method calls the standard eval() method then modifies elements of `RHS`
327 //! that correspond to algebraic constraints.
328 //!
329 //! @since New in %Cantera 4.0.
330 virtual void evalSteady(double t, span<double> LHS, span<double> RHS) {
331 throw NotImplementedError("ReactorBase::evalSteady",
332 "Not implemented for reactor type '{}'.", type());
333 }
334
335 //! Given a vector of length neq(), mark which variables should be
336 //! considered algebraic constraints
337 virtual void getConstraints(span<double> constraints) {
338 throw NotImplementedError("ReactorBase::getConstraints");
339 }
340
341 //! Initialize the reactor before solving a steady-state problem.
342 //!
343 //! This method is responsible for storing the initial value for any algebraic
344 //! constraints and returning the indices of those constraints.
345 //!
346 //! @return Indices of equations that are algebraic constraints when solving the
347 //! steady-state problem.
348 //!
349 //! @warning This method is an experimental part of the %Cantera API and may be
350 //! changed or removed without notice.
351 //! @since New in %Cantera 3.2. Renamed from `steadyConstraints` in %Cantera 4.0.
352 virtual vector<size_t> initializeSteady() {
353 throw NotImplementedError("ReactorBase::initializeSteady");
354 }
355
356 //! Set the state of the reactor to correspond to the state vector *y*.
357 virtual void updateState(span<const double> y) {
358 throw NotImplementedError("ReactorBase::updateState");
359 }
360
361 //! Add a sensitivity parameter associated with the enthalpy formation of
362 //! species *k*.
363 virtual void addSensitivitySpeciesEnthalpy(size_t k) {
364 throw NotImplementedError("ReactorBase::addSensitivitySpeciesEnthalpy");
365 }
366
367 //! Return the index in the solution vector for this reactor of the
368 //! component named *nm*.
369 virtual size_t componentIndex(const string& nm) const {
370 throw NotImplementedError("ReactorBase::componentIndex");
371 }
372
373 //! Return the name of the solution component with index *i*.
374 //! @see componentIndex()
375 virtual string componentName(size_t k) {
376 throw NotImplementedError("ReactorBase::componentName");
377 }
378
379 //! Get the upper bound on the k-th component of the local state vector.
380 virtual double upperBound(size_t k) const {
381 throw NotImplementedError("ReactorBase::upperBound");
382 }
383
384 //! Get the lower bound on the k-th component of the local state vector.
385 virtual double lowerBound(size_t k) const {
386 throw NotImplementedError("ReactorBase::lowerBound");
387 }
388
389 //! Reset physically or mathematically problematic values, such as negative species
390 //! concentrations.
391 //! @param[inout] y current state vector, to be updated; length neq()
392 virtual void resetBadValues(span<double> y) {
393 throw NotImplementedError("ReactorBase::resetBadValues");
394 }
395
396 //! Get Jacobian elements for this reactor within the full reactor network.
397 //!
398 //! Indices within `trips` are global indices within the full reactor network. The
399 //! reactor is responsible for providing all elements of the Jacobian in the rows
400 //! corresponding to its state variables, that is, all derivatives of its state
401 //! variables with respect to all state variables in the network.
402 //!
403 //! @warning This method is an experimental part of the %Cantera API and may be
404 //! changed or removed without notice.
405 virtual void getJacobianElements(SparseTriplets& trips) {};
406
407 //! Add terms proportional to derivatives with respect to this reactor's pressure.
408 //!
409 //! This method is used by connectors to apply the chain rule. For example, a
410 //! valve may know a row contribution as
411 //! `coeff * d(mdot)/d(P_in - P_out)`, while the reactor knows how its pressure
412 //! depends on its own state variables. Calling this method appends all entries
413 //! `coeff * dP/dy_j` for this reactor's state variables `y_j`.
414 //!
415 //! @param[in,out] trips Sparse Jacobian entries. Implementations append entries
416 //! using global row and column indices in the reactor network.
417 //! @param row Global row index receiving these chain-rule terms. This may be a
418 //! row belonging to another reactor when a connector creates cross-reactor
419 //! coupling.
420 //! @param coeff Multiplicative factor applied to `dP/dy_j`.
421 //! @param includeSpecies Include pressure derivatives with respect to species
422 //! state variables. These terms may be dense for mole reactors and are
423 //! controlled by preconditioner sparsity settings.
424 //! @since New in %Cantera 4.0.
425 virtual void addPressureJacobian(SparseTriplets& trips, size_t row, double coeff,
426 bool includeSpecies=true) const {};
427
428 //! Add terms proportional to derivatives with respect to this reactor's
429 //! temperature.
430 //!
431 //! This method is used by connectors to apply the chain rule. Calling this method
432 //! appends all entries `coeff * dT/dy_j` for this reactor's state variables `y_j`.
433 //!
434 //! @param[in,out] trips Sparse Jacobian entries. Implementations append entries
435 //! using global row and column indices in the reactor network.
436 //! @param row Global row index receiving these chain-rule terms. This may be a
437 //! row belonging to another reactor when a connector creates cross-reactor
438 //! coupling.
439 //! @param coeff Multiplicative factor applied to `dT/dy_j`.
440 //! @since New in %Cantera 4.0.
441 virtual void addTemperatureJacobian(SparseTriplets& trips, size_t row,
442 double coeff) const {};
443
444 //! Add terms proportional to derivatives with respect to a species mass fraction.
445 //!
446 //! Flow devices transport species according to upstream mass fractions. This
447 //! method lets a connector request chain-rule entries for `coeff * dY_k/dy_j`
448 //! while the reactor implementation handles how `Y_k` depends on its native state
449 //! variables, such as species moles.
450 //!
451 //! @param[in,out] trips Sparse Jacobian entries. Implementations append entries
452 //! using global row and column indices in the reactor network.
453 //! @param row Global row index receiving these chain-rule terms. This may be a
454 //! row belonging to another reactor when a connector creates cross-reactor
455 //! coupling.
456 //! @param k Species index in this reactor's phase.
457 //! @param coeff Multiplicative factor applied to `dY_k/dy_j`.
458 //! @since New in %Cantera 4.0.
459 virtual void addSpeciesMassFractionJacobian(SparseTriplets& trips, size_t row,
460 size_t k, double coeff) const {};
461
462 //! Add terms proportional to derivatives of this reactor's specific enthalpy.
463 //!
464 //! This method is used by flow device connectors to apply the chain rule for inlet
465 //! enthalpy contributions. It appends entries `coeff * dh/dy_j` for this reactor's
466 //! state variables `y_j`. The temperature contribution is `coeff * cp_mass`; the
467 //! optional composition contribution adds one entry per species.
468 //!
469 //! @param[in,out] trips Sparse Jacobian entries. Implementations append entries
470 //! using global row and column indices in the reactor network.
471 //! @param row Global row index receiving these chain-rule terms. This may be a
472 //! row belonging to another reactor when the inlet energy source is cross-reactor.
473 //! @param coeff Multiplicative factor applied to `dh/dy_j`.
474 //! @param includeComposition Include derivatives of specific enthalpy with respect
475 //! to species state variables. These terms may be dense and are controlled by
476 //! preconditioner sparsity settings.
477 //! @since New in %Cantera 4.0.
478 virtual void addEnthalpyJacobian(SparseTriplets& trips, size_t row, double coeff,
479 bool includeComposition=true) const {};
480
481 //! Calculate the Jacobian of a specific reactor specialization.
482 //! @warning Depending on the particular implementation, this may return an
483 //! approximate Jacobian intended only for use in forming a preconditioner for
484 //! iterative solvers.
485 //! @ingroup derivGroup
486 //!
487 //! @warning This method is an experimental part of the %Cantera
488 //! API and may be changed or removed without notice.
489 virtual Eigen::SparseMatrix<double> jacobian();
490
491 //! Use this to set the kinetics objects derivative settings
492 virtual void setDerivativeSettings(AnyMap& settings) {
493 throw NotImplementedError("ReactorBase::setDerivativeSettings");
494 }
495
496 //! Set reaction rate multipliers based on the sensitivity variables in
497 //! *params*.
498 virtual void applySensitivity(span<const double> params) {
499 throw NotImplementedError("ReactorBase::applySensitivity");
500 }
501
502 //! Reset the reaction rate multipliers
503 virtual void resetSensitivity(span<const double> params) {
504 throw NotImplementedError("ReactorBase::resetSensitivity");
505 }
506
507 //! Return a false if preconditioning is not supported or true otherwise.
508 //!
509 //! @warning This method is an experimental part of the %Cantera
510 //! API and may be changed or removed without notice.
511 //!
512 //! @since New in %Cantera 3.0
513 //!
514 virtual bool preconditionerSupported() const { return false; };
515
516 //! @}
517
518 //! Set the state of the reactor to the associated ThermoPhase object.
519 //! This method will trigger integrator reinitialization.
520 //! @deprecated To be removed after %Cantera 4.0. Use ReactorNet::reinitialize to
521 //! indicate a change in state that requires integrator reinitialization.
522 virtual void syncState();
523
524 //! Update state information needed by connected reactors, flow devices, and walls.
525 //!
526 //! Called from updateState() for normal reactor types, and from
527 //! ReactorNet::updateState for Reservoir.
528 //!
529 //! @param updatePressure Indicates whether to update #m_pressure. Should
530 //! `true` for reactors where the pressure is a dependent property,
531 //! calculated from the state, and `false` when the pressure is constant
532 //! or an independent variable.
533 virtual void updateConnected(bool updatePressure);
534
535 //! Return the residence time (s) of the contents of this reactor, based
536 //! on the outlet mass flow rates and the mass of the reactor contents.
537 double residenceTime();
538
539 //! @name Solution components
540 //!
541 //! The values returned are those after the last call to ReactorNet::advance
542 //! or ReactorNet::step.
543 //! @{
544
545 //! Returns the current volume (m^3) of the reactor.
546 double volume() const {
547 return m_vol;
548 }
549
550 //! Returns the current density (kg/m^3) of the reactor's contents.
551 double density() const;
552
553 //! Returns the current temperature (K) of the reactor's contents.
554 double temperature() const;
555
556 //! Returns the current enthalpy (J/kg) of the reactor's contents.
557 double enthalpy_mass() const {
558 return m_enthalpy;
559 }
560
561 //! Returns the current pressure (Pa) of the reactor.
562 double pressure() const {
563 return m_pressure;
564 }
565
566 //! Returns the mass (kg) of the reactor's contents.
567 double mass() const {
568 return m_mass;
569 }
570
571 //! Return the vector of species mass fractions.
572 span<const double> massFractions() const;
573
574 //! Return the mass fraction of the *k*-th species.
575 double massFraction(size_t k) const;
576
577 //! @}
578
579 //! The ReactorNet that this reactor belongs to.
581
582 //! Set the ReactorNet that this reactor belongs to.
583 void setNetwork(ReactorNet* net);
584
585 //! Get the starting offset for this reactor's state variables within the global
586 //! state vector of the ReactorNet.
587 size_t offset() const { return m_offset; }
588
589 //! Set the starting offset for this reactor's state variables within the global
590 //! state vector of the ReactorNet.
591 void setOffset(size_t offset) { m_offset = offset; }
592
593 //! Offset of the first species in the local state vector
594 size_t speciesOffset() const {
595 return m_nv - m_nsp;
596 }
597
598 //! Get scaling factors for the Jacobian matrix terms proportional to
599 //! @f$ d\dot{n}_k/dC_j @f$.
600 //!
601 //! Used to determine contribution of surface phases to the Jacobian.
602 //!
603 //! @param f_species Scaling factor for derivatives appearing in the species
604 //! equations. Equal to $1/V$.
605 //! @param f_energy Scaling factor for each species term appearing in the energy
606 //! equation.
607 virtual void getJacobianScalingFactors(double& f_species, span<double> f_energy) {
608 throw NotImplementedError("ReactorBase::getJacobianScalingFactors");
609 }
610
611 //! Add a sensitivity parameter associated with the reaction number *rxn*
612 virtual void addSensitivityReaction(size_t rxn) {
613 throw NotImplementedError("ReactorBase::addSensitivityReaction");
614 }
615
616 //! Number of sensitivity parameters associated with this reactor.
617 virtual size_t nSensParams() const {
618 return m_sensParams.size();
619 }
620
621protected:
622 explicit ReactorBase(const string& name="(none)");
623
624 //! Check whether a global Jacobian row belongs to this reactor.
625 //!
626 //! Used by Jacobian helper methods when local rows are already covered by the
627 //! reactor's own block approximation, while cross-reactor rows need explicit
628 //! connector entries.
629 bool isJacobianLocalRow(size_t row) const {
630 return row >= m_offset && row < m_offset + m_nv;
631 }
632
633 //! Number of homogeneous species in the mixture
634 size_t m_nsp = 0;
635
636 ThermoPhase* m_thermo = nullptr;
637 double m_vol = 0.0; //!< Current volume of the reactor [m^3]
638 double m_mass = 0.0; //!< Current mass of the reactor [kg]
639 double m_enthalpy = 0.0; //!< Current specific enthalpy of the reactor [J/kg]
640 double m_pressure = 0.0; //!< Current pressure in the reactor [Pa]
641 vector<double> m_state;
642 vector<FlowDevice*> m_inlet, m_outlet;
643
644 vector<WallBase*> m_wall;
645 vector<ReactorSurface*> m_surfaces;
646 vector<double> m_sdot; //!< species production rates on surfaces
647
648 //! Vector of length nWalls(), indicating whether this reactor is on the left (0)
649 //! or right (1) of each wall.
650 vector<int> m_lr;
651 string m_name; //!< Reactor name.
652 bool m_defaultNameSet = false; //!< `true` if default name has been previously set.
653 size_t m_nv = 0; //!< Number of state variables for this reactor
654
655 ReactorNet* m_net = nullptr; //!< The ReactorNet that this reactor is part of
656 size_t m_offset = 0; //!< Offset into global ReactorNet state vector
657
658 //! User-specified absolute tolerances for this reactor's state variables.
659 vector<double> m_userAtol;
660
661 //! Composite thermo/kinetics/transport handler
662 shared_ptr<Solution> m_solution;
663
664 // Data associated each sensitivity parameter
665 vector<SensitivityParameter> m_sensParams;
666};
667}
668
669#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:26
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 double upperBound(size_t k) const
Get the upper bound on the k-th component of the local state vector.
virtual void getStateDae(span< double > y, span< double > ydot)
Get the current state and derivative vector of the reactor for a DAE solver.
virtual void getJacobianScalingFactors(double &f_species, span< double > f_energy)
Get scaling factors for the Jacobian matrix terms proportional to .
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.
bool m_defaultNameSet
true if default name has been previously set.
size_t nWalls()
Return the number of Wall objects connected to this reactor.
bool getAbsoluteTolerances(span< double > atol) const
Get the absolute tolerances set for this reactor's state variables.
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.
virtual void evalDae(double t, span< const double > y, span< const double > ydot, span< double > residual)
Evaluate the reactor governing equations.
double pressure() const
Returns the current pressure (Pa) of the reactor.
virtual void addOutlet(FlowDevice &outlet)
Connect an outlet FlowDevice to this reactor.
virtual void resetBadValues(span< double > y)
Reset physically or mathematically problematic values, such as negative species concentrations.
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.
bool isJacobianLocalRow(size_t row) const
Check whether a global Jacobian row belongs to this reactor.
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.
virtual void eval(double t, span< double > LHS, span< double > RHS)
Evaluate the reactor governing equations.
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 addEnthalpyJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true) const
Add terms proportional to derivatives of this reactor's specific enthalpy.
virtual void updateState(span< const double > y)
Set the state of the reactor to correspond to the state vector y.
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.
void setAbsoluteTolerances(span< const double > atol)
Set absolute tolerances for this reactor's state variables.
double temperature() const
Returns the current temperature (K) of the reactor's contents.
virtual bool chemistryEnabled() const
Returns true if changes in the reactor composition due to chemical reactions are enabled.
virtual void evalSteady(double t, span< double > LHS, span< double > RHS)
Evaluate the governing equations with modifications for the steady-state solver.
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].
size_t m_offset
Offset into global ReactorNet state vector.
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor.
span< const double > massFractions() const
Return the vector of species mass fractions.
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.
virtual void getJacobianElements(SparseTriplets &trips)
Get Jacobian elements for this reactor within the full reactor network.
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 updateDefaultTolerances(span< double > atol, double baseAtol)
Update default absolute tolerances based on this reactor's state variables.
virtual void getState(span< double > y)
Get the current state of the reactor.
vector< double > m_userAtol
User-specified absolute tolerances for this reactor's state variables.
double m_mass
Current mass of the reactor [kg].
virtual void setEnergyEnabled(bool eflag=true)
Set the energy equation on or off.
vector< double > m_sdot
species production rates on surfaces
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.
virtual void addTemperatureJacobian(SparseTriplets &trips, size_t row, double coeff) const
Add terms proportional to derivatives with respect to this reactor's temperature.
bool hasUserTolerances() const
Returns true if user-specified absolute tolerances are set for this reactor.
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.
virtual void getConstraints(span< double > constraints)
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
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 initialize(double t0=0.0)
Initialize the reactor.
virtual void applySensitivity(span< const double > params)
Set reaction rate multipliers based on the sensitivity variables in params.
size_t nInlets()
Return the number of inlet FlowDevice objects connected to this reactor.
virtual void addPressureJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeSpecies=true) const
Add terms proportional to derivatives with respect to this reactor's pressure.
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 addSpeciesMassFractionJacobian(SparseTriplets &trips, size_t row, size_t k, double coeff) const
Add terms proportional to derivatives with respect to a species mass fraction.
void clearAbsoluteTolerances()
Clear user-specified absolute tolerances for this reactor.
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.
shared_ptr< const Solution > phase() const
Access the Solution object used to represent the contents of this reactor.
Definition ReactorBase.h:99
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 resetSensitivity(span< const double > params)
Reset the reaction rate multipliers.
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:24
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