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;
285 AnyMap settings;
286 settings["skip-nonideal"] = true;
287 setDerivativeSettings(settings);
289}
290
292{
293 integrator().setMaxSteps(nmax);
294}
295
297{
298 return integrator().maxSteps();
299}
300
301void ReactorNet::advance(double time)
302{
303 initialize();
304 m_integ->integrate(time);
305 m_time = m_integ->currentTime();
306 updateState(m_integ->solution());
307}
308
309double ReactorNet::advance(double time, bool applylimit)
310{
311 initialize();
312 if (!applylimit || !hasAdvanceLimits()) {
313 // No limit enforcement requested; integrate to requested time
314 advance(time);
315 return time;
316 }
317
318 // Enable root-based limit detection and set the base state to the current state
319 m_ybase.assign(m_nv, 0.0);
323 m_integ->setRootFunctionCount(nRootFunctions());
324
325 // Integrate toward the requested time; integrator will return early if a limit is
326 // reached (CV_ROOT_RETURN). The try/catch ensures the temporary root-finding state
327 // is cleared even when CVODE throws so subsequent calls start clean.
328 try {
329 m_integ->integrate(time);
330 } catch (...) {
331 m_limit_check_active = false;
332 m_integ->setRootFunctionCount(nRootFunctions());
333 throw;
334 }
335 m_time = m_integ->currentTime();
336
337 // Update reactor states to match the integrator solution at the time reached
338 // (which may be earlier than 'time' if a limit was triggered)
339 updateState(m_integ->solution());
340
341 // Disable limit checking after this call
342 m_limit_check_active = false;
343 m_integ->setRootFunctionCount(nRootFunctions());
344
345 // When a root event stopped integration before reaching the requested time, report
346 // the most limiting component and details about the step.
347 if (m_verbose && m_time < time) {
348 // Ensure limits are available
349 if (m_advancelimits.size() != m_nv) {
350 m_advancelimits.assign(m_nv, -1.0);
351 }
352 getAdvanceLimits(m_advancelimits);
353 auto ycurr = m_integ->solution();
354 size_t jmax = npos;
355 double max_ratio = -1.0;
356 double best_limit = 0.0;
357 for (size_t j = 0; j < m_nv; j++) {
358 double lim = m_advancelimits[j];
359 if (lim > 0.0) {
360 double delta = std::abs(ycurr[j] - m_ybase[j]);
361 double ratio = delta / lim;
362 if (ratio > max_ratio) {
363 max_ratio = ratio;
364 jmax = j;
365 best_limit = lim;
366 }
367 }
368 }
369 if (jmax != npos) {
370 double dt = m_time - m_ybase_time;
371 double y_start = m_ybase[jmax];
372 double y_end = ycurr[jmax];
373 double delta = y_end - y_start;
374 writelog(" Advance limit triggered for component {:d} (dt = {:9.4g}):"
375 " y_start = {:11.6g}, y_end = {:11.6g},"
376 " delta = {:11.6g}, limit = {:9.4g}\n",
377 jmax, dt, y_start, y_end, delta, best_limit);
378 }
379 }
380
381 // m_time is tracked via callbacks during integration
382 return m_time;
383}
384
386{
387 initialize();
388 m_time = m_integ->step(m_time + 1.0);
389 updateState(m_integ->solution());
390 return m_time;
391}
392
393void ReactorNet::solveSteady(int loglevel)
394{
395 initialize();
396 vector<double> y(neq());
397 getState(y);
398 SteadyReactorSolver solver(this, y);
400 solver.solve(loglevel);
401 solver.getState(y);
402 updateState(y);
403}
404
405Eigen::SparseMatrix<double> ReactorNet::steadyJacobian(double rdt)
406{
407 initialize();
408 vector<double> y0(neq());
409 vector<double> y1(neq());
410 getState(y0);
411 SteadyReactorSolver solver(this, y0);
412 solver.evalJacobian(y0);
413 if (rdt) {
414 solver.linearSolver()->updateTransient(rdt, solver.transientMask());
415 }
416 return std::dynamic_pointer_cast<EigenSparseJacobian>(solver.linearSolver())->jacobian();
417}
418
419void ReactorNet::getEstimate(double time, int k, span<double> yest)
420{
421 initialize();
422 auto cvode_dky = m_integ->solution();
423 std::copy(cvode_dky.begin(), cvode_dky.end(), yest.begin());
424
425 // Taylor expansion
426 double factor = 1.;
427 double deltat = time - m_time;
428 for (int n = 1; n <= k; n++) {
429 factor *= deltat / n;
430 cvode_dky = m_integ->derivative(m_time, n);
431 for (size_t j = 0; j < m_nv; j++) {
432 yest[j] += factor * cvode_dky[j];
433 }
434 }
435}
436
438{
439 if (m_integ) {
440 return m_integ->lastOrder();
441 } else {
442 return 0;
443 }
444}
445
447{
448 return (m_limit_check_active && hasAdvanceLimits()) ? 1 : 0;
449}
450
451void ReactorNet::evalRootFunctions(double t, span<const double> y, span<double> gout)
452{
453 // Default: no root detected
454 double g = 1.0;
455
457 // Ensure limits vector is current
458 if (m_advancelimits.size() != m_nv) {
459 m_advancelimits.assign(m_nv, -1.0);
460 }
461 getAdvanceLimits(m_advancelimits);
462
463 double max_ratio = 0.0;
464 for (size_t i = 0; i < m_nv; i++) {
465 double lim = m_advancelimits[i];
466 if (lim > 0.0) {
467 double delta = std::abs(y[i] - m_ybase[i]);
468 double ratio = delta / lim;
469 if (ratio > max_ratio) {
470 max_ratio = ratio;
471 }
472 }
473 }
474 g = 1.0 - max_ratio; // root at g = 0 when any component reaches its limit
475 }
476
477 gout[0] = g;
478}
479
481{
482 // ensure that reactors and components have reproducible names
484
485 for (size_t i=0; i<r.nWalls(); i++) {
486 auto& w = r.wall(i);
488 if (w.left().type() == "Reservoir") {
489 w.left().setDefaultName(m_counts);
490 }
491 if (w.right().type() == "Reservoir") {
492 w.right().setDefaultName(m_counts);
493 }
494 }
495
496 for (size_t i=0; i<r.nInlets(); i++) {
497 auto& in = r.inlet(i);
499 if (in.in().type() == "Reservoir") {
500 in.in().setDefaultName(m_counts);
501 }
502 }
503
504 for (size_t i=0; i<r.nOutlets(); i++) {
505 auto& out = r.outlet(i);
507 if (out.out().type() == "Reservoir") {
508 out.out().setDefaultName(m_counts);
509 }
510 }
511
512 for (size_t i=0; i<r.nSurfs(); i++) {
514 }
515}
516
518 if (m_integ == nullptr) {
519 throw CanteraError("ReactorNet::integrator",
520 "Integrator has not been instantiated. Add one or more reactors first.");
521 }
522 return *m_integ;
523}
524
525void ReactorNet::eval(double t, span<const double> y, span<double> ydot,
526 span<const double> p)
527{
528 m_time = t;
529 updateState(y);
530 m_LHS.assign(m_nv, 1);
531 m_RHS.assign(m_nv, 0);
532 for (auto& R : m_reactors) {
533 size_t offset = R->offset();
534 R->applySensitivity(p);
535 R->eval(t, span<double>(m_LHS.data() + offset, R->neq()),
536 span<double>(m_RHS.data() + offset, R->neq()));
537 for (size_t i = offset; i < offset + R->neq(); i++) {
538 ydot[i] = m_RHS[i] / m_LHS[i];
539 }
540 R->resetSensitivity(p);
541 }
542 checkFinite("ydot", ydot);
543}
544
545void ReactorNet::evalSteady(span<const double> y, span<double> residual)
546{
547 updateState(y);
548 m_LHS.assign(m_nv, 1);
549 m_RHS.assign(m_nv, 0);
550 for (auto& R : m_reactors) {
551 size_t offset = R->offset();
552 R->evalSteady(m_time, span<double>(m_LHS.data() + offset, R->neq()),
553 span<double>(m_RHS.data() + offset, R->neq()));
554 for (size_t i = offset; i < offset + R->neq(); i++) {
555 residual[i] = m_RHS[i] / m_LHS[i];
556 }
557 }
558 checkFinite("residual", residual);
559}
560
561void ReactorNet::evalDae(double t, span<const double> y, span<const double> ydot,
562 span<const double> p, span<double> residual)
563{
564 m_time = t;
565 updateState(y);
566 for (auto& R : m_reactors) {
567 size_t offset = R->offset();
568 R->applySensitivity(p);
569 R->evalDae(t, y.subspan(offset, R->neq()),
570 ydot.subspan(offset, R->neq()),
571 residual.subspan(offset, R->neq()));
572 R->resetSensitivity(p);
573 }
574 checkFinite("ydot", ydot);
575}
576
577void ReactorNet::getConstraints(span<double> constraints)
578{
579 for (auto& R : m_reactors) {
580 R->getConstraints(constraints.subspan(R->offset(), R->neq()));
581 }
582}
583
584double ReactorNet::sensitivity(size_t k, size_t p)
585{
586 initialize();
587 if (p >= m_sens_params.size()) {
588 throw IndexError("ReactorNet::sensitivity",
589 "m_sens_params", p, m_sens_params.size());
590 }
591 double denom = m_integ->solution(k);
592 if (denom == 0.0) {
593 denom = SmallNumber;
594 }
595 return m_integ->sensitivity(k, p) / denom;
596}
597
598void ReactorNet::evalJacobian(double t, span<double> y, span<double> ydot,
599 span<const double> p, Array2D* j)
600{
601 //evaluate the unperturbed ydot
602 eval(t, y, ydot, p);
603 for (size_t n = 0; n < m_nv; n++) {
604 // perturb x(n)
605 double ysave = y[n];
606 double dy = m_atol[n] + fabs(ysave)*m_rtol;
607 y[n] = ysave + dy;
608 dy = y[n] - ysave;
609
610 // calculate perturbed residual
611 eval(t, y, m_ydot, p);
612
613 // compute nth column of Jacobian
614 for (size_t m = 0; m < m_nv; m++) {
615 j->value(m,n) = (m_ydot[m] - ydot[m])/dy;
616 }
617 y[n] = ysave;
618 }
619}
620
621void ReactorNet::updateState(span<const double> y)
622{
623 checkFinite("y", y);
624 for (auto& R : m_reactors) {
625 R->updateState(y.subspan(R->offset(), R->neq()));
626 }
627}
628
629void ReactorNet::getDerivative(int k, span<double> dky)
630{
631 initialize();
632 auto cvode_dky = m_integ->derivative(m_time, k);
633 std::copy(cvode_dky.begin(), cvode_dky.end(), dky.begin());
634}
635
636void ReactorNet::setAdvanceLimits(span<const double> limits)
637{
638 initialize();
639 for (auto& R : m_bulkReactors) {
640 R->setAdvanceLimits(limits.subspan(R->offset(), R->neq()));
641 }
642}
643
645{
646 bool has_limit = false;
647 for (size_t n = 0; n < m_bulkReactors.size(); n++) {
648 has_limit |= m_bulkReactors[n]->hasAdvanceLimits();
649 }
650 return has_limit;
651}
652
653bool ReactorNet::getAdvanceLimits(span<double> limits) const
654{
655 bool has_limit = false;
656 for (auto& R : m_bulkReactors) {
657 has_limit |= R->getAdvanceLimits(limits.subspan(R->offset(), R->neq()));
658 }
659 return has_limit;
660}
661
662void ReactorNet::getState(span<double> y)
663{
664 for (auto& R : m_reactors) {
665 R->getState(y.subspan(R->offset(), R->neq()));
666 }
667}
668
669void ReactorNet::getStateDae(span<double> y, span<double> ydot)
670{
671 // Iterate in reverse order so that surfaces will be handled first and up-to-date
672 // values of the surface production rates of bulk species will be available when
673 // bulk reactors are processed.
674 // TODO: Replace with view::reverse once Cantera requires C++20
675 for (size_t n = m_reactors.size(); n != 0 ; n--) {
676 auto& R = m_reactors[n-1];
677 R->getStateDae(y.subspan(R->offset(), R->neq()),
678 ydot.subspan(R->offset(), R->neq()));
679 }
680}
681
682size_t ReactorNet::globalComponentIndex(const string& component, size_t reactor)
683{
684 initialize();
685 return m_reactors[reactor]->offset()
686 + m_reactors[reactor]->componentIndex(component);
687}
688
689string ReactorNet::componentName(size_t i) const
690{
691 size_t iTot = 0;
692 size_t i0 = i;
693 for (auto r : m_reactors) {
694 if (i < r->neq()) {
695 return r->name() + ": " + r->componentName(i);
696 } else {
697 i -= r->neq();
698 }
699 iTot += r->neq();
700 }
701 throw IndexError("ReactorNet::componentName", "component", i0, iTot);
702}
703
704double ReactorNet::upperBound(size_t i) const
705{
706 size_t iTot = 0;
707 size_t i0 = i;
708 for (auto r : m_reactors) {
709 if (i < r->neq()) {
710 return r->upperBound(i);
711 } else {
712 i -= r->neq();
713 }
714 iTot += r->neq();
715 }
716 throw IndexError("ReactorNet::upperBound", "upperBound", i0, iTot);
717}
718
719double ReactorNet::lowerBound(size_t i) const
720{
721 size_t iTot = 0;
722 size_t i0 = i;
723 for (auto r : m_reactors) {
724 if (i < r->neq()) {
725 return r->lowerBound(i);
726 } else {
727 i -= r->neq();
728 }
729 iTot += r->neq();
730 }
731 throw IndexError("ReactorNet::lowerBound", "lowerBound", i0, iTot);
732}
733
734void ReactorNet::resetBadValues(span<double> y) {
735 for (auto& R : m_reactors) {
736 R->resetBadValues(y.subspan(R->offset(), R->neq()));
737 }
738}
739
741 const string& name, double value, double scale)
742{
744 throw CanteraError("ReactorNet::registerSensitivityParameter",
745 "Sensitivity parameters cannot be added after the"
746 "integrator has been initialized.");
747 }
748 m_paramNames.push_back(name);
749 m_sens_params.push_back(value);
750 m_paramScales.push_back(scale);
751 return m_sens_params.size() - 1;
752}
753
755{
756 // Apply given settings to all reactors
757 for (auto& R : m_bulkReactors) {
758 R->setDerivativeSettings(settings);
759 }
760}
761
763{
764 if (m_integ) {
765 return m_integ->solverStats();
766 } else {
767 return AnyMap();
768 }
769}
770
772{
773 if (m_integ) {
774 return m_integ->linearSolverType();
775 } else {
776 return "";
777 }
778}
779
780void ReactorNet::preconditionerSolve(span<const double> rhs, span<double> output)
781{
782 if (!m_integ) {
783 throw CanteraError("ReactorNet::preconditionerSolve",
784 "Must only be called after ReactorNet is initialized.");
785 }
786 m_integ->preconditionerSolve(rhs, output);
787}
788
789void ReactorNet::preconditionerSetup(double t, span<const double> y, double gamma)
790{
791 // ensure state is up to date.
792 updateState(y);
793 // get the preconditioner
794 auto precon = std::dynamic_pointer_cast<EigenSparseJacobian>(
795 m_integ->preconditioner());
796 // Reset preconditioner
797 precon->reset();
798 // Set gamma value for M =I - gamma*J
799 precon->setGamma(gamma);
800 // Make a copy of state to adjust it for preconditioner
801 vector<double> yCopy(m_nv);
802 // Get state of reactor
803 getState(yCopy);
804 // transform state based on preconditioner rules
805 precon->stateAdjustment(yCopy);
806 // update network with adjusted state
807 updateState(yCopy);
808 // Get jacobians and give elements to preconditioners
809 vector<Eigen::Triplet<double>> trips;
810 for (auto& R : m_reactors) {
811 R->getJacobianElements(trips);
812 }
813 precon->setFromTriplets(trips);
814 // post reactor setup operations
815 precon->updatePreconditioner();
816}
817
819{
820 if (!m_integ) {
821 throw CanteraError("ReactorNet::updatePreconditioner",
822 "Must only be called after ReactorNet is initialized.");
823 }
824 auto precon = m_integ->preconditioner();
825 precon->setGamma(gamma);
826 precon->updatePreconditioner();
827}
828
830 // check for non-mole-based reactors and throw an error otherwise
831 for (auto reactor : m_bulkReactors) {
833 throw CanteraError("ReactorNet::checkPreconditionerSupported",
834 "Preconditioning is only supported for type *MoleReactor,\n"
835 "Reactor type given: '{}'.",
836 reactor->type());
837 }
838 }
839}
840
841SteadyReactorSolver::SteadyReactorSolver(ReactorNet* net, span<const double> x0)
842 : m_net(net)
843{
844 m_size = m_net->neq();
845 m_jac = newSystemJacobian("eigen-sparse-direct");
847 m_initialState.assign(x0.begin(), x0.end());
848 setInitialGuess(x0);
849 setJacobianPerturbation(m_jacobianRelPerturb, 1000 * m_net->atol(),
850 m_jacobianThreshold);
851 m_mask.assign(m_size, 1);
852 size_t start = 0;
853 for (size_t i = 0; i < net->nReactors(); i++) {
854 auto& R = net->reactor(i);
855 R.updateState(x0.subspan(start, R.neq()));
856 auto algebraic = R.initializeSteady();
857 for (auto& m : algebraic) {
858 m_mask[start + m] = 0;
859 }
860 start += R.neq();
861 }
862}
863
864void SteadyReactorSolver::eval(span<const double> x, span<double> r,
865 double rdt, int count)
866{
867 if (rdt < 0.0) {
868 rdt = m_rdt;
869 }
870 m_net->evalSteady(x, r);
871 for (size_t i = 0; i < size(); i++) {
872 r[i] -= (x[i] - m_initialState[i]) * rdt * m_mask[i];
873 }
874}
875
876void SteadyReactorSolver::initTimeInteg(double dt, span<const double> x)
877{
879 m_initialState.assign(x.begin(), x.end());
880}
881
882void SteadyReactorSolver::evalJacobian(span<const double> x0)
883{
884 m_jac->reset();
885 clock_t t0 = clock();
886 m_work1.resize(size());
887 m_work2.resize(size());
888 eval(x0, m_work1, 0.0, 0);
889 vector<double> perturbed(x0.begin(), x0.end());
890 for (size_t j = 0; j < size(); j++) {
891 // perturb x(n); preserve sign(x(n))
892 double xsave = x0[j];
893 double dx = fabs(xsave) * m_jacobianRelPerturb + m_jacobianAbsPerturb;
894 if (xsave < 0) {
895 dx = -dx;
896 }
897 perturbed[j] = xsave + dx;
898 double rdx = 1.0 / (perturbed[j] - xsave);
899
900 // calculate perturbed residual
901 eval(perturbed, m_work2, 0.0, 0);
902
903 // compute nth column of Jacobian
904 for (size_t i = 0; i < size(); i++) {
905 double delta = m_work2[i] - m_work1[i];
906 if (std::abs(delta) > m_jacobianThreshold || i == j) {
907 m_jac->setValue(i, j, delta * rdx);
908 }
909 }
910 perturbed[j] = xsave;
911 }
912 // Restore system to unperturbed state
913 m_net->updateState(x0);
914
915 m_jac->updateElapsed(double(clock() - t0) / CLOCKS_PER_SEC);
916 m_jac->incrementEvals();
917 m_jac->setAge(0);
918}
919
920double SteadyReactorSolver::weightedNorm(span<const double> step) const
921{
922 double sum = 0.0;
923 const double* x = m_state->data();
924 for (size_t i = 0; i < size(); i++) {
925 double ewt = m_net->rtol()*fabs(x[i]) + m_net->atol();
926 double f = step[i] / ewt;
927 sum += f*f;
928 }
929 return sqrt(sum / size());
930}
931
933{
934 return m_net->componentName(i);
935}
936
938{
939 return m_net->upperBound(i);
940}
941
943{
944 return m_net->lowerBound(i);
945}
946
948{
949 m_net->resetBadValues(x);
950}
951
952void SteadyReactorSolver::writeDebugInfo(const string& header_suffix,
953 const string& message, int loglevel, int attempt_counter)
954{
955 if (loglevel >= 6 && !m_state->empty()) {
956 const auto& state = *m_state;
957 writelog("Current state ({}):\n[", header_suffix);
958 for (size_t i = 0; i < state.size() - 1; i++) {
959 writelog("{}, ", state[i]);
960 }
961 writelog("{}]\n", state.back());
962 }
963 if (loglevel >= 7 && !m_xnew.empty()) {
964 writelog("Current residual ({}):\n[", header_suffix);
965 for (size_t i = 0; i < m_xnew.size() - 1; i++) {
966 writelog("{}, ", m_xnew[i]);
967 }
968 writelog("{}]\n", m_xnew.back());
969 }
970}
971
972shared_ptr<ReactorNet> newReactorNet(span<shared_ptr<ReactorBase>> reactors)
973{
974 return make_shared<ReactorNet>(reactors);
975}
976
977}
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...