Cantera 2.6.0
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
14
15using namespace std;
16
17namespace Cantera
18{
19
21 m_redo_rates(false),
22 m_surf(0),
23 m_integrator(0),
24 m_ROP_ok(false),
25 m_temp(0.0),
26 m_logtemp(0.0),
27 m_has_coverage_dependence(false),
28 m_has_electrochem_rxns(false),
29 m_has_exchange_current_density_formulation(false),
30 m_phaseExistsCheck(false),
31 m_ioFlag(0),
32 m_nDim(2)
33{
34 if (thermo != 0) {
36 }
37}
38
39InterfaceKinetics::~InterfaceKinetics()
40{
41 delete m_integrator;
42}
43
45{
47
48 for (auto& rates : m_interfaceRates) {
49 rates->resize(nTotalSpecies(), nReactions(), nPhases());
50 // @todo ensure that ReactionData are updated; calling rates->update
51 // blocks correct behavior in InterfaceKinetics::_update_rates_T
52 // and running updateROP() is premature
53 }
54}
55
57{
59 m_redo_rates = true;
60}
61
63{
64 // First task is update the electrical potentials from the Phases
68 m_rates.update_C(m_actConc.data());
69 m_redo_rates = true;
70 }
71
72 // Go find the temperature from the surface
73 doublereal T = thermo(surfacePhaseIndex()).temperature();
74 m_redo_rates = true;
75 if (T != m_temp || m_redo_rates) {
76 m_logtemp = log(T);
77
78 // Calculate the forward rate constant by calling m_rates and store it in m_rfn[]
79 m_rates.update(T, m_logtemp, m_rfn.data());
80 for (size_t n = 0; n < nPhases(); n++) {
82 }
83
84 // Use the stoichiometric manager to find deltaH for each reaction.
85 getReactionDelta(m_grt.data(), m_dH.data());
86 applyStickingCorrection(T, m_rfn.data());
87
88 // If we need to do conversions between exchange current density
89 // formulation and regular formulation (either way) do it here.
92 }
95 }
96 m_temp = T;
97 m_ROP_ok = false;
98 m_redo_rates = false;
99 }
100
101 // loop over interface MultiRate evaluators for each reaction type
102 for (auto& rates : m_interfaceRates) {
103 bool changed = rates->update(thermo(surfacePhaseIndex()), *this);
104 if (changed) {
105 rates->getRateConstants(m_rfn.data());
106 m_ROP_ok = false;
107 m_redo_rates = true;
108 }
109 }
110
111 if (!m_ROP_ok) {
112 updateKc();
113 }
114}
115
117{
118 // Store electric potentials for each phase in the array m_phi[].
119 for (size_t n = 0; n < nPhases(); n++) {
120 if (thermo(n).electricPotential() != m_phi[n]) {
122 m_redo_rates = true;
123 }
124 }
125}
126
128{
129 for (size_t n = 0; n < nPhases(); n++) {
130 const ThermoPhase* tp = m_thermo[n];
131 /*
132 * We call the getActivityConcentrations function of each ThermoPhase
133 * class that makes up this kinetics object to obtain the generalized
134 * concentrations for species within that class. This is collected in
135 * the vector m_conc. m_start[] are integer indices for that vector
136 * denoting the start of the species for each phase.
137 */
139
140 // Get regular concentrations too
141 tp->getConcentrations(m_conc.data() + m_start[n]);
142 }
143 m_ROP_ok = false;
144}
145
147{
149 copy(m_actConc.begin(), m_actConc.end(), conc);
150}
151
153{
154 fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
155
156 if (m_revindex.size() > 0) {
157 /*
158 * Get the vector of standard state electrochemical potentials for
159 * species in the Interfacial kinetics object and store it in m_mu0[]
160 * and m_mu0_Kc[]
161 */
162 updateMu0();
163 doublereal rrt = 1.0 / thermo(reactionPhaseIndex()).RT();
164
165 // compute Delta mu^0 for all reversible reactions
166 getRevReactionDelta(m_mu0_Kc.data(), m_rkcn.data());
167
168 for (size_t i = 0; i < m_revindex.size(); i++) {
169 size_t irxn = m_revindex[i];
170 if (irxn == npos || irxn >= nReactions()) {
171 throw CanteraError("InterfaceKinetics::updateKc",
172 "illegal value: irxn = {}", irxn);
173 }
174 // WARNING this may overflow HKM
175 m_rkcn[irxn] = exp(m_rkcn[irxn]*rrt);
176 }
177 for (size_t i = 0; i != m_irrev.size(); ++i) {
178 m_rkcn[ m_irrev[i] ] = 0.0;
179 }
180 }
181}
182
184{
185 // First task is update the electrical potentials from the Phases
187
188 // @todo There is significant potential to further simplify calculations
189 // once the old framework is removed
191 size_t ik = 0;
192 for (size_t n = 0; n < nPhases(); n++) {
194 for (size_t k = 0; k < thermo(n).nSpecies(); k++) {
195 m_mu0_Kc[ik] = m_mu0[ik] + Faraday * m_phi[n] * thermo(n).charge(k);
197 * thermo(n).logStandardConc(k);
198 ik++;
199 }
200 }
201}
202
204{
205 updateMu0();
206 doublereal rrt = 1.0 / thermo(reactionPhaseIndex()).RT();
207 std::fill(kc, kc + nReactions(), 0.0);
208 getReactionDelta(m_mu0_Kc.data(), kc);
209 for (size_t i = 0; i < nReactions(); i++) {
210 kc[i] = exp(-kc[i]*rrt);
211 }
212}
213
215{
216 // Calculate:
217 // - m_StandardConc[]
218 // - m_ProdStanConcReac[]
219 // - m_deltaG0[]
220 // - m_mu0[]
221
222 // First collect vectors of the standard Gibbs free energies of the
223 // species and the standard concentrations
224 // - m_mu0
225 // - m_StandardConc
226 size_t ik = 0;
227
228 for (size_t n = 0; n < nPhases(); n++) {
230 size_t nsp = thermo(n).nSpecies();
231 for (size_t k = 0; k < nsp; k++) {
233 ik++;
234 }
235 }
236
237 getReactionDelta(m_mu0.data(), m_deltaG0.data());
238
239 // Calculate the product of the standard concentrations of the reactants
240 for (size_t i = 0; i < nReactions(); i++) {
241 m_ProdStanConcReac[i] = 1.0;
242 }
243 m_reactantStoich.multiply(m_StandardConc.data(), m_ProdStanConcReac.data());
244}
245
247{
248 // Compute the electrical potential energy of each species
249 size_t ik = 0;
250 for (size_t n = 0; n < nPhases(); n++) {
251 size_t nsp = thermo(n).nSpecies();
252 for (size_t k = 0; k < nsp; k++) {
253 m_pot[ik] = Faraday * thermo(n).charge(k) * m_phi[n];
254 ik++;
255 }
256 }
257
258 // Compute the change in electrical potential energy for each reaction. This
259 // will only be non-zero if a potential difference is present.
261
262 // Modify the reaction rates. Only modify those with a non-zero activation
263 // energy. Below we decrease the activation energy below zero but in some
264 // debug modes we print out a warning message about this.
265
266 // NOTE, there is some discussion about this point. Should we decrease the
267 // activation energy below zero? I don't think this has been decided in any
268 // definitive way. The treatment below is numerically more stable, however.
269 for (size_t i = 0; i < m_beta.size(); i++) {
270 size_t irxn = m_ctrxn[i];
271
272 // Add the voltage correction to the forward reaction rate constants.
273 double eamod = m_beta[i] * deltaElectricEnergy_[irxn];
274 if (eamod != 0.0) {
275 kf[irxn] *= exp(-eamod/thermo(reactionPhaseIndex()).RT());
276 }
277 }
278}
279
281{
283 // Loop over all reactions which are defined to have a voltage transfer
284 // coefficient that affects the activity energy for the reaction
285 for (size_t i = 0; i < m_ctrxn.size(); i++) {
286 size_t irxn = m_ctrxn[i];
287
288 // Determine whether the reaction rate constant is in an exchange
289 // current density formulation format.
290 int iECDFormulation = m_ctrxn_ecdf[i];
291 if (iECDFormulation) {
292 // We need to have the straight chemical reaction rate constant to
293 // come out of this calculation.
294 double tmp = exp(- m_beta[i] * m_deltaG0[irxn]
295 / thermo(reactionPhaseIndex()).RT());
296 tmp *= 1.0 / m_ProdStanConcReac[irxn] / Faraday;
297 kfwd[irxn] *= tmp;
298 }
299 }
300}
301
303{
304 updateROP();
305 for (size_t i = 0; i < nReactions(); i++) {
306 // base rate coefficient multiplied by perturbation factor
307 kfwd[i] = m_rfn[i] * m_perturb[i];
308 }
309}
310
311void InterfaceKinetics::getRevRateConstants(doublereal* krev, bool doIrreversible)
312{
314 if (doIrreversible) {
316 for (size_t i = 0; i < nReactions(); i++) {
317 krev[i] /= m_ropnet[i];
318 }
319 } else {
320 for (size_t i = 0; i < nReactions(); i++) {
321 krev[i] *= m_rkcn[i];
322 }
323 }
324}
325
327{
328 // evaluate rate constants and equilibrium constants at temperature and phi
329 // (electric potential)
331 // get updated activities (rates updated below)
333
334 if (m_ROP_ok) {
335 return;
336 }
337
338 for (size_t i = 0; i < nReactions(); i++) {
339 // Scale the base forward rate coefficient by the perturbation factor
340 m_ropf[i] = m_rfn[i] * m_perturb[i];
341 // Multiply the scaled forward rate coefficient by the reciprocal of the
342 // equilibrium constant
343 m_ropr[i] = m_ropf[i] * m_rkcn[i];
344 }
345
346 // multiply ropf by the activity concentration reaction orders to obtain
347 // the forward rates of progress.
348 m_reactantStoich.multiply(m_actConc.data(), m_ropf.data());
349
350 // For reversible reactions, multiply ropr by the activity concentration
351 // products
352 m_revProductStoich.multiply(m_actConc.data(), m_ropr.data());
353
354 for (size_t j = 0; j != nReactions(); ++j) {
355 m_ropnet[j] = m_ropf[j] - m_ropr[j];
356 }
357
358 // For reactions involving multiple phases, we must check that the phase
359 // being consumed actually exists. This is particularly important for phases
360 // that are stoichiometric phases containing one species with a unity
361 // activity
362 if (m_phaseExistsCheck) {
363 for (size_t j = 0; j != nReactions(); ++j) {
364 if ((m_ropr[j] > m_ropf[j]) && (m_ropr[j] > 0.0)) {
365 for (size_t p = 0; p < nPhases(); p++) {
366 if (m_rxnPhaseIsProduct[j][p] && !m_phaseExists[p]) {
367 m_ropnet[j] = 0.0;
368 m_ropr[j] = m_ropf[j];
369 if (m_ropf[j] > 0.0) {
370 for (size_t rp = 0; rp < nPhases(); rp++) {
371 if (m_rxnPhaseIsReactant[j][rp] && !m_phaseExists[rp]) {
372 m_ropnet[j] = 0.0;
373 m_ropr[j] = m_ropf[j] = 0.0;
374 }
375 }
376 }
377 }
378 if (m_rxnPhaseIsReactant[j][p] && !m_phaseIsStable[p]) {
379 m_ropnet[j] = 0.0;
380 m_ropr[j] = m_ropf[j];
381 }
382 }
383 } else if ((m_ropf[j] > m_ropr[j]) && (m_ropf[j] > 0.0)) {
384 for (size_t p = 0; p < nPhases(); p++) {
385 if (m_rxnPhaseIsReactant[j][p] && !m_phaseExists[p]) {
386 m_ropnet[j] = 0.0;
387 m_ropf[j] = m_ropr[j];
388 if (m_ropf[j] > 0.0) {
389 for (size_t rp = 0; rp < nPhases(); rp++) {
390 if (m_rxnPhaseIsProduct[j][rp] && !m_phaseExists[rp]) {
391 m_ropnet[j] = 0.0;
392 m_ropf[j] = m_ropr[j] = 0.0;
393 }
394 }
395 }
396 }
397 if (m_rxnPhaseIsProduct[j][p] && !m_phaseIsStable[p]) {
398 m_ropnet[j] = 0.0;
399 m_ropf[j] = m_ropr[j];
400 }
401 }
402 }
403 }
404 }
405 m_ROP_ok = true;
406}
407
408void InterfaceKinetics::getDeltaGibbs(doublereal* deltaG)
409{
410 // Get the chemical potentials of the species in the all of the phases used
411 // in the kinetics mechanism
412 for (size_t n = 0; n < nPhases(); n++) {
413 m_thermo[n]->getChemPotentials(m_mu.data() + m_start[n]);
414 }
415
416 // Use the stoichiometric manager to find deltaG for each reaction.
417 getReactionDelta(m_mu.data(), m_deltaG.data());
418 if (deltaG != 0 && (m_deltaG.data() != deltaG)) {
419 for (size_t j = 0; j < nReactions(); ++j) {
420 deltaG[j] = m_deltaG[j];
421 }
422 }
423}
424
426{
427 // Get the chemical potentials of the species
428 for (size_t n = 0; n < nPhases(); n++) {
430 }
431
432 // Use the stoichiometric manager to find deltaG for each reaction.
433 getReactionDelta(m_grt.data(), deltaM);
434}
435
437{
438 // Get the partial molar enthalpy of all species
439 for (size_t n = 0; n < nPhases(); n++) {
441 }
442
443 // Use the stoichiometric manager to find deltaH for each reaction.
444 getReactionDelta(m_grt.data(), deltaH);
445}
446
448{
449 // Get the partial molar entropy of all species in all of the phases
450 for (size_t n = 0; n < nPhases(); n++) {
452 }
453
454 // Use the stoichiometric manager to find deltaS for each reaction.
455 getReactionDelta(m_grt.data(), deltaS);
456}
457
458void InterfaceKinetics::getDeltaSSGibbs(doublereal* deltaGSS)
459{
460 // Get the standard state chemical potentials of the species. This is the
461 // array of chemical potentials at unit activity We define these here as the
462 // chemical potentials of the pure species at the temperature and pressure
463 // of the solution.
464 for (size_t n = 0; n < nPhases(); n++) {
466 }
467
468 // Use the stoichiometric manager to find deltaG for each reaction.
469 getReactionDelta(m_mu0.data(), deltaGSS);
470}
471
473{
474 // Get the standard state enthalpies of the species. This is the array of
475 // chemical potentials at unit activity We define these here as the
476 // enthalpies of the pure species at the temperature and pressure of the
477 // solution.
478 for (size_t n = 0; n < nPhases(); n++) {
479 thermo(n).getEnthalpy_RT(m_grt.data() + m_start[n]);
480 }
481 for (size_t k = 0; k < m_kk; k++) {
483 }
484
485 // Use the stoichiometric manager to find deltaH for each reaction.
486 getReactionDelta(m_grt.data(), deltaH);
487}
488
490{
491 // Get the standard state entropy of the species. We define these here as
492 // the entropies of the pure species at the temperature and pressure of the
493 // solution.
494 for (size_t n = 0; n < nPhases(); n++) {
495 thermo(n).getEntropy_R(m_grt.data() + m_start[n]);
496 }
497 for (size_t k = 0; k < m_kk; k++) {
498 m_grt[k] *= GasConstant;
499 }
500
501 // Use the stoichiometric manager to find deltaS for each reaction.
502 getReactionDelta(m_grt.data(), deltaS);
503}
504
505bool InterfaceKinetics::addReaction(shared_ptr<Reaction> r_base, bool resize)
506{
507 if (!m_surf) {
508 init();
509 }
510
511 size_t i = nReactions();
512 bool added = Kinetics::addReaction(r_base, resize);
513 if (!added) {
514 return false;
515 }
516
517 if (r_base->reversible) {
518 m_revindex.push_back(i);
519 } else {
520 m_irrev.push_back(i);
521 }
522
523 m_rxnPhaseIsReactant.emplace_back(nPhases(), false);
524 m_rxnPhaseIsProduct.emplace_back(nPhases(), false);
525
526 for (const auto& sp : r_base->reactants) {
527 size_t k = kineticsSpeciesIndex(sp.first);
528 size_t p = speciesPhaseIndex(k);
529 m_rxnPhaseIsReactant[i][p] = true;
530 }
531 for (const auto& sp : r_base->products) {
532 size_t k = kineticsSpeciesIndex(sp.first);
533 size_t p = speciesPhaseIndex(k);
534 m_rxnPhaseIsProduct[i][p] = true;
535 }
536
537 if (!(r_base->usesLegacy())) {
538 // Set index of rate to number of reaction within kinetics
539 shared_ptr<ReactionRate> rate = r_base->rate();
540 rate->setRateIndex(nReactions() - 1);
541 rate->setContext(*r_base, *this);
542
543 // If necessary, add new interface MultiRate evaluator
544 if (m_interfaceTypes.find(rate->type()) == m_interfaceTypes.end()) {
545 m_interfaceTypes[rate->type()] = m_interfaceRates.size();
546 m_interfaceRates.push_back(rate->newMultiRate());
547 m_interfaceRates.back()->resize(m_kk, nReactions(), nPhases());
548 }
549
550 // Add reaction rate to evaluator
551 size_t index = m_interfaceTypes[rate->type()];
552 m_interfaceRates[index]->add(nReactions() - 1, *rate);
553
554 } else if (r_base->reaction_type == SURFACE_RXN) {
555 InterfaceReaction2& r = dynamic_cast<InterfaceReaction2&>(*r_base);
556 SurfaceArrhenius rate = buildSurfaceArrhenius(i, r, false);
557 m_rates.install(i, rate);
558
559 // Turn on the global flag indicating surface coverage dependence
560 if (!r.coverage_deps.empty()) {
562 }
563 ElectrochemicalReaction2* re = dynamic_cast<ElectrochemicalReaction2*>(&r);
564 if (re) {
566 m_beta.push_back(re->beta);
567 m_ctrxn.push_back(i);
568 if (re->exchange_current_density_formulation) {
570 m_ctrxn_ecdf.push_back(1);
571 } else {
572 m_ctrxn_ecdf.push_back(0);
573 }
574 }
575 } else {
576 throw NotImplementedError("InterfaceKinetics::addReaction");
577 }
578 deltaElectricEnergy_.push_back(0.0);
579 m_deltaG0.push_back(0.0);
580 m_deltaG.push_back(0.0);
581 m_ProdStanConcReac.push_back(0.0);
582
583 return true;
584}
585
586void InterfaceKinetics::modifyReaction(size_t i, shared_ptr<Reaction> r_base)
587{
588 Kinetics::modifyReaction(i, r_base);
589 if (!(r_base->usesLegacy())) {
590 shared_ptr<ReactionRate> rate = r_base->rate();
591 rate->setRateIndex(i);
592 rate->setContext(*r_base, *this);
593
594 const auto& rtype = rate->type();
595 // Ensure that interface MultiRate evaluator is available
596 if (!m_interfaceTypes.count(rtype)) {
597 throw CanteraError("InterfaceKinetics::modifyReaction",
598 "Interface evaluator not available for type '{}'.", rtype);
599 }
600 // Replace reaction rate evaluator
601 size_t index = m_interfaceTypes[rate->type()];
602 m_interfaceRates[index]->replace(i, *rate);
603
604 } else if (r_base->reaction_type == SURFACE_RXN) {
605 InterfaceReaction2& r = dynamic_cast<InterfaceReaction2&>(*r_base);
606 SurfaceArrhenius rate = buildSurfaceArrhenius(i, r, true);
607 m_rates.replace(i, rate);
608 } else {
609 throw NotImplementedError("InterfaceKinetics::modifyReaction");
610 }
611 // Invalidate cached data
612 m_redo_rates = true;
613 m_temp += 0.1;
614}
615
617 size_t i, InterfaceReaction2& r, bool replace)
618{
620 // Identify the interface phase
621 size_t iInterface = npos;
622 size_t min_dim = 4;
623 for (size_t n = 0; n < nPhases(); n++) {
624 if (thermo(n).nDim() < min_dim) {
625 iInterface = n;
626 min_dim = thermo(n).nDim();
627 }
628 }
629
630 std::string sticking_species = r.sticking_species;
631 if (sticking_species == "") {
632 // Identify the sticking species if not explicitly given
633 bool foundStick = false;
634 for (const auto& sp : r.reactants) {
635 size_t iPhase = speciesPhaseIndex(kineticsSpeciesIndex(sp.first));
636 if (iPhase != iInterface) {
637 // Non-interface species. There should be exactly one of these
638 if (foundStick) {
639 throw InputFileError("InterfaceKinetics::buildSurfaceArrhenius",
640 r.input, "Multiple non-interface species ('{}' and '{}')\n"
641 "found in sticking reaction: '{}'.\nSticking species "
642 "must be explicitly specified.",
643 sticking_species, sp.first, r.equation());
644 }
645 foundStick = true;
646 sticking_species = sp.first;
647 }
648 }
649 if (!foundStick) {
650 throw InputFileError("InterfaceKinetics::buildSurfaceArrhenius",
651 r.input, "No non-interface species found "
652 "in sticking reaction: '{}'", r.equation());
653 }
654 }
655
656 double surface_order = 0.0;
657 double multiplier = 1.0;
658 // Adjust the A-factor
659 for (const auto& sp : r.reactants) {
660 size_t iPhase = speciesPhaseIndex(kineticsSpeciesIndex(sp.first));
661 const ThermoPhase& p = thermo(iPhase);
662 size_t k = p.speciesIndex(sp.first);
663 if (sp.first == sticking_species) {
664 multiplier *= sqrt(GasConstant/(2*Pi*p.molecularWeight(k)));
665 } else {
666 // Non-sticking species. Convert from coverages used in the
667 // sticking probability expression to the concentration units
668 // used in the mass action rate expression. For surface phases,
669 // the dependence on the site density is incorporated when the
670 // rate constant is evaluated, since we don't assume that the
671 // site density is known at this time.
672 double order = getValue(r.orders, sp.first, sp.second);
673 if (&p == m_surf) {
674 multiplier *= pow(m_surf->size(k), order);
675 surface_order += order;
676 } else {
677 multiplier *= pow(p.standardConcentration(k), -order);
678 }
679 }
680 }
681
682 if (!replace) {
683 m_stickingData.emplace_back(StickData{i, surface_order, multiplier,
685 } else {
686 // Modifying an existing sticking reaction.
687 for (auto& item : m_stickingData) {
688 if (item.index == i) {
689 item.order = surface_order;
690 item.multiplier = multiplier;
691 item.use_motz_wise = r.use_motz_wise_correction;
692 break;
693 }
694 }
695 }
696 }
697
699 r.rate.temperatureExponent(),
700 r.rate.activationEnergy_R());
701
702 // Set up coverage dependencies
703 for (const auto& sp : r.coverage_deps) {
704 size_t k = thermo(reactionPhaseIndex()).speciesIndex(sp.first);
705 rate.addCoverageDependence(k, sp.second.a, sp.second.m, sp.second.E);
706 }
707 return rate;
708}
709
710void InterfaceKinetics::setIOFlag(int ioFlag)
711{
712 m_ioFlag = ioFlag;
713 if (m_integrator) {
714 m_integrator->setIOFlag(ioFlag);
715 }
716}
717
719{
721 m_phaseExists.push_back(true);
722 m_phaseIsStable.push_back(true);
723}
724
726{
727 size_t ks = reactionPhaseIndex();
728 if (ks == npos) {
729 throw CanteraError("InterfaceKinetics::init",
730 "no surface phase is present.");
731 }
732
733 // Check to see that the interface routine has a dimension of 2
734 m_surf = (SurfPhase*)&thermo(ks);
735 if (m_surf->nDim() != m_nDim) {
736 throw CanteraError("InterfaceKinetics::init",
737 "expected interface dimension = 2, but got dimension = {}",
738 m_surf->nDim());
739 }
740}
741
743{
744 size_t kOld = m_kk;
746 if (m_kk != kOld && nReactions()) {
747 throw CanteraError("InterfaceKinetics::resizeSpecies", "Cannot add"
748 " species to InterfaceKinetics after reactions have been added.");
749 }
750 m_actConc.resize(m_kk);
751 m_conc.resize(m_kk);
752 m_StandardConc.resize(m_kk, 0.0);
753 m_mu0.resize(m_kk);
754 m_mu.resize(m_kk);
755 m_mu0_Kc.resize(m_kk);
756 m_grt.resize(m_kk);
757 m_pot.resize(m_kk, 0.0);
758 m_phi.resize(nPhases(), 0.0);
759}
760
761doublereal InterfaceKinetics::electrochem_beta(size_t irxn) const
762{
763 warn_deprecated("InterfaceKinetics::electrochem_beta",
764 "This function only works for the legacy framework. "
765 "To be removed after Cantera 2.6.");
766
767 for (size_t i = 0; i < m_ctrxn.size(); i++) {
768 if (m_ctrxn[i] == irxn) {
769 return m_beta[i];
770 }
771 }
772 return 0.0;
773}
774
775void InterfaceKinetics::advanceCoverages(doublereal tstep, doublereal rtol,
776 doublereal atol, doublereal maxStepSize,
777 size_t maxSteps, size_t maxErrTestFails)
778{
779 if (m_integrator == 0) {
780 vector<InterfaceKinetics*> k{this};
782 }
783 m_integrator->setTolerances(rtol, atol);
784 m_integrator->setMaxStepSize(maxStepSize);
785 m_integrator->setMaxSteps(maxSteps);
786 m_integrator->setMaxErrTestFails(maxErrTestFails);
787 m_integrator->integrate(0.0, tstep);
788 delete m_integrator;
789 m_integrator = 0;
790}
791
793 int ifuncOverride, doublereal timeScaleOverride)
794{
795 // create our own solver object
796 if (m_integrator == 0) {
797 vector<InterfaceKinetics*> k{this};
800 }
801 m_integrator->setIOFlag(m_ioFlag);
802 // New direct method to go here
803 m_integrator->solvePseudoSteadyStateProblem(ifuncOverride, timeScaleOverride);
804}
805
806void InterfaceKinetics::setPhaseExistence(const size_t iphase, const int exists)
807{
808 checkPhaseIndex(iphase);
809 if (exists) {
810 if (!m_phaseExists[iphase]) {
813 m_phaseExists[iphase] = true;
814 }
815 m_phaseIsStable[iphase] = true;
816 } else {
817 if (m_phaseExists[iphase]) {
819 m_phaseExists[iphase] = false;
820 }
821 m_phaseIsStable[iphase] = false;
822 }
823}
824
825int InterfaceKinetics::phaseExistence(const size_t iphase) const
826{
827 checkPhaseIndex(iphase);
828 return m_phaseExists[iphase];
829}
830
831int InterfaceKinetics::phaseStability(const size_t iphase) const
832{
833 checkPhaseIndex(iphase);
834 return m_phaseIsStable[iphase];
835}
836
837void InterfaceKinetics::setPhaseStability(const size_t iphase, const int isStable)
838{
839 checkPhaseIndex(iphase);
840 if (isStable) {
841 m_phaseIsStable[iphase] = true;
842 } else {
843 m_phaseIsStable[iphase] = false;
844 }
845}
846
847void InterfaceKinetics::applyStickingCorrection(double T, double* kf)
848{
849 if (m_stickingData.empty()) {
850 return;
851 }
852
853 static const int cacheId = m_cache.getId();
854 CachedArray cached = m_cache.getArray(cacheId);
855 vector_fp& factors = cached.value;
856
857 double n0 = m_surf->siteDensity();
858 if (!cached.validate(n0)) {
859 factors.resize(m_stickingData.size());
860 for (size_t n = 0; n < m_stickingData.size(); n++) {
861 factors[n] = pow(n0, -m_stickingData[n].order);
862 }
863 }
864
865 for (size_t n = 0; n < m_stickingData.size(); n++) {
866 const StickData& item = m_stickingData[n];
867 if (item.use_motz_wise) {
868 kf[item.index] /= 1 - 0.5 * kf[item.index];
869 }
870 kf[item.index] *= factors[n] * sqrt(T) * item.multiplier;
871 }
872}
873
874}
Declarations for the implicit integration of surface site density equations (see Kinetics Managers an...
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
double activationEnergy_R() const
Return the activation energy divided by the gas constant (that is, the activation temperature) [K].
Definition: RxnRates.h:104
virtual double temperatureExponent() const
Return the temperature exponent b
Definition: Arrhenius.h:117
virtual double preExponentialFactor() const
Return the pre-exponential factor A (in m, kmol, s to powers depending on the reaction order)
Definition: Arrhenius.h:108
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
An interface reaction which involves charged species.
Definition: Reaction.h:458
doublereal beta
Forward value of the apparent Electrochemical transfer coefficient.
Definition: Reaction.h:466
Advances the surface coverages of the associated set of SurfacePhase objects in time.
void integrate(doublereal t0, doublereal t1)
Integrate from t0 to t1. The integrator is reinitialized first.
virtual void setTolerances(double rtol=1.e-7, double atol=1.e-14)
virtual void setMaxStepSize(double maxstep=0.0)
virtual void setMaxSteps(size_t maxsteps=20000)
virtual void initialize(doublereal t0=0.0)
void solvePseudoSteadyStateProblem(int ifuncOverride=-1, doublereal timeScaleOverride=1.0)
Solve for the pseudo steady-state of the surface problem.
virtual void setMaxErrTestFails(size_t maxErrTestFails=7)
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition: AnyMap.h:702
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
virtual void getDeltaSSGibbs(doublereal *deltaG)
Return the vector of values for the reaction standard state Gibbs free energy change.
std::vector< size_t > m_revindex
List of reactions numbers which are reversible reactions.
int phaseExistence(const size_t iphase) const
Gets the phase existence int for the ith phase.
SurfPhase * m_surf
Pointer to the single surface phase.
vector_fp m_beta
Electrochemical transfer coefficient for the forward direction.
virtual void getRevRateConstants(doublereal *krev, bool doIrreversible=false)
Return the reverse rate constants.
void applyVoltageKfwdCorrection(doublereal *const kfwd)
Apply modifications for the forward reaction rate for interfacial charge transfer reactions.
size_t m_nDim
Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for EdgeKinetics)
vector_fp m_deltaG0
Vector of delta G^0, the standard state Gibbs free energies for each reaction.
void setPhaseStability(const size_t iphase, const int isStable)
Set the stability of a phase in the reaction object.
vector_fp m_StandardConc
Vector of standard concentrations.
SurfaceArrhenius buildSurfaceArrhenius(size_t i, InterfaceReaction2 &r, bool replace)
Build a SurfaceArrhenius object from a Reaction, taking into account the possible sticking coefficien...
vector_fp m_deltaG
Vector of deltaG[] of reaction, the delta Gibbs free energies for each reaction.
vector_fp m_phi
Vector of phase electric potentials.
virtual void getDeltaElectrochemPotentials(doublereal *deltaM)
Return the vector of values for the reaction electrochemical free energy change.
void convertExchangeCurrentDensityFormulation(doublereal *const kfwd)
When an electrode reaction rate is optionally specified in terms of its exchange current density,...
void _update_rates_phi()
Update properties that depend on the electric potential.
doublereal m_temp
Current temperature of the data.
std::vector< std::vector< bool > > m_rxnPhaseIsReactant
Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant.
virtual void getDeltaEnthalpy(doublereal *deltaH)
Return the vector of values for the reactions change in enthalpy.
bool m_has_electrochem_rxns
Boolean flag indicating whether any reaction in the mechanism has a beta electrochemical parameter.
virtual void getFwdRateConstants(doublereal *kfwd)
Return the forward rate constants.
std::vector< unique_ptr< MultiRateBase > > m_interfaceRates
Vector of rate handlers for interface reactions.
virtual void addPhase(ThermoPhase &thermo)
Add a phase to the kinetics manager object.
std::vector< size_t > m_ctrxn
Vector of reaction indexes specifying the id of the charge transfer reactions in the mechanism.
std::vector< StickData > m_stickingData
Data for sticking reactions.
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
vector_fp m_actConc
Array of activity concentrations for each species in the kinetics object.
void setElectricPotential(int n, doublereal V)
Set the electric potential in the nth phase.
vector_int m_ctrxn_ecdf
Vector of booleans indicating whether the charge transfer reaction rate constant is described by an e...
vector_fp m_conc
Array of concentrations for each species in the kinetics mechanism.
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.
std::vector< size_t > m_irrev
Vector of irreversible reaction numbers.
vector_fp m_mu0
Vector of standard state chemical potentials for all species.
virtual void getActivityConcentrations(doublereal *const conc)
Get the vector of activity concentrations used in the kinetics object.
vector_fp m_pot
Vector of potential energies due to Voltages.
void _update_rates_C()
Update properties that depend on the species mole fractions and/or concentration,.
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
vector_fp m_mu
Vector of chemical potentials for all species.
vector_fp m_grt
Temporary work vector of length m_kk.
vector_fp deltaElectricEnergy_
Storage for the net electric energy change due to reaction.
doublereal m_logtemp
Current log of the temperature.
int phaseStability(const size_t iphase) const
Gets the phase stability int for the ith phase.
ImplicitSurfChem * m_integrator
Pointer to the Implicit surface chemistry object.
virtual void getEquilibriumConstants(doublereal *kc)
Equilibrium constant for all reactions including the voltage term.
virtual void getDeltaSSEntropy(doublereal *deltaS)
Return the vector of values for the change in the standard state entropies for each reaction.
virtual void getDeltaEntropy(doublereal *deltaS)
Return the vector of values for the reactions change in entropy.
void solvePseudoSteadyStateProblem(int ifuncOverride=-1, doublereal timeScaleOverride=1.0)
Solve for the pseudo steady-state of the surface problem.
void _update_rates_T()
Update properties that depend on temperature.
virtual void updateROP()
Internal routine that updates the Rates of Progress of the reactions.
virtual void getDeltaSSEnthalpy(doublereal *deltaH)
Return the vector of values for the change in the standard state enthalpies of reaction.
bool m_has_coverage_dependence
Boolean flag indicating whether any reaction in the mechanism has a coverage dependent forward reacti...
virtual void updateMu0()
Update the standard state chemical potentials and species equilibrium constant entries.
InterfaceKinetics(ThermoPhase *thermo=0)
Constructor.
std::vector< std::vector< bool > > m_rxnPhaseIsProduct
Vector of vector of booleans indicating whether a phase participates in a reaction as a product.
vector_fp m_mu0_Kc
Vector of standard state electrochemical potentials modified by a standard concentration term.
vector_int m_phaseIsStable
Vector of int indicating whether phases are stable or not.
void updateExchangeCurrentQuantities()
values needed to convert from exchange current density to surface reaction rate.
std::vector< bool > m_phaseExists
Vector of booleans indicating whether phases exist or not.
virtual void getDeltaGibbs(doublereal *deltaG)
Return the vector of values for the reaction Gibbs free energy change.
int m_phaseExistsCheck
Int flag to indicate that some phases in the kinetics mechanism are non-existent.
std::map< std::string, size_t > m_interfaceTypes
Rate handler mapping.
vector_fp m_ProdStanConcReac
Vector of the products of the standard concentrations of the reactants.
void advanceCoverages(doublereal 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.
doublereal electrochem_beta(size_t irxn) const
Return the charge transfer rxn Beta parameter for the ith reaction.
bool m_has_exchange_current_density_formulation
Boolean flag indicating whether any reaction in the mechanism is described by an exchange current den...
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Rate1< SurfaceArrhenius > m_rates
Templated class containing the vector of reactions for this interface.
A reaction occurring on an interface (for example, a SurfPhase or an EdgePhase)
Definition: Reaction.h:419
bool is_sticking_coefficient
Set to true if rate is a parameterization of the sticking coefficient rather than the forward rate co...
Definition: Reaction.h:442
std::map< std::string, CoverageDependency > coverage_deps
Adjustments to the Arrhenius rate expression dependent on surface species coverages.
Definition: Reaction.h:438
std::string sticking_species
For reactions with multiple non-surface species, the sticking species needs to be explicitly identifi...
Definition: Reaction.h:452
bool use_motz_wise_correction
Set to true if rate is a sticking coefficient which should be translated into a rate coefficient usin...
Definition: Reaction.h:448
size_t reactionPhaseIndex() const
Phase where the reactions occur.
Definition: Kinetics.h:218
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition: Kinetics.cpp:48
void checkPhaseIndex(size_t m) const
Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
Definition: Kinetics.cpp:75
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition: Kinetics.h:231
std::vector< size_t > m_start
m_start is a vector of integers specifying the beginning position for the species vector for the n'th...
Definition: Kinetics.h:1345
ValueCache m_cache
Cache for saved calculations within each Kinetics object.
Definition: Kinetics.h:1259
virtual void getRevReactionDelta(const double *g, double *dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition: Kinetics.cpp:450
vector_fp m_ropnet
Net rate-of-progress for each reaction.
Definition: Kinetics.h:1384
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition: Kinetics.h:1317
doublereal multiplier(size_t i) const
The current value of the multiplier for reaction i.
Definition: Kinetics.h:1206
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition: Kinetics.cpp:441
vector_fp m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition: Kinetics.h:1321
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition: Kinetics.h:171
virtual void addPhase(ThermoPhase &thermo)
Add a phase to the kinetics manager object.
Definition: Kinetics.cpp:616
vector_fp m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition: Kinetics.h:1375
vector_fp m_ropf
Forward rate-of-progress for each reaction.
Definition: Kinetics.h:1378
vector_fp m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition: Kinetics.h:1387
vector_fp m_ropr
Reverse rate-of-progress for each reaction.
Definition: Kinetics.h:1381
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition: Kinetics.cpp:660
std::vector< ThermoPhase * > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
Definition: Kinetics.h:1339
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition: Kinetics.cpp:748
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition: Kinetics.h:1306
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition: Kinetics.h:139
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:265
vector_fp m_rfn
Forward rate constant for each reaction.
Definition: Kinetics.h:1369
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition: Kinetics.h:1300
size_t surfacePhaseIndex() const
This returns the integer index of the phase which has ThermoPhase type cSurf.
Definition: Kinetics.h:206
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition: Kinetics.cpp:648
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition: Kinetics.h:243
size_t speciesPhaseIndex(size_t k) const
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of ...
Definition: Kinetics.cpp:373
An error indicating that an unimplemented function has been called.
Definition: ctexceptions.h:187
void getConcentrations(double *const c) const
Get the species concentrations (kmol/m^3).
Definition: Phase.cpp:596
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:273
doublereal 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:630
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition: Phase.h:638
doublereal molecularWeight(size_t k) const
Molecular weight of species k.
Definition: Phase.cpp:492
doublereal temperature() const
Temperature (K).
Definition: Phase.h:654
size_t speciesIndex(const std::string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition: Phase.cpp:187
Composition orders
Forward reaction order with respect to specific species.
Definition: Reaction.h:143
Composition reactants
Reactant species and stoichiometric coefficients.
Definition: Reaction.h:135
AnyMap input
Input data used for specific models.
Definition: Reaction.h:163
std::string equation() const
The chemical equation for this reaction.
Definition: Reaction.cpp:296
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition: SurfPhase.h:125
void getCoverages(doublereal *theta) const
Return a vector of surface coverages.
Definition: SurfPhase.cpp:277
virtual double size(size_t k) const
Returns the number of sites occupied by one molecule of species k.
Definition: SurfPhase.h:313
double siteDensity() const
Returns the site density.
Definition: SurfPhase.h:308
An Arrhenius rate with coverage-dependent terms.
Definition: RxnRates.h:135
void addCoverageDependence(size_t k, doublereal a, doublereal m, doublereal e)
Add a coverage dependency for species k, with exponential dependence a, power-law exponent m,...
Definition: RxnRates.cpp:86
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
Definition: ThermoPhase.h:531
void setElectricPotential(doublereal v)
Set the electric potential of this phase (V).
Definition: ThermoPhase.h:317
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
Definition: ThermoPhase.h:521
doublereal RT() const
Return the Gas Constant multiplied by the current temperature.
Definition: ThermoPhase.h:782
virtual void getEnthalpy_RT(doublereal *hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
Definition: ThermoPhase.h:590
virtual doublereal logStandardConc(size_t k=0) const
Natural logarithm of the standard concentration of the kth species.
Definition: ThermoPhase.cpp:72
virtual void getEntropy_R(doublereal *sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
Definition: ThermoPhase.h:600
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
Definition: ThermoPhase.h:404
virtual doublereal standardConcentration(size_t k=0) const
Return the standard concentration for the kth species.
Definition: ThermoPhase.h:425
void getElectrochemPotentials(doublereal *mu) const
Get the species electrochemical potentials.
Definition: ThermoPhase.cpp:93
virtual void getStandardChemPotentials(doublereal *mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
Definition: ThermoPhase.h:580
doublereal electricPotential() const
Returns the electric potential of this phase (V).
Definition: ThermoPhase.h:326
CachedArray getArray(int id)
Get a reference to a CachedValue object representing an array (vector_fp) with the given id.
Definition: ValueCache.h:171
int getId()
Get a unique id for a cached value.
Definition: ValueCache.cpp:21
const double Pi
Pi.
Definition: ct_defs.h:52
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:192
const double Faraday
Faraday constant [C/kmol].
Definition: ct_defs.h:128
const int SURFACE_RXN
A reaction occurring on a surface.
Definition: reaction_defs.h:76
void warn_deprecated(const std::string &source, const AnyBase &node, const std::string &message)
A deprecation warning for syntax in an input file.
Definition: AnyMap.cpp:1901
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition: ct_defs.h:113
const U & getValue(const std::map< T, U > &m, const T &key, const U &default_val)
Const accessor for a value in a std::map.
Definition: utilities.h:184
Values used for converting sticking coefficients into rate constants.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...