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