Cantera
2.2.1
|
A class representing a network of connected reactors. More...
#include <ReactorNet.h>
Public Member Functions | |
void | addReactor (Reactor &r) |
Add the reactor r to this reactor network. More... | |
void | addReactor (Reactor *r, bool iown=false) |
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) |
Fill the solution vector with the initial conditions at initial time t0. 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... | |
void | registerSensitivityReaction (void *reactor, size_t reactionIndex, const std::string &name, int leftright=0) |
Used by Reactor and Wall objects to register the addition of sensitivity reactions 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 (doublereal 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 (doublereal rtol, doublereal atol) |
Set the relative and absolute tolerances for the integrator. More... | |
void | setSensitivityTolerances (doublereal rtol, doublereal 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) |
Advance the state of all reactors in time. More... | |
Protected Member Functions | |
void | initialize () |
Initialize the reactor network. More... | |
Protected Attributes | |
std::vector< Reactor * > | m_reactors |
Integrator * | m_integ |
doublereal | m_time |
bool | m_init |
bool | m_integrator_init |
size_t | m_nv |
True if integrator initialization is current. More... | |
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 |
int | m_maxErrTestFails |
bool | m_verbose |
size_t | m_ntotpar |
std::vector< size_t > | m_nparams |
std::vector< std::string > | m_paramNames |
Names corresponding to each sensitivity parameter. More... | |
std::map< std::pair< void *, int >, std::map< size_t, size_t > > | m_sensOrder |
Structure used to determine the order of sensitivity parameters m_sensOrder[Reactor or Wall, leftright][reaction number] = parameter index. More... | |
std::vector< size_t > | m_sensIndex |
Mapping from the order in which sensitivity parameters were added to the ReactorNet to the order in which they occur in the integrator output. More... | |
vector_fp | m_ydot |
std::vector< bool > | m_iown |
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.
|
inline |
Set initial time.
Default = 0.0 s. Restarts integration from this time using the current mixture state as the initial condition.
Definition at line 37 of file ReactorNet.h.
References ReactorNet::time().
|
inline |
Set the maximum time step.
Definition at line 43 of file ReactorNet.h.
|
inline |
Set the maximum number of error test failures permitted by the CVODES integrator in a single time step.
Definition at line 50 of file ReactorNet.h.
|
inline |
Set the relative and absolute tolerances for the integrator.
Definition at line 56 of file ReactorNet.h.
References ReactorNet::atol(), and ReactorNet::rtol().
|
inline |
Set the relative and absolute tolerances for integrating the sensitivity equations.
Definition at line 68 of file ReactorNet.h.
References ReactorNet::atol(), and ReactorNet::rtol().
|
inline |
Current value of the simulation time.
Definition at line 79 of file ReactorNet.h.
Referenced by ReactorNet::advance(), and ReactorNet::setInitialTime().
|
inline |
Relative tolerance.
Definition at line 84 of file ReactorNet.h.
Referenced by ReactorNet::setSensitivityTolerances(), and ReactorNet::setTolerances().
|
inline |
Absolute integration tolerance.
Definition at line 89 of file ReactorNet.h.
Referenced by ReactorNet::setSensitivityTolerances(), and ReactorNet::setTolerances().
|
inline |
Relative sensitivity tolerance.
Definition at line 94 of file ReactorNet.h.
|
inline |
Absolute sensitivity tolerance.
Definition at line 99 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 113 of file ReactorNet.cpp.
References ReactorNet::initialize(), Integrator::integrate(), ReactorNet::reinitialize(), Integrator::solution(), ReactorNet::time(), and ReactorNet::updateState().
double step | ( | doublereal | time | ) |
Advance the state of all reactors in time.
Take a single timestep toward time.
Definition at line 128 of file ReactorNet.cpp.
References ReactorNet::initialize(), ReactorNet::reinitialize(), Integrator::solution(), Integrator::step(), and ReactorNet::updateState().
void addReactor | ( | Reactor & | r | ) |
Add the reactor r to this reactor network.
Definition at line 162 of file ReactorNet.cpp.
References ReactorBase::setNetwork().
void addReactor | ( | Reactor * | r, |
bool | iown = false |
||
) |
Add the reactor r to this reactor network.
Definition at line 143 of file ReactorNet.cpp.
References Cantera::int2str(), ReactorBase::name(), ReactorBase::setNetwork(), Reactor::type(), Cantera::warn_deprecated(), and Cantera::writelog().
|
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 129 of file ReactorNet.h.
Referenced by ReactorNet::globalComponentIndex(), and ReactorNet::sensitivity().
|
inline |
Returns true
if verbose logging output is enabled.
Definition at line 134 of file ReactorNet.h.
|
inline |
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition at line 140 of file ReactorNet.h.
|
inline |
Return a reference to the integrator.
Definition at line 145 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 ReactorNet::m_start.
Referenced by ReactorNet::advance(), ReactorNet::eval(), and ReactorNet::step().
|
inline |
Return the sensitivity of the k-th solution component with respect to the p-th sensitivity parameter.
The normalized sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:
\[ S_{ki} = \frac{p_i}{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).
Definition at line 166 of file ReactorNet.h.
References ReactorNet::initialize(), ReactorNet::m_sensIndex, and Integrator::solution().
Referenced by ReactorNet::sensitivity().
|
inline |
Return the sensitivity of the component named component with respect to the p-th sensitivity parameter.
The normalized sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:
\[ S_{ki} = \frac{p_i}{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).
Definition at line 180 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 183 of file ReactorNet.cpp.
References DATA_PTR, ReactorNet::eval(), and ReactorNet::m_nv.
|
inlinevirtual |
Number of equations.
Implements FuncEval.
Definition at line 197 of file ReactorNet.h.
References ReactorNet::m_nv.
Referenced by ReactorNet::initialize().
|
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 169 of file ReactorNet.cpp.
References ReactorNet::m_start, and ReactorNet::updateState().
Referenced by ReactorNet::evalJacobian().
|
virtual |
Fill the solution vector with the initial conditions at initial time t0.
Implements FuncEval.
Definition at line 217 of file ReactorNet.cpp.
References ReactorNet::m_start.
|
inlinevirtual |
Number of sensitivity parameters.
Reimplemented from FuncEval.
Definition at line 204 of file ReactorNet.h.
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 226 of file ReactorNet.cpp.
References ReactorNet::initialize(), ReactorNet::m_start, and ReactorNet::reactor().
Referenced by ReactorNet::sensitivity().
void registerSensitivityReaction | ( | void * | reactor, |
size_t | reactionIndex, | ||
const std::string & | name, | ||
int | leftright = 0 |
||
) |
Used by Reactor and Wall objects to register the addition of sensitivity reactions so that the ReactorNet can keep track of the order in which sensitivity parameters are added.
Definition at line 234 of file ReactorNet.cpp.
References ReactorNet::m_paramNames, and ReactorNet::m_sensOrder.
Referenced by Reactor::addSensitivityReaction().
|
inline |
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition at line 220 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 102 of file ReactorNet.cpp.
References ReactorNet::initialize(), and Cantera::writelog().
Referenced by ReactorNet::advance(), and ReactorNet::step().
|
inline |
Called to trigger integrator reinitialization before further integration.
Definition at line 232 of file ReactorNet.h.
Referenced by ReactorBase::syncState().
|
protected |
Initialize the reactor network.
Called automatically the first time advance or step is called.
Definition at line 41 of file ReactorNet.cpp.
References DATA_PTR, Reactor::getSensitivityOrder(), Reactor::initialize(), Integrator::initialize(), Cantera::int2str(), ReactorNet::m_nv, ReactorNet::m_sensIndex, ReactorNet::m_sensOrder, ReactorNet::m_start, Reactor::neq(), ReactorNet::neq(), Reactor::nSensParams(), Integrator::setMaxErrTestFails(), Integrator::setMaxStepSize(), Integrator::setSensitivityTolerances(), Integrator::setTolerances(), Reactor::type(), and Cantera::writelog().
Referenced by ReactorNet::advance(), ReactorNet::globalComponentIndex(), ReactorNet::reinitialize(), ReactorNet::sensitivity(), and ReactorNet::step().
|
protected |
True if integrator initialization is current.
Definition at line 248 of file ReactorNet.h.
Referenced by ReactorNet::evalJacobian(), ReactorNet::initialize(), and ReactorNet::neq().
|
protected |
m_start[n] is the starting point in the state vector for reactor n
Definition at line 251 of file ReactorNet.h.
Referenced by ReactorNet::eval(), ReactorNet::getInitialConditions(), ReactorNet::globalComponentIndex(), ReactorNet::initialize(), and ReactorNet::updateState().
|
protected |
Names corresponding to each sensitivity parameter.
Definition at line 263 of file ReactorNet.h.
Referenced by ReactorNet::registerSensitivityReaction(), and ReactorNet::sensitivityParameterName().
|
protected |
Structure used to determine the order of sensitivity parameters m_sensOrder[Reactor or Wall, leftright][reaction number] = parameter index.
Definition at line 267 of file ReactorNet.h.
Referenced by ReactorNet::initialize(), and ReactorNet::registerSensitivityReaction().
|
protected |
Mapping from the order in which sensitivity parameters were added to the ReactorNet to the order in which they occur in the integrator output.
Definition at line 272 of file ReactorNet.h.
Referenced by ReactorNet::initialize(), and ReactorNet::sensitivity().