A class representing a network of connected reactors. More...
#include <ReactorNet.h>
A class representing a network of connected reactors.
This class is used to integrate the governing equations for a network of reactors that are time dependent (Reactor, ConstPressureReactor) connected by various means, for example Wall, MassFlowController, Valve, or PressureController; or reactors dependent on a single spatial variable (FlowReactor).
Definition at line 29 of file ReactorNet.h.
Public Member Functions | |
| ReactorNet (shared_ptr< ReactorBase > reactor) | |
| Create reactor network containing single reactor. | |
| ReactorNet (span< shared_ptr< ReactorBase > > reactors) | |
| Create reactor network from multiple reactors. | |
| ReactorNet (const ReactorNet &)=delete | |
| ReactorNet & | operator= (const ReactorNet &)=delete |
| void | advance (double t) |
| Advance the state of all reactors in the independent variable (time or space). | |
| double | advance (double t, bool applylimit) |
| Advance the state of all reactors in the independent variable (time or space). | |
| double | step () |
| Advance the state of all reactors with respect to the independent variable (time or space). | |
| void | solveSteady (int loglevel=0) |
| Solve directly for the steady-state solution. | |
| Eigen::SparseMatrix< double > | steadyJacobian (double rdt=0.0) |
| Get the Jacobian used by the steady-state solver. | |
| Eigen::SparseMatrix< double > | jacobian () |
| Calculate the semi-analytical preconditioner Jacobian for the entire network. | |
| Eigen::SparseMatrix< double > | finiteDifferenceJacobian () |
| Calculate the system Jacobian for the reactor network using finite differences. | |
| ReactorBase & | reactor (int n) |
| Return a reference to the n-th reactor in this network. | |
| bool | verbose () const |
Returns true if verbose logging output is enabled. | |
| void | setVerbose (bool v=true) |
| Enable or disable verbose logging while setting up and integrating the reactor network. | |
| Integrator & | integrator () |
| Return a reference to the integrator. | |
| void | updateState (span< const double > y) |
| Update the state of all the reactors in the network to correspond to the values in the solution vector y. | |
| 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. | |
| 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 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. | |
| size_t | neq () const override |
| Number of equations. | |
| size_t | nReactors () const |
| void | eval (double t, span< const double > y, span< double > ydot, span< const double > p) override |
| Evaluate the right-hand-side ODE function. | |
| void | evalSteady (span< const double > y, span< double > residual) |
| Evaluate the governing equations adapted for the steady-state solver. | |
| 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 | getState (span< double > y) override |
| Fill in the vector y with the current state of the system. | |
| void | getStateDae (span< double > y, span< double > ydot) override |
| Fill in the vectors y and ydot with the current state of the system. | |
| virtual void | getDerivative (int k, span< double > dky) |
| Return k-th derivative at the current state of the system. | |
| void | getConstraints (span< double > constraints) override |
| Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
| size_t | nparams () const override |
| Number of sensitivity parameters. | |
| 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 the global state vector for the reactor network. | |
| string | componentName (size_t i) const |
| Return the name of the i-th component of the global state vector. | |
| double | upperBound (size_t i) const |
| Get the upper bound on the i-th component of the global state vector. | |
| double | lowerBound (size_t i) const |
| Get the lower bound on the i-th component of the global state vector. | |
| void | resetBadValues (span< double > y) |
| Reset physically or mathematically problematic values, such as negative species concentrations. | |
| 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 ReactorNet can keep track of the order in which sensitivity parameters are added. | |
| const string & | sensitivityParameterName (size_t p) const |
| The name of the p-th sensitivity parameter added to this ReactorNet. | |
| void | initialize () |
| Initialize the reactor network. | |
| void | reinitialize () |
| Reinitialize the integrator. | |
| void | setNeedsReinit () |
| Called to trigger integrator reinitialization before further integration. | |
| virtual void | setMaxSteps (int nmax) |
| Set the maximum number of internal integration steps the integrator will take before reaching the next output point. | |
| void | setAdvanceLimits (span< const double > limits) |
| Set absolute step size limits during advance. | |
| bool | hasAdvanceLimits () const |
| Check whether ReactorNet object uses advance limits. | |
| bool | getAdvanceLimits (span< double > limits) const |
| Retrieve absolute step size limits during advance. | |
| void | preconditionerSetup (double t, span< const double > y, double gamma) override |
| Evaluate the setup processes for the Jacobian preconditioner. | |
| void | preconditionerSolve (span< const double > rhs, span< double > output) override |
| Evaluate the linear system Ax=b where A is the preconditioner. | |
| AnyMap | solverStats () const |
| Get solver stats from integrator. | |
| virtual void | setDerivativeSettings (AnyMap &settings) |
| Set derivative settings of all reactors. | |
| 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. | |
Methods to set up a simulation | |
| void | setLinearSolverType (const string &linSolverType="DENSE") |
| Set the type of linear solver used in the integration. | |
| void | setPreconditioner (shared_ptr< SystemJacobian > preconditioner) |
| Set preconditioner used by the linear solver. | |
| void | setInitialTime (double time) |
| Set the initial value of the independent variable (typically time). | |
| double | getInitialTime () const |
| Get the initial value of the independent variable (typically time). | |
| double | maxTimeStep () const |
| Get the maximum integrator step. | |
| void | setMaxTimeStep (double maxstep) |
| Set the maximum integrator step. | |
| void | setMaxErrTestFails (int nmax) |
| Set the maximum number of error test failures permitted by the CVODES integrator in a single step. | |
| void | setRelativeTolerance (double rtol) |
| Set the relative tolerance for the integrator. | |
| void | setAbsoluteTolerance (double atol) |
| Set the scalar absolute tolerance for the integrator. | |
| void | clearAbsoluteTolerance () |
| Clear the user-specified scalar absolute tolerance. | |
| void | setTolerances (double rtol, double atol) |
| Set the relative and scalar absolute tolerances for the integrator. | |
| void | setSensitivityTolerances (double rtol, double atol) |
| Set the relative and absolute tolerances for integrating the sensitivity equations. | |
| double | time () |
| Current value of the simulation time [s], for reactor networks that are solved in the time domain. | |
| double | distance () |
| Current position [m] along the length of the reactor network, for reactors that are solved as a function of space. | |
| double | rtol () |
| Relative tolerance. | |
| double | atol () |
| Scalar absolute integration tolerance. | |
| double | rtolSensitivity () const |
| Relative sensitivity tolerance. | |
| double | atolSensitivity () const |
| Absolute sensitivity tolerance. | |
| string | linearSolverType () const |
| Problem type of integrator. | |
| int | maxSteps () |
| Returns the maximum number of internal integration steps the integrator will take before reaching the next output point. | |
Public Member Functions inherited from FuncEval | |
| virtual void | eval (double t, span< const double > y, span< double > ydot, span< const double > p) |
| Evaluate the right-hand-side ODE function. | |
| virtual void | evalDae (double t, span< const double > y, span< const double > ydot, span< const double > p, span< double > residual) |
| Evaluate the right-hand-side DAE function. | |
| virtual void | getConstraints (span< double > constraints) |
| Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
| int | evalNoThrow (double t, span< const double > y, span< double > ydot) |
| Evaluate the right-hand side using return code to indicate status. | |
| int | evalDaeNoThrow (double t, span< const double > y, span< const double > ydot, span< double > residual) |
| Evaluate the right-hand side using return code to indicate status. | |
| virtual void | preconditionerSetup (double t, span< const double > y, double gamma) |
| Evaluate the setup processes for the Jacobian preconditioner. | |
| virtual void | preconditionerSolve (span< const double > rhs, span< double > output) |
| Evaluate the linear system Ax=b where A is the preconditioner. | |
| virtual void | updatePreconditioner (double gamma) |
| Update the preconditioner based on already computed jacobian values. | |
| int | preconditioner_setup_nothrow (double t, span< const double > y, double gamma) |
| Preconditioner setup that doesn't throw an error but returns a CVODES flag. | |
| int | preconditioner_solve_nothrow (span< const double > rhs, span< double > output) |
| Preconditioner solve that doesn't throw an error but returns a CVODES flag. | |
| virtual size_t | nRootFunctions () const |
| Number of event/root functions exposed to the integrator. | |
| virtual void | evalRootFunctions (double t, span< const double > y, span< double > gout) |
| Evaluate the event/root functions currently in play. | |
| int | evalRootFunctionsNoThrow (double t, span< const double > y, span< double > gout) |
| Wrapper for evalRootFunctions that converts exceptions to return codes. | |
| virtual void | getState (span< double > y) |
| Fill in the vector y with the current state of the system. | |
| virtual void | getStateDae (span< double > y, span< double > ydot) |
| Fill in the vectors y and ydot with the current state of the system. | |
| virtual size_t | neq () const =0 |
| Number of equations. | |
| virtual size_t | nparams () const |
| Number of sensitivity parameters. | |
| void | suppressErrors (bool suppress) |
| Enable or disable suppression of errors when calling eval() | |
| bool | suppressErrors () const |
| Get current state of error suppression. | |
| string | getErrors () const |
| Return a string containing the text of any suppressed errors. | |
| void | clearErrors () |
| Clear any previously-stored suppressed errors. | |
Protected Member Functions | |
| virtual void | checkPreconditionerSupported () const |
| Check that preconditioning is supported by all reactors in the network. | |
| void | updateTolerances () |
| Update the integrator tolerance vector from the current scalar settings, reactor-specific user tolerances, and reactor-specific default scaling rules. | |
| void | updatePreconditioner (double gamma) override |
| Update the preconditioner based on already computed jacobian values. | |
| void | updateNames (ReactorBase &r) |
| Create reproducible names for reactors and walls/connectors. | |
| virtual void | getEstimate (double time, int k, span< double > yest) |
| Estimate a future state based on current derivatives. | |
| virtual int | lastOrder () const |
| Returns the order used for last solution step of the ODE integrator The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces. | |
Protected Attributes | |
| vector< shared_ptr< ReactorBase > > | m_reactors |
| vector< shared_ptr< Reactor > > | m_bulkReactors |
| vector< shared_ptr< ReactorBase > > | m_surfaces |
| vector< shared_ptr< ReactorBase > > | m_reservoirs |
| set< FlowDevice * > | m_flowDevices |
| set< WallBase * > | m_walls |
| map< string, int > | m_counts |
| Map used for default name generation. | |
| unique_ptr< Integrator > | m_integ |
| double | m_time = 0.0 |
| The independent variable in the system. | |
| double | m_initial_time = 0.0 |
| The initial value of the independent variable in the system. | |
| bool | m_integratorInitialized = false |
| True if the integrator has been initialized at least once. | |
| bool | m_needIntegratorInit = true |
| True if integrator needs to be (re)initialized. | |
| size_t | m_nv = 0 |
| vector< double > | m_atol |
| double | m_rtol = 1.0e-9 |
| double | m_rtolsens = 1.0e-4 |
| double | m_atols = 1.0e-15 |
| double | m_atolsens = 1.0e-6 |
| bool | m_atolUserSpecified = false |
| True if scalar absolute tolerance was user-specified. | |
| shared_ptr< SystemJacobian > | m_precon |
| string | m_linearSolverType |
| double | m_maxstep = 0.0 |
| Maximum integrator internal timestep. Default of 0.0 means infinity. | |
| bool | m_verbose = false |
| bool | m_timeIsIndependent = true |
| Indicates whether time or space is the independent variable. | |
| vector< string > | m_paramNames |
| Names corresponding to each sensitivity parameter. | |
| vector< double > | m_ydot |
| vector< double > | m_yest |
| vector< double > | m_advancelimits |
| vector< double > | m_ybase |
| Base state used for evaluating advance limits during a single advance() call when root-finding is enabled. | |
| double | m_ybase_time = 0.0 |
| Base time corresponding to m_ybase. | |
| bool | m_limit_check_active = false |
Indicates whether the advance-limit root check is active for the current call to advance(t, applylimit=true) | |
| vector< double > | m_LHS |
| m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation | |
| vector< double > | m_RHS |
Protected Attributes inherited from FuncEval | |
| bool | m_suppress_errors = false |
| vector< string > | m_errors |
| Errors occurring during function evaluations. | |
Additional Inherited Members | |
Public Attributes inherited from FuncEval | |
| vector< double > | m_sens_params |
| Values for the problem parameters for which sensitivities are computed This is the array which is perturbed and passed back as the fourth argument to eval(). | |
| vector< double > | m_paramScales |
| Scaling factors for each sensitivity parameter. | |
| ReactorNet | ( | shared_ptr< ReactorBase > | reactor | ) |
Create reactor network containing single reactor.
Definition at line 32 of file ReactorNet.cpp.
| ReactorNet | ( | span< shared_ptr< ReactorBase > > | reactors | ) |
Create reactor network from multiple reactors.
Definition at line 37 of file ReactorNet.cpp.
| ~ReactorNet | ( | ) |
Definition at line 28 of file ReactorNet.cpp.
| void setLinearSolverType | ( | const string & | linSolverType = "DENSE" | ) |
Set the type of linear solver used in the integration.
| linSolverType | type of linear solver. Default type: "DENSE" Other options include: "DIAG", "DENSE", "GMRES", "BAND" |
Definition at line 324 of file ReactorNet.cpp.
| void setPreconditioner | ( | shared_ptr< SystemJacobian > | preconditioner | ) |
Set preconditioner used by the linear solver.
| preconditioner | preconditioner object used for the linear solver |
Definition at line 330 of file ReactorNet.cpp.
| void setInitialTime | ( | double | time | ) |
Set the initial value of the independent variable (typically time).
Default = 0.0 s. Restarts integration from this value using the current mixture state as the initial condition.
Definition at line 172 of file ReactorNet.cpp.
|
inline |
Get the initial value of the independent variable (typically time).
Definition at line 63 of file ReactorNet.h.
|
inline |
Get the maximum integrator step.
Definition at line 68 of file ReactorNet.h.
| void setMaxTimeStep | ( | double | maxstep | ) |
Set the maximum integrator step.
Definition at line 179 of file ReactorNet.cpp.
| void setMaxErrTestFails | ( | int | nmax | ) |
Set the maximum number of error test failures permitted by the CVODES integrator in a single step.
Definition at line 185 of file ReactorNet.cpp.
| void setRelativeTolerance | ( | double | rtol | ) |
Set the relative tolerance for the integrator.
| rtol | Relative tolerance for all state variables; must be positive. |
Definition at line 190 of file ReactorNet.cpp.
| void setAbsoluteTolerance | ( | double | atol | ) |
Set the scalar absolute tolerance for the integrator.
| atol | Scalar absolute tolerance for state variables without reactor-specific absolute tolerances; must be positive. Superseded by absolute tolerances set for individual reactors. |
Definition at line 200 of file ReactorNet.cpp.
| void clearAbsoluteTolerance | ( | ) |
Clear the user-specified scalar absolute tolerance.
Restores the default scalar absolute tolerance and enables reactor-specific default scaling rules.
Definition at line 211 of file ReactorNet.cpp.
| void setTolerances | ( | double | rtol, |
| double | atol | ||
| ) |
Set the relative and scalar absolute tolerances for the integrator.
| rtol | If positive, set the relative tolerance for all state variables; ignored if negative. |
| atol | If positive, set the scalar absolute tolerance for all state variables; ignored if negative. Superseded by absolute tolerances set for individual reactors. |
Definition at line 218 of file ReactorNet.cpp.
| void setSensitivityTolerances | ( | double | rtol, |
| double | atol | ||
| ) |
Set the relative and absolute tolerances for integrating the sensitivity equations.
Definition at line 247 of file ReactorNet.cpp.
| double time | ( | ) |
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
Definition at line 258 of file ReactorNet.cpp.
| double distance | ( | ) |
Current position [m] along the length of the reactor network, for reactors that are solved as a function of space.
Definition at line 267 of file ReactorNet.cpp.
|
inline |
Relative tolerance.
Definition at line 126 of file ReactorNet.h.
|
inline |
Scalar absolute integration tolerance.
Definition at line 131 of file ReactorNet.h.
|
inline |
Relative sensitivity tolerance.
Definition at line 136 of file ReactorNet.h.
|
inline |
Absolute sensitivity tolerance.
Definition at line 141 of file ReactorNet.h.
| string linearSolverType | ( | ) | const |
Problem type of integrator.
Definition at line 863 of file ReactorNet.cpp.
| int maxSteps | ( | ) |
Returns the maximum number of internal integration steps the integrator will take before reaching the next output point.
Definition at line 346 of file ReactorNet.cpp.
| void advance | ( | double | t | ) |
Advance the state of all reactors in the independent variable (time or space).
Take as many internal steps as necessary to reach t.
| t | Time/distance to advance to (s or m). |
Definition at line 351 of file ReactorNet.cpp.
| double advance | ( | double | t, |
| bool | applylimit | ||
| ) |
Advance the state of all reactors in the independent variable (time or space).
Take as many internal steps as necessary towards t. If applylimit is true, the advance step will be automatically reduced if needed to stay within limits (set by setAdvanceLimit). Returns the time/distance at the end of integration.
| t | Time/distance to advance to (s or m). |
| applylimit | Limit advance step (boolean). |
Definition at line 359 of file ReactorNet.cpp.
| double step | ( | ) |
Advance the state of all reactors with respect to the independent variable (time or space).
Returns the new value of the independent variable [s or m].
Definition at line 435 of file ReactorNet.cpp.
| void solveSteady | ( | int | loglevel = 0 | ) |
Solve directly for the steady-state solution.
This approach is generally more efficient than time marching to the steady-state, but imposes a few limitations:
| loglevel | Print information about solver progress to aid in understanding cases where the solver fails to converge. Higher levels are more verbose.
|
Definition at line 443 of file ReactorNet.cpp.
| Eigen::SparseMatrix< double > steadyJacobian | ( | double | rdt = 0.0 | ) |
Get the Jacobian used by the steady-state solver.
| rdt | Reciprocal of the pseudo-timestep [1/s]. Default of 0.0 returns the steady-state Jacobian. |
Definition at line 455 of file ReactorNet.cpp.
| Eigen::SparseMatrix< double > jacobian | ( | ) |
Calculate the semi-analytical preconditioner Jacobian for the entire network.
Collects sparse Jacobian entries from each reactor's getJacobianElements() implementation and assembles them into a single network-level matrix. Reactors that do not implement getJacobianElements() contribute no entries. The matrix uses global row and column indices for the full network state vector.
Definition at line 469 of file ReactorNet.cpp.
| Eigen::SparseMatrix< double > finiteDifferenceJacobian | ( | ) |
Calculate the system Jacobian for the reactor network using finite differences.
Perturbs each element of the state vector and evaluates the network RHS to estimate the full Jacobian. Uses the integrator tolerances to scale each perturbation. This method is intended for debugging and validation of analytical Jacobian implementations, including those added via {ct}ExtensibleReactor.
Definition at line 481 of file ReactorNet.cpp.
|
inline |
Return a reference to the n-th reactor in this network.
The reactor indices are determined by the order in which the reactors were added to the reactor network.
Definition at line 245 of file ReactorNet.h.
|
inline |
Returns true if verbose logging output is enabled.
Definition at line 250 of file ReactorNet.h.
|
inline |
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition at line 256 of file ReactorNet.h.
| Integrator & integrator | ( | ) |
Return a reference to the integrator.
Only valid after adding at least one reactor to the network.
Definition at line 609 of file ReactorNet.cpp.
| void updateState | ( | span< const double > | y | ) |
Update the state of all the reactors in the network to correspond to the values in the solution vector y.
Definition at line 713 of file ReactorNet.cpp.
| 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.
The sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:
\[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} \]
For reaction sensitivities, the parameter is a multiplier on the forward rate constant (and implicitly on the reverse rate constant for reversible reactions) which has a nominal value of 1.0, and the sensitivity is nondimensional.
For species enthalpy sensitivities, the parameter is a perturbation to the molar enthalpy of formation, such that the dimensions of the sensitivity are kmol/J.
Definition at line 676 of file ReactorNet.cpp.
|
inline |
Return the sensitivity of the component named component with respect to the p-th sensitivity parameter.
The sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:
\[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} \]
For reaction sensitivities, the parameter is a multiplier on the forward rate constant (and implicitly on the reverse rate constant for reversible reactions) which has a nominal value of 1.0, and the sensitivity is nondimensional.
For species enthalpy sensitivities, the parameter is a perturbation to the molar enthalpy of formation, such that the dimensions of the sensitivity are kmol/J.
Definition at line 291 of file ReactorNet.h.
| void evalJacobian | ( | double | t, |
| span< double > | y, | ||
| span< double > | ydot, | ||
| span< const double > | p, | ||
| Array2D * | j | ||
| ) |
Evaluate the Jacobian matrix for the reactor network.
| [in] | t | Time/distance at which to evaluate the Jacobian |
| [in] | y | Global state vector at t |
| [out] | ydot | Derivative of the state vector evaluated at t, with respect to t. |
| [in] | p | sensitivity parameter vector (unused?) |
| [out] | j | Jacobian matrix, size neq() by neq(). |
Definition at line 690 of file ReactorNet.cpp.
|
inlineoverridevirtual |
|
inline |
Definition at line 313 of file ReactorNet.h.
|
overridevirtual |
Evaluate the right-hand-side ODE function.
Called by the integrator.
| [in] | t | time. |
| [in] | y | solution vector, length neq() |
| [out] | ydot | rate of change of solution vector, length neq() |
| [in] | p | sensitivity parameter vector, length nparams() |
Reimplemented from FuncEval.
Definition at line 617 of file ReactorNet.cpp.
| void evalSteady | ( | span< const double > | y, |
| span< double > | residual | ||
| ) |
Evaluate the governing equations adapted for the steady-state solver.
| [in] | y | Current state vector, length neq() |
| [out] | residual | For differential variables, this is the time derivative; for algebraic variables, this is the residual of the governing equation. |
Definition at line 637 of file ReactorNet.cpp.
|
overridevirtual |
eval coupling for IDA / DAEs
Reimplemented from FuncEval.
Definition at line 653 of file ReactorNet.cpp.
|
overridevirtual |
Fill in the vector y with the current state of the system.
Used for getting the initial state for ODE systems.
Reimplemented from FuncEval.
Definition at line 754 of file ReactorNet.cpp.
|
overridevirtual |
Fill in the vectors y and ydot with the current state of the system.
Used for getting the initial state for DAE systems.
Reimplemented from FuncEval.
Definition at line 761 of file ReactorNet.cpp.
|
virtual |
Return k-th derivative at the current state of the system.
Definition at line 721 of file ReactorNet.cpp.
|
overridevirtual |
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Reimplemented from FuncEval.
Definition at line 669 of file ReactorNet.cpp.
|
inlineoverridevirtual |
Number of sensitivity parameters.
Reimplemented from FuncEval.
Definition at line 341 of file ReactorNet.h.
| 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 the global state vector for the reactor network.
Definition at line 774 of file ReactorNet.cpp.
| string componentName | ( | size_t | i | ) | const |
Return the name of the i-th component of the global state vector.
The name returned includes both the name of the reactor and the specific component, for example ‘'reactor1: CH4’`.
Definition at line 781 of file ReactorNet.cpp.
| double upperBound | ( | size_t | i | ) | const |
Get the upper bound on the i-th component of the global state vector.
Definition at line 796 of file ReactorNet.cpp.
| double lowerBound | ( | size_t | i | ) | const |
Get the lower bound on the i-th component of the global state vector.
Definition at line 811 of file ReactorNet.cpp.
| void resetBadValues | ( | span< double > | y | ) |
Reset physically or mathematically problematic values, such as negative species concentrations.
This method is used within solveSteady() if certain errors are encountered.
| [in,out] | y | current state vector, to be updated; length neq() |
Definition at line 826 of file ReactorNet.cpp.
| 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 ReactorNet can keep track of the order in which sensitivity parameters are added.
| name | A name describing the parameter, for example the reaction string |
| value | The nominal value of the parameter |
| scale | A scaling factor to be applied to the sensitivity coefficient |
Definition at line 832 of file ReactorNet.cpp.
|
inline |
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition at line 381 of file ReactorNet.h.
| void initialize | ( | ) |
Initialize the reactor network.
Called automatically the first time advance or step is called.
Definition at line 276 of file ReactorNet.cpp.
| void reinitialize | ( | ) |
Reinitialize the integrator.
Used to solve a new problem (different initial conditions) but with the same configuration of the reactor network. Can be called manually, or automatically after calling setInitialTime or modifying a reactor's contents.
Definition at line 306 of file ReactorNet.cpp.
|
inline |
Called to trigger integrator reinitialization before further integration.
Definition at line 397 of file ReactorNet.h.
|
virtual |
Set the maximum number of internal integration steps the integrator will take before reaching the next output point.
| nmax | The maximum number of steps, setting this value to zero disables this option. |
Definition at line 341 of file ReactorNet.cpp.
| void setAdvanceLimits | ( | span< const double > | limits | ) |
Set absolute step size limits during advance.
Definition at line 728 of file ReactorNet.cpp.
| bool hasAdvanceLimits | ( | ) | const |
Check whether ReactorNet object uses advance limits.
Definition at line 736 of file ReactorNet.cpp.
| bool getAdvanceLimits | ( | span< double > | limits | ) | const |
Retrieve absolute step size limits during advance.
Definition at line 745 of file ReactorNet.cpp.
|
overridevirtual |
Evaluate the setup processes for the Jacobian preconditioner.
| [in] | t | time. |
| [in] | y | solution vector, length neq() |
| gamma | the gamma in M=I-gamma*J |
Reimplemented from FuncEval.
Definition at line 881 of file ReactorNet.cpp.
|
overridevirtual |
Evaluate the linear system Ax=b where A is the preconditioner.
| [in] | rhs | right hand side vector used in linear system |
| [out] | output | output vector for solution |
Reimplemented from FuncEval.
Definition at line 872 of file ReactorNet.cpp.
| AnyMap solverStats | ( | ) | const |
Get solver stats from integrator.
Definition at line 854 of file ReactorNet.cpp.
|
virtual |
Set derivative settings of all reactors.
| settings | the settings map propagated to all reactors and kinetics objects |
Definition at line 846 of file ReactorNet.cpp.
|
overridevirtual |
Root finding is enabled only while enforcing advance limits.
Reimplemented from FuncEval.
Definition at line 538 of file ReactorNet.cpp.
|
overridevirtual |
Evaluate the advance-limit root function used to stop integration once a limit is met.
When limits are active, this sets gout[0] to 1 - max_i(|y[i]-y_base[i]| / limit[i]) so a zero indicates a component has reached its limit; otherwise gout[0] is positive.
Reimplemented from FuncEval.
Definition at line 543 of file ReactorNet.cpp.
|
protectedvirtual |
Check that preconditioning is supported by all reactors in the network.
Definition at line 921 of file ReactorNet.cpp.
|
protected |
Update the integrator tolerance vector from the current scalar settings, reactor-specific user tolerances, and reactor-specific default scaling rules.
Definition at line 233 of file ReactorNet.cpp.
|
overrideprotectedvirtual |
Update the preconditioner based on already computed jacobian values.
Reimplemented from FuncEval.
Definition at line 910 of file ReactorNet.cpp.
|
protected |
Create reproducible names for reactors and walls/connectors.
Definition at line 572 of file ReactorNet.cpp.
|
protectedvirtual |
Estimate a future state based on current derivatives.
The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces.
Definition at line 511 of file ReactorNet.cpp.
|
protectedvirtual |
Returns the order used for last solution step of the ODE integrator The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces.
Definition at line 529 of file ReactorNet.cpp.
|
protected |
Definition at line 462 of file ReactorNet.h.
|
protected |
Definition at line 463 of file ReactorNet.h.
|
protected |
Definition at line 464 of file ReactorNet.h.
|
protected |
Definition at line 465 of file ReactorNet.h.
|
protected |
Definition at line 466 of file ReactorNet.h.
|
protected |
Definition at line 467 of file ReactorNet.h.
|
protected |
Map used for default name generation.
Definition at line 468 of file ReactorNet.h.
|
protected |
Definition at line 469 of file ReactorNet.h.
|
protected |
The independent variable in the system.
May be either time or space depending on the type of reactors in the network.
Definition at line 473 of file ReactorNet.h.
|
protected |
The initial value of the independent variable in the system.
Definition at line 476 of file ReactorNet.h.
|
protected |
True if the integrator has been initialized at least once.
Definition at line 478 of file ReactorNet.h.
|
protected |
True if integrator needs to be (re)initialized.
Definition at line 479 of file ReactorNet.h.
|
protected |
Definition at line 480 of file ReactorNet.h.
|
protected |
Definition at line 482 of file ReactorNet.h.
|
protected |
Definition at line 483 of file ReactorNet.h.
|
protected |
Definition at line 484 of file ReactorNet.h.
|
protected |
Definition at line 485 of file ReactorNet.h.
|
protected |
Definition at line 486 of file ReactorNet.h.
|
protected |
True if scalar absolute tolerance was user-specified.
Definition at line 488 of file ReactorNet.h.
|
protected |
Definition at line 489 of file ReactorNet.h.
|
protected |
Definition at line 490 of file ReactorNet.h.
|
protected |
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition at line 493 of file ReactorNet.h.
|
protected |
Definition at line 495 of file ReactorNet.h.
|
protected |
Indicates whether time or space is the independent variable.
Definition at line 498 of file ReactorNet.h.
|
protected |
Names corresponding to each sensitivity parameter.
Definition at line 501 of file ReactorNet.h.
|
protected |
Definition at line 503 of file ReactorNet.h.
|
protected |
Definition at line 504 of file ReactorNet.h.
|
protected |
Definition at line 505 of file ReactorNet.h.
|
protected |
Base state used for evaluating advance limits during a single advance() call when root-finding is enabled.
Definition at line 508 of file ReactorNet.h.
|
protected |
Base time corresponding to m_ybase.
Definition at line 510 of file ReactorNet.h.
|
protected |
Indicates whether the advance-limit root check is active for the current call to advance(t, applylimit=true)
Definition at line 513 of file ReactorNet.h.
|
protected |
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition at line 516 of file ReactorNet.h.
|
protected |
Definition at line 517 of file ReactorNet.h.