115 return *m_reactors[n];
175 doublereal* ydot, doublereal* p,
Array2D* j);
181 virtual void eval(doublereal t, doublereal* y,
182 doublereal* ydot, doublereal* p);
184 virtual void getState(doublereal* y);
268 std::vector<Reactor*> m_reactors;
269 std::unique_ptr<Integrator> m_integ;
279 doublereal m_rtol, m_rtolsens;
280 doublereal m_atols, m_atolsens;
285 int m_maxErrTestFails;
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Virtual base class for ODE right-hand-side function evaluators.
vector_fp m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
bool suppressErrors() const
Get current state of error suppression.
Abstract base class for ODE system integrators.
A class representing a network of connected reactors.
std::vector< bool > m_have_deprecated_eval
virtual size_t nparams()
Number of sensitivity parameters.
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()
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
void initialize()
Initialize the reactor network.
void advance(doublereal time)
Advance the state of all reactors in time.
bool getAdvanceLimits(double *limits)
Retrieve absolute step size limits during advance.
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
const std::string & sensitivityParameterName(size_t p)
The name of the p-th sensitivity parameter added to this ReactorNet.
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
doublereal rtol()
Relative tolerance.
doublereal m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration time-steps the integrator will take before reaching th...
std::string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
double maxTimeStep()
Get the maximum time step.
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
virtual void getState(doublereal *y)
Fill in the vector y with the current state of the system.
virtual size_t neq()
Number of equations.
void setInitialTime(double time)
Set initial time.
virtual void getDerivative(int k, double *dky)
Return k-th derivative at the current time.
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single time ste...
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 ...
doublereal atolSensitivity() const
Absolute sensitivity tolerance.
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
virtual int maxSteps()
Returns the maximum number of internal integration time-steps the integrator will take before reachin...
virtual int lastOrder()
Returns the order used for last solution step of the ODE integrator The function is intended for inte...
doublereal rtolSensitivity() const
Relative 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.
bool hasAdvanceLimits()
Check whether ReactorNet object uses advance limits.
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
void evalJacobian(doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
doublereal time()
Current value of the simulation time.
bool verbose() const
Returns true if verbose logging output is enabled.
bool m_checked_eval_deprecation
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 setMaxTimeStep(double maxstep)
Set the maximum time step.
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...
bool m_integrator_init
True if integrator initialization is current.
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
doublereal atol()
Absolute integration tolerance.
vector_fp m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
virtual void getEstimate(double time, int k, double *yest)
Estimate a future state based on current derivatives.
Class Reactor is a general-purpose class for stirred reactors.
Namespace for the Cantera kernel.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.