Cantera  3.3.0a1
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
42size_t Kinetics::checkReactionIndex(size_t i) const
43{
44 if (i < nReactions()) {
45 return i;
46 }
47 throw IndexError("Kinetics::checkReactionIndex", "reactions", i, nReactions());
48}
49
51{
52 size_t nRxn = nReactions();
53
54 // Stoichiometry managers
58
59 m_rbuf.resize(nRxn);
60
61 // products are created for positive net rate of progress
63 // reactants are destroyed for positive net rate of progress
65}
66
67size_t Kinetics::checkPhaseIndex(size_t m) const
68{
69 if (m < nPhases()) {
70 return m;
71 }
72 throw IndexError("Kinetics::checkPhaseIndex", "phase", m, nPhases());
73}
74
75size_t Kinetics::phaseIndex(const string& ph, bool raise) const
76{
77 if (m_phaseindex.find(ph) == m_phaseindex.end()) {
78 if (raise) {
79 throw CanteraError("Kinetics::phaseIndex", "Phase '{}' not found", ph);
80 }
81 return npos;
82 } else {
83 return m_phaseindex.at(ph) - 1;
84 }
85}
86
87shared_ptr<ThermoPhase> Kinetics::reactionPhase() const
88{
89 return m_thermo[0];
90}
91
92size_t Kinetics::checkSpeciesIndex(size_t k) const
93{
94 if (k < m_kk) {
95 return k;
96 }
97 throw IndexError("Kinetics::checkSpeciesIndex", "species", k, m_kk);
98}
99
101{
102 if (flag == "warn" || flag == "error" || flag == "mark-duplicate"
103 || flag == "modify-efficiency")
104 {
105 m_explicit_third_body_duplicates = flag;
106 } else {
107 throw CanteraError("Kinetics::setExplicitThirdBodyDuplicateHandling",
108 "Invalid flag '{}'", flag);
109 }
110}
111
112pair<size_t, size_t> Kinetics::checkDuplicates(bool throw_err, bool fix)
113{
114 //! Map of (key indicating participating species) to reaction numbers
115 map<size_t, vector<size_t>> participants;
116 vector<map<int, double>> net_stoich;
117 std::unordered_set<size_t> unmatched_duplicates;
118 for (size_t i = 0; i < m_reactions.size(); i++) {
119 if (m_reactions[i]->duplicate) {
120 unmatched_duplicates.insert(i);
121 }
122 }
123
124 vector<InputFileError> errs;
125 for (size_t i = 0; i < m_reactions.size(); i++) {
126 // Get data about this reaction
127 unsigned long int key = 0;
128 Reaction& R = *m_reactions[i];
129 net_stoich.emplace_back();
130 map<int, double>& net = net_stoich.back();
131 for (const auto& [name, stoich] : R.reactants) {
132 int k = static_cast<int>(kineticsSpeciesIndex(name));
133 key += k*(k+1);
134 net[-1 -k] -= stoich;
135 }
136 for (const auto& [name, stoich] : R.products) {
137 int k = static_cast<int>(kineticsSpeciesIndex(name));
138 key += k*(k+1);
139 net[1+k] += stoich;
140 }
141
142 // Compare this reaction to others with similar participants
143 vector<size_t>& related = participants[key];
144 for (size_t m : related) {
145 Reaction& other = *m_reactions[m];
146 if (R.duplicate && other.duplicate) {
147 // marked duplicates
148 unmatched_duplicates.erase(i);
149 unmatched_duplicates.erase(m);
150 continue;
151 } else if (R.type() != other.type()) {
152 continue; // different reaction types
153 } else if (R.rate() && other.rate()
154 && R.rate()->type() != other.rate()->type())
155 {
156 continue; // different rate parameterizations
157 }
158 double c = checkDuplicateStoich(net_stoich[i], net_stoich[m]);
159 if (c == 0) {
160 continue; // stoichiometries differ (not by a multiple)
161 } else if (c < 0.0 && !R.reversible && !other.reversible) {
162 continue; // irreversible reactions in opposite directions
163 } else if (R.usesThirdBody() && other.usesThirdBody()) {
164 ThirdBody& tb1 = *(R.thirdBody());
165 ThirdBody& tb2 = *(other.thirdBody());
166 bool thirdBodyOk = true;
167 for (size_t k = 0; k < nTotalSpecies(); k++) {
168 string s = kineticsSpeciesName(k);
169 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
170 // non-zero third body efficiencies for species `s` in
171 // both reactions
172 thirdBodyOk = false;
173 break;
174 }
175 }
176 if (thirdBodyOk) {
177 continue; // No overlap in third body efficiencies
178 } else if ((tb1.name() == "M") != (tb2.name() == "M")) {
179 // Exactly one of the reactions uses M as the third body
180 if (m_explicit_third_body_duplicates == "mark-duplicate") {
181 R.duplicate = true;
182 other.duplicate = true;
183 continue;
184 } else if (m_explicit_third_body_duplicates == "modify-efficiency") {
185 if (tb1.name() == "M") {
186 tb1.efficiencies[tb2.name()] = 0.0;
187 } else {
188 tb2.efficiencies[tb1.name()] = 0.0;
189 }
190 continue;
191 } else if (m_explicit_third_body_duplicates == "warn") {
192 InputFileError msg("Kinetics::checkDuplicates",
193 R.input, other.input,
194 "Undeclared duplicate third body reactions with a common "
195 "third body detected.\nAdd the field "
196 "'explicit-third-body-duplicates: mark-duplicate' or\n"
197 "'explicit-third-body-duplicates: modify-efficiency' to "
198 "the YAML phase entry\nto choose how these reactions "
199 "should be handled and suppress this warning.\n"
200 "Reaction {}: {}\nReaction {}: {}\n",
201 m+1, R.equation(), i+1, other.equation());
202 if (!fix) {
203 warn_user("Kinetics::checkDuplicates", msg.what());
204 }
205 continue;
206 } // else m_explicit_third_body_duplicates == "error"
207 }
208 }
209 if (throw_err) {
210 errs.emplace_back("Kinetics::checkDuplicates",
211 R.input, other.input,
212 "Undeclared duplicate reactions detected:\n"
213 "Reaction {}: {}\nReaction {}: {}\n",
214 m+1, R.equation(), i+1, other.equation());
215 } else if (fix) {
216 R.duplicate = true;
217 other.duplicate = true;
218 unmatched_duplicates.erase(i);
219 unmatched_duplicates.erase(m);
220 } else {
221 return {i,m};
222 }
223 }
224 participants[key].push_back(i);
225 }
226 if (unmatched_duplicates.size()) {
227 for (auto i : unmatched_duplicates) {
228 if (throw_err) {
229 errs.emplace_back("Kinetics::checkDuplicates",
230 m_reactions[i]->input,
231 "No duplicate found for declared duplicate reaction number {}"
232 " ({})", i, m_reactions[i]->equation());
233 } else if (fix) {
234 m_reactions[i]->duplicate = false;
235 } else {
236 return {i, i};
237 }
238 }
239 }
240 if (errs.empty()) {
241 return {npos, npos};
242 } else if (errs.size() == 1) {
243 throw errs[0];
244 } else {
245 fmt::memory_buffer msg;
246 for (const auto& err : errs) {
247 fmt_append(msg, "\n{}\n", err.getMessage());
248 }
249 throw CanteraError("Kinetics::checkDuplicates", to_string(msg));
250 }
251}
252
253double Kinetics::checkDuplicateStoich(map<int, double>& r1, map<int, double>& r2) const
254{
255 std::unordered_set<int> keys; // species keys (k+1 or -k-1)
256 for (auto& [speciesKey, stoich] : r1) {
257 keys.insert(speciesKey);
258 }
259 for (auto& [speciesKey, stoich] : r2) {
260 keys.insert(speciesKey);
261 }
262 int k1 = r1.begin()->first;
263 // check for duplicate written in the same direction
264 double ratio = 0.0;
265 if (r1[k1] && r2[k1]) {
266 ratio = r2[k1]/r1[k1];
267 bool different = false;
268 for (int k : keys) {
269 if ((r1[k] && !r2[k]) ||
270 (!r1[k] && r2[k]) ||
271 (r1[k] && fabs(r2[k]/r1[k] - ratio) > 1.e-8)) {
272 different = true;
273 break;
274 }
275 }
276 if (!different) {
277 return ratio;
278 }
279 }
280
281 // check for duplicate written in the reverse direction
282 if (r1[k1] == 0.0 || r2[-k1] == 0.0) {
283 return 0.0;
284 }
285 ratio = r2[-k1]/r1[k1];
286 for (int k : keys) {
287 if ((r1[k] && !r2[-k]) ||
288 (!r1[k] && r2[-k]) ||
289 (r1[k] && fabs(r2[-k]/r1[k] - ratio) > 1.e-8)) {
290 return 0.0;
291 }
292 }
293 return ratio;
294}
295
296string Kinetics::kineticsSpeciesName(size_t k) const
297{
298 for (size_t n = m_start.size()-1; n != npos; n--) {
299 if (k >= m_start[n]) {
300 return thermo(n).speciesName(k - m_start[n]);
301 }
302 }
303 throw IndexError("Kinetics::kineticsSpeciesName", "totalSpecies", k, m_kk);
304}
305
306size_t Kinetics::kineticsSpeciesIndex(const string& nm, bool raise) const
307{
308 for (size_t n = 0; n < m_thermo.size(); n++) {
309 // Check the ThermoPhase object for a match
310 size_t k = thermo(n).speciesIndex(nm, false);
311 if (k != npos) {
312 return k + m_start[n];
313 }
314 }
315 if (raise) {
316 throw CanteraError("Kinetics::kineticsSpeciesIndex",
317 "Species '{}' not found", nm);
318 }
319 return npos;
320}
321
323{
324 for (size_t n = 0; n < m_thermo.size(); n++) {
325 size_t k = thermo(n).speciesIndex(nm, false);
326 if (k != npos) {
327 return thermo(n);
328 }
329 }
330 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
331}
332
333const ThermoPhase& Kinetics::speciesPhase(const string& nm) const
334{
335 for (const auto& thermo : m_thermo) {
336 if (thermo->speciesIndex(nm, false) != npos) {
337 return *thermo;
338 }
339 }
340 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
341}
342
343size_t Kinetics::speciesPhaseIndex(size_t k) const
344{
345 for (size_t n = m_start.size()-1; n != npos; n--) {
346 if (k >= m_start[n]) {
347 return n;
348 }
349 }
350 throw CanteraError("Kinetics::speciesPhaseIndex",
351 "illegal species index: {}", k);
352}
353
354double Kinetics::reactantStoichCoeff(size_t kSpec, size_t irxn) const
355{
356 return m_reactantStoich.stoichCoeffs().coeff(kSpec, irxn);
357}
358
359double Kinetics::productStoichCoeff(size_t kSpec, size_t irxn) const
360{
361 return m_productStoich.stoichCoeffs().coeff(kSpec, irxn);
362}
363
365{
366 updateROP();
367 std::copy(m_ropf.begin(), m_ropf.end(), fwdROP);
368}
369
371{
372 updateROP();
373 std::copy(m_ropr.begin(), m_ropr.end(), revROP);
374}
375
377{
378 updateROP();
379 std::copy(m_ropnet.begin(), m_ropnet.end(), netROP);
380}
381
382void Kinetics::getReactionDelta(const double* prop, double* deltaProp) const
383{
384 fill(deltaProp, deltaProp + nReactions(), 0.0);
385 // products add
386 m_productStoich.incrementReactions(prop, deltaProp);
387 // reactants subtract
388 m_reactantStoich.decrementReactions(prop, deltaProp);
389}
390
391void Kinetics::getRevReactionDelta(const double* prop, double* deltaProp) const
392{
393 fill(deltaProp, deltaProp + nReactions(), 0.0);
394 // products add
395 m_revProductStoich.incrementReactions(prop, deltaProp);
396 // reactants subtract
397 m_reactantStoich.decrementReactions(prop, deltaProp);
398}
399
401{
402 updateROP();
403
404 // zero out the output array
405 fill(cdot, cdot + m_kk, 0.0);
406
407 // the forward direction creates product species
408 m_productStoich.incrementSpecies(m_ropf.data(), cdot);
409
410 // the reverse direction creates reactant species
411 m_reactantStoich.incrementSpecies(m_ropr.data(), cdot);
412}
413
415{
416 updateROP();
417
418 fill(ddot, ddot + m_kk, 0.0);
419 // the reverse direction destroys products in reversible reactions
420 m_revProductStoich.incrementSpecies(m_ropr.data(), ddot);
421 // the forward direction destroys reactants
422 m_reactantStoich.incrementSpecies(m_ropf.data(), ddot);
423}
424
426{
427 updateROP();
428
429 fill(net, net + m_kk, 0.0);
430 // products are created for positive net rate of progress
431 m_productStoich.incrementSpecies(m_ropnet.data(), net);
432 // reactants are destroyed for positive net rate of progress
433 m_reactantStoich.decrementSpecies(m_ropnet.data(), net);
434}
435
437{
438 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
439 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
440 // the forward direction creates product species
441 getFwdRatesOfProgress_ddT(buf.data());
442 out = m_productStoich.stoichCoeffs() * buf;
443 // the reverse direction creates reactant species
444 getRevRatesOfProgress_ddT(buf.data());
445 out += m_reactantStoich.stoichCoeffs() * buf;
446}
447
449{
450 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
451 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
452 // the forward direction creates product species
453 getFwdRatesOfProgress_ddP(buf.data());
454 out = m_productStoich.stoichCoeffs() * buf;
455 // the reverse direction creates reactant species
456 getRevRatesOfProgress_ddP(buf.data());
457 out += m_reactantStoich.stoichCoeffs() * buf;
458}
459
461{
462 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
463 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
464 // the forward direction creates product species
465 getFwdRatesOfProgress_ddC(buf.data());
466 out = m_productStoich.stoichCoeffs() * buf;
467 // the reverse direction creates reactant species
468 getRevRatesOfProgress_ddC(buf.data());
469 out += m_reactantStoich.stoichCoeffs() * buf;
470}
471
472Eigen::SparseMatrix<double> Kinetics::creationRates_ddX()
473{
474 Eigen::SparseMatrix<double> jac;
475 // the forward direction creates product species
477 // the reverse direction creates reactant species
479 return jac;
480}
481
482Eigen::SparseMatrix<double> Kinetics::creationRates_ddCi()
483{
484 Eigen::SparseMatrix<double> jac;
485 // the forward direction creates product species
487 // the reverse direction creates reactant species
489 return jac;
490}
491
493{
494 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
495 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
496 // the reverse direction destroys products in reversible reactions
497 getRevRatesOfProgress_ddT(buf.data());
498 out = m_revProductStoich.stoichCoeffs() * buf;
499 // the forward direction destroys reactants
500 getFwdRatesOfProgress_ddT(buf.data());
501 out += m_reactantStoich.stoichCoeffs() * buf;
502}
503
505{
506 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
507 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
508 // the reverse direction destroys products in reversible reactions
509 getRevRatesOfProgress_ddP(buf.data());
510 out = m_revProductStoich.stoichCoeffs() * buf;
511 // the forward direction destroys reactants
512 getFwdRatesOfProgress_ddP(buf.data());
513 out += m_reactantStoich.stoichCoeffs() * buf;
514}
515
517{
518 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
519 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
520 // the reverse direction destroys products in reversible reactions
521 getRevRatesOfProgress_ddC(buf.data());
522 out = m_revProductStoich.stoichCoeffs() * buf;
523 // the forward direction destroys reactants
524 getFwdRatesOfProgress_ddC(buf.data());
525 out += m_reactantStoich.stoichCoeffs() * buf;
526}
527
528Eigen::SparseMatrix<double> Kinetics::destructionRates_ddX()
529{
530 Eigen::SparseMatrix<double> jac;
531 // the reverse direction destroys products in reversible reactions
533 // the forward direction destroys reactants
535 return jac;
536}
537
538Eigen::SparseMatrix<double> Kinetics::destructionRates_ddCi()
539{
540 Eigen::SparseMatrix<double> jac;
541 // the reverse direction destroys products in reversible reactions
543 // the forward direction destroys reactants
545 return jac;
546}
547
549{
550 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
551 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
552 getNetRatesOfProgress_ddT(buf.data());
553 out = m_stoichMatrix * buf;
554}
555
557{
558 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
559 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
560 getNetRatesOfProgress_ddP(buf.data());
561 out = m_stoichMatrix * buf;
562}
563
565{
566 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
567 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
568 getNetRatesOfProgress_ddC(buf.data());
569 out = m_stoichMatrix * buf;
570}
571
572Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddX()
573{
575}
576
577Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddCi()
578{
580}
581
582void Kinetics::addThermo(shared_ptr<ThermoPhase> thermo)
583{
584 // the phase with lowest dimensionality is assumed to be the
585 // phase/interface at which reactions take place
586 if (thermo->nDim() <= m_mindim) {
587 if (!m_thermo.empty()) {
588 throw CanteraError("Kinetics::addThermo",
589 "The reacting (lowest dimensional) phase must be added first.");
590 }
591 m_mindim = thermo->nDim();
592 }
593
594 m_thermo.push_back(thermo);
595 m_phaseindex[m_thermo.back()->name()] = nPhases();
597}
598
599void Kinetics::setParameters(const AnyMap& phaseNode) {
600 skipUndeclaredThirdBodies(phaseNode.getBool("skip-undeclared-third-bodies", false));
602 phaseNode.getString("explicit-third-body-duplicates", "warn"));
603
604 if (phaseNode.hasKey("rate-multipliers")) {
605 const auto& defaultMultipliers = phaseNode["rate-multipliers"];
606 for (auto& [key, val] : defaultMultipliers) {
607 if (key == "default") {
608 m_defaultPerturb[-1] = val.asDouble();
609 } else {
610 m_defaultPerturb[stoi(key)] = val.asDouble();
611 }
612 }
613 }
614}
615
617{
618 AnyMap out;
619 string name = KineticsFactory::factory()->canonicalize(kineticsType());
620 if (name != "none") {
621 out["kinetics"] = name;
622 if (nReactions() == 0) {
623 out["reactions"] = "none";
624 }
626 out["skip-undeclared-third-bodies"] = true;
627 }
628 if (m_explicit_third_body_duplicates == "error") {
629 // "warn" is the default, and does not need to be added. "mark-duplicate"
630 // and "modify-efficiency" do not need to be propagated here as their
631 // effects are already applied to the corresponding reactions.
632 out["explicit-third-body-duplicates"] = "error";
633 }
634 map<double, int> multipliers;
635 for (auto m : m_perturb) {
636 multipliers[m] += 1;
637 }
638 if (static_cast<size_t>(multipliers[1.0]) != (nReactions())) {
639 int defaultCount = 0;
640 double defaultMultiplier = 1.0;
641 for (auto& [m, count] : multipliers) {
642 if (count > defaultCount) {
643 defaultCount = count;
644 defaultMultiplier = m;
645 }
646 }
647 AnyMap multiplierMap;
648 multiplierMap["default"] = defaultMultiplier;
649 for (size_t i = 0; i < nReactions(); i++) {
650 if (m_perturb[i] != defaultMultiplier) {
651 multiplierMap[to_string(i)] = m_perturb[i];
652 }
653 }
654 out["rate-multipliers"] = multiplierMap;
655 }
656 }
657 return out;
658}
659
661{
662 m_kk = 0;
663 m_start.resize(nPhases());
664
665 for (size_t i = 0; i < m_thermo.size(); i++) {
666 m_start[i] = m_kk; // global index of first species of phase i
667 m_kk += m_thermo[i]->nSpecies();
668 }
669 invalidateCache();
670}
671
672bool Kinetics::addReaction(shared_ptr<Reaction> r, bool resize)
673{
674 r->check();
675 r->validate(*this);
676
677 if (m_kk == 0) {
678 init();
679 }
681
682 // Check validity of reaction within the context of the Kinetics object
683 if (!r->checkSpecies(*this)) {
684 // Do not add reaction
685 return false;
686 }
687
688 // For reactions created outside the context of a Kinetics object, the units
689 // of the rate coefficient can't be determined in advance. Do that here.
690 if (r->rate_units.factor() == 0) {
691 r->rate()->setRateUnits(r->calculateRateCoeffUnits(*this));
692 }
693
694 size_t irxn = nReactions(); // index of the new reaction
695
696 // indices of reactant and product species within this Kinetics object
697 vector<size_t> rk, pk;
698
699 // Reactant and product stoichiometric coefficients, such that rstoich[i] is
700 // the coefficient for species rk[i]
701 vector<double> rstoich, pstoich;
702
703 for (const auto& [name, stoich] : r->reactants) {
704 rk.push_back(kineticsSpeciesIndex(name));
705 rstoich.push_back(stoich);
706 }
707
708 for (const auto& [name, stoich] : r->products) {
709 pk.push_back(kineticsSpeciesIndex(name));
710 pstoich.push_back(stoich);
711 }
712
713 // The default order for each reactant is its stoichiometric coefficient,
714 // which can be overridden by entries in the Reaction.orders map. rorder[i]
715 // is the order for species rk[i].
716 vector<double> rorder = rstoich;
717 for (const auto& [name, order] : r->orders) {
718 size_t k = kineticsSpeciesIndex(name);
719 // Find the index of species k within rk
720 auto rloc = std::find(rk.begin(), rk.end(), k);
721 if (rloc != rk.end()) {
722 rorder[rloc - rk.begin()] = order;
723 } else {
724 // If the reaction order involves a non-reactant species, add an
725 // extra term to the reactants with zero stoichiometry so that the
726 // stoichiometry manager can be used to compute the global forward
727 // reaction rate.
728 rk.push_back(k);
729 rstoich.push_back(0.0);
730 rorder.push_back(order);
731 }
732 }
733
734 m_reactantStoich.add(irxn, rk, rorder, rstoich);
735 // product orders = product stoichiometric coefficients
736 m_productStoich.add(irxn, pk, pstoich, pstoich);
737 if (r->reversible) {
738 m_revindex.push_back(irxn);
739 m_revProductStoich.add(irxn, pk, pstoich, pstoich);
740 } else {
741 m_irrev.push_back(irxn);
742 }
743
744 m_reactions.push_back(r);
745 m_rfn.push_back(0.0);
746 m_delta_gibbs0.push_back(0.0);
747 m_rkcn.push_back(0.0);
748 m_ropf.push_back(0.0);
749 m_ropr.push_back(0.0);
750 m_ropnet.push_back(0.0);
752 m_dH.push_back(0.0);
753
754 if (resize) {
756 }
757
758 for (const auto& [id, callback] : m_reactionAddedCallbacks) {
759 callback();
760 }
761
762 return true;
763}
764
765void Kinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
766{
768 shared_ptr<Reaction>& rOld = m_reactions[i];
769
770 if (rNew->rate()->type() == "electron-collision-plasma") {
771 throw CanteraError("Kinetics::modifyReaction",
772 "Type electron-collision-plasma is not supported. "
773 "Use the rate object of the reaction to modify the data.");
774 }
775
776 if (rNew->type() != rOld->type()) {
777 throw CanteraError("Kinetics::modifyReaction",
778 "Reaction types are different: {} != {}.",
779 rOld->type(), rNew->type());
780 }
781
782 if (rNew->rate()->type() != rOld->rate()->type()) {
783 throw CanteraError("Kinetics::modifyReaction",
784 "ReactionRate types are different: {} != {}.",
785 rOld->rate()->type(), rNew->rate()->type());
786 }
787
788 if (rNew->reactants != rOld->reactants) {
789 throw CanteraError("Kinetics::modifyReaction",
790 "Reactants are different: '{}' != '{}'.",
791 rOld->reactantString(), rNew->reactantString());
792 }
793
794 if (rNew->products != rOld->products) {
795 throw CanteraError("Kinetics::modifyReaction",
796 "Products are different: '{}' != '{}'.",
797 rOld->productString(), rNew->productString());
798 }
799 m_reactions[i] = rNew;
800 invalidateCache();
801}
802
803shared_ptr<Reaction> Kinetics::reaction(size_t i)
804{
806 return m_reactions[i];
807}
808
809shared_ptr<const Reaction> Kinetics::reaction(size_t i) const
810{
812 return m_reactions[i];
813}
814
815}
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
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:50
virtual void getFwdRatesOfProgress(double *fwdROP)
Return the forward rates of progress of the reactions.
Definition Kinetics.cpp:364
virtual void setParameters(const AnyMap &phaseNode)
Set kinetics-related parameters from an AnyMap phase description.
Definition Kinetics.cpp:599
size_t checkPhaseIndex(size_t m) const
Check that the specified phase index is in range.
Definition Kinetics.cpp:67
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:239
vector< shared_ptr< Reaction > > m_reactions
Vector of Reaction objects represented by this Kinetics manager.
Definition Kinetics.h:1488
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:253
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1481
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1529
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition Kinetics.cpp:803
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:391
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:100
virtual void getNetRatesOfProgress(double *netROP)
Net rates of progress.
Definition Kinetics.cpp:376
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:1506
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:1526
shared_ptr< ThermoPhase > reactionPhase() const
Return pointer to phase where the reactions occur.
Definition Kinetics.cpp:87
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1477
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:112
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:382
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1541
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1532
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition Kinetics.h:189
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:1500
map< void *, function< void()> > m_reactionAddedCallbacks
Callback functions that are invoked when the reaction is added.
Definition Kinetics.h:1562
AnyMap parameters() const
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition Kinetics.cpp:616
size_t phaseIndex(const string &ph, bool raise=true) const
Return the index of a phase among the phases participating in this kinetic mechanism.
Definition Kinetics.cpp:75
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:672
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1537
string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition Kinetics.cpp:296
virtual void getDestructionRates(double *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:414
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition Kinetics.h:1272
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1405
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:765
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1535
void skipUndeclaredThirdBodies(bool skip)
Determine behavior when adding a new reaction that contains third-body efficiencies for species not d...
Definition Kinetics.h:1334
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition Kinetics.cpp:359
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:582
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1544
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition Kinetics.h:1472
map< string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition Kinetics.h:1514
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition Kinetics.h:1517
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition Kinetics.h:1466
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1469
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:1485
virtual void getRevRatesOfProgress(double *revROP)
Return the Reverse rates of progress of the reactions.
Definition Kinetics.cpp:370
vector< size_t > m_irrev
Indices of irreversible reactions.
Definition Kinetics.h:1538
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1553
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:273
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1520
size_t checkSpeciesIndex(size_t k) const
Check that the specified species index is in range.
Definition Kinetics.cpp:92
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1463
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:660
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:251
virtual double reactantStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a reactant in reaction i.
Definition Kinetics.cpp:354
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:425
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:322
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1523
size_t checkReactionIndex(size_t m) const
Check that the specified reaction index is in range.
Definition Kinetics.cpp:42
virtual void getCreationRates(double *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:400
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:343
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:582
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:127
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:143
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:813
Eigen::SparseMatrix< double > creationRates_ddCi()
Calculate derivatives for species creation rates with respect to species concentration at constant te...
Definition Kinetics.cpp:482
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddCi()
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
Definition Kinetics.h:959
void getCreationRates_ddT(double *dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:436
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:939
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition Kinetics.cpp:528
void getCreationRates_ddC(double *dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition Kinetics.cpp:460
void getDestructionRates_ddP(double *dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:504
void getNetProductionRates_ddC(double *dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition Kinetics.cpp:564
void getDestructionRates_ddT(double *dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:492
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:762
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition Kinetics.cpp:572
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition Kinetics.cpp:472
Eigen::SparseMatrix< double > destructionRates_ddCi()
Calculate derivatives for species destruction rates with respect to species concentration at constant...
Definition Kinetics.cpp:538
void getNetProductionRates_ddT(double *dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition Kinetics.cpp:548
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:793
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:824
void getCreationRates_ddP(double *dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:448
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:897
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:835
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:577
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:849
void getNetProductionRates_ddP(double *dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition Kinetics.cpp:556
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:908
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:886
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:751
void getDestructionRates_ddC(double *dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition Kinetics.cpp:516
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition Kinetics.h:922
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:866
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:776
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:263
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...