13 ReactorNet::ReactorNet() :
14 m_integ(0), m_time(0.0), m_init(false), m_integrator_init(false),
15 m_nv(0), m_rtol(1.0e-9), m_rtolsens(1.0e-4),
16 m_atols(1.0e-15), m_atolsens(1.0e-4),
17 m_maxstep(-1.0), m_maxErrTestFails(0),
18 m_verbose(false), m_ntotpar(0)
20 m_integ = newIntegrator(
"CVODE");
26 m_integ->setProblemType(DENSE + NOJAC);
30 ReactorNet::~ReactorNet()
32 for (
size_t n = 0; n < m_reactors.size(); n++) {
46 writelog(
"Initializing reactor network.\n", m_verbose);
47 if (m_reactors.empty())
49 "no reactors in network!");
50 size_t sensParamNumber = 0;
52 for (n = 0; n < m_reactors.size(); n++) {
58 for (
size_t i = 0; i < sens_objs.size(); i++) {
59 std::map<size_t, size_t>& s =
m_sensOrder[sens_objs[i]];
60 for (std::map<size_t, size_t>::iterator iter = s.begin();
71 sprintf(buf,
"Reactor %s: %s variables.\n",
74 sprintf(buf,
" %s sensitivity params.\n",
78 if (r.
type() == FlowReactorType && m_reactors.size() > 1) {
80 "FlowReactors must be used alone.");
84 m_ydot.resize(
m_nv,0.0);
86 fill(m_atol.begin(), m_atol.end(), m_atols);
92 sprintf(buf,
"Number of equations: %s\n",
int2str(
neq()).c_str());
94 sprintf(buf,
"Maximum time step: %14.6g\n", m_maxstep);
98 m_integrator_init =
true;
105 writelog(
"Re-initializing reactor network.\n", m_verbose);
106 m_integ->reinitialize(m_time, *
this);
107 m_integrator_init =
true;
116 if (m_maxstep < 0.0) {
117 m_maxstep = time - m_time;
120 }
else if (!m_integrator_init) {
131 if (m_maxstep < 0.0) {
132 m_maxstep = time - m_time;
135 }
else if (!m_integrator_init) {
138 m_time = m_integ->
step(time);
146 "To be removed after Cantera 2.2. Use 'addReactor(Reactor&) instead'.");
149 "Ownership of Reactors by ReactorNet is deprecated.");
152 if (r->
type() >= ReactorType) {
153 m_reactors.push_back(r);
154 m_iown.push_back(iown);
158 ", since type = "+
int2str(r->
type())+
"\n", m_verbose);
165 m_reactors.push_back(&r);
166 m_iown.push_back(
false);
170 doublereal* ydot, doublereal* p)
176 for (n = 0; n < m_reactors.size(); n++) {
177 m_reactors[n]->evalEqs(t, y +
m_start[n],
178 ydot +
m_start[n], p + pstart);
179 pstart += m_nparams[n];
184 doublereal* ydot, doublereal* p,
Array2D* j)
186 doublereal ysave, dy;
191 for (
size_t n = 0; n <
m_nv; n++) {
195 dy = m_atol[n] + fabs(ysave)*m_rtol;
203 for (
size_t m = 0; m <
m_nv; m++) {
204 jac(m,n) = (m_ydot[m] - ydot[m])/dy;
212 for (
size_t n = 0; n < m_reactors.size(); n++) {
213 m_reactors[n]->updateState(y +
m_start[n]);
218 size_t leny, doublereal* y)
220 for (
size_t n = 0; n < m_reactors.size(); n++) {
235 size_t reactionIndex,
const std::string& name,
int leftright)
237 if (m_integrator_init) {
238 throw CanteraError(
"ReactorNet::registerSensitivityReaction",
239 "Sensitivity reactions cannot be added after the"
240 "integrator has been initialized.");
242 std::pair<void*, int> R = std::make_pair(reactor, leftright);
245 throw CanteraError(
"ReactorNet::registerSensitivityReaction",
246 "Attempted to register duplicate sensitivity reaction");
Backward Differentiation.
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 ...
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
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
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
Header file for class Wall.
virtual size_t neq()
Number of equations.
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
virtual void integrate(doublereal tout)
Integrate the system 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.
virtual void setTolerances(doublereal reltol, size_t n, doublereal *abstol)
Set or reset the number of equations.
size_t m_nv
True if integrator initialization is current.
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
virtual void initialize(doublereal t0, FuncEval &func)
Initialize the integrator for a new problem.
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
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.
Base class for exceptions thrown by Cantera classes.
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.
std::string name() const
Return the name of this reactor.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
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...
virtual void setSensitivityTolerances(doublereal reltol, doublereal abstol)
Set the sensitivity error tolerances.
virtual int type() const
Return a constant indicating the type of this Reactor.
virtual doublereal & solution(size_t k)
The current value of the solution of equation k.
virtual doublereal step(doublereal tout)
Integrate the system of equations.
doublereal time()
Current value of the simulation time.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
std::vector< std::pair< void *, int > > getSensitivityOrder() const
Return a vector specifying the ordering of objects to use when determining sensitivity parameter indi...
virtual size_t neq()
Number of equations (state variables) for this reactor.
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...
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
void writelog(const std::string &msg)
Write a message to the screen.
void advance(doublereal time)
Advance the state of all reactors in time.
Class Reactor is a general-purpose class for stirred reactors.