Cantera  3.1.0
Loading...
Searching...
No Matches
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
11#include "cantera/base/Array.h"
14
15#include <cstdio>
16
17namespace Cantera
18{
19
20ReactorNet::ReactorNet()
21{
22 suppressErrors(true);
23}
24
25ReactorNet::~ReactorNet()
26{
27}
28
30{
31 m_time = time;
33 m_integrator_init = false;
34}
35
36void ReactorNet::setMaxTimeStep(double maxstep)
37{
38 m_maxstep = maxstep;
40}
41
43{
45}
46
47void ReactorNet::setTolerances(double rtol, double atol)
48{
49 if (rtol >= 0.0) {
50 m_rtol = rtol;
51 }
52 if (atol >= 0.0) {
53 m_atols = atol;
54 }
55 m_init = false;
56}
57
58void ReactorNet::setSensitivityTolerances(double rtol, double atol)
59{
60 if (rtol >= 0.0) {
61 m_rtolsens = rtol;
62 }
63 if (atol >= 0.0) {
64 m_atolsens = atol;
65 }
66 m_init = false;
67}
68
71 return m_time;
72 } else {
73 throw CanteraError("ReactorNet::time", "Time is not the independent variable"
74 " for this reactor network.");
75 }
76}
77
80 return m_time;
81 } else {
82 throw CanteraError("ReactorNet::distance", "Distance is not the independent"
83 " variable for this reactor network.");
84 }
85}
86
88{
89 m_nv = 0;
90 debuglog("Initializing reactor network.\n", m_verbose);
91 if (m_reactors.empty()) {
92 throw CanteraError("ReactorNet::initialize",
93 "no reactors in network!");
94 }
95 m_start.assign(1, 0);
96 for (size_t n = 0; n < m_reactors.size(); n++) {
97 Reactor& r = *m_reactors[n];
99 size_t nv = r.neq();
100 m_nv += nv;
101 m_start.push_back(m_nv);
102
103 if (m_verbose) {
104 writelog("Reactor {:d}: {:d} variables.\n", n, nv);
105 writelog(" {:d} sensitivity params.\n", r.nSensParams());
106 }
107 if (r.type() == "FlowReactor" && m_reactors.size() > 1) {
108 throw CanteraError("ReactorNet::initialize",
109 "FlowReactors must be used alone.");
110 }
111 }
112
113 m_ydot.resize(m_nv,0.0);
114 m_yest.resize(m_nv,0.0);
115 m_advancelimits.resize(m_nv,-1.0);
116 m_atol.resize(neq());
117 fill(m_atol.begin(), m_atol.end(), m_atols);
118 m_integ->setTolerances(m_rtol, neq(), m_atol.data());
119 m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
120 if (!m_linearSolverType.empty()) {
121 m_integ->setLinearSolverType(m_linearSolverType);
122 }
123 if (m_precon) {
124 m_integ->setPreconditioner(m_precon);
125 }
126 m_integ->initialize(m_time, *this);
127 if (m_verbose) {
128 writelog("Number of equations: {:d}\n", neq());
129 writelog("Maximum time step: {:14.6g}\n", m_maxstep);
130 }
131 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
133 }
134 m_integrator_init = true;
135 m_init = true;
136}
137
139{
140 if (m_init) {
141 debuglog("Re-initializing reactor network.\n", m_verbose);
142 m_integ->reinitialize(m_time, *this);
143 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
145 }
146 m_integrator_init = true;
147 } else {
148 initialize();
149 }
150}
151
152void ReactorNet::setLinearSolverType(const string& linSolverType)
153{
154 m_linearSolverType = linSolverType;
155 m_integrator_init = false;
156}
157
158void ReactorNet::setPreconditioner(shared_ptr<PreconditionerBase> preconditioner)
159{
160 m_precon = preconditioner;
161 m_integrator_init = false;
162}
163
165{
166 integrator().setMaxSteps(nmax);
167}
168
170{
171 return integrator().maxSteps();
172}
173
174void ReactorNet::advance(double time)
175{
176 if (!m_init) {
177 initialize();
178 } else if (!m_integrator_init) {
179 reinitialize();
180 }
181 m_integ->integrate(time);
182 m_time = time;
183 updateState(m_integ->solution());
184}
185
186double ReactorNet::advance(double time, bool applylimit)
187{
188 if (!m_init) {
189 initialize();
190 } else if (!m_integrator_init) {
191 reinitialize();
192 }
193
194 if (!applylimit) {
195 // take full step
196 advance(time);
197 return time;
198 }
199
200 if (!hasAdvanceLimits()) {
201 // take full step
202 advance(time);
203 return time;
204 }
205
206 getAdvanceLimits(m_advancelimits.data());
207
208 // ensure that gradient is available
209 while (lastOrder() < 1) {
210 step();
211 }
212
213 int k = lastOrder();
214 double t = time, delta;
215 double* y = m_integ->solution();
216
217 // reduce time step if limits are exceeded
218 while (true) {
219 bool exceeded = false;
220 getEstimate(t, k, &m_yest[0]);
221 for (size_t j = 0; j < m_nv; j++) {
222 delta = abs(m_yest[j] - y[j]);
223 if ( (m_advancelimits[j] > 0.) && ( delta > m_advancelimits[j]) ) {
224 exceeded = true;
225 if (m_verbose) {
226 writelog(" Limiting global state vector component {:d} (dt = {:9.4g}):"
227 "{:11.6g} > {:9.4g}\n",
228 j, t - m_time, delta, m_advancelimits[j]);
229 }
230 }
231 }
232 if (!exceeded) {
233 break;
234 }
235 t = .5 * (m_time + t);
236 }
237 advance(t);
238 return t;
239}
240
242{
243 if (!m_init) {
244 initialize();
245 } else if (!m_integrator_init) {
246 reinitialize();
247 }
248 m_time = m_integ->step(m_time + 1.0);
249 updateState(m_integ->solution());
250 return m_time;
251}
252
253void ReactorNet::getEstimate(double time, int k, double* yest)
254{
255 if (!m_init) {
256 initialize();
257 }
258 // initialize
259 double* cvode_dky = m_integ->solution();
260 for (size_t j = 0; j < m_nv; j++) {
261 yest[j] = cvode_dky[j];
262 }
263
264 // Taylor expansion
265 double factor = 1.;
266 double deltat = time - m_time;
267 for (int n = 1; n <= k; n++) {
268 factor *= deltat / n;
269 cvode_dky = m_integ->derivative(m_time, n);
270 for (size_t j = 0; j < m_nv; j++) {
271 yest[j] += factor * cvode_dky[j];
272 }
273 }
274}
275
277{
278 if (m_integ) {
279 return m_integ->lastOrder();
280 } else {
281 return 0;
282 }
283}
284
286{
287 for (auto current : m_reactors) {
288 if (current->isOde() != r.isOde()) {
289 throw CanteraError("ReactorNet::addReactor",
290 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
291 current->type(), r.type());
292 }
293 if (current->timeIsIndependent() != r.timeIsIndependent()) {
294 throw CanteraError("ReactorNet::addReactor",
295 "Cannot mix Reactor types using time and space as independent variables"
296 "\n({} and {})", current->type(), r.type());
297 }
298 }
300 r.setNetwork(this);
301 m_reactors.push_back(&r);
302 if (!m_integ) {
303 m_integ.reset(newIntegrator(r.isOde() ? "CVODE" : "IDA"));
304 // use backward differencing, with a full Jacobian computed
305 // numerically, and use a Newton linear iterator
306 m_integ->setMethod(BDF_Method);
307 m_integ->setLinearSolverType("DENSE");
308 }
309 updateNames(r);
310}
311
313{
314 // ensure that reactors and components have reproducible names
316
317 for (size_t i=0; i<r.nWalls(); i++) {
318 auto& w = r.wall(i);
320 if (w.left().type() == "Reservoir") {
321 w.left().setDefaultName(m_counts);
322 }
323 if (w.right().type() == "Reservoir") {
324 w.right().setDefaultName(m_counts);
325 }
326 }
327
328 for (size_t i=0; i<r.nInlets(); i++) {
329 auto& in = r.inlet(i);
331 if (in.in().type() == "Reservoir") {
332 in.in().setDefaultName(m_counts);
333 }
334 }
335
336 for (size_t i=0; i<r.nOutlets(); i++) {
337 auto& out = r.outlet(i);
339 if (out.out().type() == "Reservoir") {
340 out.out().setDefaultName(m_counts);
341 }
342 }
343
344 for (size_t i=0; i<r.nSurfs(); i++) {
346 }
347}
348
350 if (m_integ == nullptr) {
351 throw CanteraError("ReactorNet::integrator",
352 "Integrator has not been instantiated. Add one or more reactors first.");
353 }
354 return *m_integ;
355}
356
357void ReactorNet::eval(double t, double* y, double* ydot, double* p)
358{
359 m_time = t;
360 updateState(y);
361 m_LHS.assign(m_nv, 1);
362 m_RHS.assign(m_nv, 0);
363 for (size_t n = 0; n < m_reactors.size(); n++) {
364 m_reactors[n]->applySensitivity(p);
365 m_reactors[n]->eval(t, m_LHS.data() + m_start[n], m_RHS.data() + m_start[n]);
366 size_t yEnd = 0;
367 if (n == m_reactors.size() - 1) {
368 yEnd = m_RHS.size();
369 } else {
370 yEnd = m_start[n + 1];
371 }
372 for (size_t i = m_start[n]; i < yEnd; i++) {
373 ydot[i] = m_RHS[i] / m_LHS[i];
374 }
375 m_reactors[n]->resetSensitivity(p);
376 }
377 checkFinite("ydot", ydot, m_nv);
378}
379
380void ReactorNet::evalDae(double t, double* y, double* ydot, double* p, double* residual)
381{
382 m_time = t;
383 updateState(y);
384 for (size_t n = 0; n < m_reactors.size(); n++) {
385 m_reactors[n]->applySensitivity(p);
386 m_reactors[n]->evalDae(t, y, ydot, residual);
387 m_reactors[n]->resetSensitivity(p);
388 }
389 checkFinite("ydot", ydot, m_nv);
390}
391
392void ReactorNet::getConstraints(double* constraints)
393{
394 for (size_t n = 0; n < m_reactors.size(); n++) {
395 m_reactors[n]->getConstraints(constraints + m_start[n]);
396 }
397}
398
399double ReactorNet::sensitivity(size_t k, size_t p)
400{
401 if (!m_init) {
402 initialize();
403 }
404 if (p >= m_sens_params.size()) {
405 throw IndexError("ReactorNet::sensitivity",
406 "m_sens_params", p, m_sens_params.size()-1);
407 }
408 double denom = m_integ->solution(k);
409 if (denom == 0.0) {
410 denom = SmallNumber;
411 }
412 return m_integ->sensitivity(k, p) / denom;
413}
414
415void ReactorNet::evalJacobian(double t, double* y, double* ydot, double* p, Array2D* j)
416{
417 //evaluate the unperturbed ydot
418 eval(t, y, ydot, p);
419 for (size_t n = 0; n < m_nv; n++) {
420 // perturb x(n)
421 double ysave = y[n];
422 double dy = m_atol[n] + fabs(ysave)*m_rtol;
423 y[n] = ysave + dy;
424 dy = y[n] - ysave;
425
426 // calculate perturbed residual
427 eval(t, y, m_ydot.data(), p);
428
429 // compute nth column of Jacobian
430 for (size_t m = 0; m < m_nv; m++) {
431 j->value(m,n) = (m_ydot[m] - ydot[m])/dy;
432 }
433 y[n] = ysave;
434 }
435}
436
438{
439 checkFinite("y", y, m_nv);
440 for (size_t n = 0; n < m_reactors.size(); n++) {
441 m_reactors[n]->updateState(y + m_start[n]);
442 }
443}
444
445void ReactorNet::getDerivative(int k, double* dky)
446{
447 if (!m_init) {
448 initialize();
449 }
450 double* cvode_dky = m_integ->derivative(m_time, k);
451 for (size_t j = 0; j < m_nv; j++) {
452 dky[j] = cvode_dky[j];
453 }
454}
455
456void ReactorNet::setAdvanceLimits(const double *limits)
457{
458 if (!m_init) {
459 initialize();
460 }
461 for (size_t n = 0; n < m_reactors.size(); n++) {
462 m_reactors[n]->setAdvanceLimits(limits + m_start[n]);
463 }
464}
465
467{
468 bool has_limit = false;
469 for (size_t n = 0; n < m_reactors.size(); n++) {
470 has_limit |= m_reactors[n]->hasAdvanceLimits();
471 }
472 return has_limit;
473}
474
475bool ReactorNet::getAdvanceLimits(double *limits) const
476{
477 bool has_limit = false;
478 for (size_t n = 0; n < m_reactors.size(); n++) {
479 has_limit |= m_reactors[n]->getAdvanceLimits(limits + m_start[n]);
480 }
481 return has_limit;
482}
483
484void ReactorNet::getState(double* y)
485{
486 for (size_t n = 0; n < m_reactors.size(); n++) {
487 m_reactors[n]->getState(y + m_start[n]);
488 }
489}
490
491void ReactorNet::getStateDae(double* y, double* ydot)
492{
493 for (size_t n = 0; n < m_reactors.size(); n++) {
494 m_reactors[n]->getStateDae(y + m_start[n], ydot + m_start[n]);
495 }
496}
497
498size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
499{
500 if (!m_init) {
501 initialize();
502 }
503 return m_start[reactor] + m_reactors[reactor]->componentIndex(component);
504}
505
506string ReactorNet::componentName(size_t i) const
507{
508 for (auto r : m_reactors) {
509 if (i < r->neq()) {
510 return r->name() + ": " + r->componentName(i);
511 } else {
512 i -= r->neq();
513 }
514 }
515 throw CanteraError("ReactorNet::componentName", "Index out of bounds");
516}
517
519 const string& name, double value, double scale)
520{
521 if (m_integrator_init) {
522 throw CanteraError("ReactorNet::registerSensitivityParameter",
523 "Sensitivity parameters cannot be added after the"
524 "integrator has been initialized.");
525 }
526 m_paramNames.push_back(name);
527 m_sens_params.push_back(value);
528 m_paramScales.push_back(scale);
529 return m_sens_params.size() - 1;
530}
531
533{
534 // Apply given settings to all reactors
535 for (size_t i = 0; i < m_reactors.size(); i++) {
536 m_reactors[i]->setDerivativeSettings(settings);
537 }
538}
539
541{
542 if (m_integ) {
543 return m_integ->solverStats();
544 } else {
545 return AnyMap();
546 }
547}
548
550{
551 if (m_integ) {
552 return m_integ->linearSolverType();
553 } else {
554 return "";
555 }
556}
557
558void ReactorNet::preconditionerSolve(double* rhs, double* output)
559{
560 if (!m_integ) {
561 throw CanteraError("ReactorNet::preconditionerSolve",
562 "Must only be called after ReactorNet is initialized.");
563 }
564 m_integ->preconditionerSolve(m_nv, rhs, output);
565}
566
567void ReactorNet::preconditionerSetup(double t, double* y, double gamma)
568{
569 // ensure state is up to date.
570 updateState(y);
571 // get the preconditioner
572 auto precon = m_integ->preconditioner();
573 // Reset preconditioner
574 precon->reset();
575 // Set gamma value for M =I - gamma*J
576 precon->setGamma(gamma);
577 // Make a copy of state to adjust it for preconditioner
578 vector<double> yCopy(m_nv);
579 // Get state of reactor
580 getState(yCopy.data());
581 // transform state based on preconditioner rules
582 precon->stateAdjustment(yCopy);
583 // update network with adjusted state
584 updateState(yCopy.data());
585 // Get jacobians and give elements to preconditioners
586 for (size_t i = 0; i < m_reactors.size(); i++) {
587 Eigen::SparseMatrix<double> rJac = m_reactors[i]->jacobian();
588 for (int k=0; k<rJac.outerSize(); ++k) {
589 for (Eigen::SparseMatrix<double>::InnerIterator it(rJac, k); it; ++it) {
590 precon->setValue(it.row() + m_start[i], it.col() + m_start[i],
591 it.value());
592 }
593 }
594 }
595 // post reactor setup operations
596 precon->setup();
597}
598
600{
601 if (!m_integ) {
602 throw CanteraError("ReactorNet::updatePreconditioner",
603 "Must only be called after ReactorNet is initialized.");
604 }
605 auto precon = m_integ->preconditioner();
606 precon->setGamma(gamma);
607 precon->updatePreconditioner();
608}
609
611 // check for non-mole-based reactors and throw an error otherwise
612 for (auto reactor : m_reactors) {
614 throw CanteraError("ReactorNet::checkPreconditionerSupported",
615 "Preconditioning is only supported for type *MoleReactor,\n"
616 "Reactor type given: '{}'.",
617 reactor->type());
618 }
619 }
620}
621
622}
Header file for class Cantera::Array2D.
Header file for class ReactorSurface.
Header file for base class WallBase.
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition Array.h:32
double & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Definition Array.h:160
Base class for exceptions thrown by Cantera classes.
bool setDefaultName(map< string, int > &counts)
Set the default name of a flow device. Returns false if it was previously set.
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
Definition FuncEval.h:184
bool suppressErrors() const
Get current state of error suppression.
Definition FuncEval.h:166
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition FuncEval.h:181
An array index is out of range.
Abstract base class for ODE system integrators.
Definition Integrator.h:44
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
Definition Integrator.h:221
virtual void setMaxSteps(int nmax)
Set the maximum number of time-steps the integrator can take before reaching the next output time.
Definition Integrator.h:239
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
Definition Integrator.h:245
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition Integrator.h:231
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
size_t nWalls()
Return the number of Wall objects connected to this reactor.
WallBase & wall(size_t n)
Return a reference to the n-th Wall connected to this reactor.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
virtual size_t nSurfs()
Return the number of surfaces in a reactor.
size_t nOutlets()
Return the number of outlet FlowDevice objects connected to this reactor.
size_t nInlets()
Return the number of inlet FlowDevice objects connected to this reactor.
ReactorSurface * surface(size_t n)
Return a reference to the n-th ReactorSurface connected to this reactor.
void setPreconditioner(shared_ptr< PreconditionerBase > preconditioner)
Set preconditioner used by the linear solver.
void setLinearSolverType(const string &linSolverType="DENSE")
Set the type of linear solver used in the integration.
void preconditionerSetup(double t, double *y, double gamma) override
Evaluate the setup processes for the Jacobian preconditioner.
double step()
Advance the state of all reactors with respect to the independent variable (time or space).
virtual int lastOrder() const
Returns the order used for last solution step of the ODE integrator The function is intended for inte...
void eval(double t, double *y, double *ydot, double *p) override
Evaluate the right-hand-side ODE function.
void initialize()
Initialize the reactor network.
void advance(double t)
Advance the state of all reactors in the independent variable (time or space).
size_t neq() const override
Number of equations.
Definition ReactorNet.h:208
vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
Definition ReactorNet.h:337
vector< double > m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition ReactorNet.h:363
double m_initial_time
The initial value of the independent variable in the system.
Definition ReactorNet.h:330
void evalJacobian(double t, double *y, double *ydot, double *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void getConstraints(double *constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
double m_time
The independent variable in the system.
Definition ReactorNet.h:327
AnyMap solverStats() const
Get solver stats from integrator.
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition ReactorNet.h:144
map< string, int > m_counts
Map used for default name generation.
Definition ReactorNet.h:322
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration steps the integrator will take before reaching the nex...
string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
void getStateDae(double *y, double *ydot) override
Fill in the vectors y and ydot with the current state of the system.
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
virtual void getDerivative(int k, double *dky)
Return k-th derivative at the current state of the system.
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single step.
size_t registerSensitivityParameter(const string &name, double value, double scale)
Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the React...
double m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition ReactorNet.h:348
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
int maxSteps()
Returns the maximum number of internal integration steps the integrator will take before reaching the...
virtual void setDerivativeSettings(AnyMap &settings)
Set derivative settings of all reactors.
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.
void updateNames(Reactor &r)
Create reproducible names for reactors and walls/connectors.
void updateState(double *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
void getState(double *y) override
Fill in the vector y with the current state of the system.
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
double rtol()
Relative tolerance.
Definition ReactorNet.h:90
size_t globalComponentIndex(const string &component, size_t reactor=0)
Return the index corresponding to the component named component in the reactor with index reactor in ...
bool m_timeIsIndependent
Indicates whether time or space is the independent variable.
Definition ReactorNet.h:353
double atol()
Absolute integration tolerance.
Definition ReactorNet.h:95
bool hasAdvanceLimits() const
Check whether ReactorNet object uses advance limits.
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
void evalDae(double t, double *y, double *ydot, double *p, double *residual) override
eval coupling for IDA / DAEs
virtual void checkPreconditionerSupported() const
Check that preconditioning is supported by all reactors in the network.
bool m_integrator_init
True if integrator initialization is current.
Definition ReactorNet.h:333
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
string linearSolverType() const
Problem type of integrator.
void updatePreconditioner(double gamma) override
Update the preconditioner based on already computed jacobian values.
void preconditionerSolve(double *rhs, double *output) override
Evaluate the linear system Ax=b where A is the preconditioner.
vector< string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition ReactorNet.h:356
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
virtual void getEstimate(double time, int k, double *yest)
Estimate a future state based on current derivatives.
bool setDefaultName(map< string, int > &counts)
Set the default name of a wall. Returns false if it was previously set.
Class Reactor is a general-purpose class for stirred reactors.
Definition Reactor.h:47
size_t neq()
Number of equations (state variables) for this reactor.
Definition Reactor.h:107
string type() const override
String indicating the reactor model implemented.
Definition Reactor.h:52
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor (including walls)
Definition Reactor.cpp:122
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:94
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
Definition Reactor.h:242
virtual bool isOde() const
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs.
Definition Reactor.h:59
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
Definition Reactor.h:65
bool setDefaultName(map< string, int > &counts)
Set the default name of a wall. Returns false if it was previously set.
Definition Wall.cpp:13
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Definition global.h:154
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:171
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Definition utilities.h:104
Integrator * newIntegrator(const string &itype)
Create new Integrator object.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
@ BDF_Method
Backward Differentiation.
Definition Integrator.h:24
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:158
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...