Cantera  3.1.0b1
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
27void Kinetics::checkReactionIndex(size_t i) const
28{
29 if (i >= nReactions()) {
30 throw IndexError("Kinetics::checkReactionIndex", "reactions", i,
31 nReactions()-1);
32 }
33}
34
36{
37 size_t nRxn = nReactions();
38
39 // Stoichiometry managers
43
44 m_rbuf.resize(nRxn);
45
46 // products are created for positive net rate of progress
48 // reactants are destroyed for positive net rate of progress
50
51 m_ready = true;
52}
53
55{
56 if (nReactions() > ii) {
57 throw ArraySizeError("Kinetics::checkReactionArraySize", ii,
58 nReactions());
59 }
60}
61
62void Kinetics::checkPhaseIndex(size_t m) const
63{
64 if (m >= nPhases()) {
65 throw IndexError("Kinetics::checkPhaseIndex", "phase", m, nPhases()-1);
66 }
67}
68
69void Kinetics::checkPhaseArraySize(size_t mm) const
70{
71 if (nPhases() > mm) {
72 throw ArraySizeError("Kinetics::checkPhaseArraySize", mm, nPhases());
73 }
74}
75
77{
78 warn_deprecated("Kinetics::reactionPhaseIndex", "The reacting phase is always "
79 "the first phase. To be removed after Cantera 3.1.");
80 return 0;
81}
82
83shared_ptr<ThermoPhase> Kinetics::reactionPhase() const
84{
85 return m_thermo[0];
86}
87
88void Kinetics::checkSpeciesIndex(size_t k) const
89{
90 if (k >= m_kk) {
91 throw IndexError("Kinetics::checkSpeciesIndex", "species", k, m_kk-1);
92 }
93}
94
96{
97 if (m_kk > kk) {
98 throw ArraySizeError("Kinetics::checkSpeciesArraySize", kk, m_kk);
99 }
100}
101
103{
104 if (flag == "warn" || flag == "error" || flag == "mark-duplicate"
105 || flag == "modify-efficiency")
106 {
107 m_explicit_third_body_duplicates = flag;
108 } else {
109 throw CanteraError("Kinetics::setExplicitThirdBodyDuplicateHandling",
110 "Invalid flag '{}'", flag);
111 }
112}
113
114pair<size_t, size_t> Kinetics::checkDuplicates(bool throw_err) const
115{
116 //! Map of (key indicating participating species) to reaction numbers
117 map<size_t, vector<size_t>> participants;
118 vector<map<int, double>> net_stoich;
119 std::unordered_set<size_t> unmatched_duplicates;
120 for (size_t i = 0; i < m_reactions.size(); i++) {
121 if (m_reactions[i]->duplicate) {
122 unmatched_duplicates.insert(i);
123 }
124 }
125
126 vector<InputFileError> errs;
127 for (size_t i = 0; i < m_reactions.size(); i++) {
128 // Get data about this reaction
129 unsigned long int key = 0;
130 Reaction& R = *m_reactions[i];
131 net_stoich.emplace_back();
132 map<int, double>& net = net_stoich.back();
133 for (const auto& [name, stoich] : R.reactants) {
134 int k = static_cast<int>(kineticsSpeciesIndex(name));
135 key += k*(k+1);
136 net[-1 -k] -= stoich;
137 }
138 for (const auto& [name, stoich] : R.products) {
139 int k = static_cast<int>(kineticsSpeciesIndex(name));
140 key += k*(k+1);
141 net[1+k] += stoich;
142 }
143
144 // Compare this reaction to others with similar participants
145 vector<size_t>& related = participants[key];
146 for (size_t m : related) {
147 Reaction& other = *m_reactions[m];
148 if (R.duplicate && other.duplicate) {
149 // marked duplicates
150 unmatched_duplicates.erase(i);
151 unmatched_duplicates.erase(m);
152 continue;
153 } else if (R.type() != other.type()) {
154 continue; // different reaction types
155 } else if (R.rate() && other.rate()
156 && R.rate()->type() != other.rate()->type())
157 {
158 continue; // different rate parameterizations
159 }
160 double c = checkDuplicateStoich(net_stoich[i], net_stoich[m]);
161 if (c == 0) {
162 continue; // stoichiometries differ (not by a multiple)
163 } else if (c < 0.0 && !R.reversible && !other.reversible) {
164 continue; // irreversible reactions in opposite directions
165 } else if (R.usesThirdBody() && other.usesThirdBody()) {
166 ThirdBody& tb1 = *(R.thirdBody());
167 ThirdBody& tb2 = *(other.thirdBody());
168 bool thirdBodyOk = true;
169 for (size_t k = 0; k < nTotalSpecies(); k++) {
170 string s = kineticsSpeciesName(k);
171 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
172 // non-zero third body efficiencies for species `s` in
173 // both reactions
174 thirdBodyOk = false;
175 break;
176 }
177 }
178 if (thirdBodyOk) {
179 continue; // No overlap in third body efficiencies
180 } else if ((tb1.name() == "M") != (tb2.name() == "M")) {
181 // Exactly one of the reactions uses M as the third body
182 if (m_explicit_third_body_duplicates == "mark-duplicate") {
183 R.duplicate = true;
184 other.duplicate = true;
185 continue;
186 } else if (m_explicit_third_body_duplicates == "modify-efficiency") {
187 if (tb1.name() == "M") {
188 tb1.efficiencies[tb2.name()] = 0.0;
189 } else {
190 tb2.efficiencies[tb1.name()] = 0.0;
191 }
192 continue;
193 } else if (m_explicit_third_body_duplicates == "warn") {
194 InputFileError msg("Kinetics::checkDuplicates",
195 R.input, other.input,
196 "Undeclared duplicate third body reactions with a common "
197 "third body detected.\nAdd the field "
198 "'explicit-third-body-duplicates: mark-duplicate' or\n"
199 "'explicit-third-body-duplicates: modify-efficiency' to "
200 "the YAML phase entry\nto choose how these reactions "
201 "should be handled and suppress this warning.\n"
202 "Reaction {}: {}\nReaction {}: {}\n",
203 m+1, R.equation(), i+1, other.equation());
204 warn_user("Kinetics::checkDuplicates", msg.what());
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 {
216 return {i,m};
217 }
218 }
219 participants[key].push_back(i);
220 }
221 if (unmatched_duplicates.size()) {
222 size_t i = *unmatched_duplicates.begin();
223 if (throw_err) {
224 errs.emplace_back("Kinetics::checkDuplicates",
225 m_reactions[i]->input,
226 "No duplicate found for declared duplicate reaction number {}"
227 " ({})", i, m_reactions[i]->equation());
228 } else {
229 return {i, i};
230 }
231 }
232 if (errs.empty()) {
233 return {npos, npos};
234 } else if (errs.size() == 1) {
235 throw errs[0];
236 } else {
237 fmt::memory_buffer msg;
238 for (const auto& err : errs) {
239 fmt_append(msg, "\n{}\n", err.getMessage());
240 }
241 throw CanteraError("Kinetics::checkDuplicates", to_string(msg));
242 }
243}
244
245double Kinetics::checkDuplicateStoich(map<int, double>& r1, map<int, double>& r2) const
246{
247 std::unordered_set<int> keys; // species keys (k+1 or -k-1)
248 for (auto& [speciesKey, stoich] : r1) {
249 keys.insert(speciesKey);
250 }
251 for (auto& [speciesKey, stoich] : r2) {
252 keys.insert(speciesKey);
253 }
254 int k1 = r1.begin()->first;
255 // check for duplicate written in the same direction
256 double ratio = 0.0;
257 if (r1[k1] && r2[k1]) {
258 ratio = r2[k1]/r1[k1];
259 bool different = false;
260 for (int k : keys) {
261 if ((r1[k] && !r2[k]) ||
262 (!r1[k] && r2[k]) ||
263 (r1[k] && fabs(r2[k]/r1[k] - ratio) > 1.e-8)) {
264 different = true;
265 break;
266 }
267 }
268 if (!different) {
269 return ratio;
270 }
271 }
272
273 // check for duplicate written in the reverse direction
274 if (r1[k1] == 0.0 || r2[-k1] == 0.0) {
275 return 0.0;
276 }
277 ratio = r2[-k1]/r1[k1];
278 for (int k : keys) {
279 if ((r1[k] && !r2[-k]) ||
280 (!r1[k] && r2[-k]) ||
281 (r1[k] && fabs(r2[-k]/r1[k] - ratio) > 1.e-8)) {
282 return 0.0;
283 }
284 }
285 return ratio;
286}
287
288string Kinetics::kineticsSpeciesName(size_t k) const
289{
290 for (size_t n = m_start.size()-1; n != npos; n--) {
291 if (k >= m_start[n]) {
292 return thermo(n).speciesName(k - m_start[n]);
293 }
294 }
295 return "<unknown>";
296}
297
298size_t Kinetics::kineticsSpeciesIndex(const string& nm) const
299{
300 for (size_t n = 0; n < m_thermo.size(); n++) {
301 // Check the ThermoPhase object for a match
302 size_t k = thermo(n).speciesIndex(nm);
303 if (k != npos) {
304 return k + m_start[n];
305 }
306 }
307 return npos;
308}
309
311{
312 for (size_t n = 0; n < m_thermo.size(); n++) {
313 size_t k = thermo(n).speciesIndex(nm);
314 if (k != npos) {
315 return thermo(n);
316 }
317 }
318 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
319}
320
321const ThermoPhase& Kinetics::speciesPhase(const string& nm) const
322{
323 for (const auto& thermo : m_thermo) {
324 if (thermo->speciesIndex(nm) != npos) {
325 return *thermo;
326 }
327 }
328 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
329}
330
331size_t Kinetics::speciesPhaseIndex(size_t k) const
332{
333 for (size_t n = m_start.size()-1; n != npos; n--) {
334 if (k >= m_start[n]) {
335 return n;
336 }
337 }
338 throw CanteraError("Kinetics::speciesPhaseIndex",
339 "illegal species index: {}", k);
340}
341
342double Kinetics::reactantStoichCoeff(size_t kSpec, size_t irxn) const
343{
344 return m_reactantStoich.stoichCoeffs().coeff(kSpec, irxn);
345}
346
347double Kinetics::productStoichCoeff(size_t kSpec, size_t irxn) const
348{
349 return m_productStoich.stoichCoeffs().coeff(kSpec, irxn);
350}
351
353{
354 updateROP();
355 std::copy(m_ropf.begin(), m_ropf.end(), fwdROP);
356}
357
359{
360 updateROP();
361 std::copy(m_ropr.begin(), m_ropr.end(), revROP);
362}
363
365{
366 updateROP();
367 std::copy(m_ropnet.begin(), m_ropnet.end(), netROP);
368}
369
370void Kinetics::getReactionDelta(const double* prop, double* deltaProp) const
371{
372 fill(deltaProp, deltaProp + nReactions(), 0.0);
373 // products add
374 m_productStoich.incrementReactions(prop, deltaProp);
375 // reactants subtract
376 m_reactantStoich.decrementReactions(prop, deltaProp);
377}
378
379void Kinetics::getRevReactionDelta(const double* prop, double* deltaProp) const
380{
381 fill(deltaProp, deltaProp + nReactions(), 0.0);
382 // products add
383 m_revProductStoich.incrementReactions(prop, deltaProp);
384 // reactants subtract
385 m_reactantStoich.decrementReactions(prop, deltaProp);
386}
387
389{
390 updateROP();
391
392 // zero out the output array
393 fill(cdot, cdot + m_kk, 0.0);
394
395 // the forward direction creates product species
396 m_productStoich.incrementSpecies(m_ropf.data(), cdot);
397
398 // the reverse direction creates reactant species
399 m_reactantStoich.incrementSpecies(m_ropr.data(), cdot);
400}
401
403{
404 updateROP();
405
406 fill(ddot, ddot + m_kk, 0.0);
407 // the reverse direction destroys products in reversible reactions
408 m_revProductStoich.incrementSpecies(m_ropr.data(), ddot);
409 // the forward direction destroys reactants
410 m_reactantStoich.incrementSpecies(m_ropf.data(), ddot);
411}
412
414{
415 updateROP();
416
417 fill(net, net + m_kk, 0.0);
418 // products are created for positive net rate of progress
419 m_productStoich.incrementSpecies(m_ropnet.data(), net);
420 // reactants are destroyed for positive net rate of progress
421 m_reactantStoich.decrementSpecies(m_ropnet.data(), net);
422}
423
425{
426 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
427 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
428 // the forward direction creates product species
429 getFwdRatesOfProgress_ddT(buf.data());
430 out = m_productStoich.stoichCoeffs() * buf;
431 // the reverse direction creates reactant species
432 getRevRatesOfProgress_ddT(buf.data());
433 out += m_reactantStoich.stoichCoeffs() * buf;
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_ddP(buf.data());
442 out = m_productStoich.stoichCoeffs() * buf;
443 // the reverse direction creates reactant species
444 getRevRatesOfProgress_ddP(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_ddC(buf.data());
454 out = m_productStoich.stoichCoeffs() * buf;
455 // the reverse direction creates reactant species
456 getRevRatesOfProgress_ddC(buf.data());
457 out += m_reactantStoich.stoichCoeffs() * buf;
458}
459
460Eigen::SparseMatrix<double> Kinetics::creationRates_ddX()
461{
462 Eigen::SparseMatrix<double> jac;
463 // the forward direction creates product species
465 // the reverse direction creates reactant species
467 return jac;
468}
469
470Eigen::SparseMatrix<double> Kinetics::creationRates_ddCi()
471{
472 Eigen::SparseMatrix<double> jac;
473 // the forward direction creates product species
475 // the reverse direction creates reactant species
477 return jac;
478}
479
481{
482 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
483 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
484 // the reverse direction destroys products in reversible reactions
485 getRevRatesOfProgress_ddT(buf.data());
486 out = m_revProductStoich.stoichCoeffs() * buf;
487 // the forward direction destroys reactants
488 getFwdRatesOfProgress_ddT(buf.data());
489 out += m_reactantStoich.stoichCoeffs() * buf;
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_ddP(buf.data());
498 out = m_revProductStoich.stoichCoeffs() * buf;
499 // the forward direction destroys reactants
500 getFwdRatesOfProgress_ddP(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_ddC(buf.data());
510 out = m_revProductStoich.stoichCoeffs() * buf;
511 // the forward direction destroys reactants
512 getFwdRatesOfProgress_ddC(buf.data());
513 out += m_reactantStoich.stoichCoeffs() * buf;
514}
515
516Eigen::SparseMatrix<double> Kinetics::destructionRates_ddX()
517{
518 Eigen::SparseMatrix<double> jac;
519 // the reverse direction destroys products in reversible reactions
521 // the forward direction destroys reactants
523 return jac;
524}
525
526Eigen::SparseMatrix<double> Kinetics::destructionRates_ddCi()
527{
528 Eigen::SparseMatrix<double> jac;
529 // the reverse direction destroys products in reversible reactions
531 // the forward direction destroys reactants
533 return jac;
534}
535
537{
538 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
539 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
540 getNetRatesOfProgress_ddT(buf.data());
541 out = m_stoichMatrix * buf;
542}
543
545{
546 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
547 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
548 getNetRatesOfProgress_ddP(buf.data());
549 out = m_stoichMatrix * buf;
550}
551
553{
554 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
555 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
556 getNetRatesOfProgress_ddC(buf.data());
557 out = m_stoichMatrix * buf;
558}
559
560Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddX()
561{
563}
564
565Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddCi()
566{
568}
569
570void Kinetics::addThermo(shared_ptr<ThermoPhase> thermo)
571{
572 // the phase with lowest dimensionality is assumed to be the
573 // phase/interface at which reactions take place
574 if (thermo->nDim() <= m_mindim) {
575 if (!m_thermo.empty()) {
576 throw CanteraError("Kinetics::addThermo",
577 "The reacting (lowest dimensional) phase must be added first.");
578 }
579 m_mindim = thermo->nDim();
580 }
581
582 m_thermo.push_back(thermo);
583 m_phaseindex[m_thermo.back()->name()] = nPhases();
585}
586
588{
589 AnyMap out;
590 string name = KineticsFactory::factory()->canonicalize(kineticsType());
591 if (name != "none") {
592 out["kinetics"] = name;
593 if (nReactions() == 0) {
594 out["reactions"] = "none";
595 }
597 out["skip-undeclared-third-bodies"] = true;
598 }
599 if (m_explicit_third_body_duplicates == "error") {
600 // "warn" is the default, and does not need to be added. "mark-duplicate"
601 // and "modify-efficiency" do not need to be propagated here as their
602 // effects are already applied to the corresponding reactions.
603 out["explicit-third-body-duplicates"] = "error";
604 }
605 }
606 return out;
607}
608
610{
611 m_kk = 0;
612 m_start.resize(nPhases());
613
614 for (size_t i = 0; i < m_thermo.size(); i++) {
615 m_start[i] = m_kk; // global index of first species of phase i
616 m_kk += m_thermo[i]->nSpecies();
617 }
618 invalidateCache();
619}
620
621bool Kinetics::addReaction(shared_ptr<Reaction> r, bool resize)
622{
623 r->check();
624 r->validate(*this);
625
626 if (m_kk == 0) {
627 init();
628 }
630
631 // Check validity of reaction within the context of the Kinetics object
632 if (!r->checkSpecies(*this)) {
633 // Do not add reaction
634 return false;
635 }
636
637 // For reactions created outside the context of a Kinetics object, the units
638 // of the rate coefficient can't be determined in advance. Do that here.
639 if (r->rate_units.factor() == 0) {
640 r->rate()->setRateUnits(r->calculateRateCoeffUnits(*this));
641 }
642
643 size_t irxn = nReactions(); // index of the new reaction
644
645 // indices of reactant and product species within this Kinetics object
646 vector<size_t> rk, pk;
647
648 // Reactant and product stoichiometric coefficients, such that rstoich[i] is
649 // the coefficient for species rk[i]
650 vector<double> rstoich, pstoich;
651
652 for (const auto& [name, stoich] : r->reactants) {
653 rk.push_back(kineticsSpeciesIndex(name));
654 rstoich.push_back(stoich);
655 }
656
657 for (const auto& [name, stoich] : r->products) {
658 pk.push_back(kineticsSpeciesIndex(name));
659 pstoich.push_back(stoich);
660 }
661
662 // The default order for each reactant is its stoichiometric coefficient,
663 // which can be overridden by entries in the Reaction.orders map. rorder[i]
664 // is the order for species rk[i].
665 vector<double> rorder = rstoich;
666 for (const auto& [name, order] : r->orders) {
667 size_t k = kineticsSpeciesIndex(name);
668 // Find the index of species k within rk
669 auto rloc = std::find(rk.begin(), rk.end(), k);
670 if (rloc != rk.end()) {
671 rorder[rloc - rk.begin()] = order;
672 } else {
673 // If the reaction order involves a non-reactant species, add an
674 // extra term to the reactants with zero stoichiometry so that the
675 // stoichiometry manager can be used to compute the global forward
676 // reaction rate.
677 rk.push_back(k);
678 rstoich.push_back(0.0);
679 rorder.push_back(order);
680 }
681 }
682
683 m_reactantStoich.add(irxn, rk, rorder, rstoich);
684 // product orders = product stoichiometric coefficients
685 m_productStoich.add(irxn, pk, pstoich, pstoich);
686 if (r->reversible) {
687 m_revProductStoich.add(irxn, pk, pstoich, pstoich);
688 }
689
690 m_reactions.push_back(r);
691 m_rfn.push_back(0.0);
692 m_delta_gibbs0.push_back(0.0);
693 m_rkcn.push_back(0.0);
694 m_ropf.push_back(0.0);
695 m_ropr.push_back(0.0);
696 m_ropnet.push_back(0.0);
697 m_perturb.push_back(1.0);
698 m_dH.push_back(0.0);
699
700 if (resize) {
702 } else {
703 m_ready = false;
704 }
705
706 return true;
707}
708
709void Kinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
710{
712 shared_ptr<Reaction>& rOld = m_reactions[i];
713 if (rNew->type() != rOld->type()) {
714 throw CanteraError("Kinetics::modifyReaction",
715 "Reaction types are different: {} != {}.",
716 rOld->type(), rNew->type());
717 }
718
719 if (rNew->rate()->type() != rOld->rate()->type()) {
720 throw CanteraError("Kinetics::modifyReaction",
721 "ReactionRate types are different: {} != {}.",
722 rOld->rate()->type(), rNew->rate()->type());
723 }
724
725 if (rNew->reactants != rOld->reactants) {
726 throw CanteraError("Kinetics::modifyReaction",
727 "Reactants are different: '{}' != '{}'.",
728 rOld->reactantString(), rNew->reactantString());
729 }
730
731 if (rNew->products != rOld->products) {
732 throw CanteraError("Kinetics::modifyReaction",
733 "Products are different: '{}' != '{}'.",
734 rOld->productString(), rNew->productString());
735 }
736 m_reactions[i] = rNew;
737 invalidateCache();
738}
739
740shared_ptr<Reaction> Kinetics::reaction(size_t i)
741{
743 return m_reactions[i];
744}
745
746shared_ptr<const Reaction> Kinetics::reaction(size_t i) const
747{
749 return m_reactions[i];
750}
751
752}
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
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
size_t reactionPhaseIndex() const
Phase where the reactions occur.
Definition Kinetics.cpp:76
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:35
void checkPhaseIndex(size_t m) const
Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
Definition Kinetics.cpp:62
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:95
virtual void getFwdRatesOfProgress(double *fwdROP)
Return the forward rates of progress of the reactions.
Definition Kinetics.cpp:352
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:242
vector< shared_ptr< Reaction > > m_reactions
Vector of Reaction objects represented by this Kinetics manager.
Definition Kinetics.h:1464
void checkSpeciesIndex(size_t k) const
Check that the specified species index is in range Throws an exception if k is greater than nSpecies(...
Definition Kinetics.cpp:88
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:245
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1461
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1505
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition Kinetics.cpp:740
bool m_ready
Boolean indicating whether Kinetics object is fully configured.
Definition Kinetics.h:1453
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:379
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:102
virtual void getNetRatesOfProgress(double *netROP)
Net rates of progress.
Definition Kinetics.cpp:364
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:1482
virtual string kineticsType() const
Identifies the Kinetics manager type.
Definition Kinetics.h:144
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1502
shared_ptr< ThermoPhase > reactionPhase() const
Return pointer to phase where the reactions occur.
Definition Kinetics.cpp:83
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1457
virtual pair< size_t, size_t > checkDuplicates(bool throw_err=true) const
Check for unmarked duplicate reactions and unmatched marked duplicates.
Definition Kinetics.cpp:114
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:370
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:69
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1514
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1508
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition Kinetics.h:184
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:1476
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:621
string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition Kinetics.cpp:288
virtual void getDestructionRates(double *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:402
AnyMap parameters()
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition Kinetics.cpp:587
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition Kinetics.h:1280
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:709
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1511
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition Kinetics.cpp:347
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:570
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1517
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition Kinetics.h:1449
map< string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition Kinetics.h:1490
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition Kinetics.h:1493
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition Kinetics.h:1443
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1446
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition Kinetics.h:152
virtual void getRevRatesOfProgress(double *revROP)
Return the Reverse rates of progress of the reactions.
Definition Kinetics.cpp:358
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1526
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:276
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1496
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1440
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:609
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:54
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:254
virtual double reactantStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a reactant in reaction i.
Definition Kinetics.cpp:342
void checkReactionIndex(size_t m) const
Check that the specified reaction index is in range Throws an exception if i is greater than nReactio...
Definition Kinetics.cpp:27
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:413
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:310
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1499
virtual void getCreationRates(double *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:388
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:331
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition Phase.h:546
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:142
size_t speciesIndex(const string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:129
Abstract base class which stores data about a reaction and its rate parameterization so that it can b...
Definition Reaction.h:25
shared_ptr< ReactionRate > rate()
Get reaction rate pointer.
Definition Reaction.h:147
bool usesThirdBody() const
Check whether reaction involves third body collider.
Definition Reaction.h:161
shared_ptr< ThirdBody > thirdBody()
Get pointer to third-body handler.
Definition Reaction.h:155
bool reversible
True if the current reaction is reversible. False otherwise.
Definition Reaction.h:126
string type() const
The type of reaction, including reaction rate information.
Definition Reaction.cpp:518
string equation() const
The chemical equation for this reaction.
Definition Reaction.cpp:357
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:192
double efficiency(const string &k) const
Get the third-body efficiency for species k
Definition Reaction.cpp:806
Composition efficiencies
Map of species to third body efficiency.
Definition Reaction.h:235
string name() const
Name of the third body collider.
Definition Reaction.h:200
void fmt_append(fmt::memory_buffer &b, const std::string &tmpl, Args... args)
Versions 6.2.0 and 6.2.1 of fmtlib do not include this define before they include windows....
Definition fmt.h:29
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:812
Eigen::SparseMatrix< double > creationRates_ddCi()
Calculate derivatives for species creation rates with respect to species concentration at constant te...
Definition Kinetics.cpp:470
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddCi()
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
Definition Kinetics.h:958
void getCreationRates_ddT(double *dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:424
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:938
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition Kinetics.cpp:516
void getCreationRates_ddC(double *dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition Kinetics.cpp:448
void getDestructionRates_ddP(double *dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:492
void getNetProductionRates_ddC(double *dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition Kinetics.cpp:552
void getDestructionRates_ddT(double *dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:480
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:761
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition Kinetics.cpp:560
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition Kinetics.cpp:460
Eigen::SparseMatrix< double > destructionRates_ddCi()
Calculate derivatives for species destruction rates with respect to species concentration at constant...
Definition Kinetics.cpp:526
void getNetProductionRates_ddT(double *dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition Kinetics.cpp:536
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:792
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:823
void getCreationRates_ddP(double *dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:436
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:896
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:834
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:565
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:848
void getNetProductionRates_ddP(double *dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition Kinetics.cpp:544
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:907
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:885
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:750
void getDestructionRates_ddC(double *dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition Kinetics.cpp:504
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition Kinetics.h:921
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:865
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:775
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
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
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...