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