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