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