Cantera  4.0.0a1
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
12#include "cantera/base/Array.h"
20
21#include <cstdio>
22#include <deque>
23
24namespace Cantera
25{
26
27ReactorNet::~ReactorNet()
28{
29}
30
31ReactorNet::ReactorNet(shared_ptr<ReactorBase> reactor)
32 : ReactorNet(span<shared_ptr<ReactorBase>>{&reactor, 1})
33{
34}
35
36ReactorNet::ReactorNet(span<shared_ptr<ReactorBase>> reactors)
37{
38 if (reactors.empty()) {
39 throw CanteraError("ReactorNet::ReactorNet", "No reactors in network!");
40 }
41
42 suppressErrors(true);
43 bool isOde = true;
44
45 // Names of Reactors and ReactorSurfaces using each Solution; should be only one
46 // reactor per Solution.
47 map<Solution*, set<string>> solutions;
48
49 // All ReactorBase objects connected to this network, starting with the ones
50 // explicitly included.
51 std::deque<shared_ptr<ReactorBase>> reactorQueue(reactors.begin(), reactors.end());
52 std::set<shared_ptr<ReactorBase>> visited;
53
54 while (!reactorQueue.empty()) {
55 auto R = reactorQueue.front();
56 reactorQueue.pop_front();
57
58 if (visited.find(R) != visited.end()) {
59 continue;
60 }
61 visited.insert(R);
62
63 if (auto bulk = std::dynamic_pointer_cast<Reactor>(R)) {
64 m_bulkReactors.push_back(bulk);
65 m_reactors.push_back(R);
66 } else if (auto surf = std::dynamic_pointer_cast<ReactorSurface>(R)) {
67 m_surfaces.push_back(surf);
68 m_reactors.push_back(R);
69 for (size_t i = 0; i < surf->nAdjacent(); i++) {
70 reactorQueue.push_back(surf->adjacent(i));
71 }
72 } else if (auto res = std::dynamic_pointer_cast<Reservoir>(R)) {
73 m_reservoirs.push_back(res);
74 }
75
76 // Discover reservoirs, flow devices, and walls connected to this reactor
77 for (size_t i = 0; i < R->nInlets(); i++) {
78 auto& inlet = R->inlet(i);
79 m_flowDevices.insert(&inlet);
80 reactorQueue.push_back(inlet.in().shared_from_this());
81 }
82 for (size_t i = 0; i < R->nOutlets(); i++) {
83 auto& outlet = R->outlet(i);
84 m_flowDevices.insert(&outlet);
85 reactorQueue.push_back(outlet.out().shared_from_this());
86 }
87 for (size_t i = 0; i < R->nWalls(); i++) {
88 auto& wall = R->wall(i);
89 m_walls.insert(&wall);
90 reactorQueue.push_back(wall.left().shared_from_this());
91 reactorQueue.push_back(wall.right().shared_from_this());
92 }
93
94 auto phase = R->phase();
95 for (size_t i = 0; i < R->nSurfs(); i++) {
96 if (R->surface(i)->phase()->adjacent(phase->name()) != phase) {
97 throw CanteraError("ReactorNet::initialize",
98 "Bulk phase '{}' used by interface '{}' must be the same object\n"
99 "as the contents of the adjacent reactor '{}'.",
100 phase->name(), R->surface(i)->name(), R->name());
101 }
102 reactorQueue.push_back(R->surface(i)->shared_from_this());
103 }
104 R->setNetwork(this);
105 updateNames(*R);
106 solutions[phase.get()].insert(R->name());
107 }
108
109 for (auto& R : m_bulkReactors) {
110 if (R->type() == "FlowReactor" && m_bulkReactors.size() > 1) {
111 throw CanteraError("ReactorNet::initialize",
112 "FlowReactors must be used alone.");
113 }
114 m_timeIsIndependent = R->timeIsIndependent();
115 R->initialize(m_time);
116 isOde = R->isOde();
117 R->setOffset(m_nv);
118 size_t nv = R->neq();
119 m_nv += nv;
120
121 for (auto current : m_bulkReactors) {
122 if (current->isOde() != R->isOde()) {
123 throw CanteraError("ReactorNet::ReactorNet",
124 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
125 current->type(), R->type());
126 }
127 if (current->timeIsIndependent() != R->timeIsIndependent()) {
128 throw CanteraError("ReactorNet::ReactorNet",
129 "Cannot mix Reactor types using time and space as independent variables"
130 "\n({} and {})", current->type(), R->type());
131 }
132 }
133 }
134
135 for (auto surf : m_surfaces) {
136 surf->setOffset(m_nv);
137 surf->initialize(m_time);
138 m_nv += surf->neq();
139 solutions[surf->phase().get()].insert(surf->name());
140 }
141
142 for (auto& [soln, reactors] : solutions) {
143 if (reactors.size() > 1) {
144 string shared;
145 for (auto r : reactors) {
146 if (r != *reactors.begin()) {
147 shared += ", ";
148 }
149 shared += fmt::format("'{}'", r);
150 }
151 throw CanteraError("ReactorNet::initialize", "The following reactors /"
152 " reactor surfaces are using the same Solution object: {}. Use"
153 " independent Solution objects or set the 'clone' argument to 'true'"
154 " when creating the Reactor or ReactorSurface objects.", shared);
155 }
156 }
157
158 m_ydot.resize(m_nv, 0.0);
159 m_yest.resize(m_nv, 0.0);
160 m_advancelimits.resize(m_nv, -1.0);
161 m_atol.resize(m_nv);
162
163 m_integ.reset(newIntegrator(isOde ? "CVODE" : "IDA"));
164 // use backward differencing, with a full Jacobian computed
165 // numerically, and use a Newton linear iterator
166 m_integ->setMethod(BDF_Method);
167 m_integ->setLinearSolverType("DENSE");
168 setTolerances(-1.0, -1.0);
169}
170
172{
173 m_time = time;
176}
177
178void ReactorNet::setMaxTimeStep(double maxstep)
179{
180 m_maxstep = maxstep;
182}
183
185{
187}
188
189void ReactorNet::setTolerances(double rtol, double atol)
190{
191 if (rtol >= 0.0) {
192 m_rtol = rtol;
193 }
194 if (atol >= 0.0) {
195 m_atols = atol;
196 }
197 fill(m_atol.begin(), m_atol.end(), m_atols);
198 m_integ->setTolerances(m_rtol, m_atol);
199}
200
201void ReactorNet::setSensitivityTolerances(double rtol, double atol)
202{
203 if (rtol >= 0.0) {
204 m_rtolsens = rtol;
205 }
206 if (atol >= 0.0) {
207 m_atolsens = atol;
208 }
209 m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
210}
211
214 return m_time;
215 } else {
216 throw CanteraError("ReactorNet::time", "Time is not the independent variable"
217 " for this reactor network.");
218 }
219}
220
222 if (!m_timeIsIndependent) {
223 return m_time;
224 } else {
225 throw CanteraError("ReactorNet::distance", "Distance is not the independent"
226 " variable for this reactor network.");
227 }
228}
229
231{
234 reinitialize();
235 }
236 return;
237 }
238 if (!m_linearSolverType.empty()) {
239 m_integ->setLinearSolverType(m_linearSolverType);
240 }
241 if (m_precon) {
242 m_integ->setPreconditioner(m_precon);
243 }
244 for (auto& reactor : m_reactors) {
246 }
247 m_integ->initialize(m_time, *this);
248 if (m_verbose) {
249 writelog("Number of equations: {:d}\n", neq());
250 writelog("Maximum time step: {:14.6g}\n", m_maxstep);
251 }
252 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
254 }
255 m_needIntegratorInit = false;
257}
258
260{
262 debuglog("Re-initializing reactor network.\n", m_verbose);
263 for (auto& reactor : m_reactors) {
265 }
266 m_integ->reinitialize(m_time, *this);
267 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
269 }
270 m_needIntegratorInit = false;
271 } else {
272 initialize();
273 }
274}
275
276void ReactorNet::setLinearSolverType(const string& linSolverType)
277{
278 m_linearSolverType = linSolverType;
280}
281
282void ReactorNet::setPreconditioner(shared_ptr<SystemJacobian> preconditioner)
283{
284 m_precon = preconditioner;
286}
287
289{
290 integrator().setMaxSteps(nmax);
291}
292
294{
295 return integrator().maxSteps();
296}
297
298void ReactorNet::advance(double time)
299{
300 initialize();
301 m_integ->integrate(time);
302 m_time = m_integ->currentTime();
303 updateState(m_integ->solution());
304}
305
306double ReactorNet::advance(double time, bool applylimit)
307{
308 initialize();
309 if (!applylimit || !hasAdvanceLimits()) {
310 // No limit enforcement requested; integrate to requested time
311 advance(time);
312 return time;
313 }
314
315 // Enable root-based limit detection and set the base state to the current state
316 m_ybase.assign(m_nv, 0.0);
320 m_integ->setRootFunctionCount(nRootFunctions());
321
322 // Integrate toward the requested time; integrator will return early if a limit is
323 // reached (CV_ROOT_RETURN). The try/catch ensures the temporary root-finding state
324 // is cleared even when CVODE throws so subsequent calls start clean.
325 try {
326 m_integ->integrate(time);
327 } catch (...) {
328 m_limit_check_active = false;
329 m_integ->setRootFunctionCount(nRootFunctions());
330 throw;
331 }
332 m_time = m_integ->currentTime();
333
334 // Update reactor states to match the integrator solution at the time reached
335 // (which may be earlier than 'time' if a limit was triggered)
336 updateState(m_integ->solution());
337
338 // Disable limit checking after this call
339 m_limit_check_active = false;
340 m_integ->setRootFunctionCount(nRootFunctions());
341
342 // When a root event stopped integration before reaching the requested time, report
343 // the most limiting component and details about the step.
344 if (m_verbose && m_time < time) {
345 // Ensure limits are available
346 if (m_advancelimits.size() != m_nv) {
347 m_advancelimits.assign(m_nv, -1.0);
348 }
349 getAdvanceLimits(m_advancelimits);
350 auto ycurr = m_integ->solution();
351 size_t jmax = npos;
352 double max_ratio = -1.0;
353 double best_limit = 0.0;
354 for (size_t j = 0; j < m_nv; j++) {
355 double lim = m_advancelimits[j];
356 if (lim > 0.0) {
357 double delta = std::abs(ycurr[j] - m_ybase[j]);
358 double ratio = delta / lim;
359 if (ratio > max_ratio) {
360 max_ratio = ratio;
361 jmax = j;
362 best_limit = lim;
363 }
364 }
365 }
366 if (jmax != npos) {
367 double dt = m_time - m_ybase_time;
368 double y_start = m_ybase[jmax];
369 double y_end = ycurr[jmax];
370 double delta = y_end - y_start;
371 writelog(" Advance limit triggered for component {:d} (dt = {:9.4g}):"
372 " y_start = {:11.6g}, y_end = {:11.6g},"
373 " delta = {:11.6g}, limit = {:9.4g}\n",
374 jmax, dt, y_start, y_end, delta, best_limit);
375 }
376 }
377
378 // m_time is tracked via callbacks during integration
379 return m_time;
380}
381
383{
384 initialize();
385 m_time = m_integ->step(m_time + 1.0);
386 updateState(m_integ->solution());
387 return m_time;
388}
389
390void ReactorNet::solveSteady(int loglevel)
391{
392 initialize();
393 vector<double> y(neq());
394 getState(y);
395 SteadyReactorSolver solver(this, y);
397 solver.solve(loglevel);
398 solver.getState(y);
399 updateState(y);
400}
401
402Eigen::SparseMatrix<double> ReactorNet::steadyJacobian(double rdt)
403{
404 initialize();
405 vector<double> y0(neq());
406 vector<double> y1(neq());
407 getState(y0);
408 SteadyReactorSolver solver(this, y0);
409 solver.evalJacobian(y0);
410 if (rdt) {
411 solver.linearSolver()->updateTransient(rdt, solver.transientMask());
412 }
413 return std::dynamic_pointer_cast<EigenSparseJacobian>(solver.linearSolver())->jacobian();
414}
415
416void ReactorNet::getEstimate(double time, int k, span<double> yest)
417{
418 initialize();
419 auto cvode_dky = m_integ->solution();
420 std::copy(cvode_dky.begin(), cvode_dky.end(), yest.begin());
421
422 // Taylor expansion
423 double factor = 1.;
424 double deltat = time - m_time;
425 for (int n = 1; n <= k; n++) {
426 factor *= deltat / n;
427 cvode_dky = m_integ->derivative(m_time, n);
428 for (size_t j = 0; j < m_nv; j++) {
429 yest[j] += factor * cvode_dky[j];
430 }
431 }
432}
433
435{
436 if (m_integ) {
437 return m_integ->lastOrder();
438 } else {
439 return 0;
440 }
441}
442
444{
445 return (m_limit_check_active && hasAdvanceLimits()) ? 1 : 0;
446}
447
448void ReactorNet::evalRootFunctions(double t, span<const double> y, span<double> gout)
449{
450 // Default: no root detected
451 double g = 1.0;
452
454 // Ensure limits vector is current
455 if (m_advancelimits.size() != m_nv) {
456 m_advancelimits.assign(m_nv, -1.0);
457 }
458 getAdvanceLimits(m_advancelimits);
459
460 double max_ratio = 0.0;
461 for (size_t i = 0; i < m_nv; i++) {
462 double lim = m_advancelimits[i];
463 if (lim > 0.0) {
464 double delta = std::abs(y[i] - m_ybase[i]);
465 double ratio = delta / lim;
466 if (ratio > max_ratio) {
467 max_ratio = ratio;
468 }
469 }
470 }
471 g = 1.0 - max_ratio; // root at g = 0 when any component reaches its limit
472 }
473
474 gout[0] = g;
475}
476
478{
479 // ensure that reactors and components have reproducible names
481
482 for (size_t i=0; i<r.nWalls(); i++) {
483 auto& w = r.wall(i);
485 if (w.left().type() == "Reservoir") {
486 w.left().setDefaultName(m_counts);
487 }
488 if (w.right().type() == "Reservoir") {
489 w.right().setDefaultName(m_counts);
490 }
491 }
492
493 for (size_t i=0; i<r.nInlets(); i++) {
494 auto& in = r.inlet(i);
496 if (in.in().type() == "Reservoir") {
497 in.in().setDefaultName(m_counts);
498 }
499 }
500
501 for (size_t i=0; i<r.nOutlets(); i++) {
502 auto& out = r.outlet(i);
504 if (out.out().type() == "Reservoir") {
505 out.out().setDefaultName(m_counts);
506 }
507 }
508
509 for (size_t i=0; i<r.nSurfs(); i++) {
511 }
512}
513
515 if (m_integ == nullptr) {
516 throw CanteraError("ReactorNet::integrator",
517 "Integrator has not been instantiated. Add one or more reactors first.");
518 }
519 return *m_integ;
520}
521
522void ReactorNet::eval(double t, span<const double> y, span<double> ydot,
523 span<const double> p)
524{
525 m_time = t;
526 updateState(y);
527 m_LHS.assign(m_nv, 1);
528 m_RHS.assign(m_nv, 0);
529 for (auto& R : m_reactors) {
530 size_t offset = R->offset();
531 R->applySensitivity(p);
532 R->eval(t, span<double>(m_LHS.data() + offset, R->neq()),
533 span<double>(m_RHS.data() + offset, R->neq()));
534 for (size_t i = offset; i < offset + R->neq(); i++) {
535 ydot[i] = m_RHS[i] / m_LHS[i];
536 }
537 R->resetSensitivity(p);
538 }
539 checkFinite("ydot", ydot);
540}
541
542void ReactorNet::evalSteady(span<const double> y, span<double> residual)
543{
544 updateState(y);
545 m_LHS.assign(m_nv, 1);
546 m_RHS.assign(m_nv, 0);
547 for (auto& R : m_reactors) {
548 size_t offset = R->offset();
549 R->evalSteady(m_time, span<double>(m_LHS.data() + offset, R->neq()),
550 span<double>(m_RHS.data() + offset, R->neq()));
551 for (size_t i = offset; i < offset + R->neq(); i++) {
552 residual[i] = m_RHS[i] / m_LHS[i];
553 }
554 }
555 checkFinite("residual", residual);
556}
557
558void ReactorNet::evalDae(double t, span<const double> y, span<const double> ydot,
559 span<const double> p, span<double> residual)
560{
561 m_time = t;
562 updateState(y);
563 for (auto& R : m_reactors) {
564 size_t offset = R->offset();
565 R->applySensitivity(p);
566 R->evalDae(t, y.subspan(offset, R->neq()),
567 ydot.subspan(offset, R->neq()),
568 residual.subspan(offset, R->neq()));
569 R->resetSensitivity(p);
570 }
571 checkFinite("ydot", ydot);
572}
573
574void ReactorNet::getConstraints(span<double> constraints)
575{
576 for (auto& R : m_reactors) {
577 R->getConstraints(constraints.subspan(R->offset(), R->neq()));
578 }
579}
580
581double ReactorNet::sensitivity(size_t k, size_t p)
582{
583 initialize();
584 if (p >= m_sens_params.size()) {
585 throw IndexError("ReactorNet::sensitivity",
586 "m_sens_params", p, m_sens_params.size());
587 }
588 double denom = m_integ->solution(k);
589 if (denom == 0.0) {
590 denom = SmallNumber;
591 }
592 return m_integ->sensitivity(k, p) / denom;
593}
594
595void ReactorNet::evalJacobian(double t, span<double> y, span<double> ydot,
596 span<const double> p, Array2D* j)
597{
598 //evaluate the unperturbed ydot
599 eval(t, y, ydot, p);
600 for (size_t n = 0; n < m_nv; n++) {
601 // perturb x(n)
602 double ysave = y[n];
603 double dy = m_atol[n] + fabs(ysave)*m_rtol;
604 y[n] = ysave + dy;
605 dy = y[n] - ysave;
606
607 // calculate perturbed residual
608 eval(t, y, m_ydot, p);
609
610 // compute nth column of Jacobian
611 for (size_t m = 0; m < m_nv; m++) {
612 j->value(m,n) = (m_ydot[m] - ydot[m])/dy;
613 }
614 y[n] = ysave;
615 }
616}
617
618void ReactorNet::updateState(span<const double> y)
619{
620 checkFinite("y", y);
621 for (auto& R : m_reactors) {
622 R->updateState(y.subspan(R->offset(), R->neq()));
623 }
624}
625
626void ReactorNet::getDerivative(int k, span<double> dky)
627{
628 initialize();
629 auto cvode_dky = m_integ->derivative(m_time, k);
630 std::copy(cvode_dky.begin(), cvode_dky.end(), dky.begin());
631}
632
633void ReactorNet::setAdvanceLimits(span<const double> limits)
634{
635 initialize();
636 for (auto& R : m_bulkReactors) {
637 R->setAdvanceLimits(limits.subspan(R->offset(), R->neq()));
638 }
639}
640
642{
643 bool has_limit = false;
644 for (size_t n = 0; n < m_bulkReactors.size(); n++) {
645 has_limit |= m_bulkReactors[n]->hasAdvanceLimits();
646 }
647 return has_limit;
648}
649
650bool ReactorNet::getAdvanceLimits(span<double> limits) const
651{
652 bool has_limit = false;
653 for (auto& R : m_bulkReactors) {
654 has_limit |= R->getAdvanceLimits(limits.subspan(R->offset(), R->neq()));
655 }
656 return has_limit;
657}
658
659void ReactorNet::getState(span<double> y)
660{
661 for (auto& R : m_reactors) {
662 R->getState(y.subspan(R->offset(), R->neq()));
663 }
664}
665
666void ReactorNet::getStateDae(span<double> y, span<double> ydot)
667{
668 // Iterate in reverse order so that surfaces will be handled first and up-to-date
669 // values of the surface production rates of bulk species will be available when
670 // bulk reactors are processed.
671 // TODO: Replace with view::reverse once Cantera requires C++20
672 for (size_t n = m_reactors.size(); n != 0 ; n--) {
673 auto& R = m_reactors[n-1];
674 R->getStateDae(y.subspan(R->offset(), R->neq()),
675 ydot.subspan(R->offset(), R->neq()));
676 }
677}
678
679size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
680{
681 initialize();
682 return m_reactors[reactor]->offset()
683 + m_reactors[reactor]->componentIndex(component);
684}
685
686string ReactorNet::componentName(size_t i) const
687{
688 size_t iTot = 0;
689 size_t i0 = i;
690 for (auto r : m_reactors) {
691 if (i < r->neq()) {
692 return r->name() + ": " + r->componentName(i);
693 } else {
694 i -= r->neq();
695 }
696 iTot += r->neq();
697 }
698 throw IndexError("ReactorNet::componentName", "component", i0, iTot);
699}
700
701double ReactorNet::upperBound(size_t i) const
702{
703 size_t iTot = 0;
704 size_t i0 = i;
705 for (auto r : m_reactors) {
706 if (i < r->neq()) {
707 return r->upperBound(i);
708 } else {
709 i -= r->neq();
710 }
711 iTot += r->neq();
712 }
713 throw IndexError("ReactorNet::upperBound", "upperBound", i0, iTot);
714}
715
716double ReactorNet::lowerBound(size_t i) const
717{
718 size_t iTot = 0;
719 size_t i0 = i;
720 for (auto r : m_reactors) {
721 if (i < r->neq()) {
722 return r->lowerBound(i);
723 } else {
724 i -= r->neq();
725 }
726 iTot += r->neq();
727 }
728 throw IndexError("ReactorNet::lowerBound", "lowerBound", i0, iTot);
729}
730
731void ReactorNet::resetBadValues(span<double> y) {
732 for (auto& R : m_reactors) {
733 R->resetBadValues(y.subspan(R->offset(), R->neq()));
734 }
735}
736
738 const string& name, double value, double scale)
739{
741 throw CanteraError("ReactorNet::registerSensitivityParameter",
742 "Sensitivity parameters cannot be added after the"
743 "integrator has been initialized.");
744 }
745 m_paramNames.push_back(name);
746 m_sens_params.push_back(value);
747 m_paramScales.push_back(scale);
748 return m_sens_params.size() - 1;
749}
750
752{
753 // Apply given settings to all reactors
754 for (auto& R : m_bulkReactors) {
755 R->setDerivativeSettings(settings);
756 }
757}
758
760{
761 if (m_integ) {
762 return m_integ->solverStats();
763 } else {
764 return AnyMap();
765 }
766}
767
769{
770 if (m_integ) {
771 return m_integ->linearSolverType();
772 } else {
773 return "";
774 }
775}
776
777void ReactorNet::preconditionerSolve(span<const double> rhs, span<double> output)
778{
779 if (!m_integ) {
780 throw CanteraError("ReactorNet::preconditionerSolve",
781 "Must only be called after ReactorNet is initialized.");
782 }
783 m_integ->preconditionerSolve(rhs, output);
784}
785
786void ReactorNet::preconditionerSetup(double t, span<const double> y, double gamma)
787{
788 // ensure state is up to date.
789 updateState(y);
790 // get the preconditioner
791 auto precon = std::dynamic_pointer_cast<EigenSparseJacobian>(
792 m_integ->preconditioner());
793 // Reset preconditioner
794 precon->reset();
795 // Set gamma value for M =I - gamma*J
796 precon->setGamma(gamma);
797 // Make a copy of state to adjust it for preconditioner
798 vector<double> yCopy(m_nv);
799 // Get state of reactor
800 getState(yCopy);
801 // transform state based on preconditioner rules
802 precon->stateAdjustment(yCopy);
803 // update network with adjusted state
804 updateState(yCopy);
805 // Get jacobians and give elements to preconditioners
806 vector<Eigen::Triplet<double>> trips;
807 for (auto& R : m_reactors) {
808 R->getJacobianElements(trips);
809 }
810 precon->setFromTriplets(trips);
811 // post reactor setup operations
812 precon->updatePreconditioner();
813}
814
816{
817 if (!m_integ) {
818 throw CanteraError("ReactorNet::updatePreconditioner",
819 "Must only be called after ReactorNet is initialized.");
820 }
821 auto precon = m_integ->preconditioner();
822 precon->setGamma(gamma);
823 precon->updatePreconditioner();
824}
825
827 // check for non-mole-based reactors and throw an error otherwise
828 for (auto reactor : m_bulkReactors) {
830 throw CanteraError("ReactorNet::checkPreconditionerSupported",
831 "Preconditioning is only supported for type *MoleReactor,\n"
832 "Reactor type given: '{}'.",
833 reactor->type());
834 }
835 }
836}
837
838SteadyReactorSolver::SteadyReactorSolver(ReactorNet* net, span<const double> x0)
839 : m_net(net)
840{
841 m_size = m_net->neq();
842 m_jac = newSystemJacobian("eigen-sparse-direct");
844 m_initialState.assign(x0.begin(), x0.end());
845 setInitialGuess(x0);
846 setJacobianPerturbation(m_jacobianRelPerturb, 1000 * m_net->atol(),
847 m_jacobianThreshold);
848 m_mask.assign(m_size, 1);
849 size_t start = 0;
850 for (size_t i = 0; i < net->nReactors(); i++) {
851 auto& R = net->reactor(i);
852 R.updateState(x0.subspan(start, R.neq()));
853 auto algebraic = R.initializeSteady();
854 for (auto& m : algebraic) {
855 m_mask[start + m] = 0;
856 }
857 start += R.neq();
858 }
859}
860
861void SteadyReactorSolver::eval(span<const double> x, span<double> r,
862 double rdt, int count)
863{
864 if (rdt < 0.0) {
865 rdt = m_rdt;
866 }
867 m_net->evalSteady(x, r);
868 for (size_t i = 0; i < size(); i++) {
869 r[i] -= (x[i] - m_initialState[i]) * rdt * m_mask[i];
870 }
871}
872
873void SteadyReactorSolver::initTimeInteg(double dt, span<const double> x)
874{
876 m_initialState.assign(x.begin(), x.end());
877}
878
879void SteadyReactorSolver::evalJacobian(span<const double> x0)
880{
881 m_jac->reset();
882 clock_t t0 = clock();
883 m_work1.resize(size());
884 m_work2.resize(size());
885 eval(x0, m_work1, 0.0, 0);
886 vector<double> perturbed(x0.begin(), x0.end());
887 for (size_t j = 0; j < size(); j++) {
888 // perturb x(n); preserve sign(x(n))
889 double xsave = x0[j];
890 double dx = fabs(xsave) * m_jacobianRelPerturb + m_jacobianAbsPerturb;
891 if (xsave < 0) {
892 dx = -dx;
893 }
894 perturbed[j] = xsave + dx;
895 double rdx = 1.0 / (perturbed[j] - xsave);
896
897 // calculate perturbed residual
898 eval(perturbed, m_work2, 0.0, 0);
899
900 // compute nth column of Jacobian
901 for (size_t i = 0; i < size(); i++) {
902 double delta = m_work2[i] - m_work1[i];
903 if (std::abs(delta) > m_jacobianThreshold || i == j) {
904 m_jac->setValue(i, j, delta * rdx);
905 }
906 }
907 perturbed[j] = xsave;
908 }
909 // Restore system to unperturbed state
910 m_net->updateState(x0);
911
912 m_jac->updateElapsed(double(clock() - t0) / CLOCKS_PER_SEC);
913 m_jac->incrementEvals();
914 m_jac->setAge(0);
915}
916
917double SteadyReactorSolver::weightedNorm(span<const double> step) const
918{
919 double sum = 0.0;
920 const double* x = m_state->data();
921 for (size_t i = 0; i < size(); i++) {
922 double ewt = m_net->rtol()*fabs(x[i]) + m_net->atol();
923 double f = step[i] / ewt;
924 sum += f*f;
925 }
926 return sqrt(sum / size());
927}
928
930{
931 return m_net->componentName(i);
932}
933
935{
936 return m_net->upperBound(i);
937}
938
940{
941 return m_net->lowerBound(i);
942}
943
945{
946 m_net->resetBadValues(x);
947}
948
949void SteadyReactorSolver::writeDebugInfo(const string& header_suffix,
950 const string& message, int loglevel, int attempt_counter)
951{
952 if (loglevel >= 6 && !m_state->empty()) {
953 const auto& state = *m_state;
954 writelog("Current state ({}):\n[", header_suffix);
955 for (size_t i = 0; i < state.size() - 1; i++) {
956 writelog("{}, ", state[i]);
957 }
958 writelog("{}]\n", state.back());
959 }
960 if (loglevel >= 7 && !m_xnew.empty()) {
961 writelog("Current residual ({}):\n[", header_suffix);
962 for (size_t i = 0; i < m_xnew.size() - 1; i++) {
963 writelog("{}, ", m_xnew[i]);
964 }
965 writelog("{}]\n", m_xnew.back());
966 }
967}
968
969shared_ptr<ReactorNet> newReactorNet(span<shared_ptr<ReactorBase>> reactors)
970{
971 return make_shared<ReactorNet>(reactors);
972}
973
974}
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:151
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:208
bool suppressErrors() const
Get current state of error suppression.
Definition FuncEval.h:190
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:205
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:229
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:247
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
Definition Integrator.h:253
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition Integrator.h:239
Base class for reactor objects.
Definition ReactorBase.h:51
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.
virtual string type() const
String indicating the reactor model implemented.
Definition ReactorBase.h:76
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
virtual void updateState(span< const double > y)
Set the state of the reactor to correspond to the state vector y.
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.
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.
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
virtual void updateConnected(bool updatePressure)
Update state information needed by connected reactors, flow devices, and walls.
A class representing a network of connected reactors.
Definition ReactorNet.h:30
virtual void getDerivative(int k, span< double > dky)
Return k-th derivative at the current state of the system.
void setLinearSolverType(const string &linSolverType="DENSE")
Set the type of linear solver used in the integration.
bool getAdvanceLimits(span< double > limits) const
Retrieve absolute step size limits during advance.
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...
double m_ybase_time
Base time corresponding to m_ybase.
Definition ReactorNet.h:444
void initialize()
Initialize the reactor network.
void evalDae(double t, span< const double > y, span< const double > ydot, span< const double > p, span< double > residual) override
eval coupling for IDA / DAEs
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:249
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition ReactorNet.h:185
vector< double > m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
Definition ReactorNet.h:450
double m_initial_time
The initial value of the independent variable in the system.
Definition ReactorNet.h:412
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void updateNames(ReactorBase &r)
Create reproducible names for reactors and walls/connectors.
double m_time
The independent variable in the system.
Definition ReactorNet.h:409
AnyMap solverStats() const
Get solver stats from integrator.
bool m_limit_check_active
Indicates whether the advance-limit root check is active for the current call to advance(t,...
Definition ReactorNet.h:447
map< string, int > m_counts
Map used for default name generation.
Definition ReactorNet.h:404
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...
bool m_integratorInitialized
True if the integrator has been initialized at least once.
Definition ReactorNet.h:414
string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
virtual void getEstimate(double time, int k, span< double > yest)
Estimate a future state based on current derivatives.
size_t nRootFunctions() const override
Root finding is enabled only while enforcing advance limits.
void evalRootFunctions(double t, span< const double > y, span< double > gout) override
Evaluate the advance-limit root function used to stop integration once a limit is met.
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
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...
void getConstraints(span< double > constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
double m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
Definition ReactorNet.h:427
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
void evalSteady(span< const double > y, span< double > residual)
Evaluate the governing equations adapted for the steady-state solver.
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.
ReactorNet(shared_ptr< ReactorBase > reactor)
Create reactor network containing single reactor.
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 evalJacobian(double t, span< double > y, span< double > ydot, span< const double > p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
void updateState(span< const double > y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
void preconditionerSolve(span< const double > rhs, span< double > output) override
Evaluate the linear system Ax=b where A is the preconditioner.
bool m_needIntegratorInit
True if integrator needs to be (re)initialized.
Definition ReactorNet.h:415
void preconditionerSetup(double t, span< const double > y, double gamma) override
Evaluate the setup processes for the Jacobian preconditioner.
double rtol()
Relative tolerance.
Definition ReactorNet.h:95
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 resetBadValues(span< double > y)
Reset physically or mathematically problematic values, such as negative species concentrations.
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:432
double atol()
Absolute integration tolerance.
Definition ReactorNet.h:100
bool hasAdvanceLimits() const
Check whether ReactorNet object uses advance limits.
vector< double > m_ybase
Base state used for evaluating advance limits during a single advance() call when root-finding is ena...
Definition ReactorNet.h:442
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
void setAdvanceLimits(span< const double > limits)
Set absolute step size limits during advance.
double lowerBound(size_t i) const
Get the lower bound on the i-th component of the global state vector.
virtual void checkPreconditionerSupported() const
Check that preconditioning is supported by all reactors in the network.
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
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 getState(span< double > y) override
Fill in the vector y with the current state of the system.
vector< string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition ReactorNet.h:435
void getStateDae(span< double > y, span< double > ydot) override
Fill in the vectors y and ydot with the current state of the system.
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
void eval(double t, span< const double > y, span< double > ydot, span< const double > p) override
Evaluate the right-hand-side ODE function.
Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
Definition ReactorNet.h:460
double weightedNorm(span< const double > step) const override
Compute the weighted norm of step.
vector< double > m_initialState
Initial value of each state variable.
Definition ReactorNet.h:479
string componentName(size_t i) const override
Get the name of the i-th component of the state vector.
void initTimeInteg(double dt, span< const double > x) override
Prepare for time stepping beginning with solution x and timestep dt.
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
void evalJacobian(span< const double > x0) override
Evaluates the Jacobian at x0 using finite differences.
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 resetBadValues(span< double > x) override
Reset values such as negative species concentrations.
void eval(span< const double > x, span< double > r, double rdt=-1.0, int count=1) override
Evaluate the residual function.
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
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.
void getState(span< double > x) const
Get the converged steady-state solution after calling solve().
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, span< const 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.
vector< int > m_mask
Transient mask.
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
void solve(int loglevel=0)
Solve the steady-state problem, taking internal timesteps as necessary until the Newton solver can co...
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: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:118
Integrator * newIntegrator(const string &itype)
Create new Integrator object.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:183
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
shared_ptr< ReactorNet > newReactorNet(span< shared_ptr< ReactorBase > > reactors)
Create a reactor network containing one or more coupled reactors.
@ 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:161
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...