Cantera  2.5.1
ReactorNet.cpp
Go to the documentation of this file.
1 //! @file ReactorNet.cpp
2 
3 // This file is part of Cantera. See License.txt in the top-level directory or
4 // at https://cantera.org/license.txt for license and copyright information.
5 
8 #include "cantera/zeroD/Wall.h"
9 
10 #include <cstdio>
11 
12 using namespace std;
13 
14 namespace Cantera
15 {
16 
17 ReactorNet::ReactorNet() :
18  m_integ(newIntegrator("CVODE")),
19  m_time(0.0), m_init(false), m_integrator_init(false),
20  m_nv(0), m_rtol(1.0e-9), m_rtolsens(1.0e-4),
21  m_atols(1.0e-15), m_atolsens(1.0e-6),
22  m_maxstep(0.0), m_maxErrTestFails(0),
23  m_verbose(false)
24 {
25  suppressErrors(true);
26 
27  // use backward differencing, with a full Jacobian computed
28  // numerically, and use a Newton linear iterator
29  m_integ->setMethod(BDF_Method);
30  m_integ->setProblemType(DENSE + NOJAC);
31 }
32 
33 void ReactorNet::setInitialTime(double time)
34 {
35  m_time = time;
36  m_integrator_init = false;
37 }
38 
39 void ReactorNet::setMaxTimeStep(double maxstep)
40 {
41  m_maxstep = maxstep;
42  m_init = false;
43 }
44 
45 void ReactorNet::setMaxErrTestFails(int nmax)
46 {
47  m_maxErrTestFails = nmax;
48  m_init = false;
49 }
50 
51 void ReactorNet::setTolerances(double rtol, double atol)
52 {
53  if (rtol >= 0.0) {
54  m_rtol = rtol;
55  }
56  if (atol >= 0.0) {
57  m_atols = atol;
58  }
59  m_init = false;
60 }
61 
62 void ReactorNet::setSensitivityTolerances(double rtol, double atol)
63 {
64  if (rtol >= 0.0) {
65  m_rtolsens = rtol;
66  }
67  if (atol >= 0.0) {
68  m_atolsens = atol;
69  }
70  m_init = false;
71 }
72 
73 void ReactorNet::initialize()
74 {
75  m_nv = 0;
76  debuglog("Initializing reactor network.\n", m_verbose);
77  if (m_reactors.empty()) {
78  throw CanteraError("ReactorNet::initialize",
79  "no reactors in network!");
80  }
81  m_start.assign(1, 0);
82  for (size_t n = 0; n < m_reactors.size(); n++) {
83  Reactor& r = *m_reactors[n];
84  r.initialize(m_time);
85  size_t nv = r.neq();
86  m_nv += nv;
87  m_start.push_back(m_nv);
88 
89  if (m_verbose) {
90  writelog("Reactor {:d}: {:d} variables.\n", n, nv);
91  writelog(" {:d} sensitivity params.\n", r.nSensParams());
92  }
93  if (r.typeStr() == "FlowReactor" && m_reactors.size() > 1) {
94  throw CanteraError("ReactorNet::initialize",
95  "FlowReactors must be used alone.");
96  }
97  }
98 
99  m_ydot.resize(m_nv,0.0);
100  m_yest.resize(m_nv,0.0);
101  m_advancelimits.resize(m_nv,-1.0);
102  m_atol.resize(neq());
103  fill(m_atol.begin(), m_atol.end(), m_atols);
104  m_integ->setTolerances(m_rtol, neq(), m_atol.data());
105  m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
106  m_integ->setMaxStepSize(m_maxstep);
107  m_integ->setMaxErrTestFails(m_maxErrTestFails);
108  if (m_verbose) {
109  writelog("Number of equations: {:d}\n", neq());
110  writelog("Maximum time step: {:14.6g}\n", m_maxstep);
111  }
112  m_integ->initialize(m_time, *this);
113  m_integrator_init = true;
114  m_init = true;
115 }
116 
117 void ReactorNet::reinitialize()
118 {
119  if (m_init) {
120  debuglog("Re-initializing reactor network.\n", m_verbose);
121  m_integ->reinitialize(m_time, *this);
122  m_integrator_init = true;
123  } else {
124  initialize();
125  }
126 }
127 
128 void ReactorNet::advance(doublereal time)
129 {
130  if (!m_init) {
131  initialize();
132  } else if (!m_integrator_init) {
133  reinitialize();
134  }
135  m_integ->integrate(time);
136  m_time = time;
137  updateState(m_integ->solution());
138 }
139 
140 double ReactorNet::advance(double time, bool applylimit)
141 {
142  if (!m_init) {
143  initialize();
144  } else if (!m_integrator_init) {
145  reinitialize();
146  }
147 
148  if (!applylimit) {
149  // take full step
150  advance(time);
151  return time;
152  }
153 
154  if (!hasAdvanceLimits()) {
155  // take full step
156  advance(time);
157  return time;
158  }
159 
160  getAdvanceLimits(m_advancelimits.data());
161 
162  // ensure that gradient is available
163  while (lastOrder() < 1) {
164  step();
165  }
166 
167  int k = lastOrder();
168  double t = time, delta;
169  double* y = m_integ->solution();
170 
171  // reduce time step if limits are exceeded
172  while (true) {
173  bool exceeded = false;
174  getEstimate(t, k, &m_yest[0]);
175  for (size_t j = 0; j < m_nv; j++) {
176  delta = abs(m_yest[j] - y[j]);
177  if ( (m_advancelimits[j] > 0.) && ( delta > m_advancelimits[j]) ) {
178  exceeded = true;
179  if (m_verbose) {
180  writelog(" Limiting global state vector component {:d} (dt = {:9.4g}):"
181  "{:11.6g} > {:9.4g}\n",
182  j, t - m_time, delta, m_advancelimits[j]);
183  }
184  }
185  }
186  if (!exceeded) {
187  break;
188  }
189  t = .5 * (m_time + t);
190  }
191  advance(t);
192  return t;
193 }
194 
195 double ReactorNet::step()
196 {
197  if (!m_init) {
198  initialize();
199  } else if (!m_integrator_init) {
200  reinitialize();
201  }
202  m_time = m_integ->step(m_time + 1.0);
203  updateState(m_integ->solution());
204  return m_time;
205 }
206 
207 void ReactorNet::getEstimate(double time, int k, double* yest)
208 {
209  // initialize
210  double* cvode_dky = m_integ->solution();
211  for (size_t j = 0; j < m_nv; j++) {
212  yest[j] = cvode_dky[j];
213  }
214 
215  // Taylor expansion
216  double factor = 1.;
217  double deltat = time - m_time;
218  for (int n = 1; n <= k; n++) {
219  factor *= deltat / n;
220  cvode_dky = m_integ->derivative(m_time, n);
221  for (size_t j = 0; j < m_nv; j++) {
222  yest[j] += factor * cvode_dky[j];
223  }
224  }
225 }
226 
227 void ReactorNet::addReactor(Reactor& r)
228 {
229  r.setNetwork(this);
230  m_reactors.push_back(&r);
231 }
232 
233 void ReactorNet::eval(doublereal t, doublereal* y,
234  doublereal* ydot, doublereal* p)
235 {
236  m_time = t; // This will be replaced at the end of the timestep
237  updateState(y);
238  for (size_t n = 0; n < m_reactors.size(); n++) {
239  m_reactors[n]->evalEqs(t, y + m_start[n], ydot + m_start[n], p);
240  }
241  checkFinite("ydot", ydot, m_nv);
242 }
243 
244 double ReactorNet::sensitivity(size_t k, size_t p)
245 {
246  if (!m_init) {
247  initialize();
248  }
249  if (p >= m_sens_params.size()) {
250  throw IndexError("ReactorNet::sensitivity",
251  "m_sens_params", p, m_sens_params.size()-1);
252  }
253  double denom = m_integ->solution(k);
254  if (denom == 0.0) {
255  denom = SmallNumber;
256  }
257  return m_integ->sensitivity(k, p) / denom;
258 }
259 
260 void ReactorNet::evalJacobian(doublereal t, doublereal* y,
261  doublereal* ydot, doublereal* p, Array2D* j)
262 {
263  //evaluate the unperturbed ydot
264  eval(t, y, ydot, p);
265  for (size_t n = 0; n < m_nv; n++) {
266  // perturb x(n)
267  double ysave = y[n];
268  double dy = m_atol[n] + fabs(ysave)*m_rtol;
269  y[n] = ysave + dy;
270  dy = y[n] - ysave;
271 
272  // calculate perturbed residual
273  eval(t, y, m_ydot.data(), p);
274 
275  // compute nth column of Jacobian
276  for (size_t m = 0; m < m_nv; m++) {
277  j->value(m,n) = (m_ydot[m] - ydot[m])/dy;
278  }
279  y[n] = ysave;
280  }
281 }
282 
283 void ReactorNet::updateState(doublereal* y)
284 {
285  checkFinite("y", y, m_nv);
286  for (size_t n = 0; n < m_reactors.size(); n++) {
287  m_reactors[n]->updateState(y + m_start[n]);
288  }
289 }
290 
291 void ReactorNet::getState(double* y)
292 {
293  for (size_t n = 0; n < m_reactors.size(); n++) {
294  m_reactors[n]->getState(y + m_start[n]);
295  }
296 }
297 
298 void ReactorNet::getDerivative(int k, double* dky)
299 {
300  double* cvode_dky = m_integ->derivative(m_time, k);
301  for (size_t j = 0; j < m_nv; j++) {
302  dky[j] = cvode_dky[j];
303  }
304 }
305 
306 void ReactorNet::setAdvanceLimits(const double *limits)
307 {
308  if (!m_init) {
309  initialize();
310  }
311  for (size_t n = 0; n < m_reactors.size(); n++) {
312  m_reactors[n]->setAdvanceLimits(limits + m_start[n]);
313  }
314 }
315 
316 bool ReactorNet::hasAdvanceLimits()
317 {
318  bool has_limit = false;
319  for (size_t n = 0; n < m_reactors.size(); n++) {
320  has_limit |= m_reactors[n]->hasAdvanceLimits();
321  }
322  return has_limit;
323 }
324 
325 bool ReactorNet::getAdvanceLimits(double *limits)
326 {
327  bool has_limit = false;
328  for (size_t n = 0; n < m_reactors.size(); n++) {
329  has_limit |= m_reactors[n]->getAdvanceLimits(limits + m_start[n]);
330  }
331  return has_limit;
332 }
333 
334 size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
335 {
336  if (!m_init) {
337  initialize();
338  }
339  return m_start[reactor] + m_reactors[reactor]->componentIndex(component);
340 }
341 
342 std::string ReactorNet::componentName(size_t i) const
343 {
344  for (auto r : m_reactors) {
345  if (i < r->neq()) {
346  return r->name() + ": " + r->componentName(i);
347  } else {
348  i -= r->neq();
349  }
350  }
351  throw CanteraError("ReactorNet::componentName", "Index out of bounds");
352 }
353 
354 size_t ReactorNet::registerSensitivityParameter(
355  const std::string& name, double value, double scale)
356 {
357  if (m_integrator_init) {
358  throw CanteraError("ReactorNet::registerSensitivityParameter",
359  "Sensitivity parameters cannot be added after the"
360  "integrator has been initialized.");
361  }
362  m_paramNames.push_back(name);
363  m_sens_params.push_back(value);
364  m_paramScales.push_back(scale);
365  return m_sens_params.size() - 1;
366 }
367 
368 }
Header file for base class WallBase.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:32
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Definition: Array.h:231
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
An array index is out of range.
Definition: ctexceptions.h:158
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
Definition: ReactorBase.cpp:95
Class Reactor is a general-purpose class for stirred reactors.
Definition: Reactor.h:38
virtual std::string typeStr() const
String indicating the reactor model implemented.
Definition: Reactor.h:42
virtual size_t neq()
Number of equations (state variables) for this reactor.
Definition: Reactor.h:97
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
Definition: Reactor.cpp:75
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
Definition: Reactor.cpp:105
void debuglog(const std::string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Definition: global.h:140
const double SmallNumber
smallest number to compare to zero.
Definition: ct_defs.h:149
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition: global.h:158
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
Definition: checkFinite.cpp:15
@ BDF_Method
Backward Differentiation.
Definition: Integrator.h:33
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition: utilities.h:135