Cantera  3.2.0a4
Loading...
Searching...
No Matches
Kinetics.cpp
Go to the documentation of this file.
1/**
2 * @file Kinetics.cpp Declarations for the base class for kinetics managers
3 * (see @ref kineticsmgr and class @link Cantera::Kinetics Kinetics @endlink).
4 *
5 * Kinetics managers calculate rates of progress of species due to
6 * homogeneous or heterogeneous kinetics.
7 */
8
9// This file is part of Cantera. See License.txt in the top-level directory or
10// at https://cantera.org/license.txt for license and copyright information.
11
18#include "cantera/base/global.h"
19#include <unordered_set>
20#include <boost/algorithm/string.hpp>
21
22using namespace std;
23
24namespace Cantera
25{
26
27shared_ptr<Kinetics> Kinetics::clone(
28 const vector<shared_ptr<ThermoPhase>>& phases) const
29{
30 vector<AnyMap> reactionDefs;
31 for (size_t i = 0; i < nReactions(); i++) {
32 reactionDefs.push_back(reaction(i)->parameters());
33 }
34 AnyMap phaseNode = parameters();
35 phaseNode["__fix-duplicate-reactions__"] = true;
36 AnyMap rootNode;
37 rootNode["reactions"] = std::move(reactionDefs);
38 rootNode.applyUnits();
39 return newKinetics(phases, phaseNode, rootNode, phases[0]->root());
40}
41
42void Kinetics::checkReactionIndex(size_t i) const
43{
44 if (i >= nReactions()) {
45 throw IndexError("Kinetics::checkReactionIndex", "reactions", i, nReactions());
46 }
47}
48
50{
51 size_t nRxn = nReactions();
52
53 // Stoichiometry managers
57
58 m_rbuf.resize(nRxn);
59
60 // products are created for positive net rate of progress
62 // reactants are destroyed for positive net rate of progress
64}
65
67{
68 if (nReactions() > ii) {
69 throw ArraySizeError("Kinetics::checkReactionArraySize", ii,
70 nReactions());
71 }
72}
73
74void Kinetics::checkPhaseIndex(size_t m) const
75{
76 if (m >= nPhases()) {
77 throw IndexError("Kinetics::checkPhaseIndex", "phase", m, nPhases());
78 }
79}
80
81void Kinetics::checkPhaseArraySize(size_t mm) const
82{
83 if (nPhases() > mm) {
84 throw ArraySizeError("Kinetics::checkPhaseArraySize", mm, nPhases());
85 }
86}
87
88shared_ptr<ThermoPhase> Kinetics::reactionPhase() const
89{
90 return m_thermo[0];
91}
92
93void Kinetics::checkSpeciesIndex(size_t k) const
94{
95 if (k >= m_kk) {
96 throw IndexError("Kinetics::checkSpeciesIndex", "species", k, m_kk);
97 }
98}
99
101{
102 if (m_kk > kk) {
103 throw ArraySizeError("Kinetics::checkSpeciesArraySize", kk, m_kk);
104 }
105}
106
108{
109 if (flag == "warn" || flag == "error" || flag == "mark-duplicate"
110 || flag == "modify-efficiency")
111 {
112 m_explicit_third_body_duplicates = flag;
113 } else {
114 throw CanteraError("Kinetics::setExplicitThirdBodyDuplicateHandling",
115 "Invalid flag '{}'", flag);
116 }
117}
118
119pair<size_t, size_t> Kinetics::checkDuplicates(bool throw_err, bool fix)
120{
121 //! Map of (key indicating participating species) to reaction numbers
122 map<size_t, vector<size_t>> participants;
123 vector<map<int, double>> net_stoich;
124 std::unordered_set<size_t> unmatched_duplicates;
125 for (size_t i = 0; i < m_reactions.size(); i++) {
126 if (m_reactions[i]->duplicate) {
127 unmatched_duplicates.insert(i);
128 }
129 }
130
131 vector<InputFileError> errs;
132 for (size_t i = 0; i < m_reactions.size(); i++) {
133 // Get data about this reaction
134 unsigned long int key = 0;
135 Reaction& R = *m_reactions[i];
136 net_stoich.emplace_back();
137 map<int, double>& net = net_stoich.back();
138 for (const auto& [name, stoich] : R.reactants) {
139 int k = static_cast<int>(kineticsSpeciesIndex(name));
140 key += k*(k+1);
141 net[-1 -k] -= stoich;
142 }
143 for (const auto& [name, stoich] : R.products) {
144 int k = static_cast<int>(kineticsSpeciesIndex(name));
145 key += k*(k+1);
146 net[1+k] += stoich;
147 }
148
149 // Compare this reaction to others with similar participants
150 vector<size_t>& related = participants[key];
151 for (size_t m : related) {
152 Reaction& other = *m_reactions[m];
153 if (R.duplicate && other.duplicate) {
154 // marked duplicates
155 unmatched_duplicates.erase(i);
156 unmatched_duplicates.erase(m);
157 continue;
158 } else if (R.type() != other.type()) {
159 continue; // different reaction types
160 } else if (R.rate() && other.rate()
161 && R.rate()->type() != other.rate()->type())
162 {
163 continue; // different rate parameterizations
164 }
165 double c = checkDuplicateStoich(net_stoich[i], net_stoich[m]);
166 if (c == 0) {
167 continue; // stoichiometries differ (not by a multiple)
168 } else if (c < 0.0 && !R.reversible && !other.reversible) {
169 continue; // irreversible reactions in opposite directions
170 } else if (R.usesThirdBody() && other.usesThirdBody()) {
171 ThirdBody& tb1 = *(R.thirdBody());
172 ThirdBody& tb2 = *(other.thirdBody());
173 bool thirdBodyOk = true;
174 for (size_t k = 0; k < nTotalSpecies(); k++) {
175 string s = kineticsSpeciesName(k);
176 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
177 // non-zero third body efficiencies for species `s` in
178 // both reactions
179 thirdBodyOk = false;
180 break;
181 }
182 }
183 if (thirdBodyOk) {
184 continue; // No overlap in third body efficiencies
185 } else if ((tb1.name() == "M") != (tb2.name() == "M")) {
186 // Exactly one of the reactions uses M as the third body
187 if (m_explicit_third_body_duplicates == "mark-duplicate") {
188 R.duplicate = true;
189 other.duplicate = true;
190 continue;
191 } else if (m_explicit_third_body_duplicates == "modify-efficiency") {
192 if (tb1.name() == "M") {
193 tb1.efficiencies[tb2.name()] = 0.0;
194 } else {
195 tb2.efficiencies[tb1.name()] = 0.0;
196 }
197 continue;
198 } else if (m_explicit_third_body_duplicates == "warn") {
199 InputFileError msg("Kinetics::checkDuplicates",
200 R.input, other.input,
201 "Undeclared duplicate third body reactions with a common "
202 "third body detected.\nAdd the field "
203 "'explicit-third-body-duplicates: mark-duplicate' or\n"
204 "'explicit-third-body-duplicates: modify-efficiency' to "
205 "the YAML phase entry\nto choose how these reactions "
206 "should be handled and suppress this warning.\n"
207 "Reaction {}: {}\nReaction {}: {}\n",
208 m+1, R.equation(), i+1, other.equation());
209 if (!fix) {
210 warn_user("Kinetics::checkDuplicates", msg.what());
211 }
212 continue;
213 } // else m_explicit_third_body_duplicates == "error"
214 }
215 }
216 if (throw_err) {
217 errs.emplace_back("Kinetics::checkDuplicates",
218 R.input, other.input,
219 "Undeclared duplicate reactions detected:\n"
220 "Reaction {}: {}\nReaction {}: {}\n",
221 m+1, R.equation(), i+1, other.equation());
222 } else if (fix) {
223 R.duplicate = true;
224 other.duplicate = true;
225 unmatched_duplicates.erase(i);
226 unmatched_duplicates.erase(m);
227 } else {
228 return {i,m};
229 }
230 }
231 participants[key].push_back(i);
232 }
233 if (unmatched_duplicates.size()) {
234 for (auto i : unmatched_duplicates) {
235 if (throw_err) {
236 errs.emplace_back("Kinetics::checkDuplicates",
237 m_reactions[i]->input,
238 "No duplicate found for declared duplicate reaction number {}"
239 " ({})", i, m_reactions[i]->equation());
240 } else if (fix) {
241 m_reactions[i]->duplicate = false;
242 } else {
243 return {i, i};
244 }
245 }
246 }
247 if (errs.empty()) {
248 return {npos, npos};
249 } else if (errs.size() == 1) {
250 throw errs[0];
251 } else {
252 fmt::memory_buffer msg;
253 for (const auto& err : errs) {
254 fmt_append(msg, "\n{}\n", err.getMessage());
255 }
256 throw CanteraError("Kinetics::checkDuplicates", to_string(msg));
257 }
258}
259
260double Kinetics::checkDuplicateStoich(map<int, double>& r1, map<int, double>& r2) const
261{
262 std::unordered_set<int> keys; // species keys (k+1 or -k-1)
263 for (auto& [speciesKey, stoich] : r1) {
264 keys.insert(speciesKey);
265 }
266 for (auto& [speciesKey, stoich] : r2) {
267 keys.insert(speciesKey);
268 }
269 int k1 = r1.begin()->first;
270 // check for duplicate written in the same direction
271 double ratio = 0.0;
272 if (r1[k1] && r2[k1]) {
273 ratio = r2[k1]/r1[k1];
274 bool different = false;
275 for (int k : keys) {
276 if ((r1[k] && !r2[k]) ||
277 (!r1[k] && r2[k]) ||
278 (r1[k] && fabs(r2[k]/r1[k] - ratio) > 1.e-8)) {
279 different = true;
280 break;
281 }
282 }
283 if (!different) {
284 return ratio;
285 }
286 }
287
288 // check for duplicate written in the reverse direction
289 if (r1[k1] == 0.0 || r2[-k1] == 0.0) {
290 return 0.0;
291 }
292 ratio = r2[-k1]/r1[k1];
293 for (int k : keys) {
294 if ((r1[k] && !r2[-k]) ||
295 (!r1[k] && r2[-k]) ||
296 (r1[k] && fabs(r2[-k]/r1[k] - ratio) > 1.e-8)) {
297 return 0.0;
298 }
299 }
300 return ratio;
301}
302
303string Kinetics::kineticsSpeciesName(size_t k) const
304{
305 for (size_t n = m_start.size()-1; n != npos; n--) {
306 if (k >= m_start[n]) {
307 return thermo(n).speciesName(k - m_start[n]);
308 }
309 }
310 return "<unknown>";
311}
312
313size_t Kinetics::kineticsSpeciesIndex(const string& nm) const
314{
315 for (size_t n = 0; n < m_thermo.size(); n++) {
316 // Check the ThermoPhase object for a match
317 size_t k = thermo(n).speciesIndex(nm);
318 if (k != npos) {
319 return k + m_start[n];
320 }
321 }
322 return npos;
323}
324
326{
327 for (size_t n = 0; n < m_thermo.size(); n++) {
328 size_t k = thermo(n).speciesIndex(nm);
329 if (k != npos) {
330 return thermo(n);
331 }
332 }
333 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
334}
335
336const ThermoPhase& Kinetics::speciesPhase(const string& nm) const
337{
338 for (const auto& thermo : m_thermo) {
339 if (thermo->speciesIndex(nm) != npos) {
340 return *thermo;
341 }
342 }
343 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
344}
345
346size_t Kinetics::speciesPhaseIndex(size_t k) const
347{
348 for (size_t n = m_start.size()-1; n != npos; n--) {
349 if (k >= m_start[n]) {
350 return n;
351 }
352 }
353 throw CanteraError("Kinetics::speciesPhaseIndex",
354 "illegal species index: {}", k);
355}
356
357double Kinetics::reactantStoichCoeff(size_t kSpec, size_t irxn) const
358{
359 return m_reactantStoich.stoichCoeffs().coeff(kSpec, irxn);
360}
361
362double Kinetics::productStoichCoeff(size_t kSpec, size_t irxn) const
363{
364 return m_productStoich.stoichCoeffs().coeff(kSpec, irxn);
365}
366
368{
369 updateROP();
370 std::copy(m_ropf.begin(), m_ropf.end(), fwdROP);
371}
372
374{
375 updateROP();
376 std::copy(m_ropr.begin(), m_ropr.end(), revROP);
377}
378
380{
381 updateROP();
382 std::copy(m_ropnet.begin(), m_ropnet.end(), netROP);
383}
384
385void Kinetics::getReactionDelta(const double* prop, double* deltaProp) const
386{
387 fill(deltaProp, deltaProp + nReactions(), 0.0);
388 // products add
389 m_productStoich.incrementReactions(prop, deltaProp);
390 // reactants subtract
391 m_reactantStoich.decrementReactions(prop, deltaProp);
392}
393
394void Kinetics::getRevReactionDelta(const double* prop, double* deltaProp) const
395{
396 fill(deltaProp, deltaProp + nReactions(), 0.0);
397 // products add
398 m_revProductStoich.incrementReactions(prop, deltaProp);
399 // reactants subtract
400 m_reactantStoich.decrementReactions(prop, deltaProp);
401}
402
404{
405 updateROP();
406
407 // zero out the output array
408 fill(cdot, cdot + m_kk, 0.0);
409
410 // the forward direction creates product species
411 m_productStoich.incrementSpecies(m_ropf.data(), cdot);
412
413 // the reverse direction creates reactant species
414 m_reactantStoich.incrementSpecies(m_ropr.data(), cdot);
415}
416
418{
419 updateROP();
420
421 fill(ddot, ddot + m_kk, 0.0);
422 // the reverse direction destroys products in reversible reactions
423 m_revProductStoich.incrementSpecies(m_ropr.data(), ddot);
424 // the forward direction destroys reactants
425 m_reactantStoich.incrementSpecies(m_ropf.data(), ddot);
426}
427
429{
430 updateROP();
431
432 fill(net, net + m_kk, 0.0);
433 // products are created for positive net rate of progress
434 m_productStoich.incrementSpecies(m_ropnet.data(), net);
435 // reactants are destroyed for positive net rate of progress
436 m_reactantStoich.decrementSpecies(m_ropnet.data(), net);
437}
438
440{
441 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
442 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
443 // the forward direction creates product species
444 getFwdRatesOfProgress_ddT(buf.data());
445 out = m_productStoich.stoichCoeffs() * buf;
446 // the reverse direction creates reactant species
447 getRevRatesOfProgress_ddT(buf.data());
448 out += m_reactantStoich.stoichCoeffs() * buf;
449}
450
452{
453 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
454 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
455 // the forward direction creates product species
456 getFwdRatesOfProgress_ddP(buf.data());
457 out = m_productStoich.stoichCoeffs() * buf;
458 // the reverse direction creates reactant species
459 getRevRatesOfProgress_ddP(buf.data());
460 out += m_reactantStoich.stoichCoeffs() * buf;
461}
462
464{
465 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
466 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
467 // the forward direction creates product species
468 getFwdRatesOfProgress_ddC(buf.data());
469 out = m_productStoich.stoichCoeffs() * buf;
470 // the reverse direction creates reactant species
471 getRevRatesOfProgress_ddC(buf.data());
472 out += m_reactantStoich.stoichCoeffs() * buf;
473}
474
475Eigen::SparseMatrix<double> Kinetics::creationRates_ddX()
476{
477 Eigen::SparseMatrix<double> jac;
478 // the forward direction creates product species
480 // the reverse direction creates reactant species
482 return jac;
483}
484
485Eigen::SparseMatrix<double> Kinetics::creationRates_ddCi()
486{
487 Eigen::SparseMatrix<double> jac;
488 // the forward direction creates product species
490 // the reverse direction creates reactant species
492 return jac;
493}
494
496{
497 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
498 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
499 // the reverse direction destroys products in reversible reactions
500 getRevRatesOfProgress_ddT(buf.data());
501 out = m_revProductStoich.stoichCoeffs() * buf;
502 // the forward direction destroys reactants
503 getFwdRatesOfProgress_ddT(buf.data());
504 out += m_reactantStoich.stoichCoeffs() * buf;
505}
506
508{
509 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
510 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
511 // the reverse direction destroys products in reversible reactions
512 getRevRatesOfProgress_ddP(buf.data());
513 out = m_revProductStoich.stoichCoeffs() * buf;
514 // the forward direction destroys reactants
515 getFwdRatesOfProgress_ddP(buf.data());
516 out += m_reactantStoich.stoichCoeffs() * buf;
517}
518
520{
521 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
522 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
523 // the reverse direction destroys products in reversible reactions
524 getRevRatesOfProgress_ddC(buf.data());
525 out = m_revProductStoich.stoichCoeffs() * buf;
526 // the forward direction destroys reactants
527 getFwdRatesOfProgress_ddC(buf.data());
528 out += m_reactantStoich.stoichCoeffs() * buf;
529}
530
531Eigen::SparseMatrix<double> Kinetics::destructionRates_ddX()
532{
533 Eigen::SparseMatrix<double> jac;
534 // the reverse direction destroys products in reversible reactions
536 // the forward direction destroys reactants
538 return jac;
539}
540
541Eigen::SparseMatrix<double> Kinetics::destructionRates_ddCi()
542{
543 Eigen::SparseMatrix<double> jac;
544 // the reverse direction destroys products in reversible reactions
546 // the forward direction destroys reactants
548 return jac;
549}
550
552{
553 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
554 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
555 getNetRatesOfProgress_ddT(buf.data());
556 out = m_stoichMatrix * buf;
557}
558
560{
561 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
562 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
563 getNetRatesOfProgress_ddP(buf.data());
564 out = m_stoichMatrix * buf;
565}
566
568{
569 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
570 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
571 getNetRatesOfProgress_ddC(buf.data());
572 out = m_stoichMatrix * buf;
573}
574
575Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddX()
576{
578}
579
580Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddCi()
581{
583}
584
585void Kinetics::addThermo(shared_ptr<ThermoPhase> thermo)
586{
587 // the phase with lowest dimensionality is assumed to be the
588 // phase/interface at which reactions take place
589 if (thermo->nDim() <= m_mindim) {
590 if (!m_thermo.empty()) {
591 throw CanteraError("Kinetics::addThermo",
592 "The reacting (lowest dimensional) phase must be added first.");
593 }
594 m_mindim = thermo->nDim();
595 }
596
597 m_thermo.push_back(thermo);
598 m_phaseindex[m_thermo.back()->name()] = nPhases();
600}
601
602void Kinetics::setParameters(const AnyMap& phaseNode) {
603 skipUndeclaredThirdBodies(phaseNode.getBool("skip-undeclared-third-bodies", false));
605 phaseNode.getString("explicit-third-body-duplicates", "warn"));
606
607 if (phaseNode.hasKey("rate-multipliers")) {
608 const auto& defaultMultipliers = phaseNode["rate-multipliers"];
609 for (auto& [key, val] : defaultMultipliers) {
610 if (key == "default") {
611 m_defaultPerturb[-1] = val.asDouble();
612 } else {
613 m_defaultPerturb[stoi(key)] = val.asDouble();
614 }
615 }
616 }
617}
618
620{
621 AnyMap out;
622 string name = KineticsFactory::factory()->canonicalize(kineticsType());
623 if (name != "none") {
624 out["kinetics"] = name;
625 if (nReactions() == 0) {
626 out["reactions"] = "none";
627 }
629 out["skip-undeclared-third-bodies"] = true;
630 }
631 if (m_explicit_third_body_duplicates == "error") {
632 // "warn" is the default, and does not need to be added. "mark-duplicate"
633 // and "modify-efficiency" do not need to be propagated here as their
634 // effects are already applied to the corresponding reactions.
635 out["explicit-third-body-duplicates"] = "error";
636 }
637 map<double, int> multipliers;
638 for (auto m : m_perturb) {
639 multipliers[m] += 1;
640 }
641 if (multipliers[1.0] != nReactions()) {
642 int defaultCount = 0;
643 double defaultMultiplier = 1.0;
644 for (auto& [m, count] : multipliers) {
645 if (count > defaultCount) {
646 defaultCount = count;
647 defaultMultiplier = m;
648 }
649 }
650 AnyMap multiplierMap;
651 multiplierMap["default"] = defaultMultiplier;
652 for (size_t i = 0; i < nReactions(); i++) {
653 if (m_perturb[i] != defaultMultiplier) {
654 multiplierMap[to_string(i)] = m_perturb[i];
655 }
656 }
657 out["rate-multipliers"] = multiplierMap;
658 }
659 }
660 return out;
661}
662
664{
665 m_kk = 0;
666 m_start.resize(nPhases());
667
668 for (size_t i = 0; i < m_thermo.size(); i++) {
669 m_start[i] = m_kk; // global index of first species of phase i
670 m_kk += m_thermo[i]->nSpecies();
671 }
672 invalidateCache();
673}
674
675bool Kinetics::addReaction(shared_ptr<Reaction> r, bool resize)
676{
677 r->check();
678 r->validate(*this);
679
680 if (m_kk == 0) {
681 init();
682 }
684
685 // Check validity of reaction within the context of the Kinetics object
686 if (!r->checkSpecies(*this)) {
687 // Do not add reaction
688 return false;
689 }
690
691 // For reactions created outside the context of a Kinetics object, the units
692 // of the rate coefficient can't be determined in advance. Do that here.
693 if (r->rate_units.factor() == 0) {
694 r->rate()->setRateUnits(r->calculateRateCoeffUnits(*this));
695 }
696
697 size_t irxn = nReactions(); // index of the new reaction
698
699 // indices of reactant and product species within this Kinetics object
700 vector<size_t> rk, pk;
701
702 // Reactant and product stoichiometric coefficients, such that rstoich[i] is
703 // the coefficient for species rk[i]
704 vector<double> rstoich, pstoich;
705
706 for (const auto& [name, stoich] : r->reactants) {
707 rk.push_back(kineticsSpeciesIndex(name));
708 rstoich.push_back(stoich);
709 }
710
711 for (const auto& [name, stoich] : r->products) {
712 pk.push_back(kineticsSpeciesIndex(name));
713 pstoich.push_back(stoich);
714 }
715
716 // The default order for each reactant is its stoichiometric coefficient,
717 // which can be overridden by entries in the Reaction.orders map. rorder[i]
718 // is the order for species rk[i].
719 vector<double> rorder = rstoich;
720 for (const auto& [name, order] : r->orders) {
721 size_t k = kineticsSpeciesIndex(name);
722 // Find the index of species k within rk
723 auto rloc = std::find(rk.begin(), rk.end(), k);
724 if (rloc != rk.end()) {
725 rorder[rloc - rk.begin()] = order;
726 } else {
727 // If the reaction order involves a non-reactant species, add an
728 // extra term to the reactants with zero stoichiometry so that the
729 // stoichiometry manager can be used to compute the global forward
730 // reaction rate.
731 rk.push_back(k);
732 rstoich.push_back(0.0);
733 rorder.push_back(order);
734 }
735 }
736
737 m_reactantStoich.add(irxn, rk, rorder, rstoich);
738 // product orders = product stoichiometric coefficients
739 m_productStoich.add(irxn, pk, pstoich, pstoich);
740 if (r->reversible) {
741 m_revindex.push_back(irxn);
742 m_revProductStoich.add(irxn, pk, pstoich, pstoich);
743 } else {
744 m_irrev.push_back(irxn);
745 }
746
747 m_reactions.push_back(r);
748 m_rfn.push_back(0.0);
749 m_delta_gibbs0.push_back(0.0);
750 m_rkcn.push_back(0.0);
751 m_ropf.push_back(0.0);
752 m_ropr.push_back(0.0);
753 m_ropnet.push_back(0.0);
755 m_dH.push_back(0.0);
756
757 if (resize) {
759 }
760
761 for (const auto& [id, callback] : m_reactionAddedCallbacks) {
762 callback();
763 }
764
765 return true;
766}
767
768void Kinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
769{
771 shared_ptr<Reaction>& rOld = m_reactions[i];
772
773 if (rNew->rate()->type() == "electron-collision-plasma") {
774 throw CanteraError("Kinetics::modifyReaction",
775 "Type electron-collision-plasma is not supported. "
776 "Use the rate object of the reaction to modify the data.");
777 }
778
779 if (rNew->type() != rOld->type()) {
780 throw CanteraError("Kinetics::modifyReaction",
781 "Reaction types are different: {} != {}.",
782 rOld->type(), rNew->type());
783 }
784
785 if (rNew->rate()->type() != rOld->rate()->type()) {
786 throw CanteraError("Kinetics::modifyReaction",
787 "ReactionRate types are different: {} != {}.",
788 rOld->rate()->type(), rNew->rate()->type());
789 }
790
791 if (rNew->reactants != rOld->reactants) {
792 throw CanteraError("Kinetics::modifyReaction",
793 "Reactants are different: '{}' != '{}'.",
794 rOld->reactantString(), rNew->reactantString());
795 }
796
797 if (rNew->products != rOld->products) {
798 throw CanteraError("Kinetics::modifyReaction",
799 "Products are different: '{}' != '{}'.",
800 rOld->productString(), rNew->productString());
801 }
802 m_reactions[i] = rNew;
803 invalidateCache();
804}
805
806shared_ptr<Reaction> Kinetics::reaction(size_t i)
807{
809 return m_reactions[i];
810}
811
812shared_ptr<const Reaction> Kinetics::reaction(size_t i) const
813{
815 return m_reactions[i];
816}
817
818}
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
Definition AnyMap.cpp:1477
void applyUnits()
Use the supplied UnitSystem to set the default units, and recursively process overrides from nodes na...
Definition AnyMap.cpp:1761
bool getBool(const string &key, bool default_) const
If key exists, return it as a bool, otherwise return default_.
Definition AnyMap.cpp:1575
const string & getString(const string &key, const string &default_) const
If key exists, return it as a string, otherwise return default_.
Definition AnyMap.cpp:1590
Array size error.
Base class for exceptions thrown by Cantera classes.
const char * what() const override
Get a description of the error.
An array index is out of range.
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition AnyMap.h:749
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:49
void checkPhaseIndex(size_t m) const
Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
Definition Kinetics.cpp:74
void checkSpeciesArraySize(size_t mm) const
Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies().
Definition Kinetics.cpp:100
virtual void getFwdRatesOfProgress(double *fwdROP)
Return the forward rates of progress of the reactions.
Definition Kinetics.cpp:367
virtual void setParameters(const AnyMap &phaseNode)
Set kinetics-related parameters from an AnyMap phase description.
Definition Kinetics.cpp:602
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:248
vector< shared_ptr< Reaction > > m_reactions
Vector of Reaction objects represented by this Kinetics manager.
Definition Kinetics.h:1502
void checkSpeciesIndex(size_t k) const
Check that the specified species index is in range Throws an exception if k is greater than nSpecies(...
Definition Kinetics.cpp:93
double checkDuplicateStoich(map< int, double > &r1, map< int, double > &r2) const
Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reac...
Definition Kinetics.cpp:260
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1495
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1543
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition Kinetics.cpp:806
virtual void getRevReactionDelta(const double *g, double *dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition Kinetics.cpp:394
void setExplicitThirdBodyDuplicateHandling(const string &flag)
Specify how to handle duplicate third body reactions where one reaction has an explicit third body an...
Definition Kinetics.cpp:107
virtual void getNetRatesOfProgress(double *netROP)
Net rates of progress.
Definition Kinetics.cpp:379
vector< size_t > m_start
m_start is a vector of integers specifying the beginning position for the species vector for the n'th...
Definition Kinetics.h:1520
virtual string kineticsType() const
Identifies the Kinetics manager type.
Definition Kinetics.h:153
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1540
shared_ptr< ThermoPhase > reactionPhase() const
Return pointer to phase where the reactions occur.
Definition Kinetics.cpp:88
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1491
virtual pair< size_t, size_t > checkDuplicates(bool throw_err=true, bool fix=false)
Check for unmarked duplicate reactions and unmatched marked duplicates.
Definition Kinetics.cpp:119
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:385
void checkPhaseArraySize(size_t mm) const
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases().
Definition Kinetics.cpp:81
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1555
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1546
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition Kinetics.h:193
vector< shared_ptr< ThermoPhase > > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
Definition Kinetics.h:1514
map< void *, function< void()> > m_reactionAddedCallbacks
Callback functions that are invoked when the reaction is added.
Definition Kinetics.h:1576
AnyMap parameters() const
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition Kinetics.cpp:619
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:675
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1551
string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition Kinetics.cpp:303
virtual void getDestructionRates(double *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:417
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition Kinetics.h:1286
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1419
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:768
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1549
void skipUndeclaredThirdBodies(bool skip)
Determine behavior when adding a new reaction that contains third-body efficiencies for species not d...
Definition Kinetics.h:1348
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition Kinetics.cpp:362
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:585
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1558
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition Kinetics.h:1486
map< string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition Kinetics.h:1528
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition Kinetics.h:1531
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition Kinetics.h:1480
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1483
shared_ptr< Kinetics > clone(const vector< shared_ptr< ThermoPhase > > &phases) const
Create a new Kinetics object with the same kinetics model and reactions as this one.
Definition Kinetics.cpp:27
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition Kinetics.h:161
map< size_t, double > m_defaultPerturb
Default values for perturbations defined in the phase definition's rate-multipliers field.
Definition Kinetics.h:1499
virtual void getRevRatesOfProgress(double *revROP)
Return the Reverse rates of progress of the reactions.
Definition Kinetics.cpp:373
vector< size_t > m_irrev
Indices of irreversible reactions.
Definition Kinetics.h:1552
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1567
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:282
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1534
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1477
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:663
void checkReactionArraySize(size_t ii) const
Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions()...
Definition Kinetics.cpp:66
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:260
virtual double reactantStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a reactant in reaction i.
Definition Kinetics.cpp:357
void checkReactionIndex(size_t m) const
Check that the specified reaction index is in range Throws an exception if i is greater than nReactio...
Definition Kinetics.cpp:42
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:428
ThermoPhase & speciesPhase(const string &nm)
This function looks up the name of a species and returns a reference to the ThermoPhase object of the...
Definition Kinetics.cpp:325
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1537
virtual void getCreationRates(double *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:403
size_t speciesPhaseIndex(size_t k) const
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of ...
Definition Kinetics.cpp:346
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:571
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:142
size_t speciesIndex(const string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:129
Abstract base class which stores data about a reaction and its rate parameterization so that it can b...
Definition Reaction.h:25
shared_ptr< ReactionRate > rate()
Get reaction rate pointer.
Definition Reaction.h:147
bool usesThirdBody() const
Check whether reaction involves third body collider.
Definition Reaction.h:161
shared_ptr< ThirdBody > thirdBody()
Get pointer to third-body handler.
Definition Reaction.h:155
bool reversible
True if the current reaction is reversible. False otherwise.
Definition Reaction.h:126
string type() const
The type of reaction, including reaction rate information.
Definition Reaction.cpp:557
string equation() const
The chemical equation for this reaction.
Definition Reaction.cpp:380
Composition products
Product species and stoichiometric coefficients.
Definition Reaction.h:114
Composition reactants
Reactant species and stoichiometric coefficients.
Definition Reaction.h:111
AnyMap input
Input data used for specific models.
Definition Reaction.h:139
bool duplicate
True if the current reaction is marked as duplicate.
Definition Reaction.h:129
void resizeCoeffs(size_t nSpc, size_t nRxn)
Resize the sparse coefficient matrix.
void add(size_t rxn, const vector< size_t > &k)
Add a single reaction to the list of reactions that this stoichiometric manager object handles.
const Eigen::SparseMatrix< double > & stoichCoeffs() const
Return matrix containing stoichiometric coefficients.
Base class for a phase with thermodynamic properties.
A class for managing third-body efficiencies, including default values.
Definition Reaction.h:207
double efficiency(const string &k) const
Get the third-body efficiency for species k
Definition Reaction.cpp:845
Composition efficiencies
Map of species to third body efficiency.
Definition Reaction.h:250
string name() const
Name of the third body collider.
Definition Reaction.h:215
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:818
Eigen::SparseMatrix< double > creationRates_ddCi()
Calculate derivatives for species creation rates with respect to species concentration at constant te...
Definition Kinetics.cpp:485
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddCi()
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
Definition Kinetics.h:964
void getCreationRates_ddT(double *dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:439
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX()
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
Definition Kinetics.h:944
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition Kinetics.cpp:531
void getCreationRates_ddC(double *dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition Kinetics.cpp:463
void getDestructionRates_ddP(double *dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:507
void getNetProductionRates_ddC(double *dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition Kinetics.cpp:567
void getDestructionRates_ddT(double *dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:495
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:767
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition Kinetics.cpp:575
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition Kinetics.cpp:475
Eigen::SparseMatrix< double > destructionRates_ddCi()
Calculate derivatives for species destruction rates with respect to species concentration at constant...
Definition Kinetics.cpp:541
void getNetProductionRates_ddT(double *dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition Kinetics.cpp:551
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:798
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:829
void getCreationRates_ddP(double *dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:451
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:902
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:840
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:580
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:854
void getNetProductionRates_ddP(double *dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition Kinetics.cpp:559
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:913
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:891
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:756
void getDestructionRates_ddC(double *dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition Kinetics.cpp:519
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition Kinetics.h:927
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:871
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:781
shared_ptr< Kinetics > newKinetics(const string &model)
Create a new Kinetics instance.
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition global.h:268
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
const U & getValue(const map< T, U > &m, const T &key, const U &default_val)
Const accessor for a value in a map.
Definition utilities.h:190
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...