Cantera  3.1.0a1
Loading...
Searching...
No Matches
InterfaceKinetics Class Reference

A kinetics manager for heterogeneous reaction mechanisms. More...

#include <InterfaceKinetics.h>

Inheritance diagram for InterfaceKinetics:
[legend]

Detailed Description

A kinetics manager for heterogeneous reaction mechanisms.

The reactions are assumed to occur at a 2D interface between two 3D phases.

There are some important additions to the behavior of the kinetics class due to the presence of multiple phases and a heterogeneous interface. If a reactant phase doesn't exists, that is, has a mole number of zero, a heterogeneous reaction can not proceed from reactants to products. Note it could perhaps proceed from products to reactants if all of the product phases exist.

In order to make the determination of whether a phase exists or not actually involves the specification of additional information to the kinetics object., which heretofore has only had access to intrinsic field information about the phases (for example, temperature, pressure, and mole fraction).

The extrinsic specification of whether a phase exists or not must be specified on top of the intrinsic calculation of the reaction rate. This class carries a set of booleans indicating whether a phase in the heterogeneous mechanism exists or not.

Additionally, the class carries a set of booleans around indicating whether a product phase is stable or not. If a phase is not thermodynamically stable, it may be the case that a particular reaction in a heterogeneous mechanism will create a product species in the unstable phase. However, other reactions in the mechanism will destruct that species. This may cause oscillations in the formation of the unstable phase from time step to time step within a ODE solver, in practice. In order to avoid this situation, a set of booleans is tracked which sets the stability of a phase. If a phase is deemed to be unstable, then species in that phase will not be allowed to be birthed by the kinetics operator. Nonexistent phases are deemed to be unstable by default, but this can be changed.

Definition at line 55 of file InterfaceKinetics.h.

Public Member Functions

 InterfaceKinetics ()=default
 Constructor.
 
void resizeReactions () override
 Finalize Kinetics object and associated objects.
 
string kineticsType () const override
 Identifies the Kinetics manager type.
 
void setElectricPotential (int n, double V)
 Set the electric potential in the nth phase.
 
void updateROP () override
 Internal routine that updates the Rates of Progress of the reactions.
 
void _update_rates_T ()
 Update properties that depend on temperature.
 
void _update_rates_phi ()
 Update properties that depend on the electric potential.
 
void _update_rates_C ()
 Update properties that depend on the species mole fractions and/or concentration,.
 
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 solvePseudoSteadyStateProblem (int ifuncOverride=-1, double timeScaleOverride=1.0)
 Solve for the pseudo steady-state of the surface problem.
 
void setIOFlag (int ioFlag)
 
virtual void updateMu0 ()
 Update the standard state chemical potentials and species equilibrium constant entries.
 
void updateKc ()
 Update the equilibrium constants and stored electrochemical potentials in molar units for all reversible reactions and for all species.
 
void setPhaseExistence (const size_t iphase, const int exists)
 Set the existence of a phase in the reaction object.
 
void setPhaseStability (const size_t iphase, const int isStable)
 Set the stability of a phase in the reaction object.
 
int phaseExistence (const size_t iphase) const
 Gets the phase existence int for the ith phase.
 
int phaseStability (const size_t iphase) const
 Gets the phase stability int for the ith phase.
 
double interfaceCurrent (const size_t iphase)
 Gets the interface current for the ith phase.
 
void setDerivativeSettings (const AnyMap &settings) override
 Set/modify derivative settings.
 
void getDerivativeSettings (AnyMap &settings) const override
 Retrieve derivative settings.
 
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi () override
 Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.
 
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi () override
 Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.
 
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi () override
 Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations.
 
Reaction Rates Of Progress
void getEquilibriumConstants (double *kc) override
 Equilibrium constant for all reactions including the voltage term.
 
void getDeltaGibbs (double *deltaG) override
 Return the vector of values for the reaction Gibbs free energy change.
 
void getDeltaElectrochemPotentials (double *deltaM) override
 Return the vector of values for the reaction electrochemical free energy change.
 
void getDeltaEnthalpy (double *deltaH) override
 Return the vector of values for the reactions change in enthalpy.
 
void getDeltaEntropy (double *deltaS) override
 Return the vector of values for the reactions change in entropy.
 
void getDeltaSSGibbs (double *deltaG) override
 Return the vector of values for the reaction standard state Gibbs free energy change.
 
void getDeltaSSEnthalpy (double *deltaH) override
 Return the vector of values for the change in the standard state enthalpies of reaction.
 
void getDeltaSSEntropy (double *deltaS) override
 Return the vector of values for the change in the standard state entropies for each reaction.
 
Reaction Mechanism Informational Query Routines
void getActivityConcentrations (double *const conc) override
 Get the vector of activity concentrations used in the kinetics object.
 
bool isReversible (size_t i) override
 True if reaction i has been declared to be reversible.
 
void getFwdRateConstants (double *kfwd) override
 Return the forward rate constants.
 
void getRevRateConstants (double *krev, bool doIrreversible=false) override
 Return the reverse rate constants.
 
Reaction Mechanism Construction
void addThermo (shared_ptr< ThermoPhase > thermo) override
 Add a thermo phase to the kinetics manager object.
 
void init () override
 Prepare the class for the addition of reactions, after all phases have been added.
 
void resizeSpecies () override
 Resize arrays with sizes that depend on the total number of species.
 
bool addReaction (shared_ptr< Reaction > r, bool resize=true) override
 Add a single reaction to the mechanism.
 
void modifyReaction (size_t i, shared_ptr< Reaction > rNew) override
 Modify the rate expression associated with a reaction.
 
void setMultiplier (size_t i, double f) override
 Set the multiplier for reaction i to f.
 
- Public Member Functions inherited from Kinetics
virtual pair< size_t, size_t > checkDuplicates (bool throw_err=true) const
 Check for unmarked duplicate reactions and unmatched marked duplicates.
 
virtual void setRoot (shared_ptr< Solution > root)
 Set root Solution holding all phase information.
 
shared_ptr< Solutionroot () const
 Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
 
 Kinetics ()=default
 Default constructor.
 
 Kinetics (const Kinetics &)=delete
 Kinetics objects are not copyable or assignable.
 
Kineticsoperator= (const Kinetics &)=delete
 
size_t nReactions () const
 Number of reactions in the reaction mechanism.
 
void checkReactionIndex (size_t m) const
 Check that the specified reaction index is in range Throws an exception if i is greater than nReactions()
 
void checkReactionArraySize (size_t ii) const
 Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions().
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range Throws an exception if k is greater than nSpecies()-1.
 
void checkSpeciesArraySize (size_t mm) const
 Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies().
 
size_t nPhases () const
 The number of phases participating in the reaction mechanism.
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases().
 
size_t phaseIndex (const string &ph) const
 Return the phase index of a phase in the list of phases defined within the object.
 
size_t reactionPhaseIndex () const
 Phase where the reactions occur.
 
shared_ptr< ThermoPhasereactionPhase () const
 Return pointer to phase where the reactions occur.
 
ThermoPhasethermo (size_t n=0)
 This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
 
const ThermoPhasethermo (size_t n=0) const
 
size_t nTotalSpecies () const
 The total number of species in all phases participating in the kinetics mechanism.
 
size_t kineticsSpeciesIndex (size_t k, size_t n) const
 The location of species k of phase n in species arrays.
 
string kineticsSpeciesName (size_t k) const
 Return the name of the kth species in the kinetics manager.
 
size_t kineticsSpeciesIndex (const string &nm) const
 This routine will look up a species number based on the input string nm.
 
ThermoPhasespeciesPhase (const string &nm)
 This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides.
 
const ThermoPhasespeciesPhase (const string &nm) const
 
ThermoPhasespeciesPhase (size_t k)
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the species' owning ThermoPhase object.
 
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 species in the kinetics manager) and returns the index of the phase owning the species.
 
virtual void getFwdRatesOfProgress (double *fwdROP)
 Return the forward rates of progress of the reactions.
 
virtual void getRevRatesOfProgress (double *revROP)
 Return the Reverse rates of progress of the reactions.
 
virtual void getNetRatesOfProgress (double *netROP)
 Net rates of progress.
 
virtual void getReactionDelta (const double *property, double *deltaProperty) const
 Change in species properties.
 
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 reversible reactions.
 
virtual void getThirdBodyConcentrations (double *concm)
 Return a vector of values of effective concentrations of third-body collision partners of any reaction.
 
virtual const vector< double > & thirdBodyConcentrations () const
 Provide direct access to current third-body concentration values.
 
virtual void getCreationRates (double *cdot)
 Species creation rates [kmol/m^3/s or kmol/m^2/s].
 
virtual void getDestructionRates (double *ddot)
 Species destruction rates [kmol/m^3/s or kmol/m^2/s].
 
virtual void getNetProductionRates (double *wdot)
 Species net production rates [kmol/m^3/s or kmol/m^2/s].
 
virtual void getFwdRateConstants_ddT (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions.
 
virtual void getFwdRateConstants_ddP (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions.
 
virtual void getFwdRateConstants_ddC (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions.
 
virtual void getFwdRatesOfProgress_ddT (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
 
virtual void getFwdRatesOfProgress_ddP (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
 
virtual void getFwdRatesOfProgress_ddC (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
 
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX ()
 Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
virtual void getRevRatesOfProgress_ddT (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
 
virtual void getRevRatesOfProgress_ddP (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
 
virtual void getRevRatesOfProgress_ddC (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
 
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX ()
 Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
virtual void getNetRatesOfProgress_ddT (double *drop)
 Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.
 
virtual void getNetRatesOfProgress_ddP (double *drop)
 Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.
 
virtual void getNetRatesOfProgress_ddC (double *drop)
 Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.
 
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX ()
 Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
void getCreationRates_ddT (double *dwdot)
 Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions.
 
void getCreationRates_ddP (double *dwdot)
 Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions.
 
void getCreationRates_ddC (double *dwdot)
 Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions.
 
Eigen::SparseMatrix< double > creationRates_ddX ()
 Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
Eigen::SparseMatrix< double > creationRates_ddCi ()
 Calculate derivatives for species creation rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
 
void getDestructionRates_ddT (double *dwdot)
 Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions.
 
void getDestructionRates_ddP (double *dwdot)
 Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions.
 
void getDestructionRates_ddC (double *dwdot)
 Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions.
 
Eigen::SparseMatrix< double > destructionRates_ddX ()
 Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
Eigen::SparseMatrix< double > destructionRates_ddCi ()
 Calculate derivatives for species destruction rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
 
void getNetProductionRates_ddT (double *dwdot)
 Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions.
 
void getNetProductionRates_ddP (double *dwdot)
 Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions.
 
void getNetProductionRates_ddC (double *dwdot)
 Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions.
 
Eigen::SparseMatrix< double > netProductionRates_ddX ()
 Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration.
 
Eigen::SparseMatrix< double > netProductionRates_ddCi ()
 Calculate derivatives for species net production rates with respect to species concentration at constant temperature, pressure, and concentration of all other species.
 
virtual double reactantStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a reactant in reaction i.
 
Eigen::SparseMatrix< double > reactantStoichCoeffs () const
 Stoichiometric coefficient matrix for reactants.
 
virtual double productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i.
 
Eigen::SparseMatrix< double > productStoichCoeffs () const
 Stoichiometric coefficient matrix for products.
 
Eigen::SparseMatrix< double > revProductStoichCoeffs () const
 Stoichiometric coefficient matrix for products of reversible reactions.
 
virtual double reactantOrder (size_t k, size_t i) const
 Reactant order of species k in reaction i.
 
virtual double productOrder (int k, int i) const
 product Order of species k in reaction i.
 
AnyMap parameters ()
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newKinetics function.
 
shared_ptr< Reactionreaction (size_t i)
 Return the Reaction object for reaction i.
 
shared_ptr< const Reactionreaction (size_t i) const
 
void skipUndeclaredSpecies (bool skip)
 Determine behavior when adding a new reaction that contains species not defined in any of the phases associated with this kinetics manager.
 
bool skipUndeclaredSpecies () const
 
void skipUndeclaredThirdBodies (bool skip)
 Determine behavior when adding a new reaction that contains third-body efficiencies for species not defined in any of the phases associated with this kinetics manager.
 
bool skipUndeclaredThirdBodies () const
 
double multiplier (size_t i) const
 The current value of the multiplier for reaction i.
 
virtual void invalidateCache ()
 

Protected Member Functions

Internal service methods
Note
These methods are for internal use, and seek to avoid code duplication while evaluating terms used for rate constants, rates of progress, and their derivatives.
void applyEquilibriumConstants (double *rop)
 Multiply rate with inverse equilibrium constant.
 
Eigen::SparseMatrix< double > calculateCompositionDerivatives (StoichManagerN &stoich, const vector< double > &in)
 Process mole fraction derivative.
 
void assertDerivativesValid (const string &name)
 Helper function ensuring that all rate derivatives can be calculated.
 
- Protected Member Functions inherited from Kinetics
virtual void updateROP ()
 
double checkDuplicateStoich (map< int, double > &r1, map< int, double > &r2) const
 Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reactions are duplicates of one another, and 0.0 otherwise.
 

Protected Attributes

vector< double > m_grt
 Temporary work vector of length m_kk.
 
vector< size_t > m_revindex
 List of reactions numbers which are reversible reactions.
 
bool m_redo_rates = false
 
vector< unique_ptr< MultiRateBase > > m_interfaceRates
 Vector of rate handlers for interface reactions.
 
map< string, size_t > m_interfaceTypes
 Rate handler mapping.
 
vector< size_t > m_irrev
 Vector of irreversible reaction numbers.
 
vector< double > m_conc
 Array of concentrations for each species in the kinetics mechanism.
 
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.
 
vector< double > m_mu
 Vector of chemical potentials for all species.
 
vector< double > m_mu0_Kc
 Vector of standard state electrochemical potentials modified by a standard concentration term.
 
vector< double > m_phi
 Vector of phase electric potentials.
 
SurfPhasem_surf = nullptr
 Pointer to the single surface phase.
 
ImplicitSurfChemm_integrator = nullptr
 Pointer to the Implicit surface chemistry object.
 
bool m_ROP_ok = false
 
double m_temp = 0.0
 Current temperature of the data.
 
int m_phaseExistsCheck = false
 Int flag to indicate that some phases in the kinetics mechanism are non-existent.
 
vector< bool > m_phaseExists
 Vector of booleans indicating whether phases exist or not.
 
vector< int > m_phaseIsStable
 Vector of int indicating whether phases are stable or not.
 
vector< vector< bool > > m_rxnPhaseIsReactant
 Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant.
 
vector< vector< bool > > m_rxnPhaseIsProduct
 Vector of vector of booleans indicating whether a phase participates in a reaction as a product.
 
int m_ioFlag = 0
 
size_t m_nDim = 2
 Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for EdgeKinetics)
 
vector< double > m_rbuf0
 Buffers for partial rop results with length nReactions()
 
vector< double > m_rbuf1
 
bool m_jac_skip_coverage_dependence = false
 A flag used to neglect rate coefficient coverage dependence in derivative formation.
 
bool m_jac_skip_electrochemistry = false
 A flag used to neglect electrochemical contributions in derivative formation.
 
double m_jac_rtol_delta = 1e-8
 Relative tolerance used in developing numerical portions of specific derivatives.
 
bool m_has_electrochemistry = false
 A flag stating if the object uses electrochemistry.
 
bool m_has_coverage_dependence = false
 A flag stating if the object has coverage dependent rates.
 
- Protected Attributes inherited from Kinetics
ValueCache m_cache
 Cache for saved calculations within each Kinetics object.
 
bool m_ready = false
 Boolean indicating whether Kinetics object is fully configured.
 
size_t m_kk = 0
 The number of species in all of the phases that participate in this kinetics mechanism.
 
vector< double > m_perturb
 Vector of perturbation factors for each reaction's rate of progress vector.
 
vector< shared_ptr< Reaction > > m_reactions
 Vector of Reaction objects represented by this Kinetics manager.
 
vector< shared_ptr< ThermoPhase > > m_thermo
 m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
 
vector< size_t > m_start
 m_start is a vector of integers specifying the beginning position for the species vector for the n'th phase in the kinetics class.
 
map< string, size_t > m_phaseindex
 Mapping of the phase name to the position of the phase within the kinetics object.
 
size_t m_mindim = 4
 number of spatial dimensions of lowest-dimensional phase.
 
vector< double > m_rfn
 Forward rate constant for each reaction.
 
vector< double > m_delta_gibbs0
 Delta G^0 for all reactions.
 
vector< double > m_rkcn
 Reciprocal of the equilibrium constant in concentration units.
 
vector< double > m_ropf
 Forward rate-of-progress for each reaction.
 
vector< double > m_ropr
 Reverse rate-of-progress for each reaction.
 
vector< double > m_ropnet
 Net rate-of-progress for each reaction.
 
vector< double > m_dH
 The enthalpy change for each reaction to calculate Blowers-Masel rates.
 
vector< double > m_rbuf
 Buffer used for storage of intermediate reaction-specific results.
 
bool m_skipUndeclaredSpecies = false
 See skipUndeclaredSpecies()
 
bool m_skipUndeclaredThirdBodies = false
 See skipUndeclaredThirdBodies()
 
bool m_hasUndeclaredThirdBodies = false
 Flag indicating whether reactions include undeclared third bodies.
 
std::weak_ptr< Solutionm_root
 reference to Solution
 
StoichManagerN m_reactantStoich
 Stoichiometry manager for the reactants for each reaction.
 
StoichManagerN m_productStoich
 Stoichiometry manager for the products for each reaction.
 
StoichManagerN m_revProductStoich
 Stoichiometry manager for the products of reversible reactions.
 
Eigen::SparseMatrix< double > m_stoichMatrix
 Net stoichiometry (products - reactants)
 

Constructor & Destructor Documentation

◆ InterfaceKinetics()

InterfaceKinetics ( )
default

Constructor.

◆ ~InterfaceKinetics()

~InterfaceKinetics ( )
override

Definition at line 17 of file InterfaceKinetics.cpp.

Member Function Documentation

◆ resizeReactions()

void resizeReactions ( )
overridevirtual

Finalize Kinetics object and associated objects.

Reimplemented from Kinetics.

Definition at line 22 of file InterfaceKinetics.cpp.

◆ kineticsType()

string kineticsType ( ) const
inlineoverridevirtual

Identifies the Kinetics manager type.

Each class derived from Kinetics should override this method to return a meaningful identifier.

Since
Starting in Cantera 3.0, the name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented from Kinetics.

Definition at line 65 of file InterfaceKinetics.h.

◆ setElectricPotential()

void setElectricPotential ( int  n,
double  V 
)

Set the electric potential in the nth phase.

Parameters
nphase Index in this kinetics object.
VElectric potential (volts)

Definition at line 38 of file InterfaceKinetics.cpp.

◆ getEquilibriumConstants()

void getEquilibriumConstants ( double *  kc)
overridevirtual

Equilibrium constant for all reactions including the voltage term.

Kc = exp(deltaG/RT)

where deltaG is the electrochemical potential difference between products minus reactants.

Reimplemented from Kinetics.

Definition at line 166 of file InterfaceKinetics.cpp.

◆ getDeltaGibbs()

void getDeltaGibbs ( double *  deltaG)
overridevirtual

Return the vector of values for the reaction Gibbs free energy change.

(virtual from Kinetics.h) These values depend upon the concentration of the solution.

units = J kmol-1

Parameters
deltaGOutput vector of deltaG's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 283 of file InterfaceKinetics.cpp.

◆ getDeltaElectrochemPotentials()

void getDeltaElectrochemPotentials ( double *  deltaM)
overridevirtual

Return the vector of values for the reaction electrochemical free energy change.

These values depend upon the concentration of the solution and the voltage of the phases

units = J kmol-1

Parameters
deltaMOutput vector of deltaM's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 300 of file InterfaceKinetics.cpp.

◆ getDeltaEnthalpy()

void getDeltaEnthalpy ( double *  deltaH)
overridevirtual

Return the vector of values for the reactions change in enthalpy.

These values depend upon the concentration of the solution.

units = J kmol-1

Parameters
deltaHOutput vector of deltaH's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 311 of file InterfaceKinetics.cpp.

◆ getDeltaEntropy()

void getDeltaEntropy ( double *  deltaS)
overridevirtual

Return the vector of values for the reactions change in entropy.

These values depend upon the concentration of the solution.

units = J kmol-1 Kelvin-1

Parameters
deltaSOutput vector of deltaS's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 322 of file InterfaceKinetics.cpp.

◆ getDeltaSSGibbs()

void getDeltaSSGibbs ( double *  deltaG)
overridevirtual

Return the vector of values for the reaction standard state Gibbs free energy change.

These values don't depend upon the concentration of the solution.

units = J kmol-1

Parameters
deltaGOutput vector of ss deltaG's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 333 of file InterfaceKinetics.cpp.

◆ getDeltaSSEnthalpy()

void getDeltaSSEnthalpy ( double *  deltaH)
overridevirtual

Return the vector of values for the change in the standard state enthalpies of reaction.

These values don't depend upon the concentration of the solution.

units = J kmol-1

Parameters
deltaHOutput vector of ss deltaH's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 347 of file InterfaceKinetics.cpp.

◆ getDeltaSSEntropy()

void getDeltaSSEntropy ( double *  deltaS)
overridevirtual

Return the vector of values for the change in the standard state entropies for each reaction.

These values don't depend upon the concentration of the solution.

units = J kmol-1 Kelvin-1

Parameters
deltaSOutput vector of ss deltaS's for reactions Length: nReactions().

Reimplemented from Kinetics.

Definition at line 364 of file InterfaceKinetics.cpp.

◆ getActivityConcentrations()

void getActivityConcentrations ( double *const  conc)
overridevirtual

Get the vector of activity concentrations used in the kinetics object.

Parameters
[out]concVector of activity concentrations. Length is equal to the number of species in the kinetics object

Reimplemented from Kinetics.

Definition at line 111 of file InterfaceKinetics.cpp.

◆ isReversible()

bool isReversible ( size_t  i)
inlineoverridevirtual

True if reaction i has been declared to be reversible.

If isReversible(i) is false, then the reverse rate of progress for reaction i is always zero.

Parameters
ireaction index

Reimplemented from Kinetics.

Definition at line 104 of file InterfaceKinetics.h.

◆ getFwdRateConstants()

void getFwdRateConstants ( double *  kfwd)
overridevirtual

Return the forward rate constants.

The computed values include all temperature-dependent and pressure-dependent contributions. By default, third-body concentrations are only considered if they are part of the reaction rate definition; for a legacy implementation that includes third-body concentrations see Cantera::use_legacy_rate_constants(). Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Parameters
kfwdOutput vector containing the forward reaction rate constants. Length: nReactions().

Reimplemented from Kinetics.

Definition at line 177 of file InterfaceKinetics.cpp.

◆ getRevRateConstants()

void getRevRateConstants ( double *  krev,
bool  doIrreversible = false 
)
overridevirtual

Return the reverse rate constants.

The computed values include all temperature-dependent and pressure-dependent contributions. By default, third-body concentrations are only considered if they are part of the reaction rate definition; for a legacy implementation that includes third-body concentrations see Cantera::use_legacy_rate_constants(). Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction. Note, this routine will return rate constants for irreversible reactions if the default for doIrreversible is overridden.

Parameters
krevOutput vector of reverse rate constants
doIrreversibleboolean indicating whether irreversible reactions should be included.

Reimplemented from Kinetics.

Definition at line 186 of file InterfaceKinetics.cpp.

◆ addThermo()

void addThermo ( shared_ptr< ThermoPhase thermo)
overridevirtual

Add a thermo phase to the kinetics manager object.

This must be done before the function init() is called or before any reactions are input. The lowest dimensional phase, where reactions occur, must be added first.

This function calls Kinetics::addThermo). It also sets the following fields:

   m_phaseExists[]
Parameters
thermoReference to the ThermoPhase to be added.

Reimplemented from Kinetics.

Definition at line 487 of file InterfaceKinetics.cpp.

◆ init()

void init ( )
overridevirtual

Prepare the class for the addition of reactions, after all phases have been added.

This method is called automatically when the first reaction is added. It needs to be called directly only in the degenerate case where there are no reactions. The base class method does nothing, but derived classes may use this to perform any initialization (allocating arrays, etc.) that requires knowing the phases.

Reimplemented from Kinetics.

Definition at line 494 of file InterfaceKinetics.cpp.

◆ resizeSpecies()

void resizeSpecies ( )
overridevirtual

Resize arrays with sizes that depend on the total number of species.

Automatically called before adding each Reaction and Phase.

Reimplemented from Kinetics.

Definition at line 501 of file InterfaceKinetics.cpp.

◆ addReaction()

bool addReaction ( shared_ptr< Reaction r,
bool  resize = true 
)
overridevirtual

Add a single reaction to the mechanism.

Derived classes should call the base class method in addition to handling their own specialized behavior.

Parameters
rPointer to the Reaction object to be added.
resizeIf true, resizeReactions is called after reaction is added.
Returns
true if the reaction is added or false if it was skipped

Reimplemented from Kinetics.

Definition at line 380 of file InterfaceKinetics.cpp.

◆ modifyReaction()

void modifyReaction ( size_t  i,
shared_ptr< Reaction rNew 
)
overridevirtual

Modify the rate expression associated with a reaction.

The stoichiometric equation, type of the reaction, reaction orders, third body efficiencies, reversibility, etc. must be unchanged.

Parameters
iIndex of the reaction to be modified
rNewReaction with the new rate expressions

Reimplemented from Kinetics.

Definition at line 446 of file InterfaceKinetics.cpp.

◆ setMultiplier()

void setMultiplier ( size_t  i,
double  f 
)
overridevirtual

Set the multiplier for reaction i to f.

Parameters
iindex of the reaction
fvalue of the multiplier.

Reimplemented from Kinetics.

Definition at line 473 of file InterfaceKinetics.cpp.

◆ updateROP()

void updateROP ( )
overridevirtual

Internal routine that updates the Rates of Progress of the reactions.

This is actually the guts of the functionality of the object

Reimplemented from Kinetics.

Definition at line 201 of file InterfaceKinetics.cpp.

◆ _update_rates_T()

void _update_rates_T ( )

Update properties that depend on temperature.

Current objects that this function updates: m_kdata->m_rfn updateKc();

Definition at line 44 of file InterfaceKinetics.cpp.

◆ _update_rates_phi()

void _update_rates_phi ( )

Update properties that depend on the electric potential.

Definition at line 81 of file InterfaceKinetics.cpp.

◆ _update_rates_C()

void _update_rates_C ( )

Update properties that depend on the species mole fractions and/or concentration,.

This method fills out the array of generalized concentrations by calling method getActivityConcentrations for each phase, which classes representing phases should overload to return the appropriate quantities.

Definition at line 92 of file InterfaceKinetics.cpp.

◆ advanceCoverages()

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.

This method carries out a time-accurate advancement of the surface coverages for a specified amount of time.

\[ \dot {\theta}_k = \dot s_k (\sigma_k / s_0) \]

Parameters
tstepTime value to advance the surface coverages
rtolThe relative tolerance for the integrator
atolThe absolute tolerance for the integrator
maxStepSizeThe maximum step-size the integrator is allowed to take. If zero, this option is disabled.
maxStepsThe maximum number of time-steps the integrator can take. If not supplied, uses the default value in CVodeIntegrator (20000).
maxErrTestFailsthe maximum permissible number of error test failures If not supplied, uses the default value in CVODES (7).

Definition at line 518 of file InterfaceKinetics.cpp.

◆ solvePseudoSteadyStateProblem()

void solvePseudoSteadyStateProblem ( int  ifuncOverride = -1,
double  timeScaleOverride = 1.0 
)

Solve for the pseudo steady-state of the surface problem.

This is the same thing as the advanceCoverages() function, but at infinite times.

Note, a direct solve is carried out under the hood here, to reduce the computational time.

Parameters
ifuncOverrideOne of the values defined in Surface Problem Solver Methods. The default is -1, which means that the program will decide.
timeScaleOverrideWhen a pseudo transient is selected this value can be used to override the default time scale for integration which is one. When SFLUX_TRANSIENT is used, this is equal to the time over which the equations are integrated. When SFLUX_INITIALIZE is used, this is equal to the time used in the initial transient algorithm, before the equation system is solved directly.

Definition at line 534 of file InterfaceKinetics.cpp.

◆ setIOFlag()

void setIOFlag ( int  ioFlag)

Definition at line 479 of file InterfaceKinetics.cpp.

◆ updateMu0()

void updateMu0 ( )
virtual

Update the standard state chemical potentials and species equilibrium constant entries.

Virtual because it is overridden when dealing with experimental open circuit voltage overrides

Definition at line 148 of file InterfaceKinetics.cpp.

◆ updateKc()

void updateKc ( )

Update the equilibrium constants and stored electrochemical potentials in molar units for all reversible reactions and for all species.

Irreversible reactions have their equilibrium constant set to zero. For reactions involving charged species the equilibrium constant is adjusted according to the electrostatic potential.

Definition at line 117 of file InterfaceKinetics.cpp.

◆ setPhaseExistence()

void setPhaseExistence ( const size_t  iphase,
const int  exists 
)

Set the existence of a phase in the reaction object.

Tell the kinetics object whether a phase in the object exists. This is actually an extrinsic specification that must be carried out on top of the intrinsic calculation of the reaction rate. The routine will also flip the IsStable boolean within the kinetics object as well.

Parameters
iphaseIndex of the phase. This is the order within the internal thermo vector object
existsBoolean indicating whether the phase exists or not

Definition at line 548 of file InterfaceKinetics.cpp.

◆ setPhaseStability()

void setPhaseStability ( const size_t  iphase,
const int  isStable 
)

Set the stability of a phase in the reaction object.

Tell the kinetics object whether a phase in the object is stable. Species in an unstable phase will not be allowed to have a positive rate of formation from this kinetics object. This is actually an extrinsic specification that must be carried out on top of the intrinsic calculation of the reaction rate.

While conceptually not needed since kinetics is consistent with thermo when taken as a whole, in practice it has found to be very useful to turn off the creation of phases which shouldn't be forming. Typically this can reduce the oscillations in phase formation and destruction which are observed.

Parameters
iphaseIndex of the phase. This is the order within the internal thermo vector object
isStableFlag indicating whether the phase is stable or not

Definition at line 579 of file InterfaceKinetics.cpp.

◆ phaseExistence()

int phaseExistence ( const size_t  iphase) const

Gets the phase existence int for the ith phase.

Parameters
iphasePhase Id
Returns
The int specifying whether the kinetics object thinks the phase exists or not. If it exists, then species in that phase can be a reactant in reactions.

Definition at line 567 of file InterfaceKinetics.cpp.

◆ phaseStability()

int phaseStability ( const size_t  iphase) const

Gets the phase stability int for the ith phase.

Parameters
iphasePhase Id
Returns
The int specifying whether the kinetics object thinks the phase is stable with nonzero mole numbers. If it stable, then the kinetics object will allow for rates of production of of species in that phase that are positive.

Definition at line 573 of file InterfaceKinetics.cpp.

◆ interfaceCurrent()

double interfaceCurrent ( const size_t  iphase)

Gets the interface current for the ith phase.

Parameters
iphasePhase Id
Returns
The double specifying the interface current. The interface current is useful when charge transfer reactions occur at an interface. It is defined here as the net positive charge entering the phase specified by the Phase Id. (Units: A/m^2 for a surface reaction, A/m for an edge reaction).

Definition at line 589 of file InterfaceKinetics.cpp.

◆ setDerivativeSettings()

void setDerivativeSettings ( const AnyMap settings)
overridevirtual

Set/modify derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented from Kinetics.

Definition at line 642 of file InterfaceKinetics.cpp.

◆ getDerivativeSettings()

void getDerivativeSettings ( AnyMap settings) const
overridevirtual

Retrieve derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented from Kinetics.

Definition at line 658 of file InterfaceKinetics.cpp.

◆ fwdRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi ( )
overridevirtual

Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.

The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.0.

Reimplemented from Kinetics.

Definition at line 606 of file InterfaceKinetics.cpp.

◆ revRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > revRatesOfProgress_ddCi ( )
overridevirtual

Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations.

The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.0.

Reimplemented from Kinetics.

Definition at line 616 of file InterfaceKinetics.cpp.

◆ netRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > netRatesOfProgress_ddCi ( )
overridevirtual

Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations.

The method returns a matrix with nReactions() rows and nTotalSpecies() columns. For a derivative with respect to \( c_i \), all other \( c_j \) are held constant.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.0.

Reimplemented from Kinetics.

Definition at line 627 of file InterfaceKinetics.cpp.

◆ applyEquilibriumConstants()

void applyEquilibriumConstants ( double *  rop)
protected

Multiply rate with inverse equilibrium constant.

Definition at line 683 of file InterfaceKinetics.cpp.

◆ calculateCompositionDerivatives()

Eigen::SparseMatrix< double > calculateCompositionDerivatives ( StoichManagerN stoich,
const vector< double > &  in 
)
protected

Process mole fraction derivative.

Parameters
stoichstoichiometry manager
inrate expression used for the derivative calculation
Returns
a sparse matrix of derivative contributions for each reaction of dimensions nTotalReactions by nTotalSpecies

Definition at line 665 of file InterfaceKinetics.cpp.

◆ assertDerivativesValid()

void assertDerivativesValid ( const string &  name)
protected

Helper function ensuring that all rate derivatives can be calculated.

Parameters
namemethod name used for error output
Exceptions
CanteraErrorif coverage dependence or electrochemical reactions are included

Definition at line 674 of file InterfaceKinetics.cpp.

Member Data Documentation

◆ m_grt

vector<double> m_grt
protected

Temporary work vector of length m_kk.

Definition at line 329 of file InterfaceKinetics.h.

◆ m_revindex

vector<size_t> m_revindex
protected

List of reactions numbers which are reversible reactions.

This is a vector of reaction numbers. Each reaction in the list is reversible. Length = number of reversible reactions

Definition at line 336 of file InterfaceKinetics.h.

◆ m_redo_rates

bool m_redo_rates = false
protected

Definition at line 338 of file InterfaceKinetics.h.

◆ m_interfaceRates

vector<unique_ptr<MultiRateBase> > m_interfaceRates
protected

Vector of rate handlers for interface reactions.

Definition at line 341 of file InterfaceKinetics.h.

◆ m_interfaceTypes

map<string, size_t> m_interfaceTypes
protected

Rate handler mapping.

Definition at line 342 of file InterfaceKinetics.h.

◆ m_irrev

vector<size_t> m_irrev
protected

Vector of irreversible reaction numbers.

vector containing the reaction numbers of irreversible reactions.

Definition at line 348 of file InterfaceKinetics.h.

◆ m_conc

vector<double> m_conc
protected

Array of concentrations for each species in the kinetics mechanism.

An array of generalized concentrations \( C_k \) that are defined such that \( a_k = C_k / C^0_k, \) where \( C^0_k \) is a standard concentration/ These generalized concentrations are used by this kinetics manager class to compute the forward and reverse rates of elementary reactions. The "units" for the concentrations of each phase depend upon the implementation of kinetics within that phase. The order of the species within the vector is based on the order of listed ThermoPhase objects in the class, and the order of the species within each ThermoPhase class.

Definition at line 362 of file InterfaceKinetics.h.

◆ m_actConc

vector<double> m_actConc
protected

Array of activity concentrations for each species in the kinetics object.

An array of activity concentrations \( Ca_k \) that are defined such that \( a_k = Ca_k / C^0_k, \) where \( C^0_k \) is a standard concentration. These activity concentrations are used by this kinetics manager class to compute the forward and reverse rates of elementary reactions. The "units" for the concentrations of each phase depend upon the implementation of kinetics within that phase. The order of the species within the vector is based on the order of listed ThermoPhase objects in the class, and the order of the species within each ThermoPhase class.

Definition at line 376 of file InterfaceKinetics.h.

◆ m_mu0

vector<double> m_mu0
protected

Vector of standard state chemical potentials for all species.

This vector contains a temporary vector of standard state chemical potentials for all of the species in the kinetics object

Length = m_kk. Units = J/kmol.

Definition at line 385 of file InterfaceKinetics.h.

◆ m_mu

vector<double> m_mu
protected

Vector of chemical potentials for all species.

This vector contains a vector of chemical potentials for all of the species in the kinetics object

Length = m_kk. Units = J/kmol.

Definition at line 394 of file InterfaceKinetics.h.

◆ m_mu0_Kc

vector<double> m_mu0_Kc
protected

Vector of standard state electrochemical potentials modified by a standard concentration term.

This vector contains a temporary vector of standard state electrochemical potentials + RTln(Cs) for all of the species in the kinetics object

In order to get the units correct for the concentration equilibrium constant, each species needs to have an RT ln(Cs) added to its contribution to the equilibrium constant Cs is the standard concentration for the species. Frequently, for solid species, Cs is equal to 1. However, for gases Cs is P/RT. Length = m_kk. Units = J/kmol.

Definition at line 408 of file InterfaceKinetics.h.

◆ m_phi

vector<double> m_phi
protected

Vector of phase electric potentials.

Temporary vector containing the potential of each phase in the kinetics object. length = number of phases. Units = Volts.

Definition at line 415 of file InterfaceKinetics.h.

◆ m_surf

SurfPhase* m_surf = nullptr
protected

Pointer to the single surface phase.

Definition at line 418 of file InterfaceKinetics.h.

◆ m_integrator

ImplicitSurfChem* m_integrator = nullptr
protected

Pointer to the Implicit surface chemistry object.

Note this object is owned by this InterfaceKinetics object. It may only be used to solve this single InterfaceKinetics object's surface problem uncoupled from other surface phases.

Definition at line 426 of file InterfaceKinetics.h.

◆ m_ROP_ok

bool m_ROP_ok = false
protected

Definition at line 428 of file InterfaceKinetics.h.

◆ m_temp

double m_temp = 0.0
protected

Current temperature of the data.

Definition at line 431 of file InterfaceKinetics.h.

◆ m_phaseExistsCheck

int m_phaseExistsCheck = false
protected

Int flag to indicate that some phases in the kinetics mechanism are non-existent.

We change the ROP vectors to make sure that non-existent phases are treated correctly in the kinetics operator. The value of this is equal to the number of phases which don't exist.

Definition at line 440 of file InterfaceKinetics.h.

◆ m_phaseExists

vector<bool> m_phaseExists
protected

Vector of booleans indicating whether phases exist or not.

Vector of booleans indicating whether a phase exists or not. We use this to set the ROP's so that unphysical things don't happen. For example, a reaction can't go in the forwards direction if a phase in which a reactant is present doesn't exist. Because InterfaceKinetics deals with intrinsic quantities only normally, nowhere else is this extrinsic concept introduced except here.

length = number of phases in the object. By default all phases exist.

Definition at line 453 of file InterfaceKinetics.h.

◆ m_phaseIsStable

vector<int> m_phaseIsStable
protected

Vector of int indicating whether phases are stable or not.

Vector of booleans indicating whether a phase is stable or not under the current conditions. We use this to set the ROP's so that unphysical things don't happen.

length = number of phases in the object. By default all phases are stable.

Definition at line 463 of file InterfaceKinetics.h.

◆ m_rxnPhaseIsReactant

vector<vector<bool> > m_rxnPhaseIsReactant
protected

Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant.

m_rxnPhaseIsReactant[j][p] indicates whether a species in phase p participates in reaction j as a reactant.

Definition at line 471 of file InterfaceKinetics.h.

◆ m_rxnPhaseIsProduct

vector<vector<bool> > m_rxnPhaseIsProduct
protected

Vector of vector of booleans indicating whether a phase participates in a reaction as a product.

m_rxnPhaseIsReactant[j][p] indicates whether a species in phase p participates in reaction j as a product.

Definition at line 479 of file InterfaceKinetics.h.

◆ m_ioFlag

int m_ioFlag = 0
protected

Definition at line 481 of file InterfaceKinetics.h.

◆ m_nDim

size_t m_nDim = 2
protected

Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for EdgeKinetics)

Definition at line 485 of file InterfaceKinetics.h.

◆ m_rbuf0

vector<double> m_rbuf0
protected

Buffers for partial rop results with length nReactions()

Definition at line 488 of file InterfaceKinetics.h.

◆ m_rbuf1

vector<double> m_rbuf1
protected

Definition at line 489 of file InterfaceKinetics.h.

◆ m_jac_skip_coverage_dependence

bool m_jac_skip_coverage_dependence = false
protected

A flag used to neglect rate coefficient coverage dependence in derivative formation.

Definition at line 493 of file InterfaceKinetics.h.

◆ m_jac_skip_electrochemistry

bool m_jac_skip_electrochemistry = false
protected

A flag used to neglect electrochemical contributions in derivative formation.

Definition at line 495 of file InterfaceKinetics.h.

◆ m_jac_rtol_delta

double m_jac_rtol_delta = 1e-8
protected

Relative tolerance used in developing numerical portions of specific derivatives.

Definition at line 497 of file InterfaceKinetics.h.

◆ m_has_electrochemistry

bool m_has_electrochemistry = false
protected

A flag stating if the object uses electrochemistry.

Definition at line 499 of file InterfaceKinetics.h.

◆ m_has_coverage_dependence

bool m_has_coverage_dependence = false
protected

A flag stating if the object has coverage dependent rates.

Definition at line 501 of file InterfaceKinetics.h.


The documentation for this class was generated from the following files: