Cantera  4.0.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
322size_t Kinetics::speciesOffset(const ThermoPhase& phase) const
323{
324 for (size_t n = 0; n < m_thermo.size(); n++) {
325 if (&phase == m_thermo[n].get()) {
326 return m_start[n];
327 }
328 }
329 throw CanteraError("Kinetics::speciesOffset",
330 "Phase '{}' not found in Kinetics object", phase.name());
331}
332
334{
335 for (size_t n = 0; n < m_thermo.size(); n++) {
336 size_t k = thermo(n).speciesIndex(nm, false);
337 if (k != npos) {
338 return thermo(n);
339 }
340 }
341 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
342}
343
344const ThermoPhase& Kinetics::speciesPhase(const string& nm) const
345{
346 for (const auto& thermo : m_thermo) {
347 if (thermo->speciesIndex(nm, false) != npos) {
348 return *thermo;
349 }
350 }
351 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
352}
353
354size_t Kinetics::speciesPhaseIndex(size_t k) const
355{
356 for (size_t n = m_start.size()-1; n != npos; n--) {
357 if (k >= m_start[n]) {
358 return n;
359 }
360 }
361 throw CanteraError("Kinetics::speciesPhaseIndex",
362 "illegal species index: {}", k);
363}
364
365double Kinetics::reactantStoichCoeff(size_t kSpec, size_t irxn) const
366{
367 return m_reactantStoich.stoichCoeffs().coeff(kSpec, irxn);
368}
369
370double Kinetics::productStoichCoeff(size_t kSpec, size_t irxn) const
371{
372 return m_productStoich.stoichCoeffs().coeff(kSpec, irxn);
373}
374
375void Kinetics::getFwdRatesOfProgress(span<double> fwdROP)
376{
377 checkArraySize("Kinetics::getFwdRatesOfProgress", fwdROP.size(), nReactions());
378 updateROP();
379 std::copy(m_ropf.begin(), m_ropf.end(), fwdROP.begin());
380}
381
382void Kinetics::getRevRatesOfProgress(span<double> revROP)
383{
384 checkArraySize("Kinetics::getRevRatesOfProgress", revROP.size(), nReactions());
385 updateROP();
386 std::copy(m_ropr.begin(), m_ropr.end(), revROP.begin());
387}
388
389void Kinetics::getNetRatesOfProgress(span<double> netROP)
390{
391 checkArraySize("Kinetics::getNetRatesOfProgress", netROP.size(), nReactions());
392 updateROP();
393 std::copy(m_ropnet.begin(), m_ropnet.end(), netROP.begin());
394}
395
396void Kinetics::getReactionDelta(span<const double> prop,
397 span<double> deltaProp) const
398{
399 checkArraySize("Kinetics::getReactionDelta", prop.size(), nTotalSpecies());
400 checkArraySize("Kinetics::getReactionDelta", deltaProp.size(), nReactions());
401 fill(deltaProp.begin(), deltaProp.end(), 0.0);
402 // products add
403 m_productStoich.incrementReactions(prop, deltaProp);
404 // reactants subtract
405 m_reactantStoich.decrementReactions(prop, deltaProp);
406}
407
408void Kinetics::getRevReactionDelta(span<const double> prop,
409 span<double> deltaProp) const
410{
411 checkArraySize("Kinetics::getRevReactionDelta", prop.size(), nTotalSpecies());
412 checkArraySize("Kinetics::getRevReactionDelta", deltaProp.size(), nReactions());
413 fill(deltaProp.begin(), deltaProp.end(), 0.0);
414 // products add
415 m_revProductStoich.incrementReactions(prop, deltaProp);
416 // reactants subtract
417 m_reactantStoich.decrementReactions(prop, deltaProp);
418}
419
420void Kinetics::getCreationRates(span<double> cdot)
421{
422 checkArraySize("Kinetics::getCreationRates", cdot.size(), nTotalSpecies());
423 updateROP();
424
425 // zero out the output array
426 fill(cdot.begin(), cdot.end(), 0.0);
427
428 // the forward direction creates product species
429 m_productStoich.incrementSpecies(m_ropf, cdot);
430
431 // the reverse direction creates reactant species
432 m_reactantStoich.incrementSpecies(m_ropr, cdot);
433}
434
435void Kinetics::getDestructionRates(span<double> ddot)
436{
437 checkArraySize("Kinetics::getDestructionRates", ddot.size(), nTotalSpecies());
438 updateROP();
439
440 fill(ddot.begin(), ddot.end(), 0.0);
441 // the reverse direction destroys products in reversible reactions
442 m_revProductStoich.incrementSpecies(m_ropr, ddot);
443 // the forward direction destroys reactants
444 m_reactantStoich.incrementSpecies(m_ropf, ddot);
445}
446
448{
449 checkArraySize("Kinetics::getNetProductionRates", net.size(), nTotalSpecies());
450 updateROP();
451
452 fill(net.begin(), net.end(), 0.0);
453 // products are created for positive net rate of progress
454 m_productStoich.incrementSpecies(m_ropnet, net);
455 // reactants are destroyed for positive net rate of progress
456 m_reactantStoich.decrementSpecies(m_ropnet, net);
457}
458
459void Kinetics::getCreationRates_ddT(span<double> dwdot)
460{
461 checkArraySize("Kinetics::getCreationRates_ddT", dwdot.size(), nTotalSpecies());
462 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
463 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
464 // the forward direction creates product species
466 out = m_productStoich.stoichCoeffs() * buf;
467 // the reverse direction creates reactant species
469 out += m_reactantStoich.stoichCoeffs() * buf;
470}
471
472void Kinetics::getCreationRates_ddP(span<double> dwdot)
473{
474 checkArraySize("Kinetics::getCreationRates_ddP", dwdot.size(), nTotalSpecies());
475 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
476 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
477 // the forward direction creates product species
479 out = m_productStoich.stoichCoeffs() * buf;
480 // the reverse direction creates reactant species
482 out += m_reactantStoich.stoichCoeffs() * buf;
483}
484
485void Kinetics::getCreationRates_ddC(span<double> dwdot)
486{
487 checkArraySize("Kinetics::getCreationRates_ddC", dwdot.size(), nTotalSpecies());
488 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
489 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
490 // the forward direction creates product species
492 out = m_productStoich.stoichCoeffs() * buf;
493 // the reverse direction creates reactant species
495 out += m_reactantStoich.stoichCoeffs() * buf;
496}
497
498Eigen::SparseMatrix<double> Kinetics::creationRates_ddX()
499{
500 Eigen::SparseMatrix<double> jac;
501 // the forward direction creates product species
503 // the reverse direction creates reactant species
505 return jac;
506}
507
508Eigen::SparseMatrix<double> Kinetics::creationRates_ddCi()
509{
510 Eigen::SparseMatrix<double> jac;
511 // the forward direction creates product species
513 // the reverse direction creates reactant species
515 return jac;
516}
517
518void Kinetics::getDestructionRates_ddT(span<double> dwdot)
519{
520 checkArraySize("Kinetics::getDestructionRates_ddT", dwdot.size(), nTotalSpecies());
521 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
522 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
523 // the reverse direction destroys products in reversible reactions
525 out = m_revProductStoich.stoichCoeffs() * buf;
526 // the forward direction destroys reactants
528 out += m_reactantStoich.stoichCoeffs() * buf;
529}
530
531void Kinetics::getDestructionRates_ddP(span<double> dwdot)
532{
533 checkArraySize("Kinetics::getDestructionRates_ddP", dwdot.size(), nTotalSpecies());
534 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
535 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
536 // the reverse direction destroys products in reversible reactions
538 out = m_revProductStoich.stoichCoeffs() * buf;
539 // the forward direction destroys reactants
541 out += m_reactantStoich.stoichCoeffs() * buf;
542}
543
544void Kinetics::getDestructionRates_ddC(span<double> dwdot)
545{
546 checkArraySize("Kinetics::getDestructionRates_ddC", dwdot.size(), nTotalSpecies());
547 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
548 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
549 // the reverse direction destroys products in reversible reactions
551 out = m_revProductStoich.stoichCoeffs() * buf;
552 // the forward direction destroys reactants
554 out += m_reactantStoich.stoichCoeffs() * buf;
555}
556
557Eigen::SparseMatrix<double> Kinetics::destructionRates_ddX()
558{
559 Eigen::SparseMatrix<double> jac;
560 // the reverse direction destroys products in reversible reactions
562 // the forward direction destroys reactants
564 return jac;
565}
566
567Eigen::SparseMatrix<double> Kinetics::destructionRates_ddCi()
568{
569 Eigen::SparseMatrix<double> jac;
570 // the reverse direction destroys products in reversible reactions
572 // the forward direction destroys reactants
574 return jac;
575}
576
578{
579 checkArraySize("Kinetics::getNetProductionRates_ddT", dwdot.size(), nTotalSpecies());
580 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
581 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
583 out = m_stoichMatrix * buf;
584}
585
587{
588 checkArraySize("Kinetics::getNetProductionRates_ddP", dwdot.size(), nTotalSpecies());
589 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
590 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
592 out = m_stoichMatrix * buf;
593}
594
596{
597 checkArraySize("Kinetics::getNetProductionRates_ddC", dwdot.size(), nTotalSpecies());
598 Eigen::Map<Eigen::VectorXd> out(dwdot.data(), m_kk);
599 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
601 out = m_stoichMatrix * buf;
602}
603
604Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddX()
605{
607}
608
609Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddCi()
610{
612}
613
614void Kinetics::addThermo(shared_ptr<ThermoPhase> thermo)
615{
616 // the phase with lowest dimensionality is assumed to be the
617 // phase/interface at which reactions take place
618 if (thermo->nDim() <= m_mindim) {
619 if (!m_thermo.empty()) {
620 throw CanteraError("Kinetics::addThermo",
621 "The reacting (lowest dimensional) phase must be added first.");
622 }
623 m_mindim = thermo->nDim();
624 }
625
626 m_thermo.push_back(thermo);
627 m_phaseindex[m_thermo.back()->name()] = nPhases();
629}
630
631void Kinetics::setParameters(const AnyMap& phaseNode) {
632 skipUndeclaredThirdBodies(phaseNode.getBool("skip-undeclared-third-bodies", false));
634 phaseNode.getString("explicit-third-body-duplicates", "warn"));
635
636 if (phaseNode.hasKey("rate-multipliers")) {
637 const auto& defaultMultipliers = phaseNode["rate-multipliers"];
638 for (auto& [key, val] : defaultMultipliers) {
639 if (key == "default") {
640 m_defaultPerturb[-1] = val.asDouble();
641 } else {
642 m_defaultPerturb[stoi(key)] = val.asDouble();
643 }
644 }
645 }
646}
647
649{
650 AnyMap out;
651 string name = KineticsFactory::factory()->canonicalize(kineticsType());
652 if (name != "none") {
653 out["kinetics"] = name;
654 if (nReactions() == 0) {
655 out["reactions"] = "none";
656 }
658 out["skip-undeclared-third-bodies"] = true;
659 }
660 if (m_explicit_third_body_duplicates == "error") {
661 // "warn" is the default, and does not need to be added. "mark-duplicate"
662 // and "modify-efficiency" do not need to be propagated here as their
663 // effects are already applied to the corresponding reactions.
664 out["explicit-third-body-duplicates"] = "error";
665 }
666 map<double, int> multipliers;
667 for (auto m : m_perturb) {
668 multipliers[m] += 1;
669 }
670 if (static_cast<size_t>(multipliers[1.0]) != (nReactions())) {
671 int defaultCount = 0;
672 double defaultMultiplier = 1.0;
673 for (auto& [m, count] : multipliers) {
674 if (count > defaultCount) {
675 defaultCount = count;
676 defaultMultiplier = m;
677 }
678 }
679 AnyMap multiplierMap;
680 multiplierMap["default"] = defaultMultiplier;
681 for (size_t i = 0; i < nReactions(); i++) {
682 if (m_perturb[i] != defaultMultiplier) {
683 multiplierMap[to_string(i)] = m_perturb[i];
684 }
685 }
686 out["rate-multipliers"] = multiplierMap;
687 }
688 }
689 return out;
690}
691
693{
694 m_kk = 0;
695 m_start.resize(nPhases());
696
697 for (size_t i = 0; i < m_thermo.size(); i++) {
698 m_start[i] = m_kk; // global index of first species of phase i
699 m_kk += m_thermo[i]->nSpecies();
700 }
701 invalidateCache();
702}
703
704bool Kinetics::addReaction(shared_ptr<Reaction> r, bool resize)
705{
706 r->check();
707 r->validate(*this);
708
709 if (m_kk == 0) {
710 init();
711 }
713
714 // Check validity of reaction within the context of the Kinetics object
715 if (!r->checkSpecies(*this)) {
716 // Do not add reaction
717 return false;
718 }
719
720 // For reactions created outside the context of a Kinetics object, the units
721 // of the rate coefficient can't be determined in advance. Do that here.
722 if (r->rate_units.factor() == 0) {
723 r->rate()->setRateUnits(r->calculateRateCoeffUnits(*this));
724 }
725
726 size_t irxn = nReactions(); // index of the new reaction
727
728 // indices of reactant and product species within this Kinetics object
729 vector<size_t> rk, pk;
730
731 // Reactant and product stoichiometric coefficients, such that rstoich[i] is
732 // the coefficient for species rk[i]
733 vector<double> rstoich, pstoich;
734
735 for (const auto& [name, stoich] : r->reactants) {
736 rk.push_back(kineticsSpeciesIndex(name));
737 rstoich.push_back(stoich);
738 }
739
740 for (const auto& [name, stoich] : r->products) {
741 pk.push_back(kineticsSpeciesIndex(name));
742 pstoich.push_back(stoich);
743 }
744
745 // The default order for each reactant is its stoichiometric coefficient,
746 // which can be overridden by entries in the Reaction.orders map. rorder[i]
747 // is the order for species rk[i].
748 vector<double> rorder = rstoich;
749 for (const auto& [name, order] : r->orders) {
750 size_t k = kineticsSpeciesIndex(name);
751 // Find the index of species k within rk
752 auto rloc = std::find(rk.begin(), rk.end(), k);
753 if (rloc != rk.end()) {
754 rorder[rloc - rk.begin()] = order;
755 } else {
756 // If the reaction order involves a non-reactant species, add an
757 // extra term to the reactants with zero stoichiometry so that the
758 // stoichiometry manager can be used to compute the global forward
759 // reaction rate.
760 rk.push_back(k);
761 rstoich.push_back(0.0);
762 rorder.push_back(order);
763 }
764 }
765
766 m_reactantStoich.add(irxn, rk, rorder, rstoich);
767 // product orders = product stoichiometric coefficients
768 m_productStoich.add(irxn, pk, pstoich, pstoich);
769 if (r->reversible) {
770 m_revindex.push_back(irxn);
771 m_revProductStoich.add(irxn, pk, pstoich, pstoich);
772 } else {
773 m_irrev.push_back(irxn);
774 }
775
776 m_reactions.push_back(r);
777 m_rfn.push_back(0.0);
778 m_delta_gibbs0.push_back(0.0);
779 m_rkcn.push_back(0.0);
780 m_ropf.push_back(0.0);
781 m_ropr.push_back(0.0);
782 m_ropnet.push_back(0.0);
784 m_dH.push_back(0.0);
785
786 if (resize) {
788 }
789
790 for (const auto& [id, callback] : m_reactionAddedCallbacks) {
791 callback();
792 }
793
794 return true;
795}
796
797void Kinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
798{
800 shared_ptr<Reaction>& rOld = m_reactions[i];
801
802 if (rNew->rate()->type() == "electron-collision-plasma") {
803 throw CanteraError("Kinetics::modifyReaction",
804 "Type electron-collision-plasma is not supported. "
805 "Use the rate object of the reaction to modify the data.");
806 }
807
808 if (rNew->type() != rOld->type()) {
809 throw CanteraError("Kinetics::modifyReaction",
810 "Reaction types are different: {} != {}.",
811 rOld->type(), rNew->type());
812 }
813
814 if (rNew->rate()->type() != rOld->rate()->type()) {
815 throw CanteraError("Kinetics::modifyReaction",
816 "ReactionRate types are different: {} != {}.",
817 rOld->rate()->type(), rNew->rate()->type());
818 }
819
820 if (rNew->reactants != rOld->reactants) {
821 throw CanteraError("Kinetics::modifyReaction",
822 "Reactants are different: '{}' != '{}'.",
823 rOld->reactantString(), rNew->reactantString());
824 }
825
826 if (rNew->products != rOld->products) {
827 throw CanteraError("Kinetics::modifyReaction",
828 "Products are different: '{}' != '{}'.",
829 rOld->productString(), rNew->productString());
830 }
831 m_reactions[i] = rNew;
832 invalidateCache();
833}
834
835shared_ptr<Reaction> Kinetics::reaction(size_t i)
836{
838 return m_reactions[i];
839}
840
841shared_ptr<const Reaction> Kinetics::reaction(size_t i) const
842{
844 return m_reactions[i];
845}
846
847}
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
shared_ptr< ThermoPhase > phase(size_t n=0) const
Return pointer to phase associated with Kinetics by index.
Definition Kinetics.h:226
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:50
virtual void getRevRatesOfProgress(span< double > revROP)
Return the Reverse rates of progress of the reactions.
Definition Kinetics.cpp:382
virtual void getNetRatesOfProgress(span< double > netROP)
Net rates of progress.
Definition Kinetics.cpp:389
virtual void setParameters(const AnyMap &phaseNode)
Set kinetics-related parameters from an AnyMap phase description.
Definition Kinetics.cpp:631
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:1491
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:1484
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1532
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition Kinetics.cpp:835
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
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:1509
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:1529
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:1480
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
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1544
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1535
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:1503
map< void *, function< void()> > m_reactionAddedCallbacks
Callback functions that are invoked when the reaction is added.
Definition Kinetics.h:1565
AnyMap parameters() const
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition Kinetics.cpp:648
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:704
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1540
string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition Kinetics.cpp:296
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition Kinetics.h:1275
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1408
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:797
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1538
void skipUndeclaredThirdBodies(bool skip)
Determine behavior when adding a new reaction that contains third-body efficiencies for species not d...
Definition Kinetics.h:1337
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition Kinetics.cpp:370
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:614
virtual void getRevReactionDelta(span< const double > g, span< 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:408
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1547
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition Kinetics.h:1475
map< string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition Kinetics.h:1517
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition Kinetics.h:1520
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition Kinetics.h:1469
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1472
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:1488
vector< size_t > m_irrev
Indices of irreversible reactions.
Definition Kinetics.h:1541
virtual void getNetProductionRates(span< double > wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:447
virtual void getFwdRatesOfProgress(span< double > fwdROP)
Return the forward rates of progress of the reactions.
Definition Kinetics.cpp:375
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1556
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
virtual void getReactionDelta(span< const double > property, span< double > deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:396
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1523
size_t checkSpeciesIndex(size_t k) const
Check that the specified species index is in range.
Definition Kinetics.cpp:92
virtual void getDestructionRates(span< double > ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:435
virtual void getCreationRates(span< double > cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:420
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1466
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:692
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:365
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:333
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1526
size_t checkReactionIndex(size_t m) const
Check that the specified reaction index is in range.
Definition Kinetics.cpp:42
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:354
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:569
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, span< const 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,...
void getNetProductionRates_ddT(span< double > dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition Kinetics.cpp:577
virtual void getFwdRatesOfProgress_ddC(span< double > drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:780
virtual void getRevRatesOfProgress_ddT(span< double > drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:828
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:817
void getCreationRates_ddT(span< double > dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:459
Eigen::SparseMatrix< double > creationRates_ddCi()
Calculate derivatives for species creation rates with respect to species concentration at constant te...
Definition Kinetics.cpp:508
void getDestructionRates_ddT(span< double > dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:518
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddCi()
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
Definition Kinetics.h:963
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:943
void getCreationRates_ddP(span< double > dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:472
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition Kinetics.cpp:557
void getNetProductionRates_ddC(span< double > dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition Kinetics.cpp:595
void getDestructionRates_ddC(span< double > dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition Kinetics.cpp:544
virtual void getRevRatesOfProgress_ddC(span< double > drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:853
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition Kinetics.cpp:604
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition Kinetics.cpp:498
Eigen::SparseMatrix< double > destructionRates_ddCi()
Calculate derivatives for species destruction rates with respect to species concentration at constant...
Definition Kinetics.cpp:567
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:797
virtual void getNetRatesOfProgress_ddT(span< double > drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:901
virtual void getRevRatesOfProgress_ddP(span< double > drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:839
virtual void getFwdRatesOfProgress_ddT(span< double > drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:755
virtual void getFwdRatesOfProgress_ddP(span< double > drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:766
virtual void getNetRatesOfProgress_ddC(span< double > drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition Kinetics.h:926
void getNetProductionRates_ddP(span< double > dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition Kinetics.cpp:586
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:609
void getCreationRates_ddC(span< double > dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition Kinetics.cpp:485
virtual void getNetRatesOfProgress_ddP(span< double > drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:912
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:890
void getDestructionRates_ddP(span< double > dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:531
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:870
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:183
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:223
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...