Cantera  4.0.0a1
Loading...
Searching...
No Matches
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"
12
13namespace Cantera
14{
15
16class Array2D;
17class Integrator;
18class SystemJacobian;
19
20//! A class representing a network of connected reactors.
21/*!
22 * This class is used to integrate the governing equations for a network of reactors
23 * that are time dependent (Reactor, ConstPressureReactor) connected by various
24 * means, for example Wall, MassFlowController, Valve, or PressureController; or
25 * reactors dependent on a single spatial variable (FlowReactor).
26 *
27 * @ingroup zerodGroup
28 */
29class ReactorNet : public FuncEval
30{
31public:
32 //! Create reactor network containing single reactor.
33 //! @since New in %Cantera 3.2.
34 ReactorNet(shared_ptr<ReactorBase> reactor);
35 //! Create reactor network from multiple reactors.
36 //! @since New in %Cantera 3.2.
37 ReactorNet(span<shared_ptr<ReactorBase>> reactors);
38 ReactorNet(const ReactorNet&) = delete;
39 ReactorNet& operator=(const ReactorNet&) = delete;
41
42 //! @name Methods to set up a simulation
43 //! @{
44
45 //! Set the type of linear solver used in the integration.
46 //! @param linSolverType type of linear solver. Default type: "DENSE"
47 //! Other options include: "DIAG", "DENSE", "GMRES", "BAND"
48 void setLinearSolverType(const string& linSolverType="DENSE");
49
50 //! Set preconditioner used by the linear solver
51 //! @param preconditioner preconditioner object used for the linear solver
52 void setPreconditioner(shared_ptr<SystemJacobian> preconditioner);
53
54 //! Set the initial value of the independent variable (typically time).
55 //! Default = 0.0 s. Restarts integration from this value using the current mixture
56 //! state as the initial condition.
57 void setInitialTime(double time);
58
59 //! Get the initial value of the independent variable (typically time).
60 /*!
61 * @since New in %Cantera 3.0.
62 */
63 double getInitialTime() const {
64 return m_initial_time;
65 }
66
67 //! Get the maximum integrator step.
68 double maxTimeStep() const {
69 return m_maxstep;
70 }
71
72 //! Set the maximum integrator step.
73 void setMaxTimeStep(double maxstep);
74
75 //! Set the maximum number of error test failures permitted by the CVODES
76 //! integrator in a single step.
77 void setMaxErrTestFails(int nmax);
78
79 //! Set the relative tolerance for the integrator.
80 //!
81 //! @param rtol Relative tolerance for all state variables; must be positive.
82 //! @since New in %Cantera 4.0.
83 void setRelativeTolerance(double rtol);
84
85 //! Set the scalar absolute tolerance for the integrator.
86 //!
87 //! @param atol Scalar absolute tolerance for state variables without
88 //! reactor-specific absolute tolerances; must be positive. Superseded by
89 //! absolute tolerances set for individual reactors.
90 //! @see ReactorBase::setAbsoluteTolerances
91 //! @since New in %Cantera 4.0.
92 void setAbsoluteTolerance(double atol);
93
94 //! Clear the user-specified scalar absolute tolerance.
95 //!
96 //! Restores the default scalar absolute tolerance and enables
97 //! reactor-specific default scaling rules.
98 //! @see ReactorBase::updateDefaultTolerances
99 //! @since New in %Cantera 4.0.
101
102 //! Set the relative and scalar absolute tolerances for the integrator.
103 //!
104 //! @param rtol If positive, set the relative tolerance for all state variables;
105 //! ignored if negative.
106 //! @param atol If positive, set the scalar absolute tolerance for all state
107 //! variables; ignored if negative. Superseded by absolute tolerances set
108 //! for individual reactors.
109 //! @deprecated To be removed after %Cantera 4.0. Use setRelativeTolerance()
110 //! and setAbsoluteTolerance() instead.
111 void setTolerances(double rtol, double atol);
112
113 //! Set the relative and absolute tolerances for integrating the
114 //! sensitivity equations.
115 void setSensitivityTolerances(double rtol, double atol);
116
117 //! Current value of the simulation time [s], for reactor networks that are solved
118 //! in the time domain.
119 double time();
120
121 //! Current position [m] along the length of the reactor network, for reactors that
122 //! are solved as a function of space.
123 double distance();
124
125 //! Relative tolerance.
126 double rtol() {
127 return m_rtol;
128 }
129
130 //! Scalar absolute integration tolerance
131 double atol() {
132 return m_atols;
133 }
134
135 //! Relative sensitivity tolerance
136 double rtolSensitivity() const {
137 return m_rtolsens;
138 }
139
140 //! Absolute sensitivity tolerance
141 double atolSensitivity() const {
142 return m_atolsens;
143 }
144
145 //! Problem type of integrator
146 string linearSolverType() const;
147
148 //! Returns the maximum number of internal integration steps the
149 //! integrator will take before reaching the next output point
150 int maxSteps();
151
152 //! @}
153
154 /**
155 * Advance the state of all reactors in the independent variable (time or space).
156 * Take as many internal steps as necessary to reach *t*.
157 * @param t Time/distance to advance to (s or m).
158 */
159 void advance(double t);
160
161 /**
162 * Advance the state of all reactors in the independent variable (time or space).
163 * Take as many internal steps as necessary towards *t*. If *applylimit* is true,
164 * the advance step will be automatically reduced if needed to stay within limits
165 * (set by setAdvanceLimit).
166 * Returns the time/distance at the end of integration.
167 * @param t Time/distance to advance to (s or m).
168 * @param applylimit Limit advance step (boolean).
169 */
170 double advance(double t, bool applylimit);
171
172 //! Advance the state of all reactors with respect to the independent variable
173 //! (time or space). Returns the new value of the independent variable [s or m].
174 double step();
175
176 //! Solve directly for the steady-state solution.
177 //!
178 //! This approach is generally more efficient than time marching to the
179 //! steady-state, but imposes a few limitations:
180 //!
181 //! - The volume of control volume reactor types (such as Reactor and
182 //! IdealGasMoleReactor) must be constant; no moving walls can be used.
183 //! - The mass of constant pressure reactor types (such as ConstPressureReactor and
184 //! IdealGasConstPressureReactor) must be constant; if flow devices are used,
185 //! inlet and outlet flows must be balanced.
186 //! - Only ideal gas reactor types can be used for when the energy equation is
187 //! disabled (fixed temperature simulations).
188 //!
189 //! @param loglevel Print information about solver progress to aid in understanding
190 //! cases where the solver fails to converge. Higher levels are more verbose.
191 //! - 0: No logging.
192 //! - 1: Basic info about each steady-state attempt and round of time stepping.
193 //! - 2: Adds details about each time step and steady-state Newton iteration.
194 //! - 3: Adds details about Newton iterations for each time step.
195 //! - 4: Adds details about state variables that are limiting steady-state
196 //! Newton step sizes.
197 //! - 5: Adds details about state variables that are limiting time-stepping
198 //! Newton step sizes.
199 //! - 6: Print current state vector after different solver stages
200 //! - 7: Print current residual vector after different solver stages
201 //!
202 //! @see SteadyStateSystem, MultiNewton
203 //! @since New in %Cantera 3.2. Support for reacting surfaces added in %Cantera 4.0.
204 void solveSteady(int loglevel=0);
205
206 //! Get the Jacobian used by the steady-state solver.
207 //!
208 //! @param rdt Reciprocal of the pseudo-timestep [1/s]. Default of 0.0 returns the
209 //! steady-state Jacobian.
210 //! @since New in %Cantera 3.2.
211 Eigen::SparseMatrix<double> steadyJacobian(double rdt=0.0);
212
213 //! Return a reference to the *n*-th reactor in this network. The reactor
214 //! indices are determined by the order in which the reactors were added
215 //! to the reactor network.
217 return *m_reactors[n];
218 }
219
220 //! Returns `true` if verbose logging output is enabled.
221 bool verbose() const {
222 return m_verbose;
223 }
224
225 //! Enable or disable verbose logging while setting up and integrating the
226 //! reactor network.
227 void setVerbose(bool v = true) {
228 m_verbose = v;
229 suppressErrors(!m_verbose);
230 }
231
232 //! Return a reference to the integrator. Only valid after adding at least one
233 //! reactor to the network.
235
236 //! Update the state of all the reactors in the network to correspond to
237 //! the values in the solution vector *y*.
238 void updateState(span<const double> y);
239
240 //! Return the sensitivity of the *k*-th solution component with respect to
241 //! the *p*-th sensitivity parameter.
242 /*!
243 * The sensitivity coefficient @f$ S_{ki} @f$ of solution variable @f$ y_k
244 * @f$ with respect to sensitivity parameter @f$ p_i @f$ is defined as:
245 *
246 * @f[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} @f]
247 *
248 * For reaction sensitivities, the parameter is a multiplier on the forward
249 * rate constant (and implicitly on the reverse rate constant for
250 * reversible reactions) which has a nominal value of 1.0, and the
251 * sensitivity is nondimensional.
252 *
253 * For species enthalpy sensitivities, the parameter is a perturbation to
254 * the molar enthalpy of formation, such that the dimensions of the
255 * sensitivity are kmol/J.
256 */
257 double sensitivity(size_t k, size_t p);
258
259 //! Return the sensitivity of the component named *component* with respect to
260 //! the *p*-th sensitivity parameter.
261 //! @copydetails ReactorNet::sensitivity(size_t, size_t)
262 double sensitivity(const string& component, size_t p, int reactor=0) {
263 size_t k = globalComponentIndex(component, reactor);
264 return sensitivity(k, p);
265 }
266
267 //! Evaluate the Jacobian matrix for the reactor network.
268 /*!
269 * @param[in] t Time/distance at which to evaluate the Jacobian
270 * @param[in] y Global state vector at *t*
271 * @param[out] ydot Derivative of the state vector evaluated at *t*, with respect
272 * to *t*.
273 * @param[in] p sensitivity parameter vector (unused?)
274 * @param[out] j Jacobian matrix, size neq() by neq().
275 */
276 void evalJacobian(double t, span<double> y,
277 span<double> ydot, span<const double> p, Array2D* j);
278
279 // overloaded methods of class FuncEval
280 size_t neq() const override {
281 return m_nv;
282 }
283
284 size_t nReactors() const {
285 return m_reactors.size();
286 }
287
288 void eval(double t, span<const double> y, span<double> ydot,
289 span<const double> p) override;
290
291 //! Evaluate the governing equations adapted for the steady-state solver.
292 //!
293 //! @param[in] y Current state vector, length neq()
294 //! @param[out] residual For differential variables, this is the time derivative;
295 //! for algebraic variables, this is the residual of the governing equation.
296 //!
297 //! @since New in %Cantera 4.0.
298 void evalSteady(span<const double> y, span<double> residual);
299
300 //! eval coupling for IDA / DAEs
301 void evalDae(double t, span<const double> y, span<const double> ydot,
302 span<const double> p, span<double> residual) override;
303
304 void getState(span<double> y) override;
305 void getStateDae(span<double> y, span<double> ydot) override;
306
307 //! Return k-th derivative at the current state of the system
308 virtual void getDerivative(int k, span<double> dky);
309
310 void getConstraints(span<double> constraints) override;
311
312 size_t nparams() const override {
313 return m_sens_params.size();
314 }
315
316 //! Return the index corresponding to the component named *component* in the
317 //! reactor with index *reactor* in the global state vector for the
318 //! reactor network.
319 size_t globalComponentIndex(const string& component, size_t reactor=0);
320
321 //! Return the name of the i-th component of the global state vector. The
322 //! name returned includes both the name of the reactor and the specific
323 //! component, for example `'reactor1: CH4'`.
324 string componentName(size_t i) const;
325
326 //! Get the upper bound on the i-th component of the global state vector.
327 double upperBound(size_t i) const;
328
329 //! Get the lower bound on the i-th component of the global state vector.
330 double lowerBound(size_t i) const;
331
332 //! Reset physically or mathematically problematic values, such as negative species
333 //! concentrations.
334 //!
335 //! This method is used within solveSteady() if certain errors are encountered.
336 //!
337 //! @param[inout] y current state vector, to be updated; length neq()
338 void resetBadValues(span<double> y);
339
340 //! Used by Reactor and Wall objects to register the addition of
341 //! sensitivity parameters so that the ReactorNet can keep track of the
342 //! order in which sensitivity parameters are added.
343 //! @param name A name describing the parameter, for example the reaction string
344 //! @param value The nominal value of the parameter
345 //! @param scale A scaling factor to be applied to the sensitivity
346 //! coefficient
347 //! @returns the index of this parameter in the vector of sensitivity
348 //! parameters (global across all reactors)
349 size_t registerSensitivityParameter(const string& name, double value, double scale);
350
351 //! The name of the p-th sensitivity parameter added to this ReactorNet.
352 const string& sensitivityParameterName(size_t p) const {
353 return m_paramNames.at(p);
354 }
355
356 //! Initialize the reactor network. Called automatically the first time
357 //! advance or step is called.
358 void initialize();
359
360 //! Reinitialize the integrator. Used to solve a new problem (different
361 //! initial conditions) but with the same configuration of the reactor
362 //! network. Can be called manually, or automatically after calling
363 //! setInitialTime or modifying a reactor's contents.
364 void reinitialize();
365
366 //! Called to trigger integrator reinitialization before further
367 //! integration.
370 }
371
372 //! Set the maximum number of internal integration steps the
373 //! integrator will take before reaching the next output point
374 //! @param nmax The maximum number of steps, setting this value
375 //! to zero disables this option.
376 virtual void setMaxSteps(int nmax);
377
378 //! Set absolute step size limits during advance
379 void setAdvanceLimits(span<const double> limits);
380
381 //! Check whether ReactorNet object uses advance limits
382 bool hasAdvanceLimits() const;
383
384 //! Retrieve absolute step size limits during advance
385 bool getAdvanceLimits(span<double> limits) const;
386
387 void preconditionerSetup(double t, span<const double> y, double gamma) override;
388
389 void preconditionerSolve(span<const double> rhs, span<double> output) override;
390
391 //! Get solver stats from integrator
392 AnyMap solverStats() const;
393
394 //! Set derivative settings of all reactors
395 //! @param settings the settings map propagated to all reactors and kinetics objects
396 virtual void setDerivativeSettings(AnyMap& settings);
397
398 //! Root finding is enabled only while enforcing advance limits
399 size_t nRootFunctions() const override;
400
401 //! Evaluate the advance-limit root function used to stop integration once a limit
402 //! is met.
403 //!
404 //! When limits are active, this sets `gout[0]` to
405 //! `1 - max_i(|y[i]-y_base[i]| / limit[i])` so a zero indicates a component has
406 //! reached its limit; otherwise `gout[0]` is positive.
407 void evalRootFunctions(double t, span<const double> y,
408 span<double> gout) override;
409
410protected:
411 //! Check that preconditioning is supported by all reactors in the network
412 virtual void checkPreconditionerSupported() const;
413
414 //! Update the integrator tolerance vector from the current scalar settings,
415 //! reactor-specific user tolerances, and reactor-specific default scaling rules.
416 void updateTolerances();
417
418 void updatePreconditioner(double gamma) override;
419
420 //! Create reproducible names for reactors and walls/connectors.
421 void updateNames(ReactorBase& r);
422
423 //! Estimate a future state based on current derivatives.
424 //! The function is intended for internal use by ReactorNet::advance
425 //! and deliberately not exposed in external interfaces.
426 virtual void getEstimate(double time, int k, span<double> yest);
427
428 //! Returns the order used for last solution step of the ODE integrator
429 //! The function is intended for internal use by ReactorNet::advance
430 //! and deliberately not exposed in external interfaces.
431 virtual int lastOrder() const;
432
433 vector<shared_ptr<ReactorBase>> m_reactors;
434 vector<shared_ptr<Reactor>> m_bulkReactors;
435 vector<shared_ptr<ReactorBase>> m_surfaces;
436 vector<shared_ptr<ReactorBase>> m_reservoirs;
437 set<FlowDevice*> m_flowDevices;
438 set<WallBase*> m_walls;
439 map<string, int> m_counts; //!< Map used for default name generation
440 unique_ptr<Integrator> m_integ;
441
442 //! The independent variable in the system. May be either time or space depending
443 //! on the type of reactors in the network.
444 double m_time = 0.0;
445
446 //! The initial value of the independent variable in the system.
447 double m_initial_time = 0.0;
448
449 bool m_integratorInitialized = false; //!< True if the integrator has been initialized at least once
450 bool m_needIntegratorInit = true; //!< True if integrator needs to be (re)initialized
451 size_t m_nv = 0;
452
453 vector<double> m_atol;
454 double m_rtol = 1.0e-9;
455 double m_rtolsens = 1.0e-4;
456 double m_atols = 1.0e-15;
457 double m_atolsens = 1.0e-6;
458 //! True if scalar absolute tolerance was user-specified
460 shared_ptr<SystemJacobian> m_precon;
461 string m_linearSolverType;
462
463 //! Maximum integrator internal timestep. Default of 0.0 means infinity.
464 double m_maxstep = 0.0;
465
466 bool m_verbose = false;
467
468 //! Indicates whether time or space is the independent variable
470
471 //! Names corresponding to each sensitivity parameter
472 vector<string> m_paramNames;
473
474 vector<double> m_ydot;
475 vector<double> m_yest;
476 vector<double> m_advancelimits;
477 //! Base state used for evaluating advance limits during a single advance()
478 //! call when root-finding is enabled
479 vector<double> m_ybase;
480 //! Base time corresponding to #m_ybase
481 double m_ybase_time = 0.0;
482 //! Indicates whether the advance-limit root check is active for the
483 //! current call to `advance(t, applylimit=true)`
485 //! m_LHS is a vector representing the coefficients on the
486 //! "left hand side" of each governing equation
487 vector<double> m_LHS;
488 vector<double> m_RHS;
489};
490
491
492//! Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
493//!
494//! @see ReactorNet::solveSteady
495//! @since New in %Cantera 3.2.
497{
498public:
499 SteadyReactorSolver(ReactorNet* net, span<const double> x0);
500 void eval(span<const double> x, span<double> r,
501 double rdt=-1.0, int count=1) override;
502 void initTimeInteg(double dt, span<const double> x) override;
503 void evalJacobian(span<const double> x0) override;
504 double weightedNorm(span<const double> step) const override;
505 string componentName(size_t i) const override;
506 double upperBound(size_t i) const override;
507 double lowerBound(size_t i) const override;
508 void resetBadValues(span<double> x) override;
509 void writeDebugInfo(const string& header_suffix, const string& message,
510 int loglevel, int attempt_counter) override;
511
512private:
513 ReactorNet* m_net = nullptr;
514
515 //! Initial value of each state variable
516 vector<double> m_initialState;
517};
518
519
520/**
521 * Create a reactor network containing one or more coupled reactors.
522 * Wall and FlowDevice objects should be installed prior to calling newReactorNet().
523 * @param[in] reactors A vector of shared pointers to the reactors to be linked
524 * together.
525 * @since New in %Cantera 3.2.
526 */
527shared_ptr<ReactorNet> newReactorNet(span<shared_ptr<ReactorBase>> reactors);
528
529}
530
531#endif
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
Virtual base class for ODE/DAE right-hand-side function evaluators.
Definition FuncEval.h:32
bool suppressErrors() const
Get current state of error suppression.
Definition FuncEval.h:190
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition FuncEval.h:205
Abstract base class for ODE system integrators.
Definition Integrator.h:44
Base class for reactor objects.
Definition ReactorBase.h:51
A class representing a network of connected reactors.
Definition ReactorNet.h:30
virtual void getDerivative(int k, span< double > dky)
Return k-th derivative at the current state of the system.
void setLinearSolverType(const string &linSolverType="DENSE")
Set the type of linear solver used in the integration.
bool getAdvanceLimits(span< double > limits) const
Retrieve absolute step size limits during advance.
double step()
Advance the state of all reactors with respect to the independent variable (time or space).
virtual int lastOrder() const
Returns the order used for last solution step of the ODE integrator The function is intended for inte...
const string & sensitivityParameterName(size_t p) const
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition ReactorNet.h:352
double m_ybase_time
Base time corresponding to m_ybase.
Definition ReactorNet.h:481
size_t nparams() const override
Number of sensitivity parameters.
Definition ReactorNet.h:312
void initialize()
Initialize the reactor network.
void evalDae(double t, span< const double > y, span< const double > ydot, span< const double > p, span< double > residual) override
eval coupling for IDA / DAEs
void advance(double t)
Advance the state of all reactors in the independent variable (time or space).
size_t neq() const override
Number of equations.
Definition ReactorNet.h:280
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition ReactorNet.h:216
vector< double > m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition ReactorNet.h:487
double m_initial_time
The initial value of the independent variable in the system.
Definition ReactorNet.h:447
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void updateNames(ReactorBase &r)
Create reproducible names for reactors and walls/connectors.
double getInitialTime() const
Get the initial value of the independent variable (typically time).
Definition ReactorNet.h:63
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
Definition ReactorNet.h:368
double m_time
The independent variable in the system.
Definition ReactorNet.h:444
AnyMap solverStats() const
Get solver stats from integrator.
void updateTolerances()
Update the integrator tolerance vector from the current scalar settings, reactor-specific user tolera...
bool m_limit_check_active
Indicates whether the advance-limit root check is active for the current call to advance(t,...
Definition ReactorNet.h:484
map< string, int > m_counts
Map used for default name generation.
Definition ReactorNet.h:439
double upperBound(size_t i) const
Get the upper bound on the i-th component of the global state vector.
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration steps the integrator will take before reaching the nex...
bool m_integratorInitialized
True if the integrator has been initialized at least once.
Definition ReactorNet.h:449
string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
virtual void getEstimate(double time, int k, span< double > yest)
Estimate a future state based on current derivatives.
size_t nRootFunctions() const override
Root finding is enabled only while enforcing advance limits.
void evalRootFunctions(double t, span< const double > y, span< double > gout) override
Evaluate the advance-limit root function used to stop integration once a limit is met.
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single step.
size_t registerSensitivityParameter(const string &name, double value, double scale)
Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the React...
void getConstraints(span< double > constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
double maxTimeStep() const
Get the maximum integrator step.
Definition ReactorNet.h:68
double m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition ReactorNet.h:464
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
void evalSteady(span< const double > y, span< double > residual)
Evaluate the governing equations adapted for the steady-state solver.
double atolSensitivity() const
Absolute sensitivity tolerance.
Definition ReactorNet.h:141
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
int maxSteps()
Returns the maximum number of internal integration steps the integrator will take before reaching the...
virtual void setDerivativeSettings(AnyMap &settings)
Set derivative settings of all reactors.
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.
void evalJacobian(double t, span< double > y, span< double > ydot, span< const double > p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
void updateState(span< const double > y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
void preconditionerSolve(span< const double > rhs, span< double > output) override
Evaluate the linear system Ax=b where A is the preconditioner.
bool m_needIntegratorInit
True if integrator needs to be (re)initialized.
Definition ReactorNet.h:450
void preconditionerSetup(double t, span< const double > y, double gamma) override
Evaluate the setup processes for the Jacobian preconditioner.
double rtol()
Relative tolerance.
Definition ReactorNet.h:126
bool verbose() const
Returns true if verbose logging output is enabled.
Definition ReactorNet.h:221
size_t globalComponentIndex(const string &component, size_t reactor=0)
Return the index corresponding to the component named component in the reactor with index reactor in ...
void resetBadValues(span< double > y)
Reset physically or mathematically problematic values, such as negative species concentrations.
void solveSteady(int loglevel=0)
Solve directly for the steady-state solution.
bool m_timeIsIndependent
Indicates whether time or space is the independent variable.
Definition ReactorNet.h:469
double atol()
Scalar absolute integration tolerance.
Definition ReactorNet.h:131
bool hasAdvanceLimits() const
Check whether ReactorNet object uses advance limits.
vector< double > m_ybase
Base state used for evaluating advance limits during a single advance() call when root-finding is ena...
Definition ReactorNet.h:479
double sensitivity(const 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:262
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
void setAdvanceLimits(span< const double > limits)
Set absolute step size limits during advance.
double rtolSensitivity() const
Relative sensitivity tolerance.
Definition ReactorNet.h:136
bool m_atolUserSpecified
True if scalar absolute tolerance was user-specified.
Definition ReactorNet.h:459
double lowerBound(size_t i) const
Get the lower bound on the i-th component of the global state vector.
virtual void checkPreconditionerSupported() const
Check that preconditioning is supported by all reactors in the network.
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition ReactorNet.h:227
Eigen::SparseMatrix< double > steadyJacobian(double rdt=0.0)
Get the Jacobian used by the steady-state solver.
string linearSolverType() const
Problem type of integrator.
void updatePreconditioner(double gamma) override
Update the preconditioner based on already computed jacobian values.
void setPreconditioner(shared_ptr< SystemJacobian > preconditioner)
Set preconditioner used by the linear solver.
void getState(span< double > y) override
Fill in the vector y with the current state of the system.
vector< string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition ReactorNet.h:472
void clearAbsoluteTolerance()
Clear the user-specified scalar absolute tolerance.
void getStateDae(span< double > y, span< double > ydot) override
Fill in the vectors y and ydot with the current state of the system.
void setTolerances(double rtol, double atol)
Set the relative and scalar absolute tolerances for the integrator.
void setRelativeTolerance(double rtol)
Set the relative tolerance for the integrator.
void eval(double t, span< const double > y, span< double > ydot, span< const double > p) override
Evaluate the right-hand-side ODE function.
void setAbsoluteTolerance(double atol)
Set the scalar absolute tolerance for the integrator.
Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
Definition ReactorNet.h:497
double weightedNorm(span< const double > step) const override
Compute the weighted norm of step.
vector< double > m_initialState
Initial value of each state variable.
Definition ReactorNet.h:516
string componentName(size_t i) const override
Get the name of the i-th component of the state vector.
void initTimeInteg(double dt, span< const double > x) override
Prepare for time stepping beginning with solution x and timestep dt.
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
void evalJacobian(span< const double > x0) override
Evaluates the Jacobian at x0 using finite differences.
void writeDebugInfo(const string &header_suffix, const string &message, int loglevel, int attempt_counter) override
Write solver debugging based on the specified log level.
void resetBadValues(span< double > x) override
Reset values such as negative species concentrations.
void eval(span< const double > x, span< double > r, double rdt=-1.0, int count=1) override
Evaluate the residual function.
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
Base class for representing a system of differential-algebraic equations and solving for its steady-s...
double rdt() const
Reciprocal of the time step.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition utilities.h:118
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
shared_ptr< ReactorNet > newReactorNet(span< shared_ptr< ReactorBase > > reactors)
Create a reactor network containing one or more coupled reactors.