Cantera 2.6.0
ReactorNet.h
Go to the documentation of this file.
1//! @file ReactorNet.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_REACTORNET_H
7#define CT_REACTORNET_H
8
9#include "Reactor.h"
11
12namespace Cantera
13{
14
15class Array2D;
16class Integrator;
17
18//! A class representing a network of connected reactors.
19/*!
20 * This class is used to integrate the time-dependent governing equations for
21 * a network of reactors (Reactor, ConstPressureReactor) connected by various
22 * means, for example Wall, MassFlowController, Valve, or PressureController.
23 *
24 * @ingroup ZeroD
25 */
26class ReactorNet : public FuncEval
27{
28public:
29 ReactorNet();
30 virtual ~ReactorNet();
31 ReactorNet(const ReactorNet&) = delete;
32 ReactorNet& operator=(const ReactorNet&) = delete;
33
34 //! @name Methods to set up a simulation.
35 //! @{
36
37 //! Set initial time. Default = 0.0 s. Restarts integration from this time
38 //! using the current mixture state as the initial condition.
39 void setInitialTime(double time);
40
41 //! Get the maximum time step.
42 double maxTimeStep() {
43 return m_maxstep;
44 }
45
46 //! Set the maximum time step.
47 void setMaxTimeStep(double maxstep);
48
49 //! Set the maximum number of error test failures permitted by the CVODES
50 //! integrator in a single time step.
51 void setMaxErrTestFails(int nmax);
52
53 //! Set the relative and absolute tolerances for the integrator.
54 void setTolerances(double rtol, double atol);
55
56 //! Set the relative and absolute tolerances for integrating the
57 //! sensitivity equations.
58 void setSensitivityTolerances(double rtol, double atol);
59
60 //! @}
61
62 //! Current value of the simulation time.
63 doublereal time() {
64 return m_time;
65 }
66
67 //! Relative tolerance.
68 doublereal rtol() {
69 return m_rtol;
70 }
71
72 //! Absolute integration tolerance
73 doublereal atol() {
74 return m_atols;
75 }
76
77 //! Relative sensitivity tolerance
78 doublereal rtolSensitivity() const {
79 return m_rtolsens;
80 }
81
82 //! Absolute sensitivity tolerance
83 doublereal atolSensitivity() const {
84 return m_atolsens;
85 }
86
87 /**
88 * Advance the state of all reactors in time. Take as many internal
89 * timesteps as necessary to reach *time*.
90 * @param time Time to advance to (s).
91 */
92 void advance(doublereal time);
93
94 /**
95 * Advance the state of all reactors in time. Take as many internal
96 * timesteps as necessary towards *time*. If *applylimit* is true,
97 * the advance step will be automatically reduced if needed to
98 * stay within limits (set by setAdvanceLimit).
99 * Returns the time at the end of integration.
100 * @param time Time to advance to (s).
101 * @param applylimit Limit advance step (boolean).
102 */
103 double advance(double time, bool applylimit);
104
105 //! Advance the state of all reactors in time.
106 double step();
107
108 //! Add the reactor *r* to this reactor network.
109 void addReactor(Reactor& r);
110
111 //! Return a reference to the *n*-th reactor in this network. The reactor
112 //! indices are determined by the order in which the reactors were added
113 //! to the reactor network.
114 Reactor& reactor(int n) {
115 return *m_reactors[n];
116 }
117
118 //! Returns `true` if verbose logging output is enabled.
119 bool verbose() const {
120 return m_verbose;
121 }
122
123 //! Enable or disable verbose logging while setting up and integrating the
124 //! reactor network.
125 void setVerbose(bool v = true) {
126 m_verbose = v;
127 suppressErrors(!m_verbose);
128 }
129
130 //! Return a reference to the integrator.
132 return *m_integ;
133 }
134
135 //! Update the state of all the reactors in the network to correspond to
136 //! the values in the solution vector *y*.
137 void updateState(doublereal* y);
138
139 //! Return the sensitivity of the *k*-th solution component with respect to
140 //! the *p*-th sensitivity parameter.
141 /*!
142 * The sensitivity coefficient \f$ S_{ki} \f$ of solution variable \f$ y_k
143 * \f$ with respect to sensitivity parameter \f$ p_i \f$ is defined as:
144 *
145 * \f[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} \f]
146 *
147 * For reaction sensitivities, the parameter is a multiplier on the forward
148 * rate constant (and implicitly on the reverse rate constant for
149 * reversible reactions) which has a nominal value of 1.0, and the
150 * sensitivity is nondimensional.
151 *
152 * For species enthalpy sensitivities, the parameter is a perturbation to
153 * the molar enthalpy of formation, such that the dimensions of the
154 * sensitivity are kmol/J.
155 */
156 double sensitivity(size_t k, size_t p);
157
158 //! Return the sensitivity of the component named *component* with respect to
159 //! the *p*-th sensitivity parameter.
160 //! @copydetails ReactorNet::sensitivity(size_t, size_t)
161 double sensitivity(const std::string& component, size_t p, int reactor=0) {
162 size_t k = globalComponentIndex(component, reactor);
163 return sensitivity(k, p);
164 }
165
166 //! Evaluate the Jacobian matrix for the reactor network.
167 /*!
168 * @param[in] t Time at which to evaluate the Jacobian
169 * @param[in] y Global state vector at time *t*
170 * @param[out] ydot Time derivative of the state vector evaluated at *t*.
171 * @param[in] p sensitivity parameter vector (unused?)
172 * @param[out] j Jacobian matrix, size neq() by neq().
173 */
174 void evalJacobian(doublereal t, doublereal* y,
175 doublereal* ydot, doublereal* p, Array2D* j);
176
177 // overloaded methods of class FuncEval
178 virtual size_t neq() {
179 return m_nv;
180 }
181 virtual void eval(doublereal t, doublereal* y,
182 doublereal* ydot, doublereal* p);
183
184 virtual void getState(doublereal* y);
185
186 //! Return k-th derivative at the current time
187 virtual void getDerivative(int k, double* dky);
188
189 virtual size_t nparams() {
190 return m_sens_params.size();
191 }
192
193 //! Return the index corresponding to the component named *component* in the
194 //! reactor with index *reactor* in the global state vector for the
195 //! reactor network.
196 size_t globalComponentIndex(const std::string& component, size_t reactor=0);
197
198 //! Return the name of the i-th component of the global state vector. The
199 //! name returned includes both the name of the reactor and the specific
200 //! component, for example `'reactor1: CH4'`.
201 std::string componentName(size_t i) const;
202
203 //! Used by Reactor and Wall objects to register the addition of
204 //! sensitivity parameters so that the ReactorNet can keep track of the
205 //! order in which sensitivity parameters are added.
206 //! @param name A name describing the parameter, for example the reaction string
207 //! @param value The nominal value of the parameter
208 //! @param scale A scaling factor to be applied to the sensitivity
209 //! coefficient
210 //! @returns the index of this parameter in the vector of sensitivity
211 //! parameters (global across all reactors)
212 size_t registerSensitivityParameter(const std::string& name, double value,
213 double scale);
214
215 //! The name of the p-th sensitivity parameter added to this ReactorNet.
216 const std::string& sensitivityParameterName(size_t p) {
217 return m_paramNames.at(p);
218 }
219
220 //! Initialize the reactor network. Called automatically the first time
221 //! advance or step is called.
222 void initialize();
223
224 //! Reinitialize the integrator. Used to solve a new problem (different
225 //! initial conditions) but with the same configuration of the reactor
226 //! network. Can be called manually, or automatically after calling
227 //! setInitialTime or modifying a reactor's contents.
228 void reinitialize();
229
230 //! Called to trigger integrator reinitialization before further
231 //! integration.
233 m_integrator_init = false;
234 }
235
236 //! Set the maximum number of internal integration time-steps the
237 //! integrator will take before reaching the next output time
238 //! @param nmax The maximum number of steps, setting this value
239 //! to zero disables this option.
240 virtual void setMaxSteps(int nmax);
241
242 //! Returns the maximum number of internal integration time-steps the
243 //! integrator will take before reaching the next output time
244 //!
245 virtual int maxSteps();
246
247 //! Set absolute step size limits during advance
248 void setAdvanceLimits(const double* limits);
249
250 //! Check whether ReactorNet object uses advance limits
251 bool hasAdvanceLimits();
252
253 //! Retrieve absolute step size limits during advance
254 bool getAdvanceLimits(double* limits);
255
256protected:
257
258 //! Estimate a future state based on current derivatives.
259 //! The function is intended for internal use by ReactorNet::advance
260 //! and deliberately not exposed in external interfaces.
261 virtual void getEstimate(double time, int k, double* yest);
262
263 //! Returns the order used for last solution step of the ODE integrator
264 //! The function is intended for internal use by ReactorNet::advance
265 //! and deliberately not exposed in external interfaces.
266 virtual int lastOrder();
267
268 std::vector<Reactor*> m_reactors;
269 std::unique_ptr<Integrator> m_integ;
270 doublereal m_time;
271 bool m_init;
272 bool m_integrator_init; //!< True if integrator initialization is current
273 size_t m_nv;
274
275 //! m_start[n] is the starting point in the state vector for reactor n
276 std::vector<size_t> m_start;
277
278 vector_fp m_atol;
279 doublereal m_rtol, m_rtolsens;
280 doublereal m_atols, m_atolsens;
281
282 //! Maximum integrator internal timestep. Default of 0.0 means infinity.
283 doublereal m_maxstep;
284
285 int m_maxErrTestFails;
286 bool m_verbose;
287
288 //! Names corresponding to each sensitivity parameter
289 std::vector<std::string> m_paramNames;
290
291 vector_fp m_ydot;
292 vector_fp m_yest;
293 vector_fp m_advancelimits;
294 //! m_LHS is a vector representing the coefficients on the
295 //! "left hand side" of each governing equation
297 vector_fp m_RHS;
298 bool m_checked_eval_deprecation; //!< @todo Remove after Cantera 2.6
299 std::vector<bool> m_have_deprecated_eval; //!< @todo Remove after Cantera 2.6
300};
301}
302
303#endif
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:30
Virtual base class for ODE right-hand-side function evaluators.
Definition: FuncEval.h:27
vector_fp m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition: FuncEval.h:81
bool suppressErrors() const
Get current state of error suppression.
Definition: FuncEval.h:71
Abstract base class for ODE system integrators.
Definition: Integrator.h:53
A class representing a network of connected reactors.
Definition: ReactorNet.h:27
std::vector< bool > m_have_deprecated_eval
Definition: ReactorNet.h:299
virtual size_t nparams()
Number of sensitivity parameters.
Definition: ReactorNet.h:189
void updateState(doublereal *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
Definition: ReactorNet.cpp:348
double step()
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:210
std::vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
Definition: ReactorNet.h:276
void initialize()
Initialize the reactor network.
Definition: ReactorNet.cpp:78
void advance(doublereal time)
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:143
bool getAdvanceLimits(double *limits)
Retrieve absolute step size limits during advance.
Definition: ReactorNet.cpp:390
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition: ReactorNet.h:289
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
Definition: ReactorNet.h:232
const std::string & sensitivityParameterName(size_t p)
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition: ReactorNet.h:216
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition: ReactorNet.h:114
doublereal rtol()
Relative tolerance.
Definition: ReactorNet.h:68
doublereal m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition: ReactorNet.h:283
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration time-steps the integrator will take before reaching th...
Definition: ReactorNet.cpp:133
std::string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
Definition: ReactorNet.cpp:407
double maxTimeStep()
Get the maximum time step.
Definition: ReactorNet.h:42
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
Definition: ReactorNet.cpp:247
virtual void getState(doublereal *y)
Fill in the vector y with the current state of the system.
Definition: ReactorNet.cpp:356
virtual size_t neq()
Number of equations.
Definition: ReactorNet.h:178
void setInitialTime(double time)
Set initial time.
Definition: ReactorNet.cpp:38
virtual void getDerivative(int k, double *dky)
Return k-th derivative at the current time.
Definition: ReactorNet.cpp:363
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single time ste...
Definition: ReactorNet.cpp:50
size_t globalComponentIndex(const std::string &component, size_t reactor=0)
Return the index corresponding to the component named component in the reactor with index reactor in ...
Definition: ReactorNet.cpp:399
doublereal atolSensitivity() const
Absolute sensitivity tolerance.
Definition: ReactorNet.h:83
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
Definition: ReactorNet.cpp:67
virtual int maxSteps()
Returns the maximum number of internal integration time-steps the integrator will take before reachin...
Definition: ReactorNet.cpp:138
virtual int lastOrder()
Returns the order used for last solution step of the ODE integrator The function is intended for inte...
Definition: ReactorNet.cpp:242
doublereal rtolSensitivity() const
Relative sensitivity tolerance.
Definition: ReactorNet.h:78
double sensitivity(size_t k, size_t p)
Return the sensitivity of the k-th solution component with respect to the p-th sensitivity parameter.
Definition: ReactorNet.cpp:309
bool hasAdvanceLimits()
Check whether ReactorNet object uses advance limits.
Definition: ReactorNet.cpp:381
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
Definition: ReactorNet.cpp:371
void evalJacobian(doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
Definition: ReactorNet.cpp:325
doublereal time()
Current value of the simulation time.
Definition: ReactorNet.h:63
bool verbose() const
Returns true if verbose logging output is enabled.
Definition: ReactorNet.h:119
bool m_checked_eval_deprecation
Definition: ReactorNet.h:298
double sensitivity(const std::string &component, size_t p, int reactor=0)
Return the sensitivity of the component named component with respect to the p-th sensitivity paramete...
Definition: ReactorNet.h:161
void setMaxTimeStep(double maxstep)
Set the maximum time step.
Definition: ReactorNet.cpp:44
size_t registerSensitivityParameter(const std::string &name, double value, double scale)
Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the React...
Definition: ReactorNet.cpp:419
bool m_integrator_init
True if integrator initialization is current.
Definition: ReactorNet.h:272
void reinitialize()
Reinitialize the integrator.
Definition: ReactorNet.cpp:122
Integrator & integrator()
Return a reference to the integrator.
Definition: ReactorNet.h:131
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition: ReactorNet.h:125
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
Definition: ReactorNet.cpp:253
doublereal atol()
Absolute integration tolerance.
Definition: ReactorNet.h:73
vector_fp m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition: ReactorNet.h:296
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
Definition: ReactorNet.cpp:56
virtual void getEstimate(double time, int k, double *yest)
Estimate a future state based on current derivatives.
Definition: ReactorNet.cpp:222
Class Reactor is a general-purpose class for stirred reactors.
Definition: Reactor.h:41
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition: utilities.h:100
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184