Cantera 2.6.0
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
InterfaceKinetics Class Reference

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

#include <InterfaceKinetics.h>

Inheritance diagram for InterfaceKinetics:
[legend]
Collaboration diagram for InterfaceKinetics:
[legend]

Classes

struct  StickData
 Values used for converting sticking coefficients into rate constants. More...
 

Public Member Functions

 InterfaceKinetics (ThermoPhase *thermo=0)
 Constructor. More...
 
virtual void resizeReactions ()
 Finalize Kinetics object and associated objects. More...
 
virtual std::string kineticsType () const
 Identifies the Kinetics manager type. More...
 
void setElectricPotential (int n, doublereal V)
 Set the electric potential in the nth phase. More...
 
virtual void updateROP ()
 Internal routine that updates the Rates of Progress of the reactions. More...
 
void _update_rates_T ()
 Update properties that depend on temperature. More...
 
void _update_rates_phi ()
 Update properties that depend on the electric potential. More...
 
void _update_rates_C ()
 Update properties that depend on the species mole fractions and/or concentration,. More...
 
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. More...
 
void solvePseudoSteadyStateProblem (int ifuncOverride=-1, doublereal timeScaleOverride=1.0)
 Solve for the pseudo steady-state of the surface problem. More...
 
void setIOFlag (int ioFlag)
 
virtual void updateMu0 ()
 Update the standard state chemical potentials and species equilibrium constant entries. More...
 
void updateKc ()
 Update the equilibrium constants and stored electrochemical potentials in molar units for all reversible reactions and for all species. More...
 
void applyVoltageKfwdCorrection (doublereal *const kfwd)
 Apply modifications for the forward reaction rate for interfacial charge transfer reactions. More...
 
void convertExchangeCurrentDensityFormulation (doublereal *const kfwd)
 When an electrode reaction rate is optionally specified in terms of its exchange current density, adjust kfwd to the standard reaction rate constant form and units. More...
 
void setPhaseExistence (const size_t iphase, const int exists)
 Set the existence of a phase in the reaction object. More...
 
void setPhaseStability (const size_t iphase, const int isStable)
 Set the stability of a phase in the reaction object. More...
 
int phaseExistence (const size_t iphase) const
 Gets the phase existence int for the ith phase. More...
 
int phaseStability (const size_t iphase) const
 Gets the phase stability int for the ith phase. More...
 
Reaction Rates Of Progress
virtual void getEquilibriumConstants (doublereal *kc)
 Equilibrium constant for all reactions including the voltage term. More...
 
void updateExchangeCurrentQuantities ()
 values needed to convert from exchange current density to surface reaction rate. More...
 
virtual void getDeltaGibbs (doublereal *deltaG)
 Return the vector of values for the reaction Gibbs free energy change. More...
 
virtual void getDeltaElectrochemPotentials (doublereal *deltaM)
 Return the vector of values for the reaction electrochemical free energy change. More...
 
virtual void getDeltaEnthalpy (doublereal *deltaH)
 Return the vector of values for the reactions change in enthalpy. More...
 
virtual void getDeltaEntropy (doublereal *deltaS)
 Return the vector of values for the reactions change in entropy. More...
 
virtual void getDeltaSSGibbs (doublereal *deltaG)
 Return the vector of values for the reaction standard state Gibbs free energy change. More...
 
virtual void getDeltaSSEnthalpy (doublereal *deltaH)
 Return the vector of values for the change in the standard state enthalpies of reaction. More...
 
virtual void getDeltaSSEntropy (doublereal *deltaS)
 Return the vector of values for the change in the standard state entropies for each reaction. More...
 
Reaction Mechanism Informational Query Routines
virtual void getActivityConcentrations (doublereal *const conc)
 Get the vector of activity concentrations used in the kinetics object. More...
 
doublereal electrochem_beta (size_t irxn) const
 Return the charge transfer rxn Beta parameter for the ith reaction. More...
 
virtual bool isReversible (size_t i)
 True if reaction i has been declared to be reversible. More...
 
virtual void getFwdRateConstants (doublereal *kfwd)
 Return the forward rate constants. More...
 
virtual void getRevRateConstants (doublereal *krev, bool doIrreversible=false)
 Return the reverse rate constants. More...
 
double effectivePreExponentialFactor (size_t irxn)
 Return effective preexponent for the specified reaction. More...
 
double effectiveActivationEnergy_R (size_t irxn)
 Return effective activation energy for the specified reaction. More...
 
double effectiveTemperatureExponent (size_t irxn)
 Return effective temperature exponent for the specified reaction. More...
 
Reaction Mechanism Construction
virtual void addPhase (ThermoPhase &thermo)
 Add a phase to the kinetics manager object. More...
 
virtual void init ()
 Prepare the class for the addition of reactions, after all phases have been added. More...
 
virtual void resizeSpecies ()
 Resize arrays with sizes that depend on the total number of species. More...
 
virtual bool addReaction (shared_ptr< Reaction > r, bool resize=true)
 Add a single reaction to the mechanism. More...
 
virtual void modifyReaction (size_t i, shared_ptr< Reaction > rNew)
 Modify the rate expression associated with a reaction. More...
 
- Public Member Functions inherited from Kinetics
virtual std::pair< size_t, size_t > checkDuplicates (bool throw_err=true) const
 Check for unmarked duplicate reactions and unmatched marked duplicates. More...
 
void selectPhase (const double *data, const ThermoPhase *phase, double *phase_data)
 
virtual void setRoot (std::shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
virtual double reactionEnthalpy (const Composition &reactants, const Composition &products)
 Calculate the reaction enthalpy of a reaction which has not necessarily been added into the Kinetics object. More...
 
 Kinetics ()
 Default constructor. More...
 
virtual ~Kinetics ()
 
 Kinetics (const Kinetics &)=delete
 Kinetics objects are not copyable or assignable. More...
 
Kineticsoperator= (const Kinetics &)=delete
 
size_t nReactions () const
 Number of reactions in the reaction mechanism. More...
 
void checkReactionIndex (size_t m) const
 Check that the specified reaction index is in range Throws an exception if i is greater than nReactions() More...
 
void checkReactionArraySize (size_t ii) const
 Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions(). More...
 
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. More...
 
void checkSpeciesArraySize (size_t mm) const
 Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). More...
 
size_t nPhases () const
 The number of phases participating in the reaction mechanism. More...
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases() More...
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). More...
 
size_t phaseIndex (const std::string &ph) const
 Return the phase index of a phase in the list of phases defined within the object. More...
 
size_t surfacePhaseIndex () const
 This returns the integer index of the phase which has ThermoPhase type cSurf. More...
 
size_t reactionPhaseIndex () const
 Phase where the reactions occur. More...
 
ThermoPhasethermo (size_t n=0)
 This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism. More...
 
const ThermoPhasethermo (size_t n=0) const
 
size_t nTotalSpecies () const
 The total number of species in all phases participating in the kinetics mechanism. More...
 
size_t kineticsSpeciesIndex (size_t k, size_t n) const
 The location of species k of phase n in species arrays. More...
 
std::string kineticsSpeciesName (size_t k) const
 Return the name of the kth species in the kinetics manager. More...
 
size_t kineticsSpeciesIndex (const std::string &nm) const
 This routine will look up a species number based on the input std::string nm. More...
 
size_t kineticsSpeciesIndex (const std::string &nm, const std::string &ph) const
 This routine will look up a species number based on the input std::string nm. More...
 
ThermoPhasespeciesPhase (const std::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. More...
 
const ThermoPhasespeciesPhase (const std::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. More...
 
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. More...
 
virtual void getFwdRatesOfProgress (doublereal *fwdROP)
 Return the forward rates of progress of the reactions. More...
 
virtual void getRevRatesOfProgress (doublereal *revROP)
 Return the Reverse rates of progress of the reactions. More...
 
virtual void getNetRatesOfProgress (doublereal *netROP)
 Net rates of progress. More...
 
virtual void getReactionDelta (const double *property, double *deltaProperty) const
 Change in species properties. More...
 
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. More...
 
virtual void getThirdBodyConcentrations (double *concm)
 Return a vector of values of effective concentrations of third-body collision partners of any reaction. More...
 
virtual const vector_fpthirdBodyConcentrations () const
 Provide direct access to current third-body concentration values. More...
 
virtual void getCreationRates (doublereal *cdot)
 Species creation rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getDestructionRates (doublereal *ddot)
 Species destruction rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getNetProductionRates (doublereal *wdot)
 Species net production rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getDerivativeSettings (AnyMap &settings) const
 Retrieve derivative settings. More...
 
virtual void setDerivativeSettings (const AnyMap &settings)
 Set/modify derivative settings. More...
 
virtual void getFwdRateConstants_ddT (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
virtual void getFwdRateConstants_ddP (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
virtual void getFwdRateConstants_ddC (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
void getCreationRates_ddT (double *dwdot)
 Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getCreationRates_ddP (double *dwdot)
 Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getCreationRates_ddC (double *dwdot)
 Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > creationRates_ddX ()
 Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getDestructionRates_ddT (double *dwdot)
 Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getDestructionRates_ddP (double *dwdot)
 Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getDestructionRates_ddC (double *dwdot)
 Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > destructionRates_ddX ()
 Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getNetProductionRates_ddT (double *dwdot)
 Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddP (double *dwdot)
 Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddC (double *dwdot)
 Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
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. More...
 
virtual double reactantStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a reactant in reaction i. More...
 
Eigen::SparseMatrix< double > reactantStoichCoeffs () const
 Stoichiometric coefficient matrix for reactants. More...
 
virtual double productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i. More...
 
Eigen::SparseMatrix< double > productStoichCoeffs () const
 Stoichiometric coefficient matrix for products. More...
 
Eigen::SparseMatrix< double > revProductStoichCoeffs () const
 Stoichiometric coefficient matrix for products of reversible reactions. More...
 
virtual doublereal reactantOrder (size_t k, size_t i) const
 Reactant order of species k in reaction i. More...
 
virtual doublereal productOrder (int k, int i) const
 product Order of species k in reaction i. More...
 
virtual int reactionType (size_t i) const
 Flag specifying the type of reaction. More...
 
virtual std::string reactionTypeStr (size_t i) const
 String specifying the type of reaction. More...
 
std::string reactionString (size_t i) const
 Return a string representing the reaction. More...
 
std::string reactantString (size_t i) const
 Returns a string containing the reactants side of the reaction equation. More...
 
std::string productString (size_t i) const
 Returns a string containing the products side of the reaction equation. More...
 
AnyMap parameters ()
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newKinetics function. More...
 
shared_ptr< Reactionreaction (size_t i)
 Return the Reaction object for reaction i. More...
 
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. More...
 
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. More...
 
bool skipUndeclaredThirdBodies () const
 
doublereal multiplier (size_t i) const
 The current value of the multiplier for reaction i. More...
 
virtual void setMultiplier (size_t i, doublereal f)
 Set the multiplier for reaction i to f. More...
 
virtual void invalidateCache ()
 

Protected Member Functions

SurfaceArrhenius buildSurfaceArrhenius (size_t i, InterfaceReaction2 &r, bool replace)
 Build a SurfaceArrhenius object from a Reaction, taking into account the possible sticking coefficient form and coverage dependencies. More...
 
void applyStickingCorrection (double T, double *kf)
 
- Protected Member Functions inherited from Kinetics
double checkDuplicateStoich (std::map< int, double > &r1, std::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. More...
 
void checkReactionBalance (const Reaction &R)
 Check that the specified reaction is balanced (same number of atoms for each element in the reactants and products). More...
 

Protected Attributes

vector_fp m_grt
 Temporary work vector of length m_kk. More...
 
std::vector< size_t > m_revindex
 List of reactions numbers which are reversible reactions. More...
 
Rate1< SurfaceArrheniusm_rates
 Templated class containing the vector of reactions for this interface. More...
 
bool m_redo_rates
 
std::vector< unique_ptr< MultiRateBase > > m_interfaceRates
 Vector of rate handlers for interface reactions. More...
 
std::map< std::string, size_t > m_interfaceTypes
 Rate handler mapping. More...
 
std::vector< size_t > m_irrev
 Vector of irreversible reaction numbers. More...
 
vector_fp m_conc
 Array of concentrations for each species in the kinetics mechanism. More...
 
vector_fp m_actConc
 Array of activity concentrations for each species in the kinetics object. More...
 
vector_fp m_mu0
 Vector of standard state chemical potentials for all species. More...
 
vector_fp m_mu
 Vector of chemical potentials for all species. More...
 
vector_fp m_mu0_Kc
 Vector of standard state electrochemical potentials modified by a standard concentration term. More...
 
vector_fp m_phi
 Vector of phase electric potentials. More...
 
vector_fp m_pot
 Vector of potential energies due to Voltages. More...
 
vector_fp deltaElectricEnergy_
 Storage for the net electric energy change due to reaction. More...
 
SurfPhasem_surf
 Pointer to the single surface phase. More...
 
ImplicitSurfChemm_integrator
 Pointer to the Implicit surface chemistry object. More...
 
vector_fp m_beta
 Electrochemical transfer coefficient for the forward direction. More...
 
std::vector< size_t > m_ctrxn
 Vector of reaction indexes specifying the id of the charge transfer reactions in the mechanism. More...
 
vector_int m_ctrxn_ecdf
 Vector of booleans indicating whether the charge transfer reaction rate constant is described by an exchange current density rate constant expression. More...
 
vector_fp m_StandardConc
 Vector of standard concentrations. More...
 
vector_fp m_deltaG0
 Vector of delta G^0, the standard state Gibbs free energies for each reaction. More...
 
vector_fp m_deltaG
 Vector of deltaG[] of reaction, the delta Gibbs free energies for each reaction. More...
 
vector_fp m_ProdStanConcReac
 Vector of the products of the standard concentrations of the reactants. More...
 
bool m_ROP_ok
 
doublereal m_temp
 Current temperature of the data. More...
 
doublereal m_logtemp
 Current log of the temperature. More...
 
bool m_has_coverage_dependence
 Boolean flag indicating whether any reaction in the mechanism has a coverage dependent forward reaction rate. More...
 
bool m_has_electrochem_rxns
 Boolean flag indicating whether any reaction in the mechanism has a beta electrochemical parameter. More...
 
bool m_has_exchange_current_density_formulation
 Boolean flag indicating whether any reaction in the mechanism is described by an exchange current density expression. More...
 
int m_phaseExistsCheck
 Int flag to indicate that some phases in the kinetics mechanism are non-existent. More...
 
std::vector< bool > m_phaseExists
 Vector of booleans indicating whether phases exist or not. More...
 
vector_int m_phaseIsStable
 Vector of int indicating whether phases are stable or not. More...
 
std::vector< std::vector< bool > > m_rxnPhaseIsReactant
 Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant. More...
 
std::vector< std::vector< bool > > m_rxnPhaseIsProduct
 Vector of vector of booleans indicating whether a phase participates in a reaction as a product. More...
 
std::vector< StickDatam_stickingData
 Data for sticking reactions. More...
 
int m_ioFlag
 
size_t m_nDim
 Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for EdgeKinetics) More...
 
- Protected Attributes inherited from Kinetics
ValueCache m_cache
 Cache for saved calculations within each Kinetics object. More...
 
bool m_ready
 Boolean indicating whether Kinetics object is fully configured. More...
 
size_t m_kk
 The number of species in all of the phases that participate in this kinetics mechanism. More...
 
vector_fp m_perturb
 Vector of perturbation factors for each reaction's rate of progress vector. More...
 
std::vector< shared_ptr< Reaction > > m_reactions
 Vector of Reaction objects represented by this Kinetics manager. More...
 
std::vector< ThermoPhase * > m_thermo
 m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator More...
 
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 phase in the kinetics class. More...
 
std::map< std::string, size_t > m_phaseindex
 Mapping of the phase name to the position of the phase within the kinetics object. More...
 
size_t m_surfphase
 Index in the list of phases of the one surface phase. More...
 
size_t m_rxnphase
 Phase Index where reactions are assumed to be taking place. More...
 
size_t m_mindim
 number of spatial dimensions of lowest-dimensional phase. More...
 
vector_fp m_rfn
 Forward rate constant for each reaction. More...
 
vector_fp m_delta_gibbs0
 Delta G^0 for all reactions. More...
 
vector_fp m_rkcn
 Reciprocal of the equilibrium constant in concentration units. More...
 
vector_fp m_ropf
 Forward rate-of-progress for each reaction. More...
 
vector_fp m_ropr
 Reverse rate-of-progress for each reaction. More...
 
vector_fp m_ropnet
 Net rate-of-progress for each reaction. More...
 
vector_fp m_dH
 The enthalpy change for each reaction to calculate Blowers-Masel rates. More...
 
vector_fp m_rbuf
 Buffer used for storage of intermediate reaction-specific results. More...
 
bool m_skipUndeclaredSpecies
 
bool m_skipUndeclaredThirdBodies
 
std::weak_ptr< Solutionm_root
 reference to Solution More...
 
StoichManagerN m_reactantStoich
 Stoichiometry manager for the reactants for each reaction. More...
 
StoichManagerN m_productStoich
 Stoichiometry manager for the products for each reaction. More...
 
StoichManagerN m_revProductStoich
 Stoichiometry manager for the products of reversible reactions. More...
 
Eigen::SparseMatrix< double > m_stoichMatrix
 Net stoichiometry (products - reactants) More...
 

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 57 of file InterfaceKinetics.h.

Constructor & Destructor Documentation

◆ InterfaceKinetics()

InterfaceKinetics ( ThermoPhase thermo = 0)

Constructor.

Parameters
thermoThe optional parameter may be used to initialize the object with one ThermoPhase object. HKM Note -> Since the interface kinetics object will probably require multiple ThermoPhase objects, this is probably not a good idea to have this parameter.

Definition at line 20 of file InterfaceKinetics.cpp.

References InterfaceKinetics::addPhase(), and Kinetics::thermo().

◆ ~InterfaceKinetics()

~InterfaceKinetics ( )
virtual

Definition at line 39 of file InterfaceKinetics.cpp.

Member Function Documentation

◆ resizeReactions()

void resizeReactions ( )
virtual

Finalize Kinetics object and associated objects.

Reimplemented from Kinetics.

Definition at line 44 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_interfaceRates, Kinetics::nPhases(), Kinetics::nReactions(), Kinetics::nTotalSpecies(), and Kinetics::resizeReactions().

◆ kineticsType()

virtual std::string kineticsType ( ) const
inlinevirtual

Identifies the Kinetics manager type.

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

Reimplemented from Kinetics.

Reimplemented in EdgeKinetics.

Definition at line 74 of file InterfaceKinetics.h.

◆ setElectricPotential()

void setElectricPotential ( int  n,
doublereal  V 
)

Set the electric potential in the nth phase.

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

Definition at line 56 of file InterfaceKinetics.cpp.

◆ getEquilibriumConstants()

void getEquilibriumConstants ( doublereal *  kc)
virtual

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 203 of file InterfaceKinetics.cpp.

References Kinetics::getReactionDelta(), InterfaceKinetics::m_mu0_Kc, Kinetics::nReactions(), Kinetics::reactionPhaseIndex(), ThermoPhase::RT(), Kinetics::thermo(), and InterfaceKinetics::updateMu0().

Referenced by InterfaceKinetics::getRevRateConstants().

◆ updateExchangeCurrentQuantities()

void updateExchangeCurrentQuantities ( )

◆ getDeltaGibbs()

void getDeltaGibbs ( doublereal *  deltaG)
virtual

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 408 of file InterfaceKinetics.cpp.

References Kinetics::getReactionDelta(), InterfaceKinetics::m_deltaG, InterfaceKinetics::m_mu, Kinetics::m_start, Kinetics::m_thermo, Kinetics::nPhases(), and Kinetics::nReactions().

◆ getDeltaElectrochemPotentials()

void getDeltaElectrochemPotentials ( doublereal *  deltaM)
virtual

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 425 of file InterfaceKinetics.cpp.

References ThermoPhase::getElectrochemPotentials(), Kinetics::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

◆ getDeltaEnthalpy()

void getDeltaEnthalpy ( doublereal *  deltaH)
virtual

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 436 of file InterfaceKinetics.cpp.

References ThermoPhase::getPartialMolarEnthalpies(), Kinetics::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

◆ getDeltaEntropy()

void getDeltaEntropy ( doublereal *  deltaS)
virtual

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 447 of file InterfaceKinetics.cpp.

References ThermoPhase::getPartialMolarEntropies(), Kinetics::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

◆ getDeltaSSGibbs()

void getDeltaSSGibbs ( doublereal *  deltaG)
virtual

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 458 of file InterfaceKinetics.cpp.

References Kinetics::getReactionDelta(), ThermoPhase::getStandardChemPotentials(), InterfaceKinetics::m_mu0, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

◆ getDeltaSSEnthalpy()

void getDeltaSSEnthalpy ( doublereal *  deltaH)
virtual

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 472 of file InterfaceKinetics.cpp.

References ThermoPhase::getEnthalpy_RT(), Kinetics::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_kk, Kinetics::m_start, Kinetics::nPhases(), Kinetics::reactionPhaseIndex(), ThermoPhase::RT(), and Kinetics::thermo().

◆ getDeltaSSEntropy()

void getDeltaSSEntropy ( doublereal *  deltaS)
virtual

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 489 of file InterfaceKinetics.cpp.

References Cantera::GasConstant, ThermoPhase::getEntropy_R(), Kinetics::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_kk, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

◆ getActivityConcentrations()

void getActivityConcentrations ( doublereal *const  conc)
virtual

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 146 of file InterfaceKinetics.cpp.

References InterfaceKinetics::_update_rates_C(), and InterfaceKinetics::m_actConc.

◆ electrochem_beta()

doublereal electrochem_beta ( size_t  irxn) const

Return the charge transfer rxn Beta parameter for the ith reaction.

Returns the beta parameter for a charge transfer reaction. This parameter is not important for non-charge transfer reactions. Note, the parameter defaults to zero. However, a value of 0.5 should be supplied for every charge transfer reaction if no information is known, as a value of 0.5 pertains to a symmetric transition state. The value can vary between 0 to 1.

Parameters
irxnReaction number in the kinetics mechanism
Returns
Beta parameter. This defaults to zero, even for charge transfer reactions.
Deprecated:
To be removed after Cantera 2.6. Parameter should be accessed from ReactionRate object instead.

Definition at line 761 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_beta, InterfaceKinetics::m_ctrxn, and Cantera::warn_deprecated().

◆ isReversible()

virtual bool isReversible ( size_t  i)
inlinevirtual

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 139 of file InterfaceKinetics.h.

References InterfaceKinetics::m_revindex.

◆ getFwdRateConstants()

void getFwdRateConstants ( doublereal *  kfwd)
virtual

Return the forward rate constants.

The computed values include all temperature-dependent, pressure-dependent, and third body contributions. 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().
Deprecated:
Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin and Glarborg, 'Chemically Reacting Flow', Wiley Interscience, 2003). For new behavior, set 'Cantera::use_legacy_rate_constants(false)'.

Reimplemented from Kinetics.

Definition at line 302 of file InterfaceKinetics.cpp.

References Kinetics::m_perturb, Kinetics::m_rfn, Kinetics::nReactions(), and InterfaceKinetics::updateROP().

Referenced by InterfaceKinetics::getRevRateConstants().

◆ getRevRateConstants()

void getRevRateConstants ( doublereal *  krev,
bool  doIrreversible = false 
)
virtual

Return the reverse rate constants.

The computed values include all temperature-dependent, pressure-dependent, and third body contributions. 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.
Deprecated:
Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin and Glarborg, 'Chemically Reacting Flow', Wiley Interscience, 2003). For new behavior, set 'Cantera::use_legacy_rate_constants(false)'.

Reimplemented from Kinetics.

Definition at line 311 of file InterfaceKinetics.cpp.

References InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::getFwdRateConstants(), Kinetics::m_rkcn, Kinetics::m_ropnet, and Kinetics::nReactions().

◆ effectivePreExponentialFactor()

double effectivePreExponentialFactor ( size_t  irxn)
inline

Return effective preexponent for the specified reaction.

Returns effective preexponent, accounting for surface coverage dependencies.

Parameters
irxnReaction number in the kinetics mechanism
Returns
Effective preexponent
Deprecated:
To be removed after Cantera 2.6.

Definition at line 162 of file InterfaceKinetics.h.

References InterfaceKinetics::m_interfaceRates, and InterfaceKinetics::m_rates.

◆ effectiveActivationEnergy_R()

double effectiveActivationEnergy_R ( size_t  irxn)
inline

Return effective activation energy for the specified reaction.

Returns effective activation energy, accounting for surface coverage dependencies.

Parameters
irxnReaction number in the kinetics mechanism
Returns
Effective activation energy divided by the gas constant
Deprecated:
To be removed after Cantera 2.6.

Definition at line 181 of file InterfaceKinetics.h.

References InterfaceKinetics::m_interfaceRates, and InterfaceKinetics::m_rates.

◆ effectiveTemperatureExponent()

double effectiveTemperatureExponent ( size_t  irxn)
inline

Return effective temperature exponent for the specified reaction.

Returns effective temperature exponent, accounting for surface coverage dependencies. Current parameterization in SurfaceArrhenius does not change this parameter with the change in surface coverages.

Parameters
irxnReaction number in the kinetics mechanism
Returns
Effective temperature exponent
Deprecated:
To be removed after Cantera 2.6.

Definition at line 201 of file InterfaceKinetics.h.

References InterfaceKinetics::m_interfaceRates, and InterfaceKinetics::m_rates.

◆ addPhase()

void addPhase ( ThermoPhase thermo)
virtual

Add a phase to the kinetics manager object.

This must be done before the function init() is called or before any reactions are input.

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

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

Reimplemented from Kinetics.

Definition at line 718 of file InterfaceKinetics.cpp.

References Kinetics::addPhase(), InterfaceKinetics::m_phaseExists, InterfaceKinetics::m_phaseIsStable, and Kinetics::thermo().

Referenced by InterfaceKinetics::InterfaceKinetics().

◆ init()

void init ( )
virtual

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 725 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_nDim, InterfaceKinetics::m_surf, Phase::nDim(), Cantera::npos, Kinetics::reactionPhaseIndex(), and Kinetics::thermo().

Referenced by InterfaceKinetics::addReaction().

◆ resizeSpecies()

void resizeSpecies ( )
virtual

◆ addReaction()

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

◆ modifyReaction()

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

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 586 of file InterfaceKinetics.cpp.

References InterfaceKinetics::buildSurfaceArrhenius(), InterfaceKinetics::m_interfaceRates, InterfaceKinetics::m_interfaceTypes, InterfaceKinetics::m_rates, InterfaceKinetics::m_temp, Kinetics::modifyReaction(), and Cantera::SURFACE_RXN.

◆ updateROP()

void updateROP ( )
virtual

◆ _update_rates_T()

void _update_rates_T ( )

Update properties that depend on temperature.

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

Definition at line 62 of file InterfaceKinetics.cpp.

References InterfaceKinetics::_update_rates_phi(), SurfPhase::getCoverages(), InterfaceKinetics::m_actConc, InterfaceKinetics::m_has_coverage_dependence, InterfaceKinetics::m_rates, and InterfaceKinetics::m_surf.

Referenced by InterfaceKinetics::updateROP().

◆ _update_rates_phi()

void _update_rates_phi ( )

Update properties that depend on the electric potential.

Definition at line 116 of file InterfaceKinetics.cpp.

References ThermoPhase::electricPotential(), InterfaceKinetics::m_phi, Kinetics::nPhases(), and Kinetics::thermo().

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::updateMu0().

◆ _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 127 of file InterfaceKinetics.cpp.

References ThermoPhase::getActivityConcentrations(), Phase::getConcentrations(), InterfaceKinetics::m_actConc, InterfaceKinetics::m_conc, Kinetics::m_start, Kinetics::m_thermo, and Kinetics::nPhases().

Referenced by InterfaceKinetics::getActivityConcentrations(), and InterfaceKinetics::updateROP().

◆ advanceCoverages()

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.

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 775 of file InterfaceKinetics.cpp.

References ImplicitSurfChem::integrate(), InterfaceKinetics::m_integrator, ImplicitSurfChem::setMaxErrTestFails(), ImplicitSurfChem::setMaxSteps(), ImplicitSurfChem::setMaxStepSize(), and ImplicitSurfChem::setTolerances().

◆ solvePseudoSteadyStateProblem()

void solvePseudoSteadyStateProblem ( int  ifuncOverride = -1,
doublereal  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 792 of file InterfaceKinetics.cpp.

References ImplicitSurfChem::initialize(), InterfaceKinetics::m_integrator, and ImplicitSurfChem::solvePseudoSteadyStateProblem().

◆ setIOFlag()

void setIOFlag ( int  ioFlag)

Definition at line 710 of file InterfaceKinetics.cpp.

◆ updateMu0()

void updateMu0 ( )
virtual

◆ 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 152 of file InterfaceKinetics.cpp.

References Kinetics::getRevReactionDelta(), InterfaceKinetics::m_irrev, InterfaceKinetics::m_mu0_Kc, InterfaceKinetics::m_revindex, Kinetics::m_rkcn, Cantera::npos, Kinetics::nReactions(), Kinetics::reactionPhaseIndex(), ThermoPhase::RT(), Kinetics::thermo(), and InterfaceKinetics::updateMu0().

◆ applyVoltageKfwdCorrection()

void applyVoltageKfwdCorrection ( doublereal *const  kfwd)

Apply modifications for the forward reaction rate for interfacial charge transfer reactions.

For reactions that transfer charge across a potential difference, the activation energies are modified by the potential difference. (see, for example, ...). This method applies this correction.

Parameters
kfwdVector of forward reaction rate constants on which to have the voltage correction applied
Deprecated:
To be removed after Cantera 2.6.

Definition at line 246 of file InterfaceKinetics.cpp.

References Phase::charge(), InterfaceKinetics::deltaElectricEnergy_, Cantera::Faraday, Kinetics::getReactionDelta(), InterfaceKinetics::m_beta, InterfaceKinetics::m_ctrxn, InterfaceKinetics::m_phi, InterfaceKinetics::m_pot, Kinetics::nPhases(), Phase::nSpecies(), Kinetics::reactionPhaseIndex(), and Kinetics::thermo().

◆ convertExchangeCurrentDensityFormulation()

void convertExchangeCurrentDensityFormulation ( doublereal *const  kfwd)

When an electrode reaction rate is optionally specified in terms of its exchange current density, adjust kfwd to the standard reaction rate constant form and units.

When the BV reaction types are used, keep the exchange current density form.

For a reaction rate constant that was given in units of Amps/m2 (exchange current density formulation with iECDFormulation == true), convert the rate to kmoles/m2/s.

For a reaction rate constant that was given in units of kmol/m2/sec when the reaction type is a Butler-Volmer form, convert it to exchange current density form (amps/m2).

Parameters
kfwdVector of forward reaction rate constants, given in either normal form or in exchange current density form.
Deprecated:
To be removed after Cantera 2.6.

Definition at line 280 of file InterfaceKinetics.cpp.

References Cantera::Faraday, InterfaceKinetics::m_beta, InterfaceKinetics::m_ctrxn, InterfaceKinetics::m_ctrxn_ecdf, InterfaceKinetics::m_deltaG0, InterfaceKinetics::m_ProdStanConcReac, Kinetics::reactionPhaseIndex(), Kinetics::thermo(), and InterfaceKinetics::updateExchangeCurrentQuantities().

◆ 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 806 of file InterfaceKinetics.cpp.

References Kinetics::checkPhaseIndex(), InterfaceKinetics::m_phaseExists, InterfaceKinetics::m_phaseExistsCheck, and InterfaceKinetics::m_phaseIsStable.

◆ 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 837 of file InterfaceKinetics.cpp.

References Kinetics::checkPhaseIndex(), and InterfaceKinetics::m_phaseIsStable.

◆ 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 825 of file InterfaceKinetics.cpp.

References Kinetics::checkPhaseIndex(), and InterfaceKinetics::m_phaseExists.

◆ 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 831 of file InterfaceKinetics.cpp.

References Kinetics::checkPhaseIndex(), and InterfaceKinetics::m_phaseIsStable.

◆ buildSurfaceArrhenius()

SurfaceArrhenius buildSurfaceArrhenius ( size_t  i,
InterfaceReaction2 r,
bool  replace 
)
protected

◆ applyStickingCorrection()

void applyStickingCorrection ( double  T,
double *  kf 
)
protected

Definition at line 847 of file InterfaceKinetics.cpp.

Member Data Documentation

◆ m_grt

vector_fp m_grt
protected

◆ m_revindex

std::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 429 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::isReversible(), and InterfaceKinetics::updateKc().

◆ m_rates

Rate1<SurfaceArrhenius> m_rates
protected

◆ m_redo_rates

bool m_redo_rates
protected

Definition at line 438 of file InterfaceKinetics.h.

◆ m_interfaceRates

std::vector<unique_ptr<MultiRateBase> > m_interfaceRates
protected

◆ m_interfaceTypes

std::map<std::string, size_t> m_interfaceTypes
protected

Rate handler mapping.

Definition at line 442 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::modifyReaction().

◆ m_irrev

std::vector<size_t> m_irrev
protected

Vector of irreversible reaction numbers.

vector containing the reaction numbers of irreversible reactions.

Definition at line 448 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::updateKc().

◆ m_conc

vector_fp 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 462 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_C(), and InterfaceKinetics::resizeSpecies().

◆ m_actConc

vector_fp 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 476 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_C(), InterfaceKinetics::_update_rates_T(), InterfaceKinetics::getActivityConcentrations(), InterfaceKinetics::resizeSpecies(), and InterfaceKinetics::updateROP().

◆ m_mu0

vector_fp 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 485 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::getDeltaSSGibbs(), InterfaceKinetics::resizeSpecies(), InterfaceKinetics::updateExchangeCurrentQuantities(), and InterfaceKinetics::updateMu0().

◆ m_mu

vector_fp 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 494 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::getDeltaGibbs(), and InterfaceKinetics::resizeSpecies().

◆ m_mu0_Kc

vector_fp 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 508 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::resizeSpecies(), InterfaceKinetics::updateKc(), and InterfaceKinetics::updateMu0().

◆ m_phi

vector_fp 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 515 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_phi(), InterfaceKinetics::applyVoltageKfwdCorrection(), InterfaceKinetics::resizeSpecies(), and InterfaceKinetics::updateMu0().

◆ m_pot

vector_fp m_pot
protected

Vector of potential energies due to Voltages.

Length is the number of species in kinetics mech. It's used to store the potential energy due to the voltage.

Definition at line 522 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::applyVoltageKfwdCorrection(), and InterfaceKinetics::resizeSpecies().

◆ deltaElectricEnergy_

vector_fp deltaElectricEnergy_
protected

Storage for the net electric energy change due to reaction.

Length is number of reactions. It's used to store the net electric potential energy change due to the reaction.

deltaElectricEnergy_[jrxn] = sum_i ( F V_i z_i nu_ij)

Deprecated:
To be removed after Cantera 2.6.

Definition at line 533 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::applyVoltageKfwdCorrection().

◆ m_surf

SurfPhase* m_surf
protected

◆ m_integrator

ImplicitSurfChem* m_integrator
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 544 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::advanceCoverages(), and InterfaceKinetics::solvePseudoSteadyStateProblem().

◆ m_beta

vector_fp m_beta
protected

Electrochemical transfer coefficient for the forward direction.

Electrochemical transfer coefficient for all reactions that have transfer reactions the reaction is given by m_ctrxn[i]

Deprecated:
To be removed after Cantera 2.6.

Definition at line 553 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::applyVoltageKfwdCorrection(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), and InterfaceKinetics::electrochem_beta().

◆ m_ctrxn

std::vector<size_t> m_ctrxn
protected

Vector of reaction indexes specifying the id of the charge transfer reactions in the mechanism.

Vector of reaction indices which involve charge transfers. This provides an index into the m_beta array.

  irxn = m_ctrxn[i]
Deprecated:
To be removed after Cantera 2.6.

Definition at line 565 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::applyVoltageKfwdCorrection(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), and InterfaceKinetics::electrochem_beta().

◆ m_ctrxn_ecdf

vector_int m_ctrxn_ecdf
protected

Vector of booleans indicating whether the charge transfer reaction rate constant is described by an exchange current density rate constant expression.

Length is equal to the number of reactions with charge transfer coefficients, m_ctrxn[]

m_ctrxn_ecdf[irxn] = 0 This means that the rate coefficient calculator will calculate the rate constant as a chemical forward rate constant, a standard format. m_ctrxn_ecdf[irxn] = 1 this means that the rate coefficient calculator will calculate the rate constant as an exchange current density rate constant expression.

Deprecated:
To be removed after Cantera 2.6.

Definition at line 579 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::convertExchangeCurrentDensityFormulation().

◆ m_StandardConc

vector_fp m_StandardConc
protected

Vector of standard concentrations.

Length number of kinetic species units depend on the definition of the standard concentration within each phase

Deprecated:
To be removed after Cantera 2.6.

Definition at line 588 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::resizeSpecies(), and InterfaceKinetics::updateExchangeCurrentQuantities().

◆ m_deltaG0

vector_fp m_deltaG0
protected

Vector of delta G^0, the standard state Gibbs free energies for each reaction.

Length is the number of reactions units are Joule kmol-1

Deprecated:
To be removed after Cantera 2.6.

Definition at line 597 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), and InterfaceKinetics::updateExchangeCurrentQuantities().

◆ m_deltaG

vector_fp m_deltaG
protected

Vector of deltaG[] of reaction, the delta Gibbs free energies for each reaction.

Length is the number of reactions units are Joule kmol-1

Deprecated:
To be removed after Cantera 2.6.

Definition at line 606 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::getDeltaGibbs().

◆ m_ProdStanConcReac

vector_fp m_ProdStanConcReac
protected

Vector of the products of the standard concentrations of the reactants.

Units vary wrt what the units of the standard concentrations are Length = number of reactions.

Deprecated:
To be removed after Cantera 2.6.

Definition at line 615 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), and InterfaceKinetics::updateExchangeCurrentQuantities().

◆ m_ROP_ok

bool m_ROP_ok
protected

Definition at line 617 of file InterfaceKinetics.h.

◆ m_temp

doublereal m_temp
protected

Current temperature of the data.

Definition at line 620 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::modifyReaction().

◆ m_logtemp

doublereal m_logtemp
protected

Current log of the temperature.

Definition at line 623 of file InterfaceKinetics.h.

◆ m_has_coverage_dependence

bool m_has_coverage_dependence
protected

Boolean flag indicating whether any reaction in the mechanism has a coverage dependent forward reaction rate.

If this is true, then the coverage dependence is multiplied into the forward reaction rates constant

Definition at line 631 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::addReaction().

◆ m_has_electrochem_rxns

bool m_has_electrochem_rxns
protected

Boolean flag indicating whether any reaction in the mechanism has a beta electrochemical parameter.

If this is true, the Butler-Volmer correction is applied to the forward reaction rate for those reactions.

fac = exp ( - beta * (delta_phi))

Deprecated:
To be removed after Cantera 2.6.

Definition at line 643 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction().

◆ m_has_exchange_current_density_formulation

bool m_has_exchange_current_density_formulation
protected

Boolean flag indicating whether any reaction in the mechanism is described by an exchange current density expression.

If this is true, the standard state Gibbs free energy of the reaction and the product of the reactant standard concentrations must be precalculated in order to calculate the rate constant.

Deprecated:
To be removed after Cantera 2.6.

Definition at line 654 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction().

◆ m_phaseExistsCheck

int m_phaseExistsCheck
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 663 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::setPhaseExistence(), and InterfaceKinetics::updateROP().

◆ m_phaseExists

std::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 676 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addPhase(), InterfaceKinetics::phaseExistence(), InterfaceKinetics::setPhaseExistence(), and InterfaceKinetics::updateROP().

◆ 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 686 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addPhase(), InterfaceKinetics::phaseStability(), InterfaceKinetics::setPhaseExistence(), InterfaceKinetics::setPhaseStability(), and InterfaceKinetics::updateROP().

◆ m_rxnPhaseIsReactant

std::vector<std::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 694 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::updateROP().

◆ m_rxnPhaseIsProduct

std::vector<std::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 702 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), and InterfaceKinetics::updateROP().

◆ m_stickingData

std::vector<StickData> m_stickingData
protected

Data for sticking reactions.

Definition at line 713 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::buildSurfaceArrhenius().

◆ m_ioFlag

int m_ioFlag
protected

Definition at line 717 of file InterfaceKinetics.h.

◆ m_nDim

size_t m_nDim
protected

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

Definition at line 721 of file InterfaceKinetics.h.

Referenced by EdgeKinetics::EdgeKinetics(), and InterfaceKinetics::init().


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