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),
30 m_integ->setProblemType(DENSE + NOJAC);
33 void ReactorNet::setInitialTime(
double time)
36 m_integrator_init =
false;
39 void ReactorNet::setMaxTimeStep(
double maxstep)
45 void ReactorNet::setMaxErrTestFails(
int nmax)
47 m_maxErrTestFails = nmax;
51 void ReactorNet::setTolerances(
double rtol,
double atol)
62 void ReactorNet::setSensitivityTolerances(
double rtol,
double atol)
73 void ReactorNet::initialize()
76 debuglog(
"Initializing reactor network.\n", m_verbose);
77 if (m_reactors.empty()) {
79 "no reactors in network!");
82 for (
size_t n = 0; n < m_reactors.size(); n++) {
87 m_start.push_back(m_nv);
90 writelog(
"Reactor {:d}: {:d} variables.\n", n, nv);
93 if (r.
typeStr() ==
"FlowReactor" && m_reactors.size() > 1) {
95 "FlowReactors must be used alone.");
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);
109 writelog(
"Number of equations: {:d}\n", neq());
110 writelog(
"Maximum time step: {:14.6g}\n", m_maxstep);
112 m_integ->initialize(m_time, *
this);
113 m_integrator_init =
true;
117 void ReactorNet::reinitialize()
120 debuglog(
"Re-initializing reactor network.\n", m_verbose);
121 m_integ->reinitialize(m_time, *
this);
122 m_integrator_init =
true;
128 void ReactorNet::advance(doublereal time)
132 }
else if (!m_integrator_init) {
135 m_integ->integrate(time);
137 updateState(m_integ->solution());
140 double ReactorNet::advance(
double time,
bool applylimit)
144 }
else if (!m_integrator_init) {
154 if (!hasAdvanceLimits()) {
160 getAdvanceLimits(m_advancelimits.data());
163 while (lastOrder() < 1) {
168 double t = time, delta;
169 double* y = m_integ->solution();
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]) ) {
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]);
189 t = .5 * (m_time + t);
195 double ReactorNet::step()
199 }
else if (!m_integrator_init) {
202 m_time = m_integ->step(m_time + 1.0);
203 updateState(m_integ->solution());
207 void ReactorNet::getEstimate(
double time,
int k,
double* yest)
210 double* cvode_dky = m_integ->solution();
211 for (
size_t j = 0; j < m_nv; j++) {
212 yest[j] = cvode_dky[j];
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];
230 m_reactors.push_back(&r);
233 void ReactorNet::eval(doublereal t, doublereal* y,
234 doublereal* ydot, doublereal* p)
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);
244 double ReactorNet::sensitivity(
size_t k,
size_t p)
249 if (p >= m_sens_params.size()) {
251 "m_sens_params", p, m_sens_params.size()-1);
253 double denom = m_integ->solution(k);
257 return m_integ->sensitivity(k, p) / denom;
260 void ReactorNet::evalJacobian(doublereal t, doublereal* y,
261 doublereal* ydot, doublereal* p,
Array2D* j)
265 for (
size_t n = 0; n < m_nv; n++) {
268 double dy = m_atol[n] + fabs(ysave)*m_rtol;
273 eval(t, y, m_ydot.data(), p);
276 for (
size_t m = 0; m < m_nv; m++) {
277 j->
value(m,n) = (m_ydot[m] - ydot[m])/dy;
283 void ReactorNet::updateState(doublereal* y)
286 for (
size_t n = 0; n < m_reactors.size(); n++) {
287 m_reactors[n]->updateState(y + m_start[n]);
291 void ReactorNet::getState(
double* y)
293 for (
size_t n = 0; n < m_reactors.size(); n++) {
294 m_reactors[n]->getState(y + m_start[n]);
298 void ReactorNet::getDerivative(
int k,
double* dky)
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];
306 void ReactorNet::setAdvanceLimits(
const double *limits)
311 for (
size_t n = 0; n < m_reactors.size(); n++) {
312 m_reactors[n]->setAdvanceLimits(limits + m_start[n]);
316 bool ReactorNet::hasAdvanceLimits()
318 bool has_limit =
false;
319 for (
size_t n = 0; n < m_reactors.size(); n++) {
320 has_limit |= m_reactors[n]->hasAdvanceLimits();
325 bool ReactorNet::getAdvanceLimits(
double *limits)
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]);
334 size_t ReactorNet::globalComponentIndex(
const string& component,
size_t reactor)
339 return m_start[reactor] + m_reactors[reactor]->componentIndex(component);
342 std::string ReactorNet::componentName(
size_t i)
const
344 for (
auto r : m_reactors) {
346 return r->name() +
": " + r->componentName(i);
351 throw CanteraError(
"ReactorNet::componentName",
"Index out of bounds");
354 size_t ReactorNet::registerSensitivityParameter(
355 const std::string& name,
double value,
double scale)
357 if (m_integrator_init) {
358 throw CanteraError(
"ReactorNet::registerSensitivityParameter",
359 "Sensitivity parameters cannot be added after the"
360 "integrator has been initialized.");
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;
Header file for base class WallBase.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Base class for exceptions thrown by Cantera classes.
An array index is out of range.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
Class Reactor is a general-purpose class for stirred reactors.
virtual std::string typeStr() const
String indicating the reactor model implemented.
virtual size_t neq()
Number of equations (state variables) for this reactor.
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
const double SmallNumber
smallest number to compare to zero.
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
Namespace for the Cantera kernel.
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
@ BDF_Method
Backward Differentiation.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.