22InterfaceKinetics::~InterfaceKinetics()
57 if (T !=
m_temp || m_redo_rates) {
59 for (
size_t n = 0; n <
nPhases(); n++) {
73 bool changed = rates->update(
thermo(0), *
this);
75 rates->getRateConstants(
m_rfn.data());
89 for (
size_t n = 0; n <
nPhases(); n++) {
99 for (
size_t n = 0; n <
nPhases(); n++) {
100 const auto& tp =
thermo(n);
132 for (
size_t i = 0; i <
m_revindex.size(); i++) {
136 "illegal value: irxn = {}", irxn);
141 for (
size_t i = 0; i !=
m_irrev.size(); ++i) {
155 for (
size_t n = 0; n <
nPhases(); n++) {
172 kc[i] = exp(-kc[i]*rrt);
188 if (doIrreversible) {
221 rates->modifyRateConstants(
m_ropf.data(),
m_ropr.data());
243 for (
size_t p = 0; p <
nPhases(); p++) {
248 for (
size_t rp = 0; rp <
nPhases(); rp++) {
262 for (
size_t p = 0; p <
nPhases(); p++) {
267 for (
size_t rp = 0; rp <
nPhases(); rp++) {
290 for (
size_t n = 0; n <
nPhases(); n++) {
296 if (deltaG != 0 && (
m_rbuf.data() != deltaG)) {
306 for (
size_t n = 0; n <
nPhases(); n++) {
317 for (
size_t n = 0; n <
nPhases(); n++) {
328 for (
size_t n = 0; n <
nPhases(); n++) {
342 for (
size_t n = 0; n <
nPhases(); n++) {
356 for (
size_t n = 0; n <
nPhases(); n++) {
359 for (
size_t k = 0; k <
m_kk; k++) {
372 for (
size_t n = 0; n <
nPhases(); n++) {
375 for (
size_t k = 0; k <
m_kk; k++) {
390 shared_ptr<ReactionRate> rate = r_base->rate();
392 rate->setContext(*r_base, *
this);
403 for (
const auto& [name, stoich] : r_base->reactants) {
408 for (
const auto& [name, stoich] : r_base->products) {
417 string rtype = rate->subType();
419 rtype = rate->type();
434 if (rate->compositionDependent()) {
439 if (r_base->usesElectrochemistry(*
this)) {
450 shared_ptr<ReactionRate> rate = r_base->rate();
451 rate->setRateIndex(i);
452 rate->setContext(*r_base, *
this);
454 string rtype = rate->subType();
456 rtype = rate->type();
462 "Interface evaluator not available for type '{}'.", rtype);
488 if (
thermo(0).nDim() > 2) {
489 throw CanteraError(
"InterfaceKinetics::init",
"no interface phase is present.");
498 throw CanteraError(
"InterfaceKinetics::resizeSpecies",
"Cannot add"
499 " species to InterfaceKinetics after reactions have been added.");
513 if (!
phase->root()) {
515 "Phase '{}' is not attached to a Solution.",
phase->name());
518 vector<shared_ptr<ReactorBase>> reservoirs;
519 for (
size_t i = 1; i <
nPhases(); i++) {
521 reservoirs.push_back(r);
528 double maxStepSize,
size_t maxSteps,
size_t maxErrTestFails)
531 vector<vector<double>> savedStates(
nPhases());
532 for (
size_t i = 1; i <
nPhases(); i++) {
550 for (
size_t i = 1; i <
nPhases(); i++) {
558 vector<vector<double>> savedStates(
nPhases());
559 for (
size_t i = 1; i <
nPhases(); i++) {
573 for (
size_t i = 1; i <
nPhases(); i++) {
621 vector<double> charges(
m_kk, 0.0);
622 vector<double> netProdRates(
m_kk, 0.0);
623 double dotProduct = 0.0;
630 dotProduct += charges[k] * netProdRates[
m_start[iphase] + k];
641 vector<double>& rop_rates =
m_rbuf0;
651 vector<double>& rop_rates =
m_rbuf0;
662 vector<double>& rop_rates =
m_rbuf0;
674 bool force = settings.
empty();
675 if (force || settings.
hasKey(
"skip-coverage-dependence")) {
679 if (force || settings.
hasKey(
"skip-electrochemistry")) {
683 if (force || settings.
hasKey(
"rtol-delta")) {
698 vector<double>& outV = m_rbuf1;
700 copy(in.begin(), in.end(), outV.begin());
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.
double getDouble(const string &key, double default_) const
If key exists, return it as a double, otherwise return default_.
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
bool empty() const
Return boolean indicating whether AnyMap is empty.
bool getBool(const string &key, bool default_) const
If key exists, return it as a bool, otherwise return default_.
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.
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
size_t checkPhaseIndex(size_t m) const
Check that the specified phase index is in range.
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
vector< double > m_ropf
Forward rate-of-progress for each reaction.
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...
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
vector< size_t > m_start
m_start is a vector of integers specifying the beginning position for the species vector for the n'th...
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
size_t nPhases() const
The number of phases participating in the reaction mechanism.
vector< shared_ptr< ThermoPhase > > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
vector< size_t > m_revindex
Indices of reversible reactions.
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
map< string, size_t > m_rateTypes
Mapping of rate handlers.
vector< double > m_ropnet
Net rate-of-progress for each reaction.
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
size_t nReactions() const
Number of reactions in the reaction mechanism.
vector< size_t > m_irrev
Indices of irreversible reactions.
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
virtual void setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
vector< double > m_rfn
Forward rate constant for each reaction.
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
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 ...
An error indicating that an unimplemented function has been called.
void getCharges(double *charges) const
Copy the vector of species charges into array charges.
virtual void restorePartialState(size_t lenstate, const double *state)
Set the internal thermodynamic state of the phase, excluding composition.
size_t nSpecies() const
Returns the number of species in the phase.
double temperature() const
Temperature (K).
virtual void savePartialState(size_t lenstate, double *state) const
Save the current thermodynamic state of the phase, excluding composition.
virtual size_t partialStateSize() const
Get the size of the partial state vector of the phase.
double charge(size_t k) const
Dimensionless electrical charge of a single molecule of species k The charge is normalized by the the...
A class representing a network of connected reactors.
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.
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].
const double GasConstant
Universal Gas Constant [J/kmol/K].
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.
const size_t npos
index returned by functions to indicate "no position"
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...