6 #ifndef CT_REACTORNET_H 7 #define CT_REACTORNET_H 97 return *m_reactors[n];
157 doublereal* ydot, doublereal* p,
Array2D* j);
163 virtual void eval(doublereal t, doublereal* y,
164 doublereal* ydot, doublereal* p);
169 virtual void getState(doublereal* y);
219 std::vector<Reactor*> m_reactors;
220 std::unique_ptr<Integrator> m_integ;
230 doublereal m_rtol, m_rtolsens;
231 doublereal m_atols, m_atolsens;
236 int m_maxErrTestFails;
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 setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
void updateState(doublereal *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
double step(doublereal time=-999)
Advance the state of all reactors in time.
std::vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
doublereal rtolSensitivity() const
Relative sensitivity tolerance.
bool verbose() const
Returns true if verbose logging output is enabled.
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.
void setMaxTimeStep(double maxstep)
Set the maximum time step.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
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 React...
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
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.
std::string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
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.
bool m_integrator_init
True if integrator initialization is current.
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
vector_fp m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
virtual void getState(doublereal *y)
Fill in the vector y with the current state of the system.
doublereal m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
bool suppressErrors() const
Get current state of error suppression.
virtual size_t nparams()
Number of sensitivity parameters.
doublereal rtol()
Relative tolerance.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
doublereal time()
Current value of the simulation time.
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
void setInitialTime(double time)
Set initial time.
Virtual base class for ODE right-hand-side function evaluators.
doublereal atolSensitivity() const
Absolute sensitivity tolerance.
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...
Namespace for the Cantera kernel.
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.