Cantera  4.0.0a1
Loading...
Searching...
No Matches
BulkKinetics.cpp
1// This file is part of Cantera. See License.txt in the top-level directory or
2// at https://cantera.org/license.txt for license and copyright information.
3
7#include "cantera/numerics/eigen_dense.h"
8
9namespace Cantera
10{
11
12BulkKinetics::BulkKinetics() {
13 setDerivativeSettings(AnyMap()); // use default settings
14}
15
16bool BulkKinetics::addReaction(shared_ptr<Reaction> r, bool resize)
17{
18 shared_ptr<ReactionRate> rate = r->rate();
19 if (rate) {
20 rate->setContext(*r, *this);
21 }
22
23 bool added = Kinetics::addReaction(r, resize);
24 if (!added) {
25 // undeclared species, etc.
26 return false;
27 }
28 double dn = 0.0;
29 for (const auto& [name, stoich] : r->products) {
30 dn += stoich;
31 }
32 for (const auto& [name, stoich] : r->reactants) {
33 dn -= stoich;
34 }
35
36 m_dn.push_back(dn);
37
38 string rtype = rate->subType();
39 if (rtype == "") {
40 rtype = rate->type();
41 }
42
43 // If necessary, add new MultiRate evaluator
44 if (m_rateTypes.find(rtype) == m_rateTypes.end()) {
45 m_rateTypes[rtype] = m_rateHandlers.size();
46 m_rateHandlers.push_back(rate->newMultiRate());
47 m_rateHandlers.back()->resize(*this);
48 }
49
50 // Set index of rate to number of reaction within kinetics
51 rate->setRateIndex(nReactions() - 1);
52
53 // Add reaction rate to evaluator
54 size_t index = m_rateTypes[rtype];
55 m_rateHandlers[index]->add(nReactions() - 1, *rate);
56
57 // Add reaction to third-body evaluator
58 if (r->thirdBody() != nullptr) {
59 addThirdBody(r);
60 }
61
62 m_concm.push_back(NAN);
63 return true;
64}
65
66void BulkKinetics::addThirdBody(shared_ptr<Reaction> r)
67{
68 map<size_t, double> efficiencies;
69 for (const auto& [name, efficiency] : r->thirdBody()->efficiencies) {
70 size_t k = kineticsSpeciesIndex(name, false);
71 if (k != npos) {
72 efficiencies[k] = efficiency;
73 } else if (!m_skipUndeclaredThirdBodies) {
74 throw CanteraError("BulkKinetics::addThirdBody", "Found third-body"
75 " efficiency for undefined species '{}' while adding reaction '{}'",
76 name, r->equation());
77 } else {
79 }
80 }
81 m_multi_concm.install(nReactions() - 1, efficiencies,
82 r->thirdBody()->default_efficiency,
83 r->thirdBody()->mass_action);
84}
85
86void BulkKinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
87{
88 // operations common to all reaction types
90
91 shared_ptr<ReactionRate> rate = rNew->rate();
92 string rtype = rate->subType();
93 if (rtype == "") {
94 rtype = rate->type();
95 }
96
97 // Ensure that MultiRate evaluator is available
98 if (m_rateTypes.find(rtype) == m_rateTypes.end()) {
99 throw CanteraError("BulkKinetics::modifyReaction",
100 "Evaluator not available for type '{}'.", rtype);
101 }
102
103 // Replace reaction rate to evaluator
104 size_t index = m_rateTypes[rtype];
105 rate->setRateIndex(i);
106 rate->setContext(*rNew, *this);
107 m_rateHandlers[index]->replace(i, *rate);
108 invalidateCache();
109}
110
112{
114 m_act_conc.resize(m_kk);
115 m_phys_conc.resize(m_kk);
116 m_grt.resize(m_kk);
117 for (auto& rates : m_rateHandlers) {
118 rates->resize(*this);
119 }
120}
121
123{
125 m_rbuf0.resize(nReactions());
126 m_rbuf1.resize(nReactions());
127 m_rbuf2.resize(nReactions());
128 m_kf0.resize(nReactions());
129 m_sbuf0.resize(nTotalSpecies());
130 m_state.resize(thermo().stateSize());
132 for (auto& rates : m_rateHandlers) {
133 rates->resize(*this);
134 // @todo ensure that ReactionData are updated; calling rates->update
135 // blocks correct behavior in update_rates_T
136 // and running updateROP() is premature
137 }
138}
139
140void BulkKinetics::setMultiplier(size_t i, double f)
141{
143 m_ROP_ok = false;
144}
145
146void BulkKinetics::invalidateCache()
147{
148 Kinetics::invalidateCache();
149 m_ROP_ok = false;
150}
151
153{
154 checkArraySize("BulkKinetics::getFwdRateConstants", kfwd.size(), nReactions());
155 updateROP();
156 copy(m_rfn.begin(), m_rfn.end(), kfwd.begin());
158 processThirdBodies(kfwd);
159 }
160}
161
163{
164 checkArraySize("BulkKinetics::getEquilibriumConstants", kc.size(), nReactions());
165 updateROP();
166
167 vector<double>& delta_gibbs0 = m_rbuf0;
168 fill(delta_gibbs0.begin(), delta_gibbs0.end(), 0.0);
169
170 // compute Delta G^0 for all reactions
171 getReactionDelta(m_grt, delta_gibbs0);
172
173 double rrt = 1.0 / thermo().RT();
174 double logStandConc = log(thermo().standardConcentration());
175 for (size_t i = 0; i < nReactions(); i++) {
176 kc[i] = exp(-delta_gibbs0[i] * rrt + m_dn[i] * logStandConc);
177 }
178}
179
180void BulkKinetics::getRevRateConstants(span<double> krev, bool doIrreversible)
181{
182 // go get the forward rate constants. -> note, we don't really care about
183 // speed or redundancy in these informational routines.
185
186 if (doIrreversible) {
188 for (size_t i = 0; i < nReactions(); i++) {
189 krev[i] /= m_rbuf0[i];
190 }
191 } else {
192 // m_rkcn[] is zero for irreversible reactions
193 for (size_t i = 0; i < nReactions(); i++) {
194 krev[i] *= m_rkcn[i];
195 }
196 }
197}
198
199void BulkKinetics::getDeltaGibbs(span<double> deltaG)
200{
201 // Get the chemical potentials for each species
202 thermo().getChemPotentials(m_sbuf0);
203 // Use the stoichiometric manager to find deltaG for each reaction.
204 getReactionDelta(m_sbuf0, deltaG);
205}
206
207void BulkKinetics::getDeltaEnthalpy(span<double> deltaH)
208{
209 // Get the partial molar enthalpy for each species
211 // Use the stoichiometric manager to find deltaH for each reaction.
212 getReactionDelta(m_sbuf0, deltaH);
213}
214
215void BulkKinetics::getDeltaEntropy(span<double> deltaS)
216{
217 // Get the partial molar entropy for each species
219 // Use the stoichiometric manager to find deltaS for each reaction.
220 getReactionDelta(m_sbuf0, deltaS);
221}
222
223void BulkKinetics::getDeltaSSGibbs(span<double> deltaG)
224{
225 // Get the standard state chemical potentials of the species. This is the
226 // array of chemical potentials at unit activity. We define these here as
227 // the chemical potentials of the pure species at the temperature and
228 // pressure of the solution.
230 // Use the stoichiometric manager to find deltaG for each reaction.
231 getReactionDelta(m_sbuf0, deltaG);
232}
233
234void BulkKinetics::getDeltaSSEnthalpy(span<double> deltaH)
235{
236 // Get the standard state enthalpies of the species.
237 thermo().getEnthalpy_RT(m_sbuf0);
238 for (size_t k = 0; k < m_kk; k++) {
239 m_sbuf0[k] *= thermo().RT();
240 }
241 // Use the stoichiometric manager to find deltaH for each reaction.
242 getReactionDelta(m_sbuf0, deltaH);
243}
244
245void BulkKinetics::getDeltaSSEntropy(span<double> deltaS)
246{
247 // Get the standard state entropy of the species. We define these here as
248 // the entropies of the pure species at the temperature and pressure of the
249 // solution.
250 thermo().getEntropy_R(m_sbuf0);
251 for (size_t k = 0; k < m_kk; k++) {
252 m_sbuf0[k] *= GasConstant;
253 }
254 // Use the stoichiometric manager to find deltaS for each reaction.
255 getReactionDelta(m_sbuf0, deltaS);
256}
257
259{
260 settings["skip-third-bodies"] = m_jac_skip_third_bodies;
261 settings["skip-falloff"] = m_jac_skip_falloff;
262 settings["skip-nonideal"] = m_jac_skip_nonideal;
263 settings["rtol-delta"] = m_jac_rtol_delta;
264}
265
267{
268 bool force = settings.empty();
269 if (force || settings.hasKey("skip-third-bodies")) {
270 m_jac_skip_third_bodies = settings.getBool("skip-third-bodies", false);
271 }
272 if (force || settings.hasKey("skip-falloff")) {
273 m_jac_skip_falloff = settings.getBool("skip-falloff", false);
274 }
275 if (force || settings.hasKey("skip-nonideal")) {
276 m_jac_skip_nonideal = settings.getBool("skip-nonideal", false);
277 }
278 if (force || settings.hasKey("rtol-delta")) {
279 m_jac_rtol_delta = settings.getDouble("rtol-delta", 1e-8);
280 }
281}
282
284{
285 assertDerivativesValid("BulkKinetics::getFwdRateConstants_ddT");
286 updateROP();
287 process_ddT(m_rfn, dkfwd);
288}
289
291{
292 assertDerivativesValid("BulkKinetics::getFwdRatesOfProgress_ddT");
293 updateROP();
294 process_ddT(m_ropf, drop);
295}
296
298{
299 assertDerivativesValid("BulkKinetics::getRevRatesOfProgress_ddT");
300 updateROP();
301 process_ddT(m_ropr, drop);
302
303 // reverse rop times scaled inverse equilibrium constant derivatives
304 copy(m_ropr.begin(), m_ropr.end(), m_rbuf2.begin());
306 Eigen::Map<Eigen::VectorXd> dRevRop(drop.data(), nReactions());
307 Eigen::Map<Eigen::VectorXd> dRevRop2(m_rbuf2.data(), nReactions());
308 dRevRop += dRevRop2;
309}
310
312{
313 assertDerivativesValid("BulkKinetics::getNetRatesOfProgress_ddT");
314 updateROP();
315 process_ddT(m_ropnet, drop);
316
317 // reverse rop times scaled inverse equilibrium constant derivatives
318 copy(m_ropr.begin(), m_ropr.end(), m_rbuf2.begin());
320 Eigen::Map<Eigen::VectorXd> dNetRop(drop.data(), nReactions());
321 Eigen::Map<Eigen::VectorXd> dNetRop2(m_rbuf2.data(), nReactions());
322 dNetRop -= dNetRop2;
323}
324
326{
327 assertDerivativesValid("BulkKinetics::getFwdRateConstants_ddP");
328 updateROP();
329 process_ddP(m_rfn, dkfwd);
330}
331
333{
334 assertDerivativesValid("BulkKinetics::getFwdRatesOfProgress_ddP");
335 updateROP();
336 process_ddP(m_ropf, drop);
337}
338
340{
341 assertDerivativesValid("BulkKinetics::getRevRatesOfProgress_ddP");
342 updateROP();
343 process_ddP(m_ropr, drop);
344}
345
347{
348 assertDerivativesValid("BulkKinetics::getNetRatesOfProgress_ddP");
349 updateROP();
350 process_ddP(m_ropnet, drop);
351}
352
354{
355 assertDerivativesValid("BulkKinetics::getFwdRateConstants_ddC");
356 updateROP();
357 process_ddC(m_reactantStoich, m_rfn, dkfwd, false);
358}
359
361{
362 assertDerivativesValid("BulkKinetics::getFwdRatesOfProgress_ddC");
363 updateROP();
365}
366
368{
369 assertDerivativesValid("BulkKinetics::getRevRatesOfProgress_ddC");
370 updateROP();
372}
373
375{
376 assertDerivativesValid("BulkKinetics::getNetRatesOfProgress_ddC");
377 updateROP();
379 Eigen::Map<Eigen::VectorXd> dNetRop(drop.data(), nReactions());
380
382 Eigen::Map<Eigen::VectorXd> dNetRop2(m_rbuf2.data(), nReactions());
383 dNetRop -= dNetRop2;
384}
385
386Eigen::SparseMatrix<double> BulkKinetics::fwdRatesOfProgress_ddX()
387{
388 assertDerivativesValid("BulkKinetics::fwdRatesOfProgress_ddX");
389
390 // forward reaction rate coefficients
391 vector<double>& rop_rates = m_rbuf0;
392 getFwdRateConstants(rop_rates);
394}
395
396Eigen::SparseMatrix<double> BulkKinetics::revRatesOfProgress_ddX()
397{
398 assertDerivativesValid("BulkKinetics::revRatesOfProgress_ddX");
399
400 // reverse reaction rate coefficients
401 vector<double>& rop_rates = m_rbuf0;
402 getFwdRateConstants(rop_rates);
403 applyEquilibriumConstants(rop_rates);
405}
406
407Eigen::SparseMatrix<double> BulkKinetics::netRatesOfProgress_ddX()
408{
409 assertDerivativesValid("BulkKinetics::netRatesOfProgress_ddX");
410
411 // forward reaction rate coefficients
412 vector<double>& rop_rates = m_rbuf0;
413 getFwdRateConstants(rop_rates);
415
416 // reverse reaction rate coefficients
417 applyEquilibriumConstants(rop_rates);
419}
420
421Eigen::SparseMatrix<double> BulkKinetics::fwdRatesOfProgress_ddCi()
422{
423 assertDerivativesValid("BulkKinetics::fwdRatesOfProgress_ddCi");
424
425 // forward reaction rate coefficients
426 vector<double>& rop_rates = m_rbuf0;
427 getFwdRateConstants(rop_rates);
428 return calculateCompositionDerivatives(m_reactantStoich, rop_rates, false);
429}
430
431Eigen::SparseMatrix<double> BulkKinetics::revRatesOfProgress_ddCi()
432{
433 assertDerivativesValid("BulkKinetics::revRatesOfProgress_ddCi");
434
435 // reverse reaction rate coefficients
436 vector<double>& rop_rates = m_rbuf0;
437 getFwdRateConstants(rop_rates);
438 applyEquilibriumConstants(rop_rates);
439 return calculateCompositionDerivatives(m_revProductStoich, rop_rates, false);
440}
441
442Eigen::SparseMatrix<double> BulkKinetics::netRatesOfProgress_ddCi()
443{
444 assertDerivativesValid("BulkKinetics::netRatesOfProgress_ddCi");
445
446 // forward reaction rate coefficients
447 vector<double>& rop_rates = m_rbuf0;
448 getFwdRateConstants(rop_rates);
449 auto jac = calculateCompositionDerivatives(m_reactantStoich, rop_rates, false);
450
451 // reverse reaction rate coefficients
452 applyEquilibriumConstants(rop_rates);
453 return jac - calculateCompositionDerivatives(m_revProductStoich, rop_rates, false);
454}
455
456void BulkKinetics::updateROP()
457{
458 static const int cacheId = m_cache.getId();
459 CachedScalar last = m_cache.getScalar(cacheId);
460 double T = thermo().temperature();
461 double rho = thermo().density();
462 int statenum = thermo().stateMFNumber();
463
464 if (last.state1 != T || last.state2 != rho) {
465 // Update properties that are independent of the composition
467 fill(m_delta_gibbs0.begin(), m_delta_gibbs0.end(), 0.0);
468 double logStandConc = log(thermo().standardConcentration());
469
470 // compute Delta G^0 for all reversible reactions
472
473 double rrt = 1.0 / thermo().RT();
474 for (size_t i = 0; i < m_revindex.size(); i++) {
475 size_t irxn = m_revindex[i];
476 m_rkcn[irxn] = std::min(
477 exp(m_delta_gibbs0[irxn] * rrt - m_dn[irxn] * logStandConc), BigNumber);
478 }
479
480 for (size_t i = 0; i != m_irrev.size(); ++i) {
481 m_rkcn[ m_irrev[i] ] = 0.0;
482 }
483 }
484
485 if (!last.validate(T, rho, statenum)) {
486 // Update terms dependent on species concentrations and temperature
489 double ctot = thermo().molarDensity();
490
491 // Third-body objects interacting with MultiRate evaluator
493 m_ROP_ok = false;
494 }
495
496 // loop over MultiRate evaluators for each reaction type
497 for (auto& rates : m_rateHandlers) {
498 bool changed = rates->update(thermo(), *this);
499 if (changed) {
500 rates->getRateConstants(m_kf0);
501 m_ROP_ok = false;
502 }
503 }
504
505 if (m_ROP_ok) {
506 // rates of progress are up-to-date only if both the thermodynamic state
507 // and m_perturb are unchanged
508 return;
509 }
510
511 // Scale the forward rate coefficient by the perturbation factor
512 for (size_t i = 0; i < nReactions(); ++i) {
513 m_rfn[i] = m_kf0[i] * m_perturb[i];
514 }
515
516 copy(m_rfn.begin(), m_rfn.end(), m_ropf.begin());
518 copy(m_ropf.begin(), m_ropf.end(), m_ropr.begin());
519
520 // for reversible reactions, multiply ropr by concentration products
522
523 for (auto& rates : m_rateHandlers) {
524 rates->modifyRateConstants(m_ropf, m_ropr);
525 }
526
527 // multiply ropf and ropr by concentration products
530
531 for (size_t j = 0; j != nReactions(); ++j) {
532 m_ropnet[j] = m_ropf[j] - m_ropr[j];
533 }
534
535 for (size_t i = 0; i < m_rfn.size(); i++) {
536 AssertFinite(m_rfn[i], "BulkKinetics::updateROP",
537 "m_rfn[{}] is not finite.", i);
538 AssertFinite(m_ropf[i], "BulkKinetics::updateROP",
539 "m_ropf[{}] is not finite.", i);
540 AssertFinite(m_ropr[i], "BulkKinetics::updateROP",
541 "m_ropr[{}] is not finite.", i);
542 }
543 m_ROP_ok = true;
544}
545
547{
548 checkArraySize("BulkKinetics::getThirdBodyConcentrations", concm.size(), nReactions());
549 updateROP();
550 std::copy(m_concm.begin(), m_concm.end(), concm.begin());
551}
552
554{
555 // reactions involving third body
556 if (!m_concm.empty()) {
558 }
559}
560
562{
563 // For reverse rates computed from thermochemistry, multiply the forward
564 // rate coefficients by the reciprocals of the equilibrium constants
565 for (size_t i = 0; i < nReactions(); ++i) {
566 rop[i] *= m_rkcn[i];
567 }
568}
569
571{
572 double T = thermo().temperature();
573 double P = thermo().pressure();
574 double rrt = 1. / thermo().RT();
575
576 vector<double>& grt = m_sbuf0;
577 vector<double>& delta_gibbs0 = m_rbuf1;
578 fill(delta_gibbs0.begin(), delta_gibbs0.end(), 0.0);
579
580 // compute perturbed Delta G^0 for all reversible reactions
581 thermo().saveState(m_state);
582 thermo().setState_TP(T * (1. + m_jac_rtol_delta), P);
584 getRevReactionDelta(grt, delta_gibbs0);
585
586 // apply scaling for derivative of inverse equilibrium constant
587 double Tinv = 1. / T;
588 double rrt_dTinv = rrt * Tinv / m_jac_rtol_delta;
589 double rrtt = rrt * Tinv;
590 for (size_t i = 0; i < m_revindex.size(); i++) {
591 size_t irxn = m_revindex[i];
592 double factor = delta_gibbs0[irxn] - m_delta_gibbs0[irxn];
593 factor *= rrt_dTinv;
594 factor += m_dn[irxn] * Tinv - m_delta_gibbs0[irxn] * rrtt;
595 drkcn[irxn] *= factor;
596 }
597
598 for (size_t i = 0; i < m_irrev.size(); ++i) {
599 drkcn[m_irrev[i]] = 0.0;
600 }
601
602 thermo().restoreState(m_state);
603}
604
605void BulkKinetics::process_ddT(span<const double> in, span<double> drop)
606{
607 // apply temperature derivative
608 copy(in.begin(), in.end(), drop.begin());
609 for (auto& rates : m_rateHandlers) {
610 rates->processRateConstants_ddT(drop, m_rfn, m_jac_rtol_delta);
611 }
612}
613
614void BulkKinetics::process_ddP(span<const double> in, span<double> drop)
615{
616 // apply pressure derivative
617 copy(in.begin(), in.end(), drop.begin());
618 for (auto& rates : m_rateHandlers) {
619 rates->processRateConstants_ddP(drop, m_rfn, m_jac_rtol_delta);
620 }
621}
622
623void BulkKinetics::process_ddC(StoichManagerN& stoich, span<const double> in,
624 span<double> drop, bool mass_action)
625{
626 Eigen::Map<Eigen::VectorXd> out(drop.data(), nReactions());
627 out.setZero();
628 double ctot_inv = 1. / thermo().molarDensity();
629
630 // derivatives due to concentrations in law of mass action
631 if (mass_action) {
632 stoich.scale(in, span<double>(out.data(), nReactions()), ctot_inv);
633 }
635 return;
636 }
637
638 // derivatives due to third-body colliders in law of mass action
639 Eigen::Map<Eigen::VectorXd> outM(m_rbuf1.data(), nReactions());
640 if (mass_action) {
641 outM.fill(0.);
642 m_multi_concm.scale(in, asSpan(outM), ctot_inv);
643 out += outM;
644 }
645
646 // derivatives due to reaction rates depending on third-body colliders
647 if (!m_jac_skip_falloff) {
648 m_multi_concm.scaleM(in, asSpan(outM), m_concm, ctot_inv);
649 for (auto& rates : m_rateHandlers) {
650 // processing step assigns zeros to entries not dependent on M
651 rates->processRateConstants_ddM(asSpan(outM), m_rfn, m_jac_rtol_delta);
652 }
653 out += outM;
654 }
655}
656
658 StoichManagerN& stoich, span<const double> in, bool ddX)
659{
660 Eigen::SparseMatrix<double> out;
661 vector<double>& scaled = m_rbuf1;
662 vector<double>& outV = m_rbuf2;
663
664 // convert from concentration to mole fraction output
665 copy(in.begin(), in.end(), scaled.begin());
666 if (ddX) {
667 double ctot = thermo().molarDensity();
668 for (size_t i = 0; i < nReactions(); ++i) {
669 scaled[i] *= ctot;
670 }
671 }
672
673 // derivatives handled by StoichManagerN
674 copy(scaled.begin(), scaled.end(), outV.begin());
675 processThirdBodies(outV);
676 out = stoich.derivatives(m_act_conc, outV);
678 return out;
679 }
680
681 // derivatives due to law of mass action
682 copy(scaled.begin(), scaled.end(), outV.begin());
683 stoich.multiply(m_act_conc, outV);
684
685 // derivatives due to reaction rates depending on third-body colliders
686 if (!m_jac_skip_falloff) {
687 for (auto& rates : m_rateHandlers) {
688 // processing step does not modify entries not dependent on M
689 rates->processRateConstants_ddM(outV, m_rfn, m_jac_rtol_delta, false);
690 }
691 }
692
693 // derivatives handled by ThirdBodyCalc
694 out += m_multi_concm.derivatives(outV);
695 return out;
696}
697
699{
700 if (!m_jac_skip_nonideal && !thermo().isIdeal()) {
701 throw NotImplementedError(name,
702 "Not supported for non-ideal ThermoPhase models.");
703 }
704}
705
706}
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
double getDouble(const string &key, double default_) const
If key exists, return it as a double, otherwise return default_.
Definition AnyMap.cpp:1580
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
Definition AnyMap.cpp:1477
bool empty() const
Return boolean indicating whether AnyMap is empty.
Definition AnyMap.cpp:1468
bool getBool(const string &key, bool default_) const
If key exists, return it as a bool, otherwise return default_.
Definition AnyMap.cpp:1575
bool addReaction(shared_ptr< Reaction > r, bool resize=true) override
Add a single reaction to the mechanism.
Eigen::SparseMatrix< double > netRatesOfProgress_ddX() override
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
void applyEquilibriumConstants(span< double > rop)
Multiply rate with inverse equilibrium constant.
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void getRevRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
void process_ddT(const span< const double > in, span< double > drop)
Process temperature derivative.
Eigen::SparseMatrix< double > calculateCompositionDerivatives(StoichManagerN &stoich, span< const double > in, bool ddX=true)
Process derivatives.
void getDeltaEnthalpy(span< double > deltaH) override
Return the vector of values for the reactions change in enthalpy.
void getEquilibriumConstants(span< double > kc) override
Return a vector of Equilibrium constants.
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX() override
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
void getFwdRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
vector< double > m_phys_conc
Physical concentrations, as calculated by ThermoPhase::getConcentrations.
void getRevRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Eigen::SparseMatrix< double > revRatesOfProgress_ddX() override
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
void getNetRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
void getNetRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
void assertDerivativesValid(const string &name)
Helper function ensuring that all rate derivatives can be calculated.
void getFwdRateConstants(span< double > kfwd) override
Return the forward rate constants.
void process_ddC(StoichManagerN &stoich, span< const double > in, span< double > drop, bool mass_action=true)
Process concentration (molar density) derivative.
void applyEquilibriumConstants_ddT(span< double > drkcn)
Multiply rate with scaled temperature derivatives of the inverse equilibrium constant.
void setMultiplier(size_t i, double f) override
Set the multiplier for reaction i to f.
vector< double > m_dn
Difference between the global reactants order and the global products order.
void process_ddP(const span< const double > in, span< double > drop)
Process pressure derivative.
void getDeltaSSGibbs(span< double > deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
void processThirdBodies(span< double > rop)
Multiply rate with third-body collider concentrations.
bool m_jac_skip_third_bodies
Derivative settings.
void resizeSpecies() override
Resize arrays with sizes that depend on the total number of species.
vector< double > m_grt
Standard chemical potentials for each species.
vector< double > m_act_conc
Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations.
void getFwdRateConstants_ddC(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to molar concentration at constant temp...
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
void getNetRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void modifyReaction(size_t i, shared_ptr< Reaction > rNew) override
Modify the rate expression associated with a reaction.
void getRevRateConstants(span< double > krev, bool doIrreversible=false) override
Return the reverse rate constants.
vector< double > m_rbuf0
Buffers for partial rop results with length nReactions()
void getDeltaSSEnthalpy(span< double > deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
void getFwdRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
void getDeltaEntropy(span< double > deltaS) override
Return the vector of values for the reactions change in entropy.
ThirdBodyCalc m_multi_concm
used with MultiRate evaluator
void getDeltaSSEntropy(span< double > deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
void resizeReactions() override
Finalize Kinetics object and associated objects.
void getFwdRateConstants_ddP(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to pressure at constant temperature,...
void getRevRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
void getThirdBodyConcentrations(span< double > concm) override
Return a vector of values of effective concentrations of third-body collision partners of any reactio...
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
void getFwdRateConstants_ddT(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to temperature at constant pressure,...
void getDeltaGibbs(span< double > deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
vector< double > m_concm
Third body concentrations.
vector< double > m_kf0
Forward rate constants without perturbation.
void setDerivativeSettings(const AnyMap &settings) override
Set/modify derivative settings.
void getFwdRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Base class for exceptions thrown by Cantera classes.
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:50
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
ValueCache m_cache
Cache for saved calculations within each Kinetics object.
Definition Kinetics.h:1438
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1488
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1536
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
Definition Kinetics.h:1459
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1533
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1484
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1539
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:704
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1544
bool m_skipUndeclaredThirdBodies
See skipUndeclaredThirdBodies()
Definition Kinetics.h:1557
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:797
map< string, size_t > m_rateTypes
Mapping of rate handlers.
Definition Kinetics.h:1460
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1542
virtual void getRevReactionDelta(span< const double > g, span< double > dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition Kinetics.cpp:408
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1476
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition Kinetics.h:161
vector< size_t > m_irrev
Indices of irreversible reactions.
Definition Kinetics.h:1545
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1560
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:273
virtual void getReactionDelta(span< const double > property, span< double > deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:396
virtual void setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
Definition Kinetics.h:1381
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1527
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1470
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:692
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:251
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1530
An error indicating that an unimplemented function has been called.
virtual double molarDensity() const
Molar density (kmol/m^3).
Definition Phase.cpp:587
double temperature() const
Temperature (K).
Definition Phase.h:585
virtual void getConcentrations(span< double > c) const
Get the species concentrations (kmol/m^3).
Definition Phase.cpp:491
virtual double density() const
Density (kg/m^3).
Definition Phase.h:610
int stateMFNumber() const
Return the State Mole Fraction Number.
Definition Phase.h:794
virtual void restoreState(span< const double > state)
Restore the state of the phase from a previously saved state vector.
Definition Phase.cpp:269
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Definition Phase.h:603
virtual void saveState(span< double > state) const
Write to array 'state' the current internal state.
Definition Phase.cpp:257
This class handles operations involving the stoichiometric coefficients on one side of a reaction (re...
Eigen::SparseMatrix< double > derivatives(span< const double > conc, span< const double > rates)
Calculate derivatives with respect to species concentrations.
void scale(span< const double > in, span< double > out, double factor) const
Scale input by reaction order and factor.
virtual void getPartialMolarEnthalpies(span< double > hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void setState_TP(double t, double p)
Set the temperature (K) and pressure (Pa)
double RT() const
Return the Gas Constant multiplied by the current temperature.
virtual void getPartialMolarEntropies(span< double > sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void getStandardChemPotentials(span< double > mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
virtual void getActivityConcentrations(span< double > c) const
This method returns an array of generalized concentrations.
virtual void getEntropy_R(span< double > sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
virtual void getEnthalpy_RT(span< double > hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
virtual void getChemPotentials(span< double > mu) const
Get the species chemical potentials. Units: J/kmol.
bool empty() const
Return boolean indicating whether ThirdBodyCalc is empty.
void install(size_t rxnNumber, const map< size_t, double > &efficiencies, double default_efficiency, bool mass_action)
Install reaction that uses third-body effects in ThirdBodyCalc manager.
void update(span< const double > conc, double ctot, span< double > concm) const
Update third-body concentrations in full vector.
Eigen::SparseMatrix< double > derivatives(span< const double > product)
Calculate derivatives with respect to species concentrations.
void multiply(span< double > output, span< const double > concm)
Multiply output with effective third-body concentration.
void scale(span< const double > in, span< double > out, double factor) const
Scale entries involving third-body collider in law of mass action by factor.
void resizeCoeffs(size_t nSpc, size_t nRxn)
Resize the sparse coefficient matrix.
void scaleM(span< const double > in, span< double > out, span< const double > concm, double factor) const
Scale entries involving third-body collider in rate expression by third-body concentration and factor...
CachedScalar getScalar(int id)
Get a reference to a CachedValue object representing a scalar (double) with the given id.
Definition ValueCache.h:161
int getId()
Get a unique id for a cached value.
#define AssertFinite(expr, procedure,...)
Throw an exception if the specified exception is not a finite number.
bool legacy_rate_constants_used()
Returns true if legacy rate constant definition is used.
Definition global.cpp:107
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:183
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
Definition eigen_dense.h:46
const double BigNumber
largest number to compare to inf.
Definition ct_defs.h:163
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
A cached property value and the state at which it was evaluated.
Definition ValueCache.h:33
double state2
Value of the second state variable for the state at which value was evaluated, for example density or...
Definition ValueCache.h:106
bool validate(double state1New)
Check whether the currently cached value is valid based on a single state variable.
Definition ValueCache.h:39
double state1
Value of the first state variable for the state at which value was evaluated, for example temperature...
Definition ValueCache.h:102