14 ReactorNet::ReactorNet() : Cantera::FuncEval(), m_nr(0), m_nreactors(0),
15 m_integ(0), m_time(0.0), m_init(false),
16 m_nv(0), m_rtol(1.0e-9), m_rtolsens(1.0e-4),
17 m_atols(1.0e-15), m_atolsens(1.0e-4),
18 m_maxstep(-1.0), m_maxErrTestFails(0),
19 m_verbose(false), m_ntotpar(0)
24 m_integ = newIntegrator(
"CVODE");
30 m_integ->setProblemType(DENSE + NOJAC);
34 ReactorNet::~ReactorNet()
36 for (
size_t n = 0; n < m_nr; n++) {
54 writelog(
"Initializing reactor network.\n", m_verbose);
57 "no reactors in network!");
58 size_t sensParamNumber = 0;
59 for (n = 0; n < m_nr; n++) {
60 if (m_r[n]->type() >= ReactorType) {
61 m_r[n]->initialize(m_time);
63 m_reactors.push_back(r);
68 for (
size_t i = 0; i < sens_objs.size(); i++) {
69 std::map<size_t, size_t>& s =
m_sensOrder[sens_objs[i]];
70 for (std::map<size_t, size_t>::iterator iter = s.begin();
81 sprintf(buf,
"Reactor %s: %s variables.\n",
84 sprintf(buf,
" %s sensitivity params.\n",
88 if (m_r[n]->type() == FlowReactorType && m_nr > 1) {
90 "FlowReactors must be used alone.");
95 m_connect.resize(m_nr*m_nr,0);
96 m_ydot.resize(m_nv,0.0);
97 size_t i, j, nin, nout, nw;
99 for (i = 0; i < m_nr; i++) {
101 for (j = 0; j < m_nr; j++) {
107 for (n = 0; n < nin; n++) {
113 for (n = 0; n < nout; n++) {
119 for (n = 0; n < nw; n++) {
132 m_atol.resize(
neq());
133 fill(m_atol.begin(), m_atol.end(), m_atols);
139 sprintf(buf,
"Number of equations: %s\n",
int2str(
neq()).c_str());
141 sprintf(buf,
"Maximum time step: %14.6g\n", m_maxstep);
151 if (m_maxstep < 0.0) {
152 m_maxstep = time - m_time;
164 if (m_maxstep < 0.0) {
165 m_maxstep = time - m_time;
169 m_time = m_integ->
step(time);
177 if (r->
type() >= ReactorType) {
179 m_iown.push_back(iown);
184 ", since type = "+
int2str(r->
type())+
"\n", m_verbose);
189 doublereal* ydot, doublereal* p)
196 for (n = 0; n < m_nreactors; n++) {
197 m_reactors[n]->evalEqs(t, y + start,
198 ydot + start, p + pstart);
200 pstart += m_nparams[n];
205 doublereal* ydot, doublereal* p,
Array2D* j)
207 doublereal ysave, dy;
215 for (
size_t n = 0; n < m_nv; n++) {
219 dy = m_atol[n] + fabs(ysave)*m_rtol;
227 for (
size_t m = 0; m < m_nv; m++) {
228 jac(m,n) = (m_ydot[m] - ydot[m])/dy;
233 std::cerr << err.
what() << std::endl;
234 error(
"Terminating execution.");
241 for (
size_t n = 0; n < m_nreactors; n++) {
242 m_reactors[n]->updateState(y + start);
248 size_t leny, doublereal* y)
251 for (
size_t n = 0; n < m_nreactors; n++) {
252 m_reactors[n]->getInitialConditions(t0, m_size[n], y + start);
264 for (n = 0; n <
reactor; n++) {
267 return start + m_reactors[n]->componentIndex(species);
271 size_t reactionIndex,
const std::string& name,
int leftright)
273 std::pair<void*, int> R = std::make_pair(reactor, leftright);
276 throw CanteraError(
"ReactorNet::registerSensitivityReaction",
277 "Attempted to register duplicate sensitivity reaction");
Backward Differentiation.
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
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.
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
const ReactorBase & right()
Return a reference to the Reactor or Reservoir to the right of the wall.
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...
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
Header file for class Wall.
const ReactorBase & out() const
Return a const reference to the downstream reactor.
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 ...
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
virtual void integrate(doublereal tout)
Integrate the system of equations.
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 nInlets()
Return the number of inlet FlowDevice objects connected to this reactor.
virtual void initialize(doublereal t0, FuncEval &func)
Initialize the integrator for a new problem.
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.
const char * what() const
Get a description of the error.
void addReactor(ReactorBase *r, bool iown=false)
Add the reactor r to this reactor network.
size_t nWalls()
Return the number of Wall objects connected to this reactor.
Wall & wall(size_t n)
Return a reference to the n-th Wall connected to this reactor.
void error(const std::string &msg)
Write an error message and quit.
ReactorBase & in() const
Return a reference to the upstream reactor.
Base class for exceptions thrown by Cantera classes.
Base class for stirred reactors.
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
std::string name() const
Return the name of this reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to 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 int type() const
Return a constant indicating the type of this Reactor.
virtual void setSensitivityTolerances(doublereal reltol, doublereal abstol)
Set the sensitivity error tolerances.
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...
size_t nOutlets()
Return the number of outlet FlowDevice objects connected to this reactor.
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.
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
ReactorBase & left() const
Return a reference to the Reactor or Reservoir to the left of the wall.
void advance(doublereal time)
Advance the state of all reactors in time.
Class Reactor is a general-purpose class for stirred reactors.