Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ReactorNet.cpp
Go to the documentation of this file.
1 //! @file ReactorNet.cpp
4 #include "cantera/zeroD/Wall.h"
5 
6 #include <cstdio>
7 
8 using namespace std;
9 
10 namespace Cantera
11 {
12 
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)
19 {
20  m_integ = newIntegrator("CVODE");
21 
22  // use backward differencing, with a full Jacobian computed
23  // numerically, and use a Newton linear iterator
24 
25  m_integ->setMethod(BDF_Method);
26  m_integ->setProblemType(DENSE + NOJAC);
27  m_integ->setIterator(Newton_Iter);
28 }
29 
30 ReactorNet::~ReactorNet()
31 {
32  for (size_t n = 0; n < m_reactors.size(); n++) {
33  if (m_iown[n]) {
34  delete m_reactors[n];
35  }
36  m_reactors[n] = 0;
37  }
38  delete m_integ;
39 }
40 
42 {
43  size_t n, nv;
44  char buf[100];
45  m_nv = 0;
46  writelog("Initializing reactor network.\n", m_verbose);
47  if (m_reactors.empty())
48  throw CanteraError("ReactorNet::initialize",
49  "no reactors in network!");
50  size_t sensParamNumber = 0;
51  m_start.assign(1, 0);
52  for (n = 0; n < m_reactors.size(); n++) {
53  Reactor& r = *m_reactors[n];
54  r.initialize(m_time);
55  nv = r.neq();
56  m_nparams.push_back(r.nSensParams());
57  std::vector<std::pair<void*, int> > sens_objs = r.getSensitivityOrder();
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();
61  iter != s.end();
62  ++iter) {
63  m_sensIndex.resize(std::max(iter->second + 1, m_sensIndex.size()));
64  m_sensIndex[iter->second] = sensParamNumber++;
65  }
66  }
67  m_nv += nv;
68  m_start.push_back(m_nv);
69 
70  if (m_verbose) {
71  sprintf(buf,"Reactor %s: %s variables.\n",
72  int2str(n).c_str(), int2str(nv).c_str());
73  writelog(buf);
74  sprintf(buf," %s sensitivity params.\n",
75  int2str(r.nSensParams()).c_str());
76  writelog(buf);
77  }
78  if (r.type() == FlowReactorType && m_reactors.size() > 1) {
79  throw CanteraError("ReactorNet::initialize",
80  "FlowReactors must be used alone.");
81  }
82  }
83 
84  m_ydot.resize(m_nv,0.0);
85  m_atol.resize(neq());
86  fill(m_atol.begin(), m_atol.end(), m_atols);
87  m_integ->setTolerances(m_rtol, neq(), DATA_PTR(m_atol));
88  m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
89  m_integ->setMaxStepSize(m_maxstep);
90  m_integ->setMaxErrTestFails(m_maxErrTestFails);
91  if (m_verbose) {
92  sprintf(buf, "Number of equations: %s\n", int2str(neq()).c_str());
93  writelog(buf);
94  sprintf(buf, "Maximum time step: %14.6g\n", m_maxstep);
95  writelog(buf);
96  }
97  m_integ->initialize(m_time, *this);
98  m_integrator_init = true;
99  m_init = true;
100 }
101 
103 {
104  if (m_init) {
105  writelog("Re-initializing reactor network.\n", m_verbose);
106  m_integ->reinitialize(m_time, *this);
107  m_integrator_init = true;
108  } else {
109  initialize();
110  }
111 }
112 
113 void ReactorNet::advance(doublereal time)
114 {
115  if (!m_init) {
116  if (m_maxstep < 0.0) {
117  m_maxstep = time - m_time;
118  }
119  initialize();
120  } else if (!m_integrator_init) {
121  reinitialize();
122  }
123  m_integ->integrate(time);
124  m_time = time;
125  updateState(m_integ->solution());
126 }
127 
128 double ReactorNet::step(doublereal time)
129 {
130  if (!m_init) {
131  if (m_maxstep < 0.0) {
132  m_maxstep = time - m_time;
133  }
134  initialize();
135  } else if (!m_integrator_init) {
136  reinitialize();
137  }
138  m_time = m_integ->step(time);
139  updateState(m_integ->solution());
140  return m_time;
141 }
142 
143 void ReactorNet::addReactor(Reactor* r, bool iown)
144 {
145  warn_deprecated("ReactorNet::addReactor(Reactor*)",
146  "To be removed after Cantera 2.2. Use 'addReactor(Reactor&) instead'.");
147  if (iown) {
148  warn_deprecated("ReactorNet::addReactor",
149  "Ownership of Reactors by ReactorNet is deprecated.");
150  }
151  r->setNetwork(this);
152  if (r->type() >= ReactorType) {
153  m_reactors.push_back(r);
154  m_iown.push_back(iown);
155  writelog("Adding reactor "+r->name()+"\n", m_verbose);
156  } else {
157  writelog("Not adding reactor "+r->name()+
158  ", since type = "+int2str(r->type())+"\n", m_verbose);
159  }
160 }
161 
163 {
164  r.setNetwork(this);
165  m_reactors.push_back(&r);
166  m_iown.push_back(false);
167 }
168 
169 void ReactorNet::eval(doublereal t, doublereal* y,
170  doublereal* ydot, doublereal* p)
171 {
172  size_t n;
173  size_t pstart = 0;
174 
175  updateState(y);
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];
180  }
181 }
182 
183 void ReactorNet::evalJacobian(doublereal t, doublereal* y,
184  doublereal* ydot, doublereal* p, Array2D* j)
185 {
186  doublereal ysave, dy;
187  Array2D& jac = *j;
188 
189  //evaluate the unperturbed ydot
190  eval(t, y, ydot, p);
191  for (size_t n = 0; n < m_nv; n++) {
192 
193  // perturb x(n)
194  ysave = y[n];
195  dy = m_atol[n] + fabs(ysave)*m_rtol;
196  y[n] = ysave + dy;
197  dy = y[n] - ysave;
198 
199  // calculate perturbed residual
200  eval(t, y, DATA_PTR(m_ydot), p);
201 
202  // compute nth column of Jacobian
203  for (size_t m = 0; m < m_nv; m++) {
204  jac(m,n) = (m_ydot[m] - ydot[m])/dy;
205  }
206  y[n] = ysave;
207  }
208 }
209 
210 void ReactorNet::updateState(doublereal* y)
211 {
212  for (size_t n = 0; n < m_reactors.size(); n++) {
213  m_reactors[n]->updateState(y + m_start[n]);
214  }
215 }
216 
218  size_t leny, doublereal* y)
219 {
220  for (size_t n = 0; n < m_reactors.size(); n++) {
221  m_reactors[n]->getInitialConditions(t0, m_start[n+1]-m_start[n],
222  y + m_start[n]);
223  }
224 }
225 
226 size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
227 {
228  if (!m_init) {
229  initialize();
230  }
231  return m_start[reactor] + m_reactors[reactor]->componentIndex(component);
232 }
233 
235  size_t reactionIndex, const std::string& name, int leftright)
236 {
237  if (m_integrator_init) {
238  throw CanteraError("ReactorNet::registerSensitivityReaction",
239  "Sensitivity reactions cannot be added after the"
240  "integrator has been initialized.");
241  }
242  std::pair<void*, int> R = std::make_pair(reactor, leftright);
243  if (m_sensOrder.count(R) &&
244  m_sensOrder[R].count(reactionIndex)) {
245  throw CanteraError("ReactorNet::registerSensitivityReaction",
246  "Attempted to register duplicate sensitivity reaction");
247  }
248  m_paramNames.push_back(name);
249  m_sensOrder[R][reactionIndex] = m_ntotpar;
250  m_ntotpar++;
251 }
252 
253 }
Backward Differentiation.
Definition: Integrator.h:32
std::vector< size_t > m_sensIndex
Mapping from the order in which sensitivity parameters were added to the ReactorNet to the order in w...
Definition: ReactorNet.h:272
double step(doublereal time)
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:128
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 ...
Definition: ReactorNet.cpp:226
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
Definition: ReactorNet.cpp:162
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition: ReactorNet.h:263
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:39
void updateState(doublereal *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
Definition: ReactorNet.cpp:210
std::vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
Definition: ReactorNet.h:251
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
Definition: Reactor.cpp:107
Header file for class Wall.
virtual size_t neq()
Number of equations.
Definition: ReactorNet.h:197
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
Definition: Integrator.h:175
virtual void integrate(doublereal tout)
Integrate the system of equations.
Definition: Integrator.h:121
void reinitialize()
Reinitialize the integrator.
Definition: ReactorNet.cpp:102
void initialize()
Initialize the reactor network.
Definition: ReactorNet.cpp:41
void evalJacobian(doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
Definition: ReactorNet.cpp:183
virtual void setTolerances(doublereal reltol, size_t n, doublereal *abstol)
Set or reset the number of equations.
Definition: Integrator.h:72
size_t m_nv
True if integrator initialization is current.
Definition: ReactorNet.h:248
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:78
virtual void initialize(doublereal t0, FuncEval &func)
Initialize the integrator for a new problem.
Definition: Integrator.h:108
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
Definition: Reactor.cpp:68
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:29
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Fill the solution vector with the initial conditions at initial time t0.
Definition: ReactorNet.cpp:217
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
Definition: ReactorNet.cpp:169
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition: ReactorNet.h:129
std::string name() const
Return the name of this reactor.
Definition: ReactorBase.h:42
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
Definition: ReactorBase.cpp:82
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...
Definition: ReactorNet.cpp:234
virtual void setSensitivityTolerances(doublereal reltol, doublereal abstol)
Set the sensitivity error tolerances.
Definition: Integrator.h:91
virtual int type() const
Return a constant indicating the type of this Reactor.
Definition: Reactor.h:44
virtual doublereal & solution(size_t k)
The current value of the solution of equation k.
Definition: Integrator.h:136
Newton Iteration.
Definition: Integrator.h:42
virtual doublereal step(doublereal tout)
Integrate the system of equations.
Definition: Integrator.h:130
doublereal time()
Current value of the simulation time.
Definition: ReactorNet.h:79
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
std::vector< std::pair< void *, int > > getSensitivityOrder() const
Return a vector specifying the ordering of objects to use when determining sensitivity parameter indi...
Definition: Reactor.cpp:341
virtual size_t neq()
Number of equations (state variables) for this reactor.
Definition: Reactor.h:93
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...
Definition: ReactorNet.h:267
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition: Integrator.h:185
void writelog(const std::string &msg)
Write a message to the screen.
Definition: global.cpp:33
void advance(doublereal time)
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:113
Class Reactor is a general-purpose class for stirred reactors.
Definition: Reactor.h:39