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(nTotalSpecies(), nReactions(), nPhases());
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 }
62
63 // Use the stoichiometric manager to find deltaH for each reaction.
64 getReactionDelta(m_grt.data(), m_dH.data());
65
66 m_temp = T;
67 m_ROP_ok = false;
68 m_redo_rates = false;
69 }
70
71 // loop over interface MultiRate evaluators for each reaction type
72 for (auto& rates : m_rateHandlers) {
73 bool changed = rates->update(thermo(0), *this);
74 if (changed) {
75 rates->getRateConstants(m_rfn.data());
76 m_ROP_ok = false;
77 m_redo_rates = true;
78 }
79 }
80
81 if (!m_ROP_ok) {
82 updateKc();
83 }
84}
85
87{
88 // Store electric potentials for each phase in the array m_phi[].
89 for (size_t n = 0; n < nPhases(); n++) {
90 if (thermo(n).electricPotential() != m_phi[n]) {
92 m_redo_rates = true;
93 }
94 }
95}
96
98{
99 for (size_t n = 0; n < nPhases(); n++) {
100 const auto& tp = thermo(n);
101 /*
102 * We call the getActivityConcentrations function of each ThermoPhase
103 * class that makes up this kinetics object to obtain the generalized
104 * concentrations for species within that class. This is collected in
105 * the vector m_conc. m_start[] are integer indices for that vector
106 * denoting the start of the species for each phase.
107 */
108 tp.getActivityConcentrations(m_actConc.data() + m_start[n]);
109
110 // Get regular concentrations too
111 tp.getConcentrations(m_conc.data() + m_start[n]);
112 }
113 m_ROP_ok = false;
114}
115
117{
118 fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
119
120 if (m_revindex.size() > 0) {
121 /*
122 * Get the vector of standard state electrochemical potentials for
123 * species in the Interfacial kinetics object and store it in m_mu0[]
124 * and m_mu0_Kc[]
125 */
126 updateMu0();
127 double rrt = 1.0 / thermo(0).RT();
128
129 // compute Delta mu^0 for all reversible reactions
130 getRevReactionDelta(m_mu0_Kc.data(), m_rkcn.data());
131
132 for (size_t i = 0; i < m_revindex.size(); i++) {
133 size_t irxn = m_revindex[i];
134 if (irxn == npos || irxn >= nReactions()) {
135 throw CanteraError("InterfaceKinetics::updateKc",
136 "illegal value: irxn = {}", irxn);
137 }
138 // WARNING this may overflow HKM
139 m_rkcn[irxn] = exp(m_rkcn[irxn]*rrt);
140 }
141 for (size_t i = 0; i != m_irrev.size(); ++i) {
142 m_rkcn[ m_irrev[i] ] = 0.0;
143 }
144 }
145}
146
148{
149 // First task is update the electrical potentials from the Phases
151
152 // @todo There is significant potential to further simplify calculations
153 // once the old framework is removed
154 size_t ik = 0;
155 for (size_t n = 0; n < nPhases(); n++) {
157 for (size_t k = 0; k < thermo(n).nSpecies(); k++) {
158 m_mu0_Kc[ik] = m_mu0[ik] + Faraday * m_phi[n] * thermo(n).charge(k);
159 m_mu0_Kc[ik] -= thermo(0).RT() * thermo(n).logStandardConc(k);
160 ik++;
161 }
162 }
163}
164
166{
167 updateMu0();
168 double rrt = 1.0 / thermo(0).RT();
169 std::fill(kc, kc + nReactions(), 0.0);
170 getReactionDelta(m_mu0_Kc.data(), kc);
171 for (size_t i = 0; i < nReactions(); i++) {
172 kc[i] = exp(-kc[i]*rrt);
173 }
174}
175
177{
178 updateROP();
179 for (size_t i = 0; i < nReactions(); i++) {
180 // base rate coefficient multiplied by perturbation factor
181 kfwd[i] = m_rfn[i] * m_perturb[i];
182 }
183}
184
185void InterfaceKinetics::getRevRateConstants(double* krev, bool doIrreversible)
186{
188 if (doIrreversible) {
190 for (size_t i = 0; i < nReactions(); i++) {
191 krev[i] /= m_ropnet[i];
192 }
193 } else {
194 for (size_t i = 0; i < nReactions(); i++) {
195 krev[i] *= m_rkcn[i];
196 }
197 }
198}
199
201{
202 // evaluate rate constants and equilibrium constants at temperature and phi
203 // (electric potential)
205 // get updated activities (rates updated below)
207
208 if (m_ROP_ok) {
209 return;
210 }
211
212 for (size_t i = 0; i < nReactions(); i++) {
213 // Scale the base forward rate coefficient by the perturbation factor
214 m_ropf[i] = m_rfn[i] * m_perturb[i];
215 // Multiply the scaled forward rate coefficient by the reciprocal of the
216 // equilibrium constant
217 m_ropr[i] = m_ropf[i] * m_rkcn[i];
218 }
219
220 for (auto& rates : m_rateHandlers) {
221 rates->modifyRateConstants(m_ropf.data(), m_ropr.data());
222 }
223
224 // multiply ropf by the activity concentration reaction orders to obtain
225 // the forward rates of progress.
226 m_reactantStoich.multiply(m_actConc.data(), m_ropf.data());
227
228 // For reversible reactions, multiply ropr by the activity concentration
229 // products
230 m_revProductStoich.multiply(m_actConc.data(), m_ropr.data());
231
232 for (size_t j = 0; j != nReactions(); ++j) {
233 m_ropnet[j] = m_ropf[j] - m_ropr[j];
234 }
235
236 // For reactions involving multiple phases, we must check that the phase
237 // being consumed actually exists. This is particularly important for phases
238 // that are stoichiometric phases containing one species with a unity
239 // activity
240 if (m_phaseExistsCheck) {
241 for (size_t j = 0; j != nReactions(); ++j) {
242 if ((m_ropr[j] > m_ropf[j]) && (m_ropr[j] > 0.0)) {
243 for (size_t p = 0; p < nPhases(); p++) {
244 if (m_rxnPhaseIsProduct[j][p] && !m_phaseExists[p]) {
245 m_ropnet[j] = 0.0;
246 m_ropr[j] = m_ropf[j];
247 if (m_ropf[j] > 0.0) {
248 for (size_t rp = 0; rp < nPhases(); rp++) {
249 if (m_rxnPhaseIsReactant[j][rp] && !m_phaseExists[rp]) {
250 m_ropnet[j] = 0.0;
251 m_ropr[j] = m_ropf[j] = 0.0;
252 }
253 }
254 }
255 }
256 if (m_rxnPhaseIsReactant[j][p] && !m_phaseIsStable[p]) {
257 m_ropnet[j] = 0.0;
258 m_ropr[j] = m_ropf[j];
259 }
260 }
261 } else if ((m_ropf[j] > m_ropr[j]) && (m_ropf[j] > 0.0)) {
262 for (size_t p = 0; p < nPhases(); p++) {
263 if (m_rxnPhaseIsReactant[j][p] && !m_phaseExists[p]) {
264 m_ropnet[j] = 0.0;
265 m_ropf[j] = m_ropr[j];
266 if (m_ropf[j] > 0.0) {
267 for (size_t rp = 0; rp < nPhases(); rp++) {
268 if (m_rxnPhaseIsProduct[j][rp] && !m_phaseExists[rp]) {
269 m_ropnet[j] = 0.0;
270 m_ropf[j] = m_ropr[j] = 0.0;
271 }
272 }
273 }
274 }
275 if (m_rxnPhaseIsProduct[j][p] && !m_phaseIsStable[p]) {
276 m_ropnet[j] = 0.0;
277 m_ropf[j] = m_ropr[j];
278 }
279 }
280 }
281 }
282 }
283 m_ROP_ok = true;
284}
285
287{
288 // Get the chemical potentials of the species in the all of the phases used
289 // in the kinetics mechanism
290 for (size_t n = 0; n < nPhases(); n++) {
291 m_thermo[n]->getChemPotentials(m_mu.data() + m_start[n]);
292 }
293
294 // Use the stoichiometric manager to find deltaG for each reaction.
295 getReactionDelta(m_mu.data(), m_rbuf.data());
296 if (deltaG != 0 && (m_rbuf.data() != deltaG)) {
297 for (size_t j = 0; j < nReactions(); ++j) {
298 deltaG[j] = m_rbuf[j];
299 }
300 }
301}
302
304{
305 // Get the chemical potentials of the species
306 for (size_t n = 0; n < nPhases(); n++) {
308 }
309
310 // Use the stoichiometric manager to find deltaG for each reaction.
311 getReactionDelta(m_grt.data(), deltaM);
312}
313
315{
316 // Get the partial molar enthalpy of all species
317 for (size_t n = 0; n < nPhases(); n++) {
319 }
320
321 // Use the stoichiometric manager to find deltaH for each reaction.
322 getReactionDelta(m_grt.data(), deltaH);
323}
324
326{
327 // Get the partial molar entropy of all species in all of the phases
328 for (size_t n = 0; n < nPhases(); n++) {
330 }
331
332 // Use the stoichiometric manager to find deltaS for each reaction.
333 getReactionDelta(m_grt.data(), deltaS);
334}
335
337{
338 // Get the standard state chemical potentials of the species. This is the
339 // array of chemical potentials at unit activity We define these here as the
340 // chemical potentials of the pure species at the temperature and pressure
341 // of the solution.
342 for (size_t n = 0; n < nPhases(); n++) {
344 }
345
346 // Use the stoichiometric manager to find deltaG for each reaction.
347 getReactionDelta(m_mu0.data(), deltaGSS);
348}
349
351{
352 // Get the standard state enthalpies of the species. This is the array of
353 // chemical potentials at unit activity We define these here as the
354 // enthalpies of the pure species at the temperature and pressure of the
355 // solution.
356 for (size_t n = 0; n < nPhases(); n++) {
357 thermo(n).getEnthalpy_RT(m_grt.data() + m_start[n]);
358 }
359 for (size_t k = 0; k < m_kk; k++) {
360 m_grt[k] *= thermo(0).RT();
361 }
362
363 // Use the stoichiometric manager to find deltaH for each reaction.
364 getReactionDelta(m_grt.data(), deltaH);
365}
366
368{
369 // Get the standard state entropy of the species. We define these here as
370 // the entropies of the pure species at the temperature and pressure of the
371 // solution.
372 for (size_t n = 0; n < nPhases(); n++) {
373 thermo(n).getEntropy_R(m_grt.data() + m_start[n]);
374 }
375 for (size_t k = 0; k < m_kk; k++) {
376 m_grt[k] *= GasConstant;
377 }
378
379 // Use the stoichiometric manager to find deltaS for each reaction.
380 getReactionDelta(m_grt.data(), deltaS);
381}
382
383bool InterfaceKinetics::addReaction(shared_ptr<Reaction> r_base, bool resize)
384{
385 if (!m_surf) {
386 init();
387 }
388
389 size_t i = nReactions();
390 shared_ptr<ReactionRate> rate = r_base->rate();
391 if (rate) {
392 rate->setContext(*r_base, *this);
393 }
394
395 bool added = Kinetics::addReaction(r_base, resize);
396 if (!added) {
397 return false;
398 }
399
400 m_rxnPhaseIsReactant.emplace_back(nPhases(), false);
401 m_rxnPhaseIsProduct.emplace_back(nPhases(), false);
402
403 for (const auto& [name, stoich] : r_base->reactants) {
404 size_t k = kineticsSpeciesIndex(name);
405 size_t p = speciesPhaseIndex(k);
406 m_rxnPhaseIsReactant[i][p] = true;
407 }
408 for (const auto& [name, stoich] : r_base->products) {
409 size_t k = kineticsSpeciesIndex(name);
410 size_t p = speciesPhaseIndex(k);
411 m_rxnPhaseIsProduct[i][p] = true;
412 }
413
414 // Set index of rate to number of reaction within kinetics
415 rate->setRateIndex(nReactions() - 1);
416
417 string rtype = rate->subType();
418 if (rtype == "") {
419 rtype = rate->type();
420 }
421
422 // If necessary, add new interface MultiRate evaluator
423 if (m_rateTypes.find(rtype) == m_rateTypes.end()) {
424 m_rateTypes[rtype] = m_rateHandlers.size();
425 m_rateHandlers.push_back(rate->newMultiRate());
426 m_rateHandlers.back()->resize(m_kk, nReactions(), nPhases());
427 }
428
429 // Add reaction rate to evaluator
430 size_t index = m_rateTypes[rtype];
431 m_rateHandlers[index]->add(nReactions() - 1, *rate);
432
433 // Set flag for coverage dependence to true
434 if (rate->compositionDependent()) {
436 }
437
438 // Set flag for electrochemistry to true
439 if (r_base->usesElectrochemistry(*this)) {
441 }
442
443 return true;
444}
445
446void InterfaceKinetics::modifyReaction(size_t i, shared_ptr<Reaction> r_base)
447{
448 Kinetics::modifyReaction(i, r_base);
449
450 shared_ptr<ReactionRate> rate = r_base->rate();
451 rate->setRateIndex(i);
452 rate->setContext(*r_base, *this);
453
454 string rtype = rate->subType();
455 if (rtype == "") {
456 rtype = rate->type();
457 }
458
459 // Ensure that interface MultiRate evaluator is available
460 if (!m_rateTypes.count(rtype)) {
461 throw CanteraError("InterfaceKinetics::modifyReaction",
462 "Interface evaluator not available for type '{}'.", rtype);
463 }
464 // Replace reaction rate evaluator
465 size_t index = m_rateTypes[rate->type()];
466 m_rateHandlers[index]->replace(i, *rate);
467
468 // Invalidate cached data
469 m_redo_rates = true;
470 m_temp += 0.1;
471}
472
473void InterfaceKinetics::setMultiplier(size_t i, double f)
474{
476 m_ROP_ok = false;
477}
478
479void InterfaceKinetics::addThermo(shared_ptr<ThermoPhase> thermo)
480{
482 m_phaseExists.push_back(true);
483 m_phaseIsStable.push_back(true);
484}
485
487{
488 if (thermo(0).nDim() > 2) {
489 throw CanteraError("InterfaceKinetics::init", "no interface phase is present.");
490 }
491}
492
494{
495 size_t kOld = m_kk;
497 if (m_kk != kOld && nReactions()) {
498 throw CanteraError("InterfaceKinetics::resizeSpecies", "Cannot add"
499 " species to InterfaceKinetics after reactions have been added.");
500 }
501 m_actConc.resize(m_kk);
502 m_conc.resize(m_kk);
503 m_mu0.resize(m_kk);
504 m_mu.resize(m_kk);
505 m_mu0_Kc.resize(m_kk);
506 m_grt.resize(m_kk);
507 m_phi.resize(nPhases(), 0.0);
508}
509
511{
512 for (auto& phase : m_thermo) {
513 if (!phase->root()) {
514 throw CanteraError("InterfaceKinetics::buildNetwork",
515 "Phase '{}' is not attached to a Solution.", phase->name());
516 }
517 }
518 vector<shared_ptr<ReactorBase>> reservoirs;
519 for (size_t i = 1; i < nPhases(); i++) {
520 auto r = newReservoir(thermo(i).root(), false);
521 reservoirs.push_back(r);
522 }
523 auto rsurf = newReactorSurface(thermo(0).root(), reservoirs, false);
524 m_integrator = new ReactorNet(rsurf);
525}
526
527void InterfaceKinetics::advanceCoverages(double tstep, double rtol, double atol,
528 double maxStepSize, size_t maxSteps, size_t maxErrTestFails)
529{
530 // Stash the state of adjacent phases, and set their T and P to match the surface
531 vector<vector<double>> savedStates(nPhases());
532 for (size_t i = 1; i < nPhases(); i++) {
533 savedStates[i].resize(thermo(i).partialStateSize());
534 thermo(i).savePartialState(savedStates[i].size(), savedStates[i].data());
535 thermo(i).setState_TP(thermo(0).temperature(), thermo(0).pressure());
536 }
537
538 if (!m_integrator) {
539 buildNetwork();
540 }
541
542 m_integrator->setTolerances(rtol, atol);
543 m_integrator->setMaxTimeStep(maxStepSize);
544 m_integrator->setMaxSteps(maxSteps);
545 m_integrator->setMaxErrTestFails(maxErrTestFails);
547 m_integrator->advance(tstep);
548
549 // Restore adjacent phases to their original states
550 for (size_t i = 1; i < nPhases(); i++) {
551 thermo(i).restorePartialState(savedStates[i].size(), savedStates[i].data());
552 }
553}
554
556{
557 // Stash the state of adjacent phases, and set their T and P to match the surface
558 vector<vector<double>> savedStates(nPhases());
559 for (size_t i = 1; i < nPhases(); i++) {
560 savedStates[i].resize(thermo(i).partialStateSize());
561 thermo(i).savePartialState(savedStates[i].size(), savedStates[i].data());
562 thermo(i).setState_TP(thermo(0).temperature(), thermo(0).pressure());
563 }
564
565 if (!m_integrator) {
566 buildNetwork();
567 }
568
569 m_integrator->setVerbose(loglevel != 0);
570 m_integrator->solveSteady(loglevel);
571
572 // Restore adjacent phases to their original states
573 for (size_t i = 1; i < nPhases(); i++) {
574 thermo(i).restorePartialState(savedStates[i].size(), savedStates[i].data());
575 }
576}
577
578void InterfaceKinetics::setPhaseExistence(const size_t iphase, const int exists)
579{
580 checkPhaseIndex(iphase);
581 if (exists) {
582 if (!m_phaseExists[iphase]) {
585 m_phaseExists[iphase] = true;
586 }
587 m_phaseIsStable[iphase] = true;
588 } else {
589 if (m_phaseExists[iphase]) {
591 m_phaseExists[iphase] = false;
592 }
593 m_phaseIsStable[iphase] = false;
594 }
595}
596
597int InterfaceKinetics::phaseExistence(const size_t iphase) const
598{
599 checkPhaseIndex(iphase);
600 return m_phaseExists[iphase];
601}
602
603int InterfaceKinetics::phaseStability(const size_t iphase) const
604{
605 checkPhaseIndex(iphase);
606 return m_phaseIsStable[iphase];
607}
608
609void InterfaceKinetics::setPhaseStability(const size_t iphase, const int isStable)
610{
611 checkPhaseIndex(iphase);
612 if (isStable) {
613 m_phaseIsStable[iphase] = true;
614 } else {
615 m_phaseIsStable[iphase] = false;
616 }
617}
618
619double InterfaceKinetics::interfaceCurrent(const size_t iphase)
620{
621 vector<double> charges(m_kk, 0.0);
622 vector<double> netProdRates(m_kk, 0.0);
623 double dotProduct = 0.0;
624
625 thermo(iphase).getCharges(charges.data());
626 getNetProductionRates(netProdRates.data());
627
628 for (size_t k = 0; k < thermo(iphase).nSpecies(); k++)
629 {
630 dotProduct += charges[k] * netProdRates[m_start[iphase] + k];
631 }
632
633 return dotProduct * Faraday;
634}
635
637{
638 // check derivatives are valid
639 assertDerivativesValid("InterfaceKinetics::fwdRatesOfProgress_ddCi");
640 // forward reaction rate coefficients
641 vector<double>& rop_rates = m_rbuf0;
642 getFwdRateConstants(rop_rates.data());
644}
645
647{
648 // check derivatives are valid
649 assertDerivativesValid("InterfaceKinetics::revRatesOfProgress_ddCi");
650 // reverse reaction rate coefficients
651 vector<double>& rop_rates = m_rbuf0;
652 getFwdRateConstants(rop_rates.data());
653 applyEquilibriumConstants(rop_rates.data());
655}
656
658{
659 // check derivatives are valid
660 assertDerivativesValid("InterfaceKinetics::netRatesOfProgress_ddCi");
661 // forward reaction rate coefficients
662 vector<double>& rop_rates = m_rbuf0;
663 getFwdRateConstants(rop_rates.data());
664 Eigen::SparseMatrix<double> jac = calculateCompositionDerivatives(m_reactantStoich,
665 rop_rates);
666
667 // reverse reaction rate coefficients
668 applyEquilibriumConstants(rop_rates.data());
670}
671
673{
674 bool force = settings.empty();
675 if (force || settings.hasKey("skip-coverage-dependence")) {
676 m_jac_skip_coverage_dependence = settings.getBool("skip-coverage-dependence",
677 false);
678 }
679 if (force || settings.hasKey("skip-electrochemistry")) {
680 m_jac_skip_electrochemistry = settings.getBool("skip-electrochemistry",
681 false);
682 }
683 if (force || settings.hasKey("rtol-delta")) {
684 m_jac_rtol_delta = settings.getDouble("rtol-delta", 1e-8);
685 }
686}
687
689{
690 settings["skip-coverage-dependence"] = m_jac_skip_coverage_dependence;
691 settings["skip-electrochemistry"] = m_jac_skip_electrochemistry;
692 settings["rtol-delta"] = m_jac_rtol_delta;
693}
694
696 StoichManagerN& stoich, const vector<double>& in)
697{
698 vector<double>& outV = m_rbuf1;
699 // derivatives handled by StoichManagerN
700 copy(in.begin(), in.end(), outV.begin());
701 return stoich.derivatives(m_actConc.data(), outV.data());
702}
703
705{
707 throw NotImplementedError(name, "Coverage-dependent reactions not supported.");
709 throw NotImplementedError(name, "Electrochemical reactions not supported.");
710 }
711}
712
714{
715 // For reverse rates computed from thermochemistry, multiply the forward
716 // rate coefficients by the reciprocals of the equilibrium constants
717 for (size_t i = 0; i < nReactions(); ++i) {
718 rop[i] *= m_rkcn[i];
719 }
720}
721
722}
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.
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...
void getFwdRateConstants(double *kfwd) override
Return the forward rate constants.
vector< bool > m_phaseExists
Vector of booleans indicating whether phases exist or not.
void getDeltaSSGibbs(double *deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
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 solvePseudoSteadyStateProblem(int loglevel=0)
Solve for the pseudo steady-state of the surface problem.
void getDeltaGibbs(double *deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
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.
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 getDeltaSSEntropy(double *deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
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 getDeltaSSEnthalpy(double *deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
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.
void getRevRateConstants(double *krev, bool doIrreversible=false) override
Return the reverse rate constants.
void getEquilibriumConstants(double *kc) override
Equilibrium constant for all reactions including the voltage term.
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 applyEquilibriumConstants(double *rop)
Multiply rate with inverse equilibrium constant.
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
void getDeltaEnthalpy(double *deltaH) override
Return the vector of values for the reactions change in enthalpy.
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.
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.
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.
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 getDeltaElectrochemPotentials(double *deltaM) override
Return the vector of values for the reaction electrochemical free energy change.
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 > calculateCompositionDerivatives(StoichManagerN &stoich, const vector< double > &in)
Process mole fraction derivative.
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
void getDeltaEntropy(double *deltaS) override
Return the vector of values for the reactions change in entropy.
vector< vector< bool > > m_rxnPhaseIsProduct
Vector of vector of booleans indicating whether a phase participates in a reaction as a product.
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:1483
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1531
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:402
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
Definition Kinetics.h:1454
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:1508
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1528
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1479
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:393
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1543
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1534
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:1502
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:683
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1539
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1407
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:776
map< string, size_t > m_rateTypes
Mapping of rate handlers.
Definition Kinetics.h:1455
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1537
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:593
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1546
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1471
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:1540
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 setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
Definition Kinetics.h:1376
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1522
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1465
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:671
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:251
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:436
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.
void getCharges(double *charges) const
Copy the vector of species charges into array charges.
Definition Phase.cpp:430
virtual void restorePartialState(size_t lenstate, const double *state)
Set the internal thermodynamic state of the phase, excluding composition.
Definition Phase.cpp:240
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:245
double temperature() const
Temperature (K).
Definition Phase.h:598
virtual void savePartialState(size_t lenstate, double *state) const
Save the current thermodynamic state of the phase, excluding composition.
Definition Phase.cpp:225
virtual size_t partialStateSize() const
Get the size of the partial state vector of the phase.
Definition Phase.h:316
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:574
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(const double *conc, const double *rates)
Calculate derivatives with respect to species concentrations.
virtual void getPartialMolarEnthalpies(double *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
double electricPotential() const
Returns the electric potential of this phase (V).
virtual void getEntropy_R(double *sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
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.
void getElectrochemPotentials(double *mu) const
Get the species electrochemical potentials.
void setElectricPotential(double v)
Set the electric potential of this phase (V).
virtual void getStandardChemPotentials(double *mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
virtual void getEnthalpy_RT(double *hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
virtual void getPartialMolarEntropies(double *sbar) const
Returns an array of partial molar entropies of the species in the solution.
const double Faraday
Faraday constant [C/kmol].
Definition ct_defs.h:133
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:122
shared_ptr< Reservoir > newReservoir(shared_ptr< Solution > phase, bool clone, const string &name)
Create a Reservoir object with the specified contents.
shared_ptr< ReactorSurface > newReactorSurface(shared_ptr< Solution > phase, const vector< shared_ptr< ReactorBase > > &reactors, bool clone, const string &name)
Create a ReactorSurface object with the specified contents and adjacent reactors participating in sur...
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:182
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...