Cantera
2.3.0
|
A class representing a network of connected reactors. More...
#include <ReactorNet.h>
Public Member Functions | |
ReactorNet (const ReactorNet &)=delete | |
void | addReactor (Reactor &r) |
Add the reactor r to this reactor network. More... | |
Reactor & | reactor (int n) |
Return a reference to the n-th reactor in this network. More... | |
bool | verbose () const |
Returns true if verbose logging output is enabled. More... | |
void | setVerbose (bool v=true) |
Enable or disable verbose logging while setting up and integrating the reactor network. More... | |
Integrator & | integrator () |
Return a reference to the integrator. More... | |
void | updateState (doublereal *y) |
Update the state of all the reactors in the network to correspond to the values in the solution vector y. More... | |
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. More... | |
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 parameter. More... | |
void | evalJacobian (doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j) |
Evaluate the Jacobian matrix for the reactor network. More... | |
virtual size_t | neq () |
Number of equations. More... | |
virtual void | eval (doublereal t, doublereal *y, doublereal *ydot, doublereal *p) |
Evaluate the right-hand-side function. More... | |
virtual void | getInitialConditions (doublereal t0, size_t leny, doublereal *y) |
virtual void | getState (doublereal *y) |
Fill in the vector y with the current state of the system. More... | |
virtual size_t | nparams () |
Number of sensitivity parameters. More... | |
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 the global state vector for the reactor network. More... | |
std::string | componentName (size_t i) const |
Return the name of the i-th component of the global state vector. More... | |
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 ReactorNet can keep track of the order in which sensitivity parameters are added. More... | |
const std::string & | sensitivityParameterName (size_t p) |
The name of the p-th sensitivity parameter added to this ReactorNet. More... | |
void | reinitialize () |
Reinitialize the integrator. More... | |
void | setNeedsReinit () |
Called to trigger integrator reinitialization before further integration. More... | |
Methods to set up a simulation. | |
void | setInitialTime (double time) |
Set initial time. More... | |
void | setMaxTimeStep (double maxstep) |
Set the maximum time step. More... | |
void | setMaxErrTestFails (int nmax) |
Set the maximum number of error test failures permitted by the CVODES integrator in a single time step. More... | |
void | setTolerances (double rtol, double atol) |
Set the relative and absolute tolerances for the integrator. More... | |
void | setSensitivityTolerances (double rtol, double atol) |
Set the relative and absolute tolerances for integrating the sensitivity equations. More... | |
doublereal | time () |
Current value of the simulation time. More... | |
doublereal | rtol () |
Relative tolerance. More... | |
doublereal | atol () |
Absolute integration tolerance. More... | |
doublereal | rtolSensitivity () const |
Relative sensitivity tolerance. More... | |
doublereal | atolSensitivity () const |
Absolute sensitivity tolerance. More... | |
void | advance (doublereal time) |
Advance the state of all reactors in time. More... | |
double | step (doublereal time=-999) |
Advance the state of all reactors in time. More... | |
Public Member Functions inherited from FuncEval | |
int | eval_nothrow (double t, double *y, double *ydot) |
Evaluate the right-hand side using return code to indicate status. More... | |
void | suppressErrors (bool suppress) |
Enable or disable suppression of errors when calling eval() More... | |
bool | suppressErrors () const |
Get current state of error suppression. More... | |
std::string | getErrors () const |
Return a string containing the text of any suppressed errors. More... | |
void | clearErrors () |
Clear any previously-stored suppressed errors. More... | |
Protected Member Functions | |
void | initialize () |
Initialize the reactor network. More... | |
Protected Attributes | |
std::vector< Reactor * > | m_reactors |
std::unique_ptr< Integrator > | m_integ |
doublereal | m_time |
bool | m_init |
bool | m_integrator_init |
True if integrator initialization is current. More... | |
size_t | m_nv |
std::vector< size_t > | m_start |
m_start[n] is the starting point in the state vector for reactor n More... | |
vector_fp | m_atol |
doublereal | m_rtol |
doublereal | m_rtolsens |
doublereal | m_atols |
doublereal | m_atolsens |
doublereal | m_maxstep |
Maximum integrator internal timestep. Default of 0.0 means infinity. More... | |
int | m_maxErrTestFails |
bool | m_verbose |
std::vector< std::string > | m_paramNames |
Names corresponding to each sensitivity parameter. More... | |
vector_fp | m_ydot |
Protected Attributes inherited from FuncEval | |
bool | m_suppress_errors |
std::vector< std::string > | m_errors |
Errors occuring during function evaluations. More... | |
Additional Inherited Members | |
Public Attributes inherited from FuncEval | |
vector_fp | 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(). More... | |
vector_fp | m_paramScales |
Scaling factors for each sensitivity parameter. More... | |
A class representing a network of connected reactors.
This class is used to integrate the time-dependent governing equations for a network of reactors (Reactor, ConstPressureReactor) connected by various means, e.g. Wall, MassFlowController, Valve, PressureController.
Definition at line 23 of file ReactorNet.h.
void setInitialTime | ( | double | time | ) |
Set initial time.
Default = 0.0 s. Restarts integration from this time using the current mixture state as the initial condition.
Definition at line 34 of file ReactorNet.cpp.
void setMaxTimeStep | ( | double | maxstep | ) |
Set the maximum time step.
Definition at line 40 of file ReactorNet.cpp.
void setMaxErrTestFails | ( | int | nmax | ) |
Set the maximum number of error test failures permitted by the CVODES integrator in a single time step.
Definition at line 46 of file ReactorNet.cpp.
void setTolerances | ( | double | rtol, |
double | atol | ||
) |
Set the relative and absolute tolerances for the integrator.
Definition at line 52 of file ReactorNet.cpp.
void setSensitivityTolerances | ( | double | rtol, |
double | atol | ||
) |
Set the relative and absolute tolerances for integrating the sensitivity equations.
Definition at line 63 of file ReactorNet.cpp.
|
inline |
Current value of the simulation time.
Definition at line 52 of file ReactorNet.h.
|
inline |
Relative tolerance.
Definition at line 57 of file ReactorNet.h.
|
inline |
Absolute integration tolerance.
Definition at line 62 of file ReactorNet.h.
|
inline |
Relative sensitivity tolerance.
Definition at line 67 of file ReactorNet.h.
|
inline |
Absolute sensitivity tolerance.
Definition at line 72 of file ReactorNet.h.
void advance | ( | doublereal | time | ) |
Advance the state of all reactors in time.
Take as many internal timesteps as necessary to reach time.
time | Time to advance to (s). |
Definition at line 127 of file ReactorNet.cpp.
double step | ( | doublereal | time = -999 | ) |
Advance the state of all reactors in time.
Definition at line 139 of file ReactorNet.cpp.
References Cantera::warn_deprecated().
void addReactor | ( | Reactor & | r | ) |
Add the reactor r to this reactor network.
Definition at line 155 of file ReactorNet.cpp.
References ReactorBase::setNetwork().
|
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 96 of file ReactorNet.h.
Referenced by ReactorNet::sensitivity().
|
inline |
Returns true
if verbose logging output is enabled.
Definition at line 101 of file ReactorNet.h.
|
inline |
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition at line 107 of file ReactorNet.h.
References FuncEval::suppressErrors().
|
inline |
Return a reference to the integrator.
Definition at line 113 of file ReactorNet.h.
void updateState | ( | doublereal * | y | ) |
Update the state of all the reactors in the network to correspond to the values in the solution vector y.
Definition at line 210 of file ReactorNet.cpp.
References Cantera::checkFinite().
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 171 of file ReactorNet.cpp.
References Cantera::SmallNumber.
Referenced by ReactorNet::sensitivity().
|
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 143 of file ReactorNet.h.
References ReactorNet::globalComponentIndex(), ReactorNet::reactor(), and ReactorNet::sensitivity().
void evalJacobian | ( | doublereal | t, |
doublereal * | y, | ||
doublereal * | ydot, | ||
doublereal * | p, | ||
Array2D * | j | ||
) |
Evaluate the Jacobian matrix for the reactor network.
[in] | t | Time at which to evaluate the Jacobian |
[in] | y | Global state vector at time t |
[out] | ydot | Time derivative of the state vector evaluated at t. |
[in] | p | sensitivity parameter vector (unused?) |
[out] | j | Jacobian matrix, size neq() by neq(). |
Definition at line 187 of file ReactorNet.cpp.
|
inlinevirtual |
|
virtual |
Evaluate the right-hand-side 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() |
Implements FuncEval.
Definition at line 161 of file ReactorNet.cpp.
References Cantera::checkFinite().
|
virtual |
Reimplemented from FuncEval.
Definition at line 218 of file ReactorNet.cpp.
References Cantera::warn_deprecated().
|
virtual |
Fill in the vector y with the current state of the system.
Reimplemented from FuncEval.
Definition at line 225 of file ReactorNet.cpp.
|
inlinevirtual |
Number of sensitivity parameters.
Reimplemented from FuncEval.
Definition at line 171 of file ReactorNet.h.
References FuncEval::m_sens_params.
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 the global state vector for the reactor network.
Definition at line 232 of file ReactorNet.cpp.
Referenced by ReactorNet::sensitivity().
std::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, e.g. ‘'reactor1: CH4’`.
Definition at line 240 of file ReactorNet.cpp.
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 ReactorNet can keep track of the order in which sensitivity parameters are added.
name | A name describing the parameter, e.g. the reaction string |
value | The nominal value of the parameter |
scale | A scaling factor to be applied to the sensitivity coefficient |
Definition at line 252 of file ReactorNet.cpp.
References Cantera::scale().
Referenced by Reactor::addSensitivityReaction(), and Reactor::addSensitivitySpeciesEnthalpy().
|
inline |
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition at line 198 of file ReactorNet.h.
References ReactorNet::m_paramNames.
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 116 of file ReactorNet.cpp.
References Cantera::debuglog().
|
inline |
Called to trigger integrator reinitialization before further integration.
Definition at line 210 of file ReactorNet.h.
References ReactorNet::m_integrator_init.
Referenced by ReactorBase::syncState().
|
protected |
Initialize the reactor network.
Called automatically the first time advance or step is called.
Definition at line 74 of file ReactorNet.cpp.
|
protected |
True if integrator initialization is current.
Definition at line 223 of file ReactorNet.h.
Referenced by ReactorNet::setNeedsReinit().
|
protected |
m_start[n] is the starting point in the state vector for reactor n
Definition at line 227 of file ReactorNet.h.
|
protected |
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition at line 234 of file ReactorNet.h.
|
protected |
Names corresponding to each sensitivity parameter.
Definition at line 240 of file ReactorNet.h.
Referenced by ReactorNet::sensitivityParameterName().