Cantera  4.0.0a1
Loading...
Searching...
No Matches
InterfaceKinetics.cpp
Go to the documentation of this file.
1/**
2 * @file InterfaceKinetics.cpp
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
16
17namespace Cantera
18{
19
20// Constructor / destructor definitions required due to forward-declared unique_ptr
21// members.
22InterfaceKinetics::~InterfaceKinetics()
23{
24 delete m_integrator;
25}
26
28{
30
31 // resize buffer
32 m_rbuf0.resize(nReactions());
33 m_rbuf1.resize(nReactions());
34
35 for (auto& rates : m_rateHandlers) {
36 rates->resize(*this);
37 // @todo ensure that ReactionData are updated; calling rates->update
38 // blocks correct behavior in InterfaceKinetics::_update_rates_T
39 // and running updateROP() is premature
40 }
41}
42
44{
46 m_redo_rates = true;
47}
48
50{
51 // First task is update the electrical potentials from the Phases
53
54 // Go find the temperature from the surface
55 double T = thermo(0).temperature();
56 m_redo_rates = true;
57 if (T != m_temp || m_redo_rates) {
58 // Calculate the forward rate constant by calling m_rates and store it in m_rfn[]
59 for (size_t n = 0; n < nPhases(); n++) {
61 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
62 }
63
64 // Use the stoichiometric manager to find deltaH for each reaction.
66
67 m_temp = T;
68 m_ROP_ok = false;
69 m_redo_rates = false;
70 }
71
72 // loop over interface MultiRate evaluators for each reaction type
73 for (auto& rates : m_rateHandlers) {
74 bool changed = rates->update(thermo(0), *this);
75 if (changed) {
76 rates->getRateConstants(m_rfn);
77 m_ROP_ok = false;
78 m_redo_rates = true;
79 }
80 }
81
82 if (!m_ROP_ok) {
83 updateKc();
84 }
85}
86
88{
89 // Store electric potentials for each phase in the array m_phi[].
90 for (size_t n = 0; n < nPhases(); n++) {
91 if (thermo(n).electricPotential() != m_phi[n]) {
93 m_redo_rates = true;
94 }
95 }
96}
97
99{
100 for (size_t n = 0; n < nPhases(); n++) {
101 const auto& tp = thermo(n);
102 /*
103 * We call the getActivityConcentrations function of each ThermoPhase
104 * class that makes up this kinetics object to obtain the generalized
105 * concentrations for species within that class. This is collected in
106 * the vector m_conc. m_start[] are integer indices for that vector
107 * denoting the start of the species for each phase.
108 */
109 tp.getActivityConcentrations(
110 span<double>(m_actConc).subspan(m_start[n], tp.nSpecies()));
111
112 // Get regular concentrations too
113 tp.getConcentrations(span<double>(m_conc).subspan(m_start[n], tp.nSpecies()));
114 }
115 m_ROP_ok = false;
116}
117
119{
120 fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
121
122 if (m_revindex.size() > 0) {
123 /*
124 * Get the vector of standard state electrochemical potentials for
125 * species in the Interfacial kinetics object and store it in m_mu0[]
126 * and m_mu0_Kc[]
127 */
128 updateMu0();
129 double rrt = 1.0 / thermo(0).RT();
130
131 // compute Delta mu^0 for all reversible reactions
133
134 for (size_t i = 0; i < m_revindex.size(); i++) {
135 size_t irxn = m_revindex[i];
136 if (irxn == npos || irxn >= nReactions()) {
137 throw CanteraError("InterfaceKinetics::updateKc",
138 "illegal value: irxn = {}", irxn);
139 }
140 // WARNING this may overflow HKM
141 m_rkcn[irxn] = exp(m_rkcn[irxn]*rrt);
142 }
143 for (size_t i = 0; i != m_irrev.size(); ++i) {
144 m_rkcn[ m_irrev[i] ] = 0.0;
145 }
146 }
147}
148
150{
151 // First task is update the electrical potentials from the Phases
153
154 // @todo There is significant potential to further simplify calculations
155 // once the old framework is removed
156 size_t ik = 0;
157 for (size_t n = 0; n < nPhases(); n++) {
159 span<double>(m_mu0).subspan(m_start[n], thermo(n).nSpecies()));
160 for (size_t k = 0; k < thermo(n).nSpecies(); k++) {
161 m_mu0_Kc[ik] = m_mu0[ik] + Faraday * m_phi[n] * thermo(n).charge(k);
162 m_mu0_Kc[ik] -= thermo(0).RT() * thermo(n).logStandardConc(k);
163 ik++;
164 }
165 }
166}
167
169{
170 updateMu0();
171 double rrt = 1.0 / thermo(0).RT();
172 std::fill(kc.begin(), kc.end(), 0.0);
174 for (size_t i = 0; i < nReactions(); i++) {
175 kc[i] = exp(-kc[i]*rrt);
176 }
177}
178
180{
181 checkArraySize("InterfaceKinetics::getFwdRateConstants", kfwd.size(), nReactions());
182 updateROP();
183 for (size_t i = 0; i < nReactions(); i++) {
184 // base rate coefficient multiplied by perturbation factor
185 kfwd[i] = m_rfn[i] * m_perturb[i];
186 }
187}
188
189void InterfaceKinetics::getRevRateConstants(span<double> krev, bool doIrreversible)
190{
192 if (doIrreversible) {
194 for (size_t i = 0; i < nReactions(); i++) {
195 krev[i] /= m_ropnet[i];
196 }
197 } else {
198 for (size_t i = 0; i < nReactions(); i++) {
199 krev[i] *= m_rkcn[i];
200 }
201 }
202}
203
205{
206 // evaluate rate constants and equilibrium constants at temperature and phi
207 // (electric potential)
209 // get updated activities (rates updated below)
211
212 if (m_ROP_ok) {
213 return;
214 }
215
216 for (size_t i = 0; i < nReactions(); i++) {
217 // Scale the base forward rate coefficient by the perturbation factor
218 m_ropf[i] = m_rfn[i] * m_perturb[i];
219 // Multiply the scaled forward rate coefficient by the reciprocal of the
220 // equilibrium constant
221 m_ropr[i] = m_ropf[i] * m_rkcn[i];
222 }
223
224 for (auto& rates : m_rateHandlers) {
225 rates->modifyRateConstants(m_ropf, m_ropr);
226 }
227
228 // multiply ropf by the activity concentration reaction orders to obtain
229 // the forward rates of progress.
231
232 // For reversible reactions, multiply ropr by the activity concentration
233 // products
235
236 for (size_t j = 0; j != nReactions(); ++j) {
237 m_ropnet[j] = m_ropf[j] - m_ropr[j];
238 }
239
240 // For reactions involving multiple phases, we must check that the phase
241 // being consumed actually exists. This is particularly important for phases
242 // that are stoichiometric phases containing one species with a unity
243 // activity
244 if (m_phaseExistsCheck) {
245 for (size_t j = 0; j != nReactions(); ++j) {
246 if ((m_ropr[j] > m_ropf[j]) && (m_ropr[j] > 0.0)) {
247 for (size_t p = 0; p < nPhases(); p++) {
248 if (m_rxnPhaseIsProduct[j][p] && !m_phaseExists[p]) {
249 m_ropnet[j] = 0.0;
250 m_ropr[j] = m_ropf[j];
251 if (m_ropf[j] > 0.0) {
252 for (size_t rp = 0; rp < nPhases(); rp++) {
253 if (m_rxnPhaseIsReactant[j][rp] && !m_phaseExists[rp]) {
254 m_ropnet[j] = 0.0;
255 m_ropr[j] = m_ropf[j] = 0.0;
256 }
257 }
258 }
259 }
260 if (m_rxnPhaseIsReactant[j][p] && !m_phaseIsStable[p]) {
261 m_ropnet[j] = 0.0;
262 m_ropr[j] = m_ropf[j];
263 }
264 }
265 } else if ((m_ropf[j] > m_ropr[j]) && (m_ropf[j] > 0.0)) {
266 for (size_t p = 0; p < nPhases(); p++) {
267 if (m_rxnPhaseIsReactant[j][p] && !m_phaseExists[p]) {
268 m_ropnet[j] = 0.0;
269 m_ropf[j] = m_ropr[j];
270 if (m_ropf[j] > 0.0) {
271 for (size_t rp = 0; rp < nPhases(); rp++) {
272 if (m_rxnPhaseIsProduct[j][rp] && !m_phaseExists[rp]) {
273 m_ropnet[j] = 0.0;
274 m_ropf[j] = m_ropr[j] = 0.0;
275 }
276 }
277 }
278 }
279 if (m_rxnPhaseIsProduct[j][p] && !m_phaseIsStable[p]) {
280 m_ropnet[j] = 0.0;
281 m_ropf[j] = m_ropr[j];
282 }
283 }
284 }
285 }
286 }
287 m_ROP_ok = true;
288}
289
290void InterfaceKinetics::getDeltaGibbs(span<double> deltaG)
291{
292 // Get the chemical potentials of the species in the all of the phases used
293 // in the kinetics mechanism
294 for (size_t n = 0; n < nPhases(); n++) {
295 m_thermo[n]->getChemPotentials(
296 span<double>(m_mu).subspan(m_start[n], thermo(n).nSpecies()));
297 }
298
299 // Use the stoichiometric manager to find deltaG for each reaction.
301 if (deltaG.data() != m_rbuf.data()) {
302 copy(m_rbuf.begin(), m_rbuf.end(), deltaG.begin());
303 }
304}
305
307{
308 // Get the chemical potentials of the species
309 for (size_t n = 0; n < nPhases(); n++) {
311 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
312 }
313
314 // Use the stoichiometric manager to find deltaG for each reaction.
315 getReactionDelta(m_grt, deltaM);
316}
317
318void InterfaceKinetics::getDeltaEnthalpy(span<double> deltaH)
319{
320 // Get the partial molar enthalpy of all species
321 for (size_t n = 0; n < nPhases(); n++) {
323 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
324 }
325
326 // Use the stoichiometric manager to find deltaH for each reaction.
327 getReactionDelta(m_grt, deltaH);
328}
329
330void InterfaceKinetics::getDeltaEntropy(span<double> deltaS)
331{
332 // Get the partial molar entropy of all species in all of the phases
333 for (size_t n = 0; n < nPhases(); n++) {
335 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
336 }
337
338 // Use the stoichiometric manager to find deltaS for each reaction.
339 getReactionDelta(m_grt, deltaS);
340}
341
342void InterfaceKinetics::getDeltaSSGibbs(span<double> deltaGSS)
343{
344 // Get the standard state chemical potentials of the species. This is the
345 // array of chemical potentials at unit activity We define these here as the
346 // chemical potentials of the pure species at the temperature and pressure
347 // of the solution.
348 for (size_t n = 0; n < nPhases(); n++) {
350 span<double>(m_mu0).subspan(m_start[n], thermo(n).nSpecies()));
351 }
352
353 // Use the stoichiometric manager to find deltaG for each reaction.
354 getReactionDelta(m_mu0, deltaGSS);
355}
356
358{
359 // Get the standard state enthalpies of the species. This is the array of
360 // chemical potentials at unit activity We define these here as the
361 // enthalpies of the pure species at the temperature and pressure of the
362 // solution.
363 for (size_t n = 0; n < nPhases(); n++) {
365 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
366 }
367 for (size_t k = 0; k < m_kk; k++) {
368 m_grt[k] *= thermo(0).RT();
369 }
370
371 // Use the stoichiometric manager to find deltaH for each reaction.
372 getReactionDelta(m_grt, deltaH);
373}
374
376{
377 // Get the standard state entropy of the species. We define these here as
378 // the entropies of the pure species at the temperature and pressure of the
379 // solution.
380 for (size_t n = 0; n < nPhases(); n++) {
382 span<double>(m_grt).subspan(m_start[n], thermo(n).nSpecies()));
383 }
384 for (size_t k = 0; k < m_kk; k++) {
385 m_grt[k] *= GasConstant;
386 }
387
388 // Use the stoichiometric manager to find deltaS for each reaction.
389 getReactionDelta(m_grt, deltaS);
390}
391
392bool InterfaceKinetics::addReaction(shared_ptr<Reaction> r_base, bool resize)
393{
394 if (!m_surf) {
395 init();
396 }
397
398 size_t i = nReactions();
399 shared_ptr<ReactionRate> rate = r_base->rate();
400 if (rate) {
401 rate->setContext(*r_base, *this);
402 }
403
404 bool added = Kinetics::addReaction(r_base, resize);
405 if (!added) {
406 return false;
407 }
408
409 m_rxnPhaseIsReactant.emplace_back(nPhases(), false);
410 m_rxnPhaseIsProduct.emplace_back(nPhases(), false);
411
412 for (const auto& [name, stoich] : r_base->reactants) {
413 size_t k = kineticsSpeciesIndex(name);
414 size_t p = speciesPhaseIndex(k);
415 m_rxnPhaseIsReactant[i][p] = true;
416 }
417 for (const auto& [name, stoich] : r_base->products) {
418 size_t k = kineticsSpeciesIndex(name);
419 size_t p = speciesPhaseIndex(k);
420 m_rxnPhaseIsProduct[i][p] = true;
421 }
422
423 // Set index of rate to number of reaction within kinetics
424 rate->setRateIndex(nReactions() - 1);
425
426 string rtype = rate->subType();
427 if (rtype == "") {
428 rtype = rate->type();
429 }
430
431 // If necessary, add new interface MultiRate evaluator
432 if (m_rateTypes.find(rtype) == m_rateTypes.end()) {
433 m_rateTypes[rtype] = m_rateHandlers.size();
434 m_rateHandlers.push_back(rate->newMultiRate());
435 m_rateHandlers.back()->resize(*this);
436 }
437
438 // Add reaction rate to evaluator
439 size_t index = m_rateTypes[rtype];
440 m_rateHandlers[index]->add(nReactions() - 1, *rate);
441
442 // Set flag for coverage dependence to true
443 if (rate->compositionDependent()) {
445 }
446
447 // Set flag for electrochemistry to true
448 if (r_base->usesElectrochemistry(*this)) {
450 }
451
452 return true;
453}
454
455void InterfaceKinetics::modifyReaction(size_t i, shared_ptr<Reaction> r_base)
456{
457 Kinetics::modifyReaction(i, r_base);
458
459 shared_ptr<ReactionRate> rate = r_base->rate();
460 rate->setRateIndex(i);
461 rate->setContext(*r_base, *this);
462
463 string rtype = rate->subType();
464 if (rtype == "") {
465 rtype = rate->type();
466 }
467
468 // Ensure that interface MultiRate evaluator is available
469 if (!m_rateTypes.count(rtype)) {
470 throw CanteraError("InterfaceKinetics::modifyReaction",
471 "Interface evaluator not available for type '{}'.", rtype);
472 }
473 // Replace reaction rate evaluator
474 size_t index = m_rateTypes[rate->type()];
475 m_rateHandlers[index]->replace(i, *rate);
476
477 // Invalidate cached data
478 m_redo_rates = true;
479 m_temp += 0.1;
480}
481
482void InterfaceKinetics::setMultiplier(size_t i, double f)
483{
485 m_ROP_ok = false;
486}
487
488void InterfaceKinetics::addThermo(shared_ptr<ThermoPhase> thermo)
489{
491 m_phaseExists.push_back(true);
492 m_phaseIsStable.push_back(true);
493}
494
496{
497 if (thermo(0).nDim() > 2) {
498 throw CanteraError("InterfaceKinetics::init", "no interface phase is present.");
499 }
500}
501
503{
504 size_t kOld = m_kk;
506 if (m_kk != kOld && nReactions()) {
507 throw CanteraError("InterfaceKinetics::resizeSpecies", "Cannot add"
508 " species to InterfaceKinetics after reactions have been added.");
509 }
510 m_actConc.resize(m_kk);
511 m_conc.resize(m_kk);
512 m_mu0.resize(m_kk);
513 m_mu.resize(m_kk);
514 m_mu0_Kc.resize(m_kk);
515 m_grt.resize(m_kk);
516 m_phi.resize(nPhases(), 0.0);
517}
518
520{
521 for (auto& phase : m_thermo) {
522 if (!phase->root()) {
523 throw CanteraError("InterfaceKinetics::buildNetwork",
524 "Phase '{}' is not attached to a Solution.", phase->name());
525 }
526 }
527 vector<shared_ptr<ReactorBase>> reservoirs;
528 for (size_t i = 1; i < nPhases(); i++) {
529 auto r = newReservoir(thermo(i).root(), false);
530 reservoirs.push_back(r);
531 }
532 auto rsurf = newReactorSurface(thermo(0).root(), reservoirs, false);
533 m_integrator = new ReactorNet(rsurf);
534}
535
536void InterfaceKinetics::advanceCoverages(double tstep, double rtol, double atol,
537 double maxStepSize, size_t maxSteps, size_t maxErrTestFails)
538{
539 // Stash the state of adjacent phases, and set their T and P to match the surface
540 vector<vector<double>> savedStates(nPhases());
541 for (size_t i = 1; i < nPhases(); i++) {
542 savedStates[i].resize(thermo(i).partialStateSize());
543 thermo(i).savePartialState(savedStates[i]);
544 thermo(i).setState_TP(thermo(0).temperature(), thermo(0).pressure());
545 }
546
547 if (!m_integrator) {
548 buildNetwork();
549 }
550
551 m_integrator->setTolerances(rtol, atol);
552 m_integrator->setMaxTimeStep(maxStepSize);
553 m_integrator->setMaxSteps(maxSteps);
554 m_integrator->setMaxErrTestFails(maxErrTestFails);
556 m_integrator->advance(tstep);
557
558 // Restore adjacent phases to their original states
559 for (size_t i = 1; i < nPhases(); i++) {
560 thermo(i).restorePartialState(savedStates[i]);
561 }
562}
563
565{
566 // Stash the state of adjacent phases, and set their T and P to match the surface
567 vector<vector<double>> savedStates(nPhases());
568 for (size_t i = 1; i < nPhases(); i++) {
569 savedStates[i].resize(thermo(i).partialStateSize());
570 thermo(i).savePartialState(savedStates[i]);
571 thermo(i).setState_TP(thermo(0).temperature(), thermo(0).pressure());
572 }
573
574 if (!m_integrator) {
575 buildNetwork();
576 }
577
578 m_integrator->setVerbose(loglevel != 0);
579 m_integrator->solveSteady(loglevel);
580
581 // Restore adjacent phases to their original states
582 for (size_t i = 1; i < nPhases(); i++) {
583 thermo(i).restorePartialState(savedStates[i]);
584 }
585}
586
587void InterfaceKinetics::setPhaseExistence(const size_t iphase, const int exists)
588{
589 checkPhaseIndex(iphase);
590 if (exists) {
591 if (!m_phaseExists[iphase]) {
594 m_phaseExists[iphase] = true;
595 }
596 m_phaseIsStable[iphase] = true;
597 } else {
598 if (m_phaseExists[iphase]) {
600 m_phaseExists[iphase] = false;
601 }
602 m_phaseIsStable[iphase] = false;
603 }
604}
605
606int InterfaceKinetics::phaseExistence(const size_t iphase) const
607{
608 checkPhaseIndex(iphase);
609 return m_phaseExists[iphase];
610}
611
612int InterfaceKinetics::phaseStability(const size_t iphase) const
613{
614 checkPhaseIndex(iphase);
615 return m_phaseIsStable[iphase];
616}
617
618void InterfaceKinetics::setPhaseStability(const size_t iphase, const int isStable)
619{
620 checkPhaseIndex(iphase);
621 if (isStable) {
622 m_phaseIsStable[iphase] = true;
623 } else {
624 m_phaseIsStable[iphase] = false;
625 }
626}
627
628double InterfaceKinetics::interfaceCurrent(const size_t iphase)
629{
630 vector<double> charges(m_kk, 0.0);
631 vector<double> netProdRates(m_kk, 0.0);
632 double dotProduct = 0.0;
633
634 thermo(iphase).getCharges(charges);
635 getNetProductionRates(netProdRates);
636
637 for (size_t k = 0; k < thermo(iphase).nSpecies(); k++)
638 {
639 dotProduct += charges[k] * netProdRates[m_start[iphase] + k];
640 }
641
642 return dotProduct * Faraday;
643}
644
646{
647 // check derivatives are valid
648 assertDerivativesValid("InterfaceKinetics::fwdRatesOfProgress_ddCi");
649 // forward reaction rate coefficients
650 vector<double>& rop_rates = m_rbuf0;
651 getFwdRateConstants(rop_rates);
653}
654
656{
657 // check derivatives are valid
658 assertDerivativesValid("InterfaceKinetics::revRatesOfProgress_ddCi");
659 // reverse reaction rate coefficients
660 vector<double>& rop_rates = m_rbuf0;
661 getFwdRateConstants(rop_rates);
662 applyEquilibriumConstants(rop_rates);
664}
665
667{
668 // check derivatives are valid
669 assertDerivativesValid("InterfaceKinetics::netRatesOfProgress_ddCi");
670 // forward reaction rate coefficients
671 vector<double>& rop_rates = m_rbuf0;
672 getFwdRateConstants(rop_rates);
673 Eigen::SparseMatrix<double> jac = calculateCompositionDerivatives(m_reactantStoich,
674 rop_rates);
675
676 // reverse reaction rate coefficients
677 applyEquilibriumConstants(rop_rates);
679}
680
682{
683 bool force = settings.empty();
684 if (force || settings.hasKey("skip-coverage-dependence")) {
685 m_jac_skip_coverage_dependence = settings.getBool("skip-coverage-dependence",
686 false);
687 }
688 if (force || settings.hasKey("skip-electrochemistry")) {
689 m_jac_skip_electrochemistry = settings.getBool("skip-electrochemistry",
690 false);
691 }
692 if (force || settings.hasKey("rtol-delta")) {
693 m_jac_rtol_delta = settings.getDouble("rtol-delta", 1e-8);
694 }
695}
696
698{
699 settings["skip-coverage-dependence"] = m_jac_skip_coverage_dependence;
700 settings["skip-electrochemistry"] = m_jac_skip_electrochemistry;
701 settings["rtol-delta"] = m_jac_rtol_delta;
702}
703
705 StoichManagerN& stoich, span<const double> in)
706{
707 vector<double>& outV = m_rbuf1;
708 // derivatives handled by StoichManagerN
709 copy(in.begin(), in.end(), outV.begin());
710 return stoich.derivatives(m_actConc, outV);
711}
712
714{
716 throw NotImplementedError(name, "Coverage-dependent reactions not supported.");
718 throw NotImplementedError(name, "Electrochemical reactions not supported.");
719 }
720}
721
723{
724 // For reverse rates computed from thermochemistry, multiply the forward
725 // rate coefficients by the reciprocals of the equilibrium constants
726 for (size_t i = 0; i < nReactions(); ++i) {
727 rop[i] *= m_rkcn[i];
728 }
729}
730
731}
Header file for class ReactorSurface.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
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
Base class for exceptions thrown by Cantera classes.
bool addReaction(shared_ptr< Reaction > r, bool resize=true) override
Add a single reaction to the mechanism.
int phaseExistence(const size_t iphase) const
Gets the phase existence int for the ith phase.
double interfaceCurrent(const size_t iphase)
Gets the interface current for the ith phase.
SurfPhase * m_surf
Pointer to the single surface phase.
vector< vector< bool > > m_rxnPhaseIsReactant
Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant.
void getDeltaElectrochemPotentials(span< double > deltaM) override
Return the vector of values for the reaction electrochemical free energy change.
void applyEquilibriumConstants(span< double > rop)
Multiply rate with inverse equilibrium constant.
vector< int > m_phaseIsStable
Vector of int indicating whether phases are stable or not.
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
vector< bool > m_phaseExists
Vector of booleans indicating whether phases exist or not.
double m_temp
Current temperature of the data.
void setPhaseStability(const size_t iphase, const int isStable)
Set the stability of a phase in the reaction object.
void getDeltaEnthalpy(span< double > deltaH) override
Return the vector of values for the reactions change in enthalpy.
void getEquilibriumConstants(span< double > kc) override
Equilibrium constant for all reactions including the voltage term.
void solvePseudoSteadyStateProblem(int loglevel=0)
Solve for the pseudo steady-state of the surface problem.
bool m_jac_skip_coverage_dependence
A flag used to neglect rate coefficient coverage dependence in derivative formation.
void _update_rates_phi()
Update properties that depend on the electric potential.
vector< double > m_phi
Vector of phase electric potentials.
void assertDerivativesValid(const string &name)
Helper function ensuring that all rate derivatives can be calculated.
double m_jac_rtol_delta
Relative tolerance used in developing numerical portions of specific derivatives.
void getFwdRateConstants(span< double > kfwd) override
Return the forward rate constants.
ReactorNet * m_integrator
Pointer to the Implicit surface chemistry object.
void buildNetwork()
Create a ReactorNet where only the surface species are active, for use with advanceCoverages() and so...
void setMultiplier(size_t i, double f) override
Set the multiplier for reaction i to f.
void advanceCoverages(double tstep, double rtol=1.e-7, double atol=1.e-14, double maxStepSize=0, size_t maxSteps=20000, size_t maxErrTestFails=7)
Advance the surface coverages in time.
void getDeltaSSGibbs(span< double > deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
void resizeSpecies() override
Resize arrays with sizes that depend on the total number of species.
vector< double > m_grt
Temporary work vector of length m_kk.
bool m_has_electrochemistry
A flag stating if the object uses electrochemistry.
void updateKc()
Update the equilibrium constants and stored electrochemical potentials in molar units for all reversi...
void setPhaseExistence(const size_t iphase, const int exists)
Set the existence of a phase in the reaction object.
void updateROP() override
Internal routine that updates the Rates of Progress of the reactions.
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
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.
bool m_jac_skip_electrochemistry
A flag used to neglect electrochemical contributions in derivative formation.
void _update_rates_C()
Update properties that depend on the species mole fractions and/or concentration,.
vector< double > m_rbuf0
Buffers for partial rop results with length nReactions()
int phaseStability(const size_t iphase) const
Gets the phase stability int for the ith phase.
void getDeltaSSEnthalpy(span< double > deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
Eigen::SparseMatrix< double > calculateCompositionDerivatives(StoichManagerN &stoich, span< const double > in)
Process mole fraction derivative.
vector< double > m_mu0_Kc
Vector of standard state electrochemical potentials modified by a standard concentration term.
void _update_rates_T()
Update properties that depend on temperature.
void getDeltaEntropy(span< double > deltaS) override
Return the vector of values for the reactions change in entropy.
void getDeltaSSEntropy(span< double > deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
bool m_has_coverage_dependence
A flag stating if the object has coverage dependent rates.
virtual void updateMu0()
Update the standard state chemical potentials and species equilibrium constant entries.
void addThermo(shared_ptr< ThermoPhase > thermo) override
Add a thermo phase to the kinetics manager object.
vector< double > m_mu
Vector of chemical potentials for all species.
vector< double > m_actConc
Array of activity concentrations for each species in the kinetics object.
vector< double > m_mu0
Vector of standard state chemical potentials for all species.
void setElectricPotential(int n, double V)
Set the electric potential in the nth phase.
void init() override
Prepare the class for the addition of reactions, after all phases have been added.
void resizeReactions() override
Finalize Kinetics object and associated objects.
int m_phaseExistsCheck
Int flag to indicate that some phases in the kinetics mechanism are non-existent.
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
vector< vector< bool > > m_rxnPhaseIsProduct
Vector of vector of booleans indicating whether a phase participates in a reaction as a product.
void getDeltaGibbs(span< double > deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
void setDerivativeSettings(const AnyMap &settings) override
Set/modify derivative settings.
vector< double > m_conc
Array of concentrations for each species in the kinetics mechanism.
shared_ptr< ThermoPhase > phase(size_t n=0) const
Return pointer to phase associated with Kinetics by index.
Definition Kinetics.h:226
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:50
size_t checkPhaseIndex(size_t m) const
Check that the specified phase index is in range.
Definition Kinetics.cpp:67
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:239
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1484
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1532
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
Definition Kinetics.h:1455
vector< size_t > m_start
m_start is a vector of integers specifying the beginning position for the species vector for the n'th...
Definition Kinetics.h:1509
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1529
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1480
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1544
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1535
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition Kinetics.h:189
vector< shared_ptr< ThermoPhase > > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
Definition Kinetics.h:1503
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:704
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1540
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1408
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:797
map< string, size_t > m_rateTypes
Mapping of rate handlers.
Definition Kinetics.h:1456
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1538
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:614
virtual void getRevReactionDelta(span< const double > g, span< double > dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition Kinetics.cpp:408
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1547
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1472
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:1541
virtual void getNetProductionRates(span< double > wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:447
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:1377
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1523
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1466
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:692
size_t speciesPhaseIndex(size_t k) const
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of ...
Definition Kinetics.cpp:354
An error indicating that an unimplemented function has been called.
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:246
double temperature() const
Temperature (K).
Definition Phase.h:585
virtual void restorePartialState(span< const double > state)
Set the internal thermodynamic state of the phase, excluding composition.
Definition Phase.cpp:237
void getCharges(span< double > charges) const
Copy the vector of species charges into array charges.
Definition Phase.cpp:411
virtual size_t partialStateSize() const
Get the size of the partial state vector of the phase.
Definition Phase.h:317
virtual void savePartialState(span< double > state) const
Save the current thermodynamic state of the phase, excluding composition.
Definition Phase.cpp:225
double charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
Definition Phase.h:561
A class representing a network of connected reactors.
Definition ReactorNet.h:30
void advance(double t)
Advance the state of all reactors in the independent variable (time or space).
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration steps the integrator will take before reaching the nex...
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single step.
void solveSteady(int loglevel=0)
Solve directly for the steady-state solution.
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition ReactorNet.h:196
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
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.
double electricPotential() const
Returns the electric potential of this phase (V).
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)
virtual double logStandardConc(size_t k=0) const
Natural logarithm of the standard concentration of the kth species.
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...
void setElectricPotential(double v)
Set the electric potential of this phase (V).
void getElectrochemPotentials(span< double > mu) const
Get the species electrochemical potentials.
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...
const double Faraday
Faraday constant [C/kmol].
Definition ct_defs.h:134
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
shared_ptr< ReactorSurface > newReactorSurface(shared_ptr< Solution > phase, span< shared_ptr< ReactorBase > > reactors, bool clone, const string &name)
Create a ReactorSurface object with the specified contents and adjacent reactors participating in sur...
shared_ptr< Reservoir > newReservoir(shared_ptr< Solution > phase, bool clone, const string &name)
Create a Reservoir object with the specified contents.
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
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...