6 #ifndef CT_REACTORNET_H
7 #define CT_REACTORNET_H
39 m_integrator_init =
false;
51 m_maxErrTestFails = nmax;
130 return *m_reactors[n];
194 doublereal* ydot, doublereal* p,
Array2D* j);
200 virtual void eval(doublereal t, doublereal* y,
201 doublereal* ydot, doublereal* p);
217 const std::string& name,
int leftright=0);
233 m_integrator_init =
false;
243 std::vector<Reactor*> m_reactors;
247 bool m_integrator_init;
254 doublereal m_rtol, m_rtolsens;
255 doublereal m_atols, m_atolsens;
256 doublereal m_maxstep;
257 int m_maxErrTestFails;
260 std::vector<size_t> m_nparams;
267 std::map<std::pair<void*, int>, std::map<size_t, size_t> >
m_sensOrder;
276 std::vector<bool> m_iown;
std::vector< size_t > m_sensIndex
Mapping from the order in which sensitivity parameters were added to the ReactorNet to the order in w...
double step(doublereal time)
Advance the state of all reactors in time.
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 ...
A class representing a network of connected reactors.
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
void updateState(doublereal *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
std::vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
void setTolerances(doublereal rtol, doublereal atol)
Set the relative and absolute tolerances for the integrator.
doublereal atolSensitivity() const
Absolute sensitivity tolerance.
doublereal rtolSensitivity() const
Relative sensitivity tolerance.
virtual size_t neq()
Number of equations.
void reinitialize()
Reinitialize the integrator.
void initialize()
Initialize the reactor network.
void evalJacobian(doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
Integrator & integrator()
Return a reference to the integrator.
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
size_t m_nv
True if integrator initialization is current.
void setMaxTimeStep(double maxstep)
Set the maximum time step.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Fill the solution vector with the initial conditions at initial time t0.
Header file for class Cantera::Array2D.
double sensitivity(const std::string &component, size_t p, int reactor=0)
Return the sensitivity of the component named component with respect to the p-th sensitivity paramete...
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
doublereal atol()
Absolute integration tolerance.
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single time ste...
Abstract base class for ODE system integrators.
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
void setSensitivityTolerances(doublereal rtol, doublereal atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
virtual size_t nparams()
Number of sensitivity parameters.
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 Reacto...
doublereal rtol()
Relative tolerance.
bool verbose() const
Returns true if verbose logging output is enabled.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
virtual doublereal & solution(size_t k)
The current value of the solution of equation k.
doublereal time()
Current value of the simulation time.
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...
void setInitialTime(doublereal time)
Set initial time.
Virtual base class for ODE right-hand-side function evaluators.
An array index is out of range.
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...
const std::string & sensitivityParameterName(size_t p)
The name of the p-th sensitivity parameter added to this ReactorNet.
void advance(doublereal time)
Advance the state of all reactors in time.
Class Reactor is a general-purpose class for stirred reactors.