6 #ifndef CT_REACTORNET_H
7 #define CT_REACTORNET_H
51 m_maxErrTestFails = nmax;
171 doublereal* ydot, doublereal* p,
Array2D* j);
177 virtual void eval(doublereal t, doublereal* y,
178 doublereal* ydot, doublereal* p);
194 const std::string& name,
int leftright=0);
201 void connect(
size_t i,
size_t j) {
202 m_connect[j*m_nr + i] = 1;
203 m_connect[i*m_nr + j] = 1;
206 bool connected(
size_t i,
size_t j) {
207 return (m_connect[m_nr*i + j] == 1);
217 std::vector<ReactorBase*> m_r;
218 std::vector<Reactor*> m_reactors;
225 std::vector<size_t> m_size;
227 doublereal m_rtol, m_rtolsens;
228 doublereal m_atols, m_atolsens;
229 doublereal m_maxstep;
230 int m_maxErrTestFails;
233 std::vector<size_t> m_nparams;
240 std::map<std::pair<void*, int>, std::map<size_t, size_t> >
m_sensOrder;
250 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...
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
double step(doublereal time)
Advance the state of all reactors in time.
A class representing a network of connected reactors.
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...
double sensitivity(const std::string &species, size_t p, int reactor=0)
Return the sensitivity of the species named species with respect to the p-th sensitivity parameter...
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.
size_t globalComponentIndex(const std::string &species, size_t reactor=0)
Return the index corresponding to the species named species in the reactor with index reactor in the ...
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 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.
std::vector< int > vector_int
Vector of ints.
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
doublereal atol()
Absolute integration tolerance.
void addReactor(ReactorBase *r, bool iown=false)
Add the reactor r to this reactor network.
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.
Base class for stirred reactors.
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
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.
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.