Cantera  3.2.0a4
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"
19
20#include <cstdio>
21
22namespace Cantera
23{
24
25ReactorNet::ReactorNet()
26{
27 suppressErrors(true);
28}
29
30ReactorNet::ReactorNet(shared_ptr<ReactorBase> reactor)
31{
32 suppressErrors(true);
34}
35
36ReactorNet::ReactorNet(vector<shared_ptr<ReactorBase>>& reactors)
37{
38 suppressErrors(true);
39 for (auto& reactor : reactors) {
41 }
42}
43
44ReactorNet::~ReactorNet()
45{
46}
47
49{
50 m_time = time;
52 m_integrator_init = false;
53}
54
55void ReactorNet::setMaxTimeStep(double maxstep)
56{
57 m_maxstep = maxstep;
59}
60
62{
64}
65
66void ReactorNet::setTolerances(double rtol, double atol)
67{
68 if (rtol >= 0.0) {
69 m_rtol = rtol;
70 }
71 if (atol >= 0.0) {
72 m_atols = atol;
73 }
74 m_init = false;
75}
76
77void ReactorNet::setSensitivityTolerances(double rtol, double atol)
78{
79 if (rtol >= 0.0) {
80 m_rtolsens = rtol;
81 }
82 if (atol >= 0.0) {
83 m_atolsens = atol;
84 }
85 m_init = false;
86}
87
90 return m_time;
91 } else {
92 throw CanteraError("ReactorNet::time", "Time is not the independent variable"
93 " for this reactor network.");
94 }
95}
96
99 return m_time;
100 } else {
101 throw CanteraError("ReactorNet::distance", "Distance is not the independent"
102 " variable for this reactor network.");
103 }
104}
105
107{
108 m_nv = 0;
109 debuglog("Initializing reactor network.\n", m_verbose);
110 if (m_reactors.empty()) {
111 throw CanteraError("ReactorNet::initialize",
112 "no reactors in network!");
113 }
114 // Names of Reactors and ReactorSurfaces using each Solution; should be only one
115 map<Solution*, vector<string>> solutions;
116 // Unique ReactorSurface objects. Can be attached to multiple Reactor objects
117 set<ReactorBase*> surfaces;
118 m_start.assign(1, 0);
119 for (size_t n = 0; n < m_reactors.size(); n++) {
120 Reactor& r = *m_reactors[n];
121 shared_ptr<Solution> bulk = r.contents4();
123 size_t nv = r.neq();
124 m_nv += nv;
125 m_start.push_back(m_nv);
126
127 if (m_verbose) {
128 writelog("Reactor {:d}: {:d} variables.\n", n, nv);
129 writelog(" {:d} sensitivity params.\n", r.nSensParams());
130 }
131 if (r.type() == "FlowReactor" && m_reactors.size() > 1) {
132 throw CanteraError("ReactorNet::initialize",
133 "FlowReactors must be used alone.");
134 }
135 solutions[bulk.get()].push_back(r.name());
136 for (size_t i = 0; i < r.nSurfs(); i++) {
137 if (r.surface(i)->contents4()->adjacent(bulk->name()) != bulk) {
138 throw CanteraError("ReactorNet::initialize",
139 "Bulk phase '{}' used by interface '{}' must be the same object\n"
140 "as the contents of the adjacent reactor '{}'.",
141 bulk->name(), r.surface(i)->name(), r.name());
142 }
143 surfaces.insert(r.surface(i));
144 }
145 }
146 for (auto surf : surfaces) {
147 solutions[surf->contents4().get()].push_back(surf->name());
148 }
149 for (auto& [soln, reactors] : solutions) {
150 if (reactors.size() > 1) {
151 string shared;
152 for (size_t i = 0; i < reactors.size() - 1; i++) {
153 shared += fmt::format("'{}', ", reactors[i]);
154 }
155 shared += fmt::format("'{}'", reactors.back());
156 warn_user("ReactorNet::initialize", "The following reactors / reactor"
157 " surfaces are using the same Solution object: {}. Use independent"
158 " Solution objects or set the 'clone' argument to 'true' when creating"
159 " the Reactor or ReactorSurface objects. Shared Solution objects within"
160 " a single ReactorNet will be an error after Cantera 3.2.", shared);
161 }
162 }
163
164 m_ydot.resize(m_nv,0.0);
165 m_yest.resize(m_nv,0.0);
166 m_advancelimits.resize(m_nv,-1.0);
167 m_atol.resize(neq());
168 fill(m_atol.begin(), m_atol.end(), m_atols);
169 m_integ->setTolerances(m_rtol, neq(), m_atol.data());
170 m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
171 if (!m_linearSolverType.empty()) {
172 m_integ->setLinearSolverType(m_linearSolverType);
173 }
174 if (m_precon) {
175 m_integ->setPreconditioner(m_precon);
176 }
177 m_integ->initialize(m_time, *this);
178 if (m_verbose) {
179 writelog("Number of equations: {:d}\n", neq());
180 writelog("Maximum time step: {:14.6g}\n", m_maxstep);
181 }
182 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
184 }
185 m_integrator_init = true;
186 m_init = true;
187}
188
190{
191 if (m_init) {
192 debuglog("Re-initializing reactor network.\n", m_verbose);
193 m_integ->reinitialize(m_time, *this);
194 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
196 }
197 m_integrator_init = true;
198 } else {
199 initialize();
200 }
201}
202
203void ReactorNet::setLinearSolverType(const string& linSolverType)
204{
205 m_linearSolverType = linSolverType;
206 m_integrator_init = false;
207}
208
209void ReactorNet::setPreconditioner(shared_ptr<SystemJacobian> preconditioner)
210{
211 m_precon = preconditioner;
212 m_integrator_init = false;
213}
214
216{
217 integrator().setMaxSteps(nmax);
218}
219
221{
222 return integrator().maxSteps();
223}
224
225void ReactorNet::advance(double time)
226{
227 if (!m_init) {
228 initialize();
229 } else if (!m_integrator_init) {
230 reinitialize();
231 }
232 m_integ->integrate(time);
233 m_time = time;
234 updateState(m_integ->solution());
235}
236
237double ReactorNet::advance(double time, bool applylimit)
238{
239 if (!m_init) {
240 initialize();
241 } else if (!m_integrator_init) {
242 reinitialize();
243 }
244
245 if (!applylimit) {
246 // take full step
247 advance(time);
248 return time;
249 }
250
251 if (!hasAdvanceLimits()) {
252 // take full step
253 advance(time);
254 return time;
255 }
256
257 getAdvanceLimits(m_advancelimits.data());
258
259 // ensure that gradient is available
260 while (lastOrder() < 1) {
261 step();
262 }
263
264 int k = lastOrder();
265 double t = time, delta;
266 double* y = m_integ->solution();
267
268 // reduce time step if limits are exceeded
269 while (true) {
270 bool exceeded = false;
271 getEstimate(t, k, &m_yest[0]);
272 for (size_t j = 0; j < m_nv; j++) {
273 delta = abs(m_yest[j] - y[j]);
274 if ( (m_advancelimits[j] > 0.) && ( delta > m_advancelimits[j]) ) {
275 exceeded = true;
276 if (m_verbose) {
277 writelog(" Limiting global state vector component {:d} (dt = {:9.4g}):"
278 "{:11.6g} > {:9.4g}\n",
279 j, t - m_time, delta, m_advancelimits[j]);
280 }
281 }
282 }
283 if (!exceeded) {
284 break;
285 }
286 t = .5 * (m_time + t);
287 }
288 advance(t);
289 return t;
290}
291
293{
294 if (!m_init) {
295 initialize();
296 } else if (!m_integrator_init) {
297 reinitialize();
298 }
299 m_time = m_integ->step(m_time + 1.0);
300 updateState(m_integ->solution());
301 return m_time;
302}
303
304void ReactorNet::solveSteady(int loglevel)
305{
306 if (!m_init) {
307 initialize();
308 } else if (!m_integrator_init) {
309 reinitialize();
310 }
311 vector<double> y(neq());
312 getState(y.data());
313 SteadyReactorSolver solver(this, y.data());
315 solver.solve(loglevel);
316 solver.getState(y.data());
317 updateState(y.data());
318}
319
320Eigen::SparseMatrix<double> ReactorNet::steadyJacobian(double rdt)
321{
322 if (!m_init) {
323 initialize();
324 } else if (!m_integrator_init) {
325 reinitialize();
326 }
327 vector<double> y0(neq());
328 vector<double> y1(neq());
329 getState(y0.data());
330 SteadyReactorSolver solver(this, y0.data());
331 solver.evalJacobian(y0.data());
332 if (rdt) {
333 solver.linearSolver()->updateTransient(rdt, solver.transientMask().data());
334 }
335 return std::dynamic_pointer_cast<EigenSparseJacobian>(solver.linearSolver())->jacobian();
336}
337
338void ReactorNet::getEstimate(double time, int k, double* yest)
339{
340 if (!m_init) {
341 initialize();
342 }
343 // initialize
344 double* cvode_dky = m_integ->solution();
345 for (size_t j = 0; j < m_nv; j++) {
346 yest[j] = cvode_dky[j];
347 }
348
349 // Taylor expansion
350 double factor = 1.;
351 double deltat = time - m_time;
352 for (int n = 1; n <= k; n++) {
353 factor *= deltat / n;
354 cvode_dky = m_integ->derivative(m_time, n);
355 for (size_t j = 0; j < m_nv; j++) {
356 yest[j] += factor * cvode_dky[j];
357 }
358 }
359}
360
362{
363 if (m_integ) {
364 return m_integ->lastOrder();
365 } else {
366 return 0;
367 }
368}
369
371{
372 warn_deprecated("ReactorNet::addReactor",
373 "To be removed after Cantera 3.2. Replaceable by reactor net "
374 "instantiation with contents.");
375 for (auto current : m_reactors) {
376 if (current->isOde() != r.isOde()) {
377 throw CanteraError("ReactorNet::addReactor",
378 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
379 current->type(), r.type());
380 }
381 if (current->timeIsIndependent() != r.timeIsIndependent()) {
382 throw CanteraError("ReactorNet::addReactor",
383 "Cannot mix Reactor types using time and space as independent variables"
384 "\n({} and {})", current->type(), r.type());
385 }
386 }
388 r.setNetwork(this);
389 m_reactors.push_back(&r);
390 if (!m_integ) {
391 m_integ.reset(newIntegrator(r.isOde() ? "CVODE" : "IDA"));
392 // use backward differencing, with a full Jacobian computed
393 // numerically, and use a Newton linear iterator
394 m_integ->setMethod(BDF_Method);
395 m_integ->setLinearSolverType("DENSE");
396 }
397 updateNames(r);
398}
399
400void ReactorNet::addReactor(shared_ptr<ReactorBase> reactor)
401{
402 auto r = std::dynamic_pointer_cast<Reactor>(reactor);
403 if (!r) {
404 throw CanteraError("ReactorNet::addReactor",
405 "Reactor with type '{}' cannot be added to network.",
406 reactor->type());
407 }
408
409 for (auto current : m_reactors) {
410 if (current->isOde() != r->isOde()) {
411 throw CanteraError("ReactorNet::addReactor",
412 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
413 current->type(), r->type());
414 }
415 if (current->timeIsIndependent() != r->timeIsIndependent()) {
416 throw CanteraError("ReactorNet::addReactor",
417 "Cannot mix Reactor types using time and space as independent variables"
418 "\n({} and {})", current->type(), r->type());
419 }
420 }
421 m_timeIsIndependent = r->timeIsIndependent();
422 r->setNetwork(this);
423 m_reactors.push_back(r.get());
424 if (!m_integ) {
425 m_integ.reset(newIntegrator(r->isOde() ? "CVODE" : "IDA"));
426 // use backward differencing, with a full Jacobian computed
427 // numerically, and use a Newton linear iterator
428 m_integ->setMethod(BDF_Method);
429 m_integ->setLinearSolverType("DENSE");
430 }
431 updateNames(*r);
432}
433
435{
436 // ensure that reactors and components have reproducible names
438
439 for (size_t i=0; i<r.nWalls(); i++) {
440 auto& w = r.wall(i);
442 if (w.left().type() == "Reservoir") {
443 w.left().setDefaultName(m_counts);
444 }
445 if (w.right().type() == "Reservoir") {
446 w.right().setDefaultName(m_counts);
447 }
448 }
449
450 for (size_t i=0; i<r.nInlets(); i++) {
451 auto& in = r.inlet(i);
453 if (in.in().type() == "Reservoir") {
454 in.in().setDefaultName(m_counts);
455 }
456 }
457
458 for (size_t i=0; i<r.nOutlets(); i++) {
459 auto& out = r.outlet(i);
461 if (out.out().type() == "Reservoir") {
462 out.out().setDefaultName(m_counts);
463 }
464 }
465
466 for (size_t i=0; i<r.nSurfs(); i++) {
468 }
469}
470
472 if (m_integ == nullptr) {
473 throw CanteraError("ReactorNet::integrator",
474 "Integrator has not been instantiated. Add one or more reactors first.");
475 }
476 return *m_integ;
477}
478
479void ReactorNet::eval(double t, double* y, double* ydot, double* p)
480{
481 m_time = t;
482 updateState(y);
483 m_LHS.assign(m_nv, 1);
484 m_RHS.assign(m_nv, 0);
485 for (size_t n = 0; n < m_reactors.size(); n++) {
486 m_reactors[n]->applySensitivity(p);
487 m_reactors[n]->eval(t, m_LHS.data() + m_start[n], m_RHS.data() + m_start[n]);
488 size_t yEnd = 0;
489 if (n == m_reactors.size() - 1) {
490 yEnd = m_RHS.size();
491 } else {
492 yEnd = m_start[n + 1];
493 }
494 for (size_t i = m_start[n]; i < yEnd; i++) {
495 ydot[i] = m_RHS[i] / m_LHS[i];
496 }
497 m_reactors[n]->resetSensitivity(p);
498 }
499 checkFinite("ydot", ydot, m_nv);
500}
501
502void ReactorNet::evalDae(double t, double* y, double* ydot, double* p, double* residual)
503{
504 m_time = t;
505 updateState(y);
506 for (size_t n = 0; n < m_reactors.size(); n++) {
507 m_reactors[n]->applySensitivity(p);
508 m_reactors[n]->evalDae(t, y, ydot, residual);
509 m_reactors[n]->resetSensitivity(p);
510 }
511 checkFinite("ydot", ydot, m_nv);
512}
513
514void ReactorNet::getConstraints(double* constraints)
515{
516 for (size_t n = 0; n < m_reactors.size(); n++) {
517 m_reactors[n]->getConstraints(constraints + m_start[n]);
518 }
519}
520
521double ReactorNet::sensitivity(size_t k, size_t p)
522{
523 if (!m_init) {
524 initialize();
525 }
526 if (p >= m_sens_params.size()) {
527 throw IndexError("ReactorNet::sensitivity",
528 "m_sens_params", p, m_sens_params.size());
529 }
530 double denom = m_integ->solution(k);
531 if (denom == 0.0) {
532 denom = SmallNumber;
533 }
534 return m_integ->sensitivity(k, p) / denom;
535}
536
537void ReactorNet::evalJacobian(double t, double* y, double* ydot, double* p, Array2D* j)
538{
539 //evaluate the unperturbed ydot
540 eval(t, y, ydot, p);
541 for (size_t n = 0; n < m_nv; n++) {
542 // perturb x(n)
543 double ysave = y[n];
544 double dy = m_atol[n] + fabs(ysave)*m_rtol;
545 y[n] = ysave + dy;
546 dy = y[n] - ysave;
547
548 // calculate perturbed residual
549 eval(t, y, m_ydot.data(), p);
550
551 // compute nth column of Jacobian
552 for (size_t m = 0; m < m_nv; m++) {
553 j->value(m,n) = (m_ydot[m] - ydot[m])/dy;
554 }
555 y[n] = ysave;
556 }
557}
558
560{
561 checkFinite("y", y, m_nv);
562 for (size_t n = 0; n < m_reactors.size(); n++) {
563 m_reactors[n]->updateState(y + m_start[n]);
564 }
565}
566
567void ReactorNet::getDerivative(int k, double* dky)
568{
569 if (!m_init) {
570 initialize();
571 }
572 double* cvode_dky = m_integ->derivative(m_time, k);
573 for (size_t j = 0; j < m_nv; j++) {
574 dky[j] = cvode_dky[j];
575 }
576}
577
578void ReactorNet::setAdvanceLimits(const double *limits)
579{
580 if (!m_init) {
581 initialize();
582 }
583 for (size_t n = 0; n < m_reactors.size(); n++) {
584 m_reactors[n]->setAdvanceLimits(limits + m_start[n]);
585 }
586}
587
589{
590 bool has_limit = false;
591 for (size_t n = 0; n < m_reactors.size(); n++) {
592 has_limit |= m_reactors[n]->hasAdvanceLimits();
593 }
594 return has_limit;
595}
596
597bool ReactorNet::getAdvanceLimits(double *limits) const
598{
599 bool has_limit = false;
600 for (size_t n = 0; n < m_reactors.size(); n++) {
601 has_limit |= m_reactors[n]->getAdvanceLimits(limits + m_start[n]);
602 }
603 return has_limit;
604}
605
606void ReactorNet::getState(double* y)
607{
608 for (size_t n = 0; n < m_reactors.size(); n++) {
609 m_reactors[n]->getState(y + m_start[n]);
610 }
611}
612
613void ReactorNet::getStateDae(double* y, double* ydot)
614{
615 for (size_t n = 0; n < m_reactors.size(); n++) {
616 m_reactors[n]->getStateDae(y + m_start[n], ydot + m_start[n]);
617 }
618}
619
620size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
621{
622 if (!m_init) {
623 initialize();
624 }
625 return m_start[reactor] + m_reactors[reactor]->componentIndex(component);
626}
627
628string ReactorNet::componentName(size_t i) const
629{
630 for (auto r : m_reactors) {
631 if (i < r->neq()) {
632 return r->name() + ": " + r->componentName(i);
633 } else {
634 i -= r->neq();
635 }
636 }
637 throw CanteraError("ReactorNet::componentName", "Index out of bounds");
638}
639
640double ReactorNet::upperBound(size_t i) const
641{
642 size_t i0 = i;
643 for (auto r : m_reactors) {
644 if (i < r->neq()) {
645 return r->upperBound(i);
646 } else {
647 i -= r->neq();
648 }
649 }
650 throw CanteraError("ReactorNet::upperBound", "Index {} out of bounds", i0);
651}
652
653double ReactorNet::lowerBound(size_t i) const
654{
655 size_t i0 = i;
656 for (auto r : m_reactors) {
657 if (i < r->neq()) {
658 return r->lowerBound(i);
659 } else {
660 i -= r->neq();
661 }
662 }
663 throw CanteraError("ReactorNet::lowerBound", "Index {} out of bounds", i0);
664}
665
667 size_t i = 0;
668 for (auto r : m_reactors) {
669 r->resetBadValues(y + m_start[i++]);
670 }
671}
672
674 const string& name, double value, double scale)
675{
676 if (m_integrator_init) {
677 throw CanteraError("ReactorNet::registerSensitivityParameter",
678 "Sensitivity parameters cannot be added after the"
679 "integrator has been initialized.");
680 }
681 m_paramNames.push_back(name);
682 m_sens_params.push_back(value);
683 m_paramScales.push_back(scale);
684 return m_sens_params.size() - 1;
685}
686
688{
689 // Apply given settings to all reactors
690 for (size_t i = 0; i < m_reactors.size(); i++) {
691 m_reactors[i]->setDerivativeSettings(settings);
692 }
693}
694
696{
697 if (m_integ) {
698 return m_integ->solverStats();
699 } else {
700 return AnyMap();
701 }
702}
703
705{
706 if (m_integ) {
707 return m_integ->linearSolverType();
708 } else {
709 return "";
710 }
711}
712
713void ReactorNet::preconditionerSolve(double* rhs, double* output)
714{
715 if (!m_integ) {
716 throw CanteraError("ReactorNet::preconditionerSolve",
717 "Must only be called after ReactorNet is initialized.");
718 }
719 m_integ->preconditionerSolve(m_nv, rhs, output);
720}
721
722void ReactorNet::preconditionerSetup(double t, double* y, double gamma)
723{
724 // ensure state is up to date.
725 updateState(y);
726 // get the preconditioner
727 auto precon = m_integ->preconditioner();
728 // Reset preconditioner
729 precon->reset();
730 // Set gamma value for M =I - gamma*J
731 precon->setGamma(gamma);
732 // Make a copy of state to adjust it for preconditioner
733 vector<double> yCopy(m_nv);
734 // Get state of reactor
735 getState(yCopy.data());
736 // transform state based on preconditioner rules
737 precon->stateAdjustment(yCopy);
738 // update network with adjusted state
739 updateState(yCopy.data());
740 // Get jacobians and give elements to preconditioners
741 for (size_t i = 0; i < m_reactors.size(); i++) {
742 Eigen::SparseMatrix<double> rJac = m_reactors[i]->jacobian();
743 for (int k=0; k<rJac.outerSize(); ++k) {
744 for (Eigen::SparseMatrix<double>::InnerIterator it(rJac, k); it; ++it) {
745 precon->setValue(it.row() + m_start[i], it.col() + m_start[i],
746 it.value());
747 }
748 }
749 }
750 // post reactor setup operations
751 precon->updatePreconditioner();
752}
753
755{
756 if (!m_integ) {
757 throw CanteraError("ReactorNet::updatePreconditioner",
758 "Must only be called after ReactorNet is initialized.");
759 }
760 auto precon = m_integ->preconditioner();
761 precon->setGamma(gamma);
762 precon->updatePreconditioner();
763}
764
766 // check for non-mole-based reactors and throw an error otherwise
767 for (auto reactor : m_reactors) {
769 throw CanteraError("ReactorNet::checkPreconditionerSupported",
770 "Preconditioning is only supported for type *MoleReactor,\n"
771 "Reactor type given: '{}'.",
772 reactor->type());
773 }
774 }
775}
776
777SteadyReactorSolver::SteadyReactorSolver(ReactorNet* net, double* x0)
778 : m_net(net)
779{
780 m_size = m_net->neq();
781 m_jac = newSystemJacobian("eigen-sparse-direct");
783 m_initialState.assign(x0, x0 + m_size);
784 setInitialGuess(x0);
785 m_mask.assign(m_size, 1);
786 size_t start = 0;
787 for (size_t i = 0; i < net->nReactors(); i++) {
788 auto& R = net->reactor(i);
789 for (auto& m : R.steadyConstraints()) {
790 m_algebraic.push_back(start + m);
791 }
792 start += R.neq();
793 }
794 for (auto& n : m_algebraic) {
795 m_mask[n] = 0;
796 }
797}
798
799void SteadyReactorSolver::eval(double* x, double* r, double rdt, int count)
800{
801 if (rdt < 0.0) {
802 rdt = m_rdt;
803 }
804 vector<double> xv(x, x + size());
805 m_net->eval(0.0, x, r, nullptr);
806 for (size_t i = 0; i < size(); i++) {
807 r[i] -= (x[i] - m_initialState[i]) * rdt;
808 }
809 // Hold algebraic constraints fixed
810 for (auto& n : m_algebraic) {
811 r[n] = x[n] - m_initialState[n];
812 }
813}
814
815void SteadyReactorSolver::initTimeInteg(double dt, double* x)
816{
818 m_initialState.assign(x, x + size());
819}
820
822{
823 m_jac->reset();
824 clock_t t0 = clock();
825 m_work1.resize(size());
826 m_work2.resize(size());
827 eval(x0, m_work1.data(), 0.0, 0);
828 for (size_t j = 0; j < size(); j++) {
829 // perturb x(n); preserve sign(x(n))
830 double xsave = x0[j];
831 double dx = fabs(xsave) * m_jacobianRelPerturb + m_jacobianAbsPerturb;
832 if (xsave < 0) {
833 dx = -dx;
834 }
835 x0[j] = xsave + dx;
836 double rdx = 1.0 / (x0[j] - xsave);
837
838 // calculate perturbed residual
839 eval(x0, m_work2.data(), 0.0, 0);
840
841 // compute nth column of Jacobian
842 for (size_t i = 0; i < size(); i++) {
843 double delta = m_work2[i] - m_work1[i];
844 if (std::abs(delta) > m_jacobianThreshold || i == j) {
845 m_jac->setValue(i, j, delta * rdx);
846 }
847 }
848 x0[j] = xsave;
849 }
850
851 m_jac->updateElapsed(double(clock() - t0) / CLOCKS_PER_SEC);
852 m_jac->incrementEvals();
853 m_jac->setAge(0);
854}
855
856double SteadyReactorSolver::weightedNorm(const double* step) const
857{
858 double sum = 0.0;
859 const double* x = m_state->data();
860 for (size_t i = 0; i < size(); i++) {
861 double ewt = m_net->rtol()*x[i] + m_net->atol();
862 double f = step[i] / ewt;
863 sum += f*f;
864 }
865 return sqrt(sum / size());
866}
867
869{
870 return m_net->componentName(i);
871}
872
874{
875 return m_net->upperBound(i);
876}
877
879{
880 return m_net->lowerBound(i);
881}
882
884{
885 m_net->resetBadValues(x);
886}
887
888void SteadyReactorSolver::writeDebugInfo(const string& header_suffix,
889 const string& message, int loglevel, int attempt_counter)
890{
891 if (loglevel >= 6 && !m_state->empty()) {
892 const auto& state = *m_state;
893 writelog("Current state ({}):\n[", header_suffix);
894 for (size_t i = 0; i < state.size() - 1; i++) {
895 writelog("{}, ", state[i]);
896 }
897 writelog("{}]\n", state.back());
898 }
899 if (loglevel >= 7 && !m_xnew.empty()) {
900 writelog("Current residual ({}):\n[", header_suffix);
901 for (size_t i = 0; i < m_xnew.size() - 1; i++) {
902 writelog("{}, ", m_xnew[i]);
903 }
904 writelog("{}]\n", m_xnew.back());
905 }
906}
907
908shared_ptr<ReactorNet> newReactorNet(vector<shared_ptr<ReactorBase>>& reactors)
909{
910 return make_shared<ReactorNet>(reactors);
911}
912
913}
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.
void setDefaultName(map< string, int > &counts)
Set the default name of a connector. 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.
virtual size_t nSurfs() const
Return the number of surfaces in a reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
shared_ptr< Solution > contents4()
Access the Solution object used to represent the contents of this 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.
string name() const
Return the name of this reactor.
Definition ReactorBase.h:82
A class representing a network of connected reactors.
Definition ReactorNet.h:30
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:258
vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
Definition ReactorNet.h:405
vector< double > m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition ReactorNet.h:431
double m_initial_time
The initial value of the independent variable in the system.
Definition ReactorNet.h:398
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:395
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:194
map< string, int > m_counts
Map used for default name generation.
Definition ReactorNet.h:390
double upperBound(size_t i) const
Get the upper bound on the i-th component of the global state vector.
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:416
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:96
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 ...
void solveSteady(int loglevel=0)
Solve directly for the steady-state solution.
bool m_timeIsIndependent
Indicates whether time or space is the independent variable.
Definition ReactorNet.h:421
double atol()
Absolute integration tolerance.
Definition ReactorNet.h:101
bool hasAdvanceLimits() const
Check whether ReactorNet object uses advance limits.
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
double lowerBound(size_t i) const
Get the lower bound on the i-th component of the global state vector.
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:401
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.
Eigen::SparseMatrix< double > steadyJacobian(double rdt=0.0)
Get the Jacobian used by the steady-state solver.
string linearSolverType() const
Problem type of integrator.
void updatePreconditioner(double gamma) override
Update the preconditioner based on already computed jacobian values.
void setPreconditioner(shared_ptr< SystemJacobian > preconditioner)
Set preconditioner used by the linear solver.
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:424
void resetBadValues(double *y)
Reset physically or mathematically problematic values, such as negative species concentrations.
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.
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:100
size_t nSensParams() const override
Number of sensitivity parameters associated with this reactor.
Definition Reactor.cpp:136
string type() const override
String indicating the reactor model implemented.
Definition Reactor.h:55
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:108
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
Definition Reactor.h:250
virtual bool isOde() const
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs.
Definition Reactor.h:62
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
Definition Reactor.h:68
Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
Definition ReactorNet.h:441
double weightedNorm(const double *step) const override
Compute the weighted norm of step.
vector< double > m_initialState
Initial value of each state variable.
Definition ReactorNet.h:459
string componentName(size_t i) const override
Get the name of the i-th component of the state vector.
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
vector< size_t > m_algebraic
Indices of variables that are held constant in the time-stepping mode of the steady-state solver.
Definition ReactorNet.h:463
void evalJacobian(double *x0) override
Evaluates the Jacobian at x0 using finite differences.
void resetBadValues(double *x) override
Reset values such as negative species concentrations.
void writeDebugInfo(const string &header_suffix, const string &message, int loglevel, int attempt_counter) override
Write solver debugging based on the specified log level.
void eval(double *x, double *r, double rdt=-1.0, int count=1) override
Evaluate the residual function.
void initTimeInteg(double dt, double *x) override
Prepare for time stepping beginning with solution x and timestep dt.
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
int solve(double *x0, double *x1, int loglevel)
Solve , where is the residual function.
virtual void resize()
Call to set the size of internal data structures after first defining the system or if the problem si...
shared_ptr< SystemJacobian > linearSolver() const
Get the the linear solver being used to hold the Jacobian matrix and solve linear systems as part of ...
vector< double > m_xnew
Work array used to hold the residual or the new solution.
double m_jacobianAbsPerturb
Absolute perturbation of each component in finite difference Jacobian.
size_t size() const
Total solution vector length;.
vector< int > & transientMask()
Access the vector indicating which equations contain a transient term.
double rdt() const
Reciprocal of the time step.
virtual void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
double m_rdt
Reciprocal of time step.
double m_jacobianThreshold
Threshold for ignoring small elements in Jacobian.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
shared_ptr< vector< double > > m_state
Solution vector.
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
void getState(double *x) const
Get the converged steady-state solution after calling solve().
double m_jacobianRelPerturb
Relative perturbation of each component in finite difference Jacobian.
vector< double > m_work1
Work arrays used during Jacobian evaluation.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
Definition global.h:159
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:176
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.
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:268
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
shared_ptr< SystemJacobian > newSystemJacobian(const string &type)
Create a SystemJacobian object of the specified type.
const double SmallNumber
smallest number to compare to zero.
Definition ct_defs.h:158
shared_ptr< ReactorNet > newReactorNet(vector< shared_ptr< ReactorBase > > &reactors)
Create a reactor network containing one or more coupled reactors.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Definition AnyMap.cpp:1997
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...