Cantera 2.6.0
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/join.hpp>
21
22using namespace std;
23
24namespace Cantera
25{
27 m_ready(false),
28 m_kk(0),
29 m_thermo(0),
30 m_surfphase(npos),
31 m_rxnphase(npos),
32 m_mindim(4),
33 m_skipUndeclaredSpecies(false),
34 m_skipUndeclaredThirdBodies(false)
35{
36}
37
38Kinetics::~Kinetics() {}
39
40void Kinetics::checkReactionIndex(size_t i) const
41{
42 if (i >= nReactions()) {
43 throw IndexError("Kinetics::checkReactionIndex", "reactions", i,
44 nReactions()-1);
45 }
46}
47
49{
50 size_t nRxn = nReactions();
51
52 // Stoichiometry managers
56
57 m_rbuf.resize(nRxn);
58
59 // products are created for positive net rate of progress
61 // reactants are destroyed for positive net rate of progress
63
64 m_ready = true;
65}
66
68{
69 if (nReactions() > ii) {
70 throw ArraySizeError("Kinetics::checkReactionArraySize", ii,
71 nReactions());
72 }
73}
74
75void Kinetics::checkPhaseIndex(size_t m) const
76{
77 if (m >= nPhases()) {
78 throw IndexError("Kinetics::checkPhaseIndex", "phase", m, nPhases()-1);
79 }
80}
81
82void Kinetics::checkPhaseArraySize(size_t mm) const
83{
84 if (nPhases() > mm) {
85 throw ArraySizeError("Kinetics::checkPhaseArraySize", mm, nPhases());
86 }
87}
88
89void Kinetics::checkSpeciesIndex(size_t k) const
90{
91 if (k >= m_kk) {
92 throw IndexError("Kinetics::checkSpeciesIndex", "species", k, m_kk-1);
93 }
94}
95
97{
98 if (m_kk > kk) {
99 throw ArraySizeError("Kinetics::checkSpeciesArraySize", kk, m_kk);
100 }
101}
102
103std::pair<size_t, size_t> Kinetics::checkDuplicates(bool throw_err) const
104{
105 //! Map of (key indicating participating species) to reaction numbers
106 std::map<size_t, std::vector<size_t> > participants;
107 std::vector<std::map<int, double> > net_stoich;
108 std::unordered_set<size_t> unmatched_duplicates;
109 for (size_t i = 0; i < m_reactions.size(); i++) {
110 if (m_reactions[i]->duplicate) {
111 unmatched_duplicates.insert(i);
112 }
113 }
114
115 for (size_t i = 0; i < m_reactions.size(); i++) {
116 // Get data about this reaction
117 unsigned long int key = 0;
118 Reaction& R = *m_reactions[i];
119 net_stoich.emplace_back();
120 std::map<int, double>& net = net_stoich.back();
121 for (const auto& sp : R.reactants) {
122 int k = static_cast<int>(kineticsSpeciesIndex(sp.first));
123 key += k*(k+1);
124 net[-1 -k] -= sp.second;
125 }
126 for (const auto& sp : R.products) {
127 int k = static_cast<int>(kineticsSpeciesIndex(sp.first));
128 key += k*(k+1);
129 net[1+k] += sp.second;
130 }
131
132 // Compare this reaction to others with similar participants
133 vector<size_t>& related = participants[key];
134 for (size_t m : related) {
135 Reaction& other = *m_reactions[m];
136 if (R.duplicate && other.duplicate) {
137 // marked duplicates
138 unmatched_duplicates.erase(i);
139 unmatched_duplicates.erase(m);
140 continue;
141 } else if (R.type() != other.type()) {
142 continue; // different reaction types
143 } else if (R.rate() && other.rate()
144 && R.rate()->type() != other.rate()->type())
145 {
146 continue; // different rate parameterizations
147 }
148 doublereal c = checkDuplicateStoich(net_stoich[i], net_stoich[m]);
149 if (c == 0) {
150 continue; // stoichiometries differ (not by a multiple)
151 } else if (c < 0.0 && !R.reversible && !other.reversible) {
152 continue; // irreversible reactions in opposite directions
153 } else if (R.type() == "falloff-legacy" ||
154 R.type() == "chemically-activated-legacy") {
155 ThirdBody& tb1 = dynamic_cast<FalloffReaction2&>(R).third_body;
156 ThirdBody& tb2 = dynamic_cast<FalloffReaction2&>(other).third_body;
157 bool thirdBodyOk = true;
158 for (size_t k = 0; k < nTotalSpecies(); k++) {
159 string s = kineticsSpeciesName(k);
160 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
161 // non-zero third body efficiencies for species `s` in
162 // both reactions
163 thirdBodyOk = false;
164 break;
165 }
166 }
167 if (thirdBodyOk) {
168 continue; // No overlap in third body efficiencies
169 }
170 } else if (R.type() == "falloff" || R.type() == "chemically-activated") {
171 auto tb1 = dynamic_cast<FalloffReaction3&>(R).thirdBody();
172 auto tb2 = dynamic_cast<FalloffReaction3&>(other).thirdBody();
173 bool thirdBodyOk = true;
174 for (size_t k = 0; k < nTotalSpecies(); k++) {
175 string s = kineticsSpeciesName(k);
176 if (tb1->efficiency(s) * tb2->efficiency(s) != 0.0) {
177 // non-zero third body efficiencies for species `s` in
178 // both reactions
179 thirdBodyOk = false;
180 break;
181 }
182 }
183 if (thirdBodyOk) {
184 continue; // No overlap in third body efficiencies
185 }
186 } else if (R.type() == "three-body") {
187 ThirdBody& tb1 = *(dynamic_cast<ThreeBodyReaction3&>(R).thirdBody());
188 ThirdBody& tb2 = *(dynamic_cast<ThreeBodyReaction3&>(other).thirdBody());
189 bool thirdBodyOk = true;
190 for (size_t k = 0; k < nTotalSpecies(); k++) {
191 string s = kineticsSpeciesName(k);
192 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
193 // non-zero third body efficiencies for species `s` in
194 // both reactions
195 thirdBodyOk = false;
196 break;
197 }
198 }
199 if (thirdBodyOk) {
200 continue; // No overlap in third body efficiencies
201 }
202 } else if (R.type() == "three-body-legacy") {
203 ThirdBody& tb1 = dynamic_cast<ThreeBodyReaction2&>(R).third_body;
204 ThirdBody& tb2 = dynamic_cast<ThreeBodyReaction2&>(other).third_body;
205 bool thirdBodyOk = true;
206 for (size_t k = 0; k < nTotalSpecies(); k++) {
207 string s = kineticsSpeciesName(k);
208 if (tb1.efficiency(s) * tb2.efficiency(s) != 0.0) {
209 // non-zero third body efficiencies for species `s` in
210 // both reactions
211 thirdBodyOk = false;
212 break;
213 }
214 }
215 if (thirdBodyOk) {
216 continue; // No overlap in third body efficiencies
217 }
218 }
219 if (throw_err) {
220 throw InputFileError("Kinetics::checkDuplicates",
221 R.input, other.input,
222 "Undeclared duplicate reactions detected:\n"
223 "Reaction {}: {}\nReaction {}: {}\n",
224 i+1, other.equation(), m+1, R.equation());
225 } else {
226 return {i,m};
227 }
228 }
229 participants[key].push_back(i);
230 }
231 if (unmatched_duplicates.size()) {
232 size_t i = *unmatched_duplicates.begin();
233 if (throw_err) {
234 throw InputFileError("Kinetics::checkDuplicates",
235 m_reactions[i]->input,
236 "No duplicate found for declared duplicate reaction number {}"
237 " ({})", i, m_reactions[i]->equation());
238 } else {
239 return {i, i};
240 }
241 }
242 return {npos, npos};
243}
244
245double Kinetics::checkDuplicateStoich(std::map<int, double>& r1,
246 std::map<int, double>& r2) const
247{
248 std::unordered_set<int> keys; // species keys (k+1 or -k-1)
249 for (auto& r : r1) {
250 keys.insert(r.first);
251 }
252 for (auto& r : r2) {
253 keys.insert(r.first);
254 }
255 int k1 = r1.begin()->first;
256 // check for duplicate written in the same direction
257 doublereal ratio = 0.0;
258 if (r1[k1] && r2[k1]) {
259 ratio = r2[k1]/r1[k1];
260 bool different = false;
261 for (int k : keys) {
262 if ((r1[k] && !r2[k]) ||
263 (!r1[k] && r2[k]) ||
264 (r1[k] && fabs(r2[k]/r1[k] - ratio) > 1.e-8)) {
265 different = true;
266 break;
267 }
268 }
269 if (!different) {
270 return ratio;
271 }
272 }
273
274 // check for duplicate written in the reverse direction
275 if (r1[k1] == 0.0 || r2[-k1] == 0.0) {
276 return 0.0;
277 }
278 ratio = r2[-k1]/r1[k1];
279 for (int k : keys) {
280 if ((r1[k] && !r2[-k]) ||
281 (!r1[k] && r2[-k]) ||
282 (r1[k] && fabs(r2[-k]/r1[k] - ratio) > 1.e-8)) {
283 return 0.0;
284 }
285 }
286 return ratio;
287}
288
290{
291 R.checkBalance(*this);
292 warn_deprecated("Kinetics::checkReactionBalance",
293 "To be removed after Cantera 2.6. Replacable by Reaction::checkBalance.");
294}
295
296void Kinetics::selectPhase(const double* data, const ThermoPhase* phase,
297 double* phase_data)
298{
299 for (size_t n = 0; n < nPhases(); n++) {
300 if (phase == m_thermo[n]) {
301 size_t nsp = phase->nSpecies();
302 copy(data + m_start[n],
303 data + m_start[n] + nsp, phase_data);
304 return;
305 }
306 }
307 throw CanteraError("Kinetics::selectPhase", "Phase not found.");
308}
309
310string Kinetics::kineticsSpeciesName(size_t k) const
311{
312 for (size_t n = m_start.size()-1; n != npos; n--) {
313 if (k >= m_start[n]) {
314 return thermo(n).speciesName(k - m_start[n]);
315 }
316 }
317 return "<unknown>";
318}
319
320size_t Kinetics::kineticsSpeciesIndex(const std::string& nm) const
321{
322 for (size_t n = 0; n < m_thermo.size(); n++) {
323 // Check the ThermoPhase object for a match
324 size_t k = thermo(n).speciesIndex(nm);
325 if (k != npos) {
326 return k + m_start[n];
327 }
328 }
329 return npos;
330}
331
332size_t Kinetics::kineticsSpeciesIndex(const std::string& nm,
333 const std::string& ph) const
334{
335 if (ph == "<any>") {
336 return kineticsSpeciesIndex(nm);
337 }
338
339 for (size_t n = 0; n < m_thermo.size(); n++) {
340 string id = thermo(n).name();
341 if (ph == id) {
342 size_t k = thermo(n).speciesIndex(nm);
343 if (k == npos) {
344 return npos;
345 }
346 return k + m_start[n];
347 }
348 }
349 return npos;
350}
351
352ThermoPhase& Kinetics::speciesPhase(const std::string& nm)
353{
354 for (size_t n = 0; n < m_thermo.size(); n++) {
355 size_t k = thermo(n).speciesIndex(nm);
356 if (k != npos) {
357 return thermo(n);
358 }
359 }
360 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
361}
362
363const ThermoPhase& Kinetics::speciesPhase(const std::string& nm) const
364{
365 for (const auto thermo : m_thermo) {
366 if (thermo->speciesIndex(nm) != npos) {
367 return *thermo;
368 }
369 }
370 throw CanteraError("Kinetics::speciesPhase", "unknown species '{}'", nm);
371}
372
373size_t Kinetics::speciesPhaseIndex(size_t k) const
374{
375 for (size_t n = m_start.size()-1; n != npos; n--) {
376 if (k >= m_start[n]) {
377 return n;
378 }
379 }
380 throw CanteraError("Kinetics::speciesPhaseIndex",
381 "illegal species index: {}", k);
382}
383
384double Kinetics::reactantStoichCoeff(size_t kSpec, size_t irxn) const
385{
386 return m_reactantStoich.stoichCoeffs().coeff(kSpec, irxn);
387}
388
389double Kinetics::productStoichCoeff(size_t kSpec, size_t irxn) const
390{
391 return m_productStoich.stoichCoeffs().coeff(kSpec, irxn);
392}
393
394int Kinetics::reactionType(size_t i) const {
395 warn_deprecated("Kinetics::reactionType",
396 "To be changed after Cantera 2.6. "
397 "Return string instead of magic number; use "
398 "Kinetics::reactionTypeStr during transition.");
399 return m_reactions[i]->reaction_type;
400}
401
402std::string Kinetics::reactionTypeStr(size_t i) const {
403 return m_reactions[i]->type();
404}
405
406std::string Kinetics::reactionString(size_t i) const
407{
408 return m_reactions[i]->equation();
409}
410
411//! Returns a string containing the reactants side of the reaction equation.
412std::string Kinetics::reactantString(size_t i) const
413{
414 return m_reactions[i]->reactantString();
415}
416
417//! Returns a string containing the products side of the reaction equation.
418std::string Kinetics::productString(size_t i) const
419{
420 return m_reactions[i]->productString();
421}
422
423void Kinetics::getFwdRatesOfProgress(doublereal* fwdROP)
424{
425 updateROP();
426 std::copy(m_ropf.begin(), m_ropf.end(), fwdROP);
427}
428
429void Kinetics::getRevRatesOfProgress(doublereal* revROP)
430{
431 updateROP();
432 std::copy(m_ropr.begin(), m_ropr.end(), revROP);
433}
434
435void Kinetics::getNetRatesOfProgress(doublereal* netROP)
436{
437 updateROP();
438 std::copy(m_ropnet.begin(), m_ropnet.end(), netROP);
439}
440
441void Kinetics::getReactionDelta(const double* prop, double* deltaProp) const
442{
443 fill(deltaProp, deltaProp + nReactions(), 0.0);
444 // products add
445 m_productStoich.incrementReactions(prop, deltaProp);
446 // reactants subtract
447 m_reactantStoich.decrementReactions(prop, deltaProp);
448}
449
450void Kinetics::getRevReactionDelta(const double* prop, double* deltaProp) const
451{
452 fill(deltaProp, deltaProp + nReactions(), 0.0);
453 // products add
454 m_revProductStoich.incrementReactions(prop, deltaProp);
455 // reactants subtract
456 m_reactantStoich.decrementReactions(prop, deltaProp);
457}
458
460{
461 updateROP();
462
463 // zero out the output array
464 fill(cdot, cdot + m_kk, 0.0);
465
466 // the forward direction creates product species
467 m_productStoich.incrementSpecies(m_ropf.data(), cdot);
468
469 // the reverse direction creates reactant species
470 m_reactantStoich.incrementSpecies(m_ropr.data(), cdot);
471}
472
473void Kinetics::getDestructionRates(doublereal* ddot)
474{
475 updateROP();
476
477 fill(ddot, ddot + m_kk, 0.0);
478 // the reverse direction destroys products in reversible reactions
479 m_revProductStoich.incrementSpecies(m_ropr.data(), ddot);
480 // the forward direction destroys reactants
481 m_reactantStoich.incrementSpecies(m_ropf.data(), ddot);
482}
483
485{
486 updateROP();
487
488 fill(net, net + m_kk, 0.0);
489 // products are created for positive net rate of progress
490 m_productStoich.incrementSpecies(m_ropnet.data(), net);
491 // reactants are destroyed for positive net rate of progress
492 m_reactantStoich.decrementSpecies(m_ropnet.data(), net);
493}
494
496{
497 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
498 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
499 // the forward direction creates product species
500 getFwdRatesOfProgress_ddT(buf.data());
501 out = m_productStoich.stoichCoeffs() * buf;
502 // the reverse direction creates reactant species
503 getRevRatesOfProgress_ddT(buf.data());
504 out += m_reactantStoich.stoichCoeffs() * buf;
505}
506
508{
509 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
510 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
511 // the forward direction creates product species
512 getFwdRatesOfProgress_ddP(buf.data());
513 out = m_productStoich.stoichCoeffs() * buf;
514 // the reverse direction creates reactant species
515 getRevRatesOfProgress_ddP(buf.data());
516 out += m_reactantStoich.stoichCoeffs() * buf;
517}
518
520{
521 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
522 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
523 // the forward direction creates product species
524 getFwdRatesOfProgress_ddC(buf.data());
525 out = m_productStoich.stoichCoeffs() * buf;
526 // the reverse direction creates reactant species
527 getRevRatesOfProgress_ddC(buf.data());
528 out += m_reactantStoich.stoichCoeffs() * buf;
529}
530
531Eigen::SparseMatrix<double> Kinetics::creationRates_ddX()
532{
533 Eigen::SparseMatrix<double> jac;
534 // the forward direction creates product species
536 // the reverse direction creates reactant species
538 return jac;
539}
540
542{
543 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
544 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
545 // the reverse direction destroys products in reversible reactions
546 getRevRatesOfProgress_ddT(buf.data());
547 out = m_revProductStoich.stoichCoeffs() * buf;
548 // the forward direction destroys reactants
549 getFwdRatesOfProgress_ddT(buf.data());
550 out += m_reactantStoich.stoichCoeffs() * buf;
551}
552
554{
555 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
556 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
557 // the reverse direction destroys products in reversible reactions
558 getRevRatesOfProgress_ddP(buf.data());
559 out = m_revProductStoich.stoichCoeffs() * buf;
560 // the forward direction destroys reactants
561 getFwdRatesOfProgress_ddP(buf.data());
562 out += m_reactantStoich.stoichCoeffs() * buf;
563}
564
566{
567 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
568 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
569 // the reverse direction destroys products in reversible reactions
570 getRevRatesOfProgress_ddC(buf.data());
571 out = m_revProductStoich.stoichCoeffs() * buf;
572 // the forward direction destroys reactants
573 getFwdRatesOfProgress_ddC(buf.data());
574 out += m_reactantStoich.stoichCoeffs() * buf;
575}
576
577Eigen::SparseMatrix<double> Kinetics::destructionRates_ddX()
578{
579 Eigen::SparseMatrix<double> jac;
580 // the reverse direction destroys products in reversible reactions
582 // the forward direction destroys reactants
584 return jac;
585}
586
588{
589 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
590 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
591 getNetRatesOfProgress_ddT(buf.data());
592 out = m_stoichMatrix * buf;
593}
594
596{
597 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
598 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
599 getNetRatesOfProgress_ddP(buf.data());
600 out = m_stoichMatrix * buf;
601}
602
604{
605 Eigen::Map<Eigen::VectorXd> out(dwdot, m_kk);
606 Eigen::Map<Eigen::VectorXd> buf(m_rbuf.data(), nReactions());
607 getNetRatesOfProgress_ddC(buf.data());
608 out = m_stoichMatrix * buf;
609}
610
611Eigen::SparseMatrix<double> Kinetics::netProductionRates_ddX()
612{
614}
615
617{
618 // the phase with lowest dimensionality is assumed to be the
619 // phase/interface at which reactions take place
620 if (thermo.nDim() <= m_mindim) {
621 m_mindim = thermo.nDim();
623 }
624
625 // there should only be one surface phase
626 if (thermo.type() == kineticsType()) {
629 }
630 m_thermo.push_back(&thermo);
631 m_phaseindex[m_thermo.back()->name()] = nPhases();
633}
634
636{
637 AnyMap out;
638 string name = KineticsFactory::factory()->canonicalize(kineticsType());
639 if (name != "none") {
640 out["kinetics"] = name;
641 if (nReactions() == 0) {
642 out["reactions"] = "none";
643 }
644 }
645 return out;
646}
647
649{
650 m_kk = 0;
651 m_start.resize(nPhases());
652
653 for (size_t i = 0; i < m_thermo.size(); i++) {
654 m_start[i] = m_kk; // global index of first species of phase i
655 m_kk += m_thermo[i]->nSpecies();
656 }
657 invalidateCache();
658}
659
660bool Kinetics::addReaction(shared_ptr<Reaction> r, bool resize)
661{
662 r->validate();
663 r->validate(*this);
664
665 if (m_kk == 0) {
666 init();
667 }
669
670 // Check validity of reaction within the context of the Kinetics object
671 if (!r->checkSpecies(*this)) {
672 // Do not add reaction
673 return false;
674 }
675
676 // For reactions created outside the context of a Kinetics object, the units
677 // of the rate coefficient can't be determined in advance. Do that here.
678 if (r->rate_units.factor() == 0) {
679 r->calculateRateCoeffUnits(*this);
680 }
681
682 size_t irxn = nReactions(); // index of the new reaction
683
684 // indices of reactant and product species within this Kinetics object
685 std::vector<size_t> rk, pk;
686
687 // Reactant and product stoichiometric coefficients, such that rstoich[i] is
688 // the coefficient for species rk[i]
689 vector_fp rstoich, pstoich;
690
691 for (const auto& sp : r->reactants) {
692 rk.push_back(kineticsSpeciesIndex(sp.first));
693 rstoich.push_back(sp.second);
694 }
695
696 for (const auto& sp : r->products) {
697 pk.push_back(kineticsSpeciesIndex(sp.first));
698 pstoich.push_back(sp.second);
699 }
700
701 // The default order for each reactant is its stoichiometric coefficient,
702 // which can be overridden by entries in the Reaction.orders map. rorder[i]
703 // is the order for species rk[i].
704 vector_fp rorder = rstoich;
705 for (const auto& sp : r->orders) {
706 size_t k = kineticsSpeciesIndex(sp.first);
707 // Find the index of species k within rk
708 auto rloc = std::find(rk.begin(), rk.end(), k);
709 if (rloc != rk.end()) {
710 rorder[rloc - rk.begin()] = sp.second;
711 } else {
712 // If the reaction order involves a non-reactant species, add an
713 // extra term to the reactants with zero stoichiometry so that the
714 // stoichiometry manager can be used to compute the global forward
715 // reaction rate.
716 rk.push_back(k);
717 rstoich.push_back(0.0);
718 rorder.push_back(sp.second);
719 }
720 }
721
722 m_reactantStoich.add(irxn, rk, rorder, rstoich);
723 // product orders = product stoichiometric coefficients
724 m_productStoich.add(irxn, pk, pstoich, pstoich);
725 if (r->reversible) {
726 m_revProductStoich.add(irxn, pk, pstoich, pstoich);
727 }
728
729 m_reactions.push_back(r);
730 m_rfn.push_back(0.0);
731 m_delta_gibbs0.push_back(0.0);
732 m_rkcn.push_back(0.0);
733 m_ropf.push_back(0.0);
734 m_ropr.push_back(0.0);
735 m_ropnet.push_back(0.0);
736 m_perturb.push_back(1.0);
737 m_dH.push_back(0.0);
738
739 if (resize) {
741 } else {
742 m_ready = false;
743 }
744
745 return true;
746}
747
748void Kinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
749{
751 shared_ptr<Reaction>& rOld = m_reactions[i];
752 if (rNew->type() != rOld->type()) {
753 throw CanteraError("Kinetics::modifyReaction",
754 "Reaction types are different: {} != {}.",
755 rOld->type(), rNew->type());
756 }
757
758 if (!(rNew->usesLegacy())) {
759 if (rNew->rate()->type() != rOld->rate()->type()) {
760 throw CanteraError("Kinetics::modifyReaction",
761 "ReactionRate types are different: {} != {}.",
762 rOld->rate()->type(), rNew->rate()->type());
763 }
764 }
765
766 if (rNew->reactants != rOld->reactants) {
767 throw CanteraError("Kinetics::modifyReaction",
768 "Reactants are different: '{}' != '{}'.",
769 rOld->reactantString(), rNew->reactantString());
770 }
771
772 if (rNew->products != rOld->products) {
773 throw CanteraError("Kinetics::modifyReaction",
774 "Products are different: '{}' != '{}'.",
775 rOld->productString(), rNew->productString());
776 }
777 m_reactions[i] = rNew;
778 invalidateCache();
779}
780
781shared_ptr<Reaction> Kinetics::reaction(size_t i)
782{
784 return m_reactions[i];
785}
786
787shared_ptr<const Reaction> Kinetics::reaction(size_t i) const
788{
790 return m_reactions[i];
791}
792
793double Kinetics::reactionEnthalpy(const Composition& reactants, const Composition& products)
794{
796 for (size_t n = 0; n < nPhases(); n++) {
798 }
799 double rxn_deltaH = 0;
800 for (const auto& product : products) {
801 size_t k = kineticsSpeciesIndex(product.first);
802 double stoich = product.second;
803 rxn_deltaH += hk[k] * stoich;
804 }
805 for (const auto& reactant : reactants) {
806 size_t k = kineticsSpeciesIndex(reactant.first);
807 double stoich = reactant.second;
808 rxn_deltaH -= hk[k] * stoich;
809 }
810 return rxn_deltaH;
811}
812
813}
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:399
Array size error.
Definition: ctexceptions.h:128
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
A reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Definition: Reaction.h:297
A falloff reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Definition: Reaction.h:504
An array index is out of range.
Definition: ctexceptions.h:158
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition: AnyMap.h:702
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition: Kinetics.cpp:48
vector_fp m_delta_gibbs0
Delta G^0 for all reactions.
Definition: Kinetics.h:1372
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:75
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:96
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition: Kinetics.h:231
std::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:1345
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:89
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition: Kinetics.cpp:781
bool m_ready
Boolean indicating whether Kinetics object is fully configured.
Definition: Kinetics.h:1313
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:450
vector_fp m_ropnet
Net rate-of-progress for each reaction.
Definition: Kinetics.h:1384
void selectPhase(const double *data, const ThermoPhase *phase, double *phase_data)
Definition: Kinetics.cpp:296
void getCreationRates_ddT(double *dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition: Kinetics.cpp:495
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:810
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition: Kinetics.cpp:577
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition: Kinetics.h:1317
ThermoPhase & speciesPhase(const std::string &nm)
This function looks up the name of a species and returns a reference to the ThermoPhase object of the...
Definition: Kinetics.cpp:352
size_t m_rxnphase
Phase Index where reactions are assumed to be taking place.
Definition: Kinetics.h:1363
void getCreationRates_ddC(double *dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition: Kinetics.cpp:519
void getDestructionRates_ddP(double *dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition: Kinetics.cpp:553
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition: Kinetics.cpp:441
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:82
vector_fp m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition: Kinetics.h:1321
void getNetProductionRates_ddC(double *dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition: Kinetics.cpp:603
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition: Kinetics.h:171
void getDestructionRates_ddT(double *dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition: Kinetics.cpp:541
std::string productString(size_t i) const
Returns a string containing the products side of the reaction equation.
Definition: Kinetics.cpp:418
virtual void addPhase(ThermoPhase &thermo)
Add a phase to the kinetics manager object.
Definition: Kinetics.cpp:616
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition: Kinetics.h:673
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition: Kinetics.cpp:611
vector_fp m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition: Kinetics.h:1375
vector_fp m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition: Kinetics.h:1390
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition: Kinetics.cpp:531
vector_fp m_ropf
Forward rate-of-progress for each reaction.
Definition: Kinetics.h:1378
vector_fp m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition: Kinetics.h:1387
void getNetProductionRates_ddT(double *dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition: Kinetics.cpp:587
vector_fp m_ropr
Reverse rate-of-progress for each reaction.
Definition: Kinetics.h:1381
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition: Kinetics.h:704
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition: Kinetics.h:715
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition: Kinetics.cpp:660
std::string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition: Kinetics.cpp:310
std::vector< ThermoPhase * > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
Definition: Kinetics.h:1339
AnyMap parameters()
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition: Kinetics.cpp:635
virtual void getNetProductionRates(doublereal *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.cpp:484
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition: Kinetics.h:1125
void getCreationRates_ddP(double *dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition: Kinetics.cpp:507
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition: Kinetics.h:768
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition: Kinetics.cpp:748
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition: Kinetics.cpp:389
std::map< std::string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition: Kinetics.h:1353
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition: Kinetics.h:1309
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition: Kinetics.h:726
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition: Kinetics.h:1366
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition: Kinetics.h:1303
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition: Kinetics.h:1306
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition: Kinetics.h:740
virtual double reactionEnthalpy(const Composition &reactants, const Composition &products)
Calculate the reaction enthalpy of a reaction which has not necessarily been added into the Kinetics ...
Definition: Kinetics.cpp:793
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition: Kinetics.h:139
virtual void getDestructionRates(doublereal *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.cpp:473
virtual int reactionType(size_t i) const
Flag specifying the type of reaction.
Definition: Kinetics.cpp:394
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:265
size_t m_surfphase
Index in the list of phases of the one surface phase.
Definition: Kinetics.h:1356
virtual void getFwdRatesOfProgress(doublereal *fwdROP)
Return the forward rates of progress of the reactions.
Definition: Kinetics.cpp:423
void getNetProductionRates_ddP(double *dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition: Kinetics.cpp:595
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition: Kinetics.h:779
vector_fp m_rfn
Forward rate constant for each reaction.
Definition: Kinetics.h:1369
virtual std::string reactionTypeStr(size_t i) const
String specifying the type of reaction.
Definition: Kinetics.cpp:402
virtual std::string kineticsType() const
Identifies the Kinetics manager type.
Definition: Kinetics.h:131
virtual std::pair< size_t, size_t > checkDuplicates(bool throw_err=true) const
Check for unmarked duplicate reactions and unmatched marked duplicates.
Definition: Kinetics.cpp:103
void checkReactionBalance(const Reaction &R)
Check that the specified reaction is balanced (same number of atoms for each element in the reactants...
Definition: Kinetics.cpp:289
std::string reactantString(size_t i) const
Returns a string containing the reactants side of the reaction equation.
Definition: Kinetics.cpp:412
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition: Kinetics.h:662
void getDestructionRates_ddC(double *dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition: Kinetics.cpp:565
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition: Kinetics.h:1300
virtual void getNetRatesOfProgress(doublereal *netROP)
Net rates of progress.
Definition: Kinetics.cpp:435
std::string reactionString(size_t i) const
Return a string representing the reaction.
Definition: Kinetics.cpp:406
Kinetics()
Default constructor.
Definition: Kinetics.cpp:26
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition: Kinetics.h:793
virtual void getCreationRates(doublereal *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.cpp:459
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition: Kinetics.h:757
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition: Kinetics.cpp:648
std::vector< shared_ptr< Reaction > > m_reactions
Vector of Reaction objects represented by this Kinetics manager.
Definition: Kinetics.h:1324
void checkReactionArraySize(size_t ii) const
Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions()...
Definition: Kinetics.cpp:67
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition: Kinetics.h:243
virtual double reactantStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a reactant in reaction i.
Definition: Kinetics.cpp:384
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:40
virtual void getRevRatesOfProgress(doublereal *revROP)
Return the Reverse rates of progress of the reactions.
Definition: Kinetics.cpp:429
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition: Kinetics.h:687
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:373
std::string name() const
Return the name of the phase.
Definition: Phase.cpp:70
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:273
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition: Phase.h:638
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:200
size_t speciesIndex(const std::string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition: Phase.cpp:187
Abstract base class which stores data about a reaction and its rate parameterization so that it can b...
Definition: Reaction.h:33
void checkBalance(const Kinetics &kin) const
Check that the specified reaction is balanced (same number of atoms for each element in the reactants...
Definition: Reaction.cpp:409
shared_ptr< ReactionRate > rate()
Get reaction rate pointer.
Definition: Reaction.h:171
bool reversible
True if the current reaction is reversible. False otherwise.
Definition: Reaction.h:150
virtual std::string type() const
The type of reaction.
Definition: Reaction.cpp:310
Composition products
Product species and stoichiometric coefficients.
Definition: Reaction.h:138
Composition reactants
Reactant species and stoichiometric coefficients.
Definition: Reaction.h:135
AnyMap input
Input data used for specific models.
Definition: Reaction.h:163
bool duplicate
True if the current reaction is marked as duplicate.
Definition: Reaction.h:153
std::string equation() const
The chemical equation for this reaction.
Definition: Reaction.cpp:296
void add(size_t rxn, const std::vector< size_t > &k)
Add a single reaction to the list of reactions that this stoichiometric manager object handles.
void resizeCoeffs(size_t nSpc, size_t nRxn)
Resize the sparse coefficient matrix.
const Eigen::SparseMatrix< double > & stoichCoeffs() const
Return matrix containing stoichiometric coefficients.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
Definition: ThermoPhase.h:521
virtual std::string type() const
String indicating the thermodynamic model implemented.
Definition: ThermoPhase.h:113
A class for managing third-body efficiencies, including default values.
Definition: Reaction.h:238
double efficiency(const std::string &k) const
Get the third-body efficiency for species k
Definition: Reaction.cpp:632
A reaction with a non-reacting third body "M" that acts to add or remove energy from the reacting spe...
Definition: Reaction.h:269
A reaction with a non-reacting third body "M" that acts to add or remove energy from the reacting spe...
Definition: Reaction.h:475
This file contains definitions for utility functions and text for modules, inputfiles,...
double checkDuplicateStoich(std::map< int, double > &r1, std::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
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:192
std::map< std::string, double > Composition
Map from string names to doubles.
Definition: ct_defs.h:180
void warn_deprecated(const std::string &source, const AnyBase &node, const std::string &message)
A deprecation warning for syntax in an input file.
Definition: AnyMap.cpp:1901
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...