Cantera  2.1.2
Public 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]

Public Member Functions

 InterfaceKinetics (thermo_t *thermo=0)
 Constructor. More...
 
virtual ~InterfaceKinetics ()
 Destructor. More...
 
 InterfaceKinetics (const InterfaceKinetics &right)
 Copy Constructor for the Kinetics object. More...
 
InterfaceKineticsoperator= (const InterfaceKinetics &right)
 Assignment operator. More...
 
virtual KineticsduplMyselfAsKinetics (const std::vector< thermo_t * > &tpVector) const
 Duplication routine for objects which inherit from Kinetics. More...
 
virtual int type () const
 Identifies the kinetics manager type. More...
 
void setElectricPotential (int n, doublereal V)
 Set the electric potential in the nth phase. More...
 
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)
 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)
 
void checkPartialEquil ()
 
size_t reactionNumber () const
 
void addElementaryReaction (ReactionData &r)
 
void addGlobalReaction (const ReactionData &r)
 
void installReagents (const ReactionData &r)
 
void updateKc ()
 Update the equilibrium constants in molar units for all reversible reactions. More...
 
void registerReaction (size_t rxnNumber, int type, size_t loc)
 Write values into m_index. More...
 
void applyButlerVolmerCorrection (doublereal *const kf)
 Apply corrections for interfacial charge transfer reactions. More...
 
void applyExchangeCurrentDensityFormulation (doublereal *const kfwd)
 When an electrode reaction rate is optionally specified in terms of its exchange current density, extra vectors need to be precalculated. 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 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 getEquilibriumConstants (doublereal *kc)
 Return a vector of Equilibrium constants. More...
 
void getExchangeCurrentQuantities ()
 
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...
 
Species Production Rates
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 *net)
 Species net production rates [kmol/m^3/s or kmol/m^2/s]. More...
 
Reaction Mechanism Informational Query Routines
virtual doublereal reactantStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a reactant in reaction i. More...
 
virtual doublereal productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i. More...
 
virtual int reactionType (size_t i) const
 Flag specifying the type of reaction. More...
 
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 std::string reactionString (size_t i) const
 Return a string representing the reaction. 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...
 
virtual void getActivationEnergies (doublereal *E)
 Return the activation energies in Kelvin. More...
 
Reaction Mechanism Construction
virtual void addPhase (thermo_t &thermo)
 Add a phase to the kinetics manager object. More...
 
virtual void init ()
 Prepare the class for the addition of reactions. More...
 
virtual void addReaction (ReactionData &r)
 Add a single reaction to the mechanism. More...
 
virtual void finalize ()
 Finish adding reactions and prepare for use. More...
 
virtual bool ready () const
 Returns true if the kinetics manager has been properly initialized and finalized. More...
 
- Public Member Functions inherited from Kinetics
void incrementRxnCount ()
 Increment the number of reactions in the mechanism by one. More...
 
void selectPhase (const doublereal *data, const thermo_t *phase, doublereal *phase_data)
 
 Kinetics ()
 Default constructor. More...
 
virtual ~Kinetics ()
 Destructor. More...
 
 Kinetics (const Kinetics &)
 Copy Constructor for the Kinetics object. More...
 
Kineticsoperator= (const Kinetics &right)
 Assignment operator. More...
 
virtual void assignShallowPointers (const std::vector< thermo_t * > &tpVector)
 Reassign the pointers within the Kinetics object. More...
 
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)
 Return the phase index of a phase in the list of phases defined within the object. More...
 
size_t surfacePhaseIndex ()
 This returns the integer index of the phase which has ThermoPhase type cSurf. More...
 
size_t reactionPhaseIndex ()
 Phase where the reactions occur. More...
 
thermo_tthermo (size_t n=0)
 This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism. More...
 
const thermo_tthermo (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...
 
thermo_tspeciesPhase (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...
 
thermo_tspeciesPhase (size_t k)
 This function takes as an argument the kineticsSpecies index (i.e., 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)
 This function takes as an argument the kineticsSpecies index (i.e., 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 getReactionDelta (const doublereal *property, doublereal *deltaProperty)
 Change in species properties. 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 const std::vector
< size_t > & 
reactants (size_t i) const
 Returns a read-only reference to the vector of reactant index numbers for reaction i. More...
 
virtual const std::vector
< size_t > & 
products (size_t i) const
 Returns a read-only reference to the vector of product index numbers for reaction i. More...
 
virtual const std::vector
< grouplist_t > & 
reactantGroups (size_t i)
 
virtual const std::vector
< grouplist_t > & 
productGroups (size_t i)
 
doublereal multiplier (size_t i) const
 The current value of the multiplier for reaction i. More...
 
void setMultiplier (size_t i, doublereal f)
 Set the multiplier for reaction i to f. 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::map< size_t, std::pair
< int, size_t > > 
m_index
 Vector of information about reactions in the mechanism. More...
 
std::vector< size_t > m_irrev
 Vector of irreversible reaction numbers. More...
 
ReactionStoichMgr m_rxnstoich
 Stoichiometric manager for the reaction mechanism. More...
 
size_t m_nirrev
 Number of irreversible reactions in the mechanism. More...
 
size_t m_nrev
 Number of reversible reactions in the mechanism. More...
 
std::vector< std::map< size_t,
doublereal > > 
m_rrxn
 m_rrxn is a vector of maps, containing the reactant stoichiometric coefficient information More...
 
std::vector< std::map< size_t,
doublereal > > 
m_prxn
 m_prxn is a vector of maps, containing the reactant stoichiometric coefficient information More...
 
std::vector< std::string > m_rxneqn
 String expression for each rxn. More...
 
vector_fp m_conc
 an array of generalized concentrations for each species More...
 
vector_fp m_mu0
 Vector of standard state chemical potentials. 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 m_rwork
 Vector temporary. More...
 
vector_fp m_E
 Vector of raw activation energies for the reactions. More...
 
SurfPhasem_surf
 Pointer to the single surface phase. More...
 
ImplicitSurfChemm_integrator
 Pointer to the Implicit surface chemistry object. More...
 
vector_fp m_beta
 
std::vector< size_t > m_ctrxn
 Vector of reaction indexes specifying the id of the current transfer reactions in the mechanism. More...
 
vector_int m_ctrxn_ecdf
 Vector of booleans indicating whether the charge transfer reaction may be described by an exchange current density expression. More...
 
vector_fp m_StandardConc
 
vector_fp m_deltaG0
 
vector_fp m_ProdStanConcReac
 
doublereal m_logp0
 
doublereal m_logc0
 
vector_fp m_ropf
 
vector_fp m_ropr
 
vector_fp m_ropnet
 
bool m_ROP_ok
 
doublereal m_temp
 Current temperature of the data. More...
 
doublereal m_logtemp
 Current log of the temperature. More...
 
vector_fp m_rfn
 
vector_fp m_rkcn
 
bool m_finalized
 boolean indicating whether mechanism has been finalized 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...
 
std::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...
 
int m_ioFlag
 
- Protected Attributes inherited from Kinetics
size_t m_ii
 Number of reactions in the mechanism. 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< std::vector
< size_t > > 
m_reactants
 This is a vector of vectors containing the reactants for each reaction. More...
 
std::vector< std::vector
< size_t > > 
m_products
 This is a vector of vectors containing the products for each reaction. More...
 
std::vector< thermo_t * > 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 id, i.e., the id attribute in the xml phase element 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...
 

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, i.e., 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 (i.e., 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 63 of file InterfaceKinetics.h.

Constructor & Destructor Documentation

InterfaceKinetics ( thermo_t 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().

Referenced by InterfaceKinetics::duplMyselfAsKinetics().

~InterfaceKinetics ( )
virtual

Destructor.

Definition at line 54 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_integrator.

Copy Constructor for the Kinetics object.

Definition at line 59 of file InterfaceKinetics.cpp.

References InterfaceKinetics::operator=().

Member Function Documentation

InterfaceKinetics & operator= ( const InterfaceKinetics right)
Kinetics * duplMyselfAsKinetics ( const std::vector< thermo_t * > &  tpVector) const
virtual

Duplication routine for objects which inherit from Kinetics.

This function can be used to duplicate objects derived from Kinetics even if the application only has a pointer to Kinetics to work with.

These routines are basically wrappers around the derived copy constructor.

Parameters
tpVectorVector of pointers to ThermoPhase objects. this is the m_thermo vector within this object

Reimplemented from Kinetics.

Reimplemented in EdgeKinetics.

Definition at line 161 of file InterfaceKinetics.cpp.

References Kinetics::assignShallowPointers(), and InterfaceKinetics::InterfaceKinetics().

int type ( ) const
virtual

Identifies the kinetics manager type.

Each class derived from Kinetics should overload this method to return a unique integer. Standard values are defined in file mix_defs.h.

Reimplemented from Kinetics.

Reimplemented in EdgeKinetics.

Definition at line 156 of file InterfaceKinetics.cpp.

Referenced by InterfaceKinetics::registerReaction().

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

References ThermoPhase::setElectricPotential(), and Kinetics::thermo().

void getFwdRatesOfProgress ( doublereal *  fwdROP)
virtual

Return the forward rates of progress of the reactions.

Forward rates of progress. Return the forward rates of progress in array fwdROP, which must be dimensioned at least as large as the total number of reactions.

Parameters
fwdROPOutput vector containing forward rates of progress of the reactions. Length: m_ii.

Reimplemented from Kinetics.

Definition at line 310 of file InterfaceKinetics.cpp.

References InterfaceKinetics::updateROP().

void getRevRatesOfProgress ( doublereal *  revROP)
virtual

Return the Reverse rates of progress of the reactions.

Return the reverse rates of progress in array revROP, which must be dimensioned at least as large as the total number of reactions.

Parameters
revROPOutput vector containing reverse rates of progress of the reactions. Length: m_ii.

Reimplemented from Kinetics.

Definition at line 316 of file InterfaceKinetics.cpp.

References InterfaceKinetics::updateROP().

void getNetRatesOfProgress ( doublereal *  netROP)
virtual

Net rates of progress.

Return the net (forward - reverse) rates of progress in array netROP, which must be dimensioned at least as large as the total number of reactions.

Parameters
netROPOutput vector of the net ROP. Length: m_ii.

Reimplemented from Kinetics.

Definition at line 322 of file InterfaceKinetics.cpp.

References InterfaceKinetics::updateROP().

void getEquilibriumConstants ( doublereal *  kc)
virtual

Return a vector of Equilibrium constants.

Return the equilibrium constants of the reactions in concentration units in array kc, which must be dimensioned at least as large as the total number of reactions.

Parameters
kcOutput vector containing the equilibrium constants. Length: m_ii.

Reimplemented from Kinetics.

Definition at line 328 of file InterfaceKinetics.cpp.

References Phase::charge(), DATA_PTR, Cantera::GasConstant, ReactionStoichMgr::getReactionDelta(), ThermoPhase::getStandardChemPotentials(), ThermoPhase::logStandardConc(), Kinetics::m_ii, InterfaceKinetics::m_mu0, InterfaceKinetics::m_phi, InterfaceKinetics::m_rxnstoich, Kinetics::m_start, Kinetics::nPhases(), Phase::nSpecies(), Phase::temperature(), and Kinetics::thermo().

Referenced by InterfaceKinetics::getRevRateConstants().

void getDeltaGibbs ( doublereal *  deltaG)
virtual

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

These values depend upon the concentration of the solution.

units = J kmol-1

Parameters
deltaGOutput vector of deltaG's for reactions Length: m_ii.

Reimplemented from Kinetics.

Definition at line 609 of file InterfaceKinetics.cpp.

References DATA_PTR, ThermoPhase::getChemPotentials(), ReactionStoichMgr::getReactionDelta(), InterfaceKinetics::m_grt, InterfaceKinetics::m_rxnstoich, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

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: m_ii.

Reimplemented from Kinetics.

Definition at line 629 of file InterfaceKinetics.cpp.

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

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: m_ii.

Reimplemented from Kinetics.

Definition at line 646 of file InterfaceKinetics.cpp.

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

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: m_ii.

Reimplemented from Kinetics.

Definition at line 662 of file InterfaceKinetics.cpp.

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

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: m_ii.

Reimplemented from Kinetics.

Definition at line 678 of file InterfaceKinetics.cpp.

References DATA_PTR, ReactionStoichMgr::getReactionDelta(), ThermoPhase::getStandardChemPotentials(), InterfaceKinetics::m_grt, InterfaceKinetics::m_rxnstoich, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

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: m_ii.

Reimplemented from Kinetics.

Definition at line 696 of file InterfaceKinetics.cpp.

References DATA_PTR, Cantera::GasConstant, ThermoPhase::getEnthalpy_RT(), ReactionStoichMgr::getReactionDelta(), InterfaceKinetics::m_grt, Kinetics::m_kk, InterfaceKinetics::m_rxnstoich, Kinetics::m_start, Kinetics::nPhases(), Phase::temperature(), and Kinetics::thermo().

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: m_ii.

Reimplemented from Kinetics.

Definition at line 718 of file InterfaceKinetics.cpp.

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

void getCreationRates ( doublereal *  cdot)
virtual

Species creation rates [kmol/m^3/s or kmol/m^2/s].

Return the species creation rates in array cdot, which must be dimensioned at least as large as the total number of species in all phases.

See Also
nTotalSpecies.
Parameters
cdotOutput vector of creation rates. Length: m_kk.

Reimplemented from Kinetics.

Definition at line 381 of file InterfaceKinetics.cpp.

References ReactionStoichMgr::getCreationRates(), Kinetics::m_kk, InterfaceKinetics::m_rxnstoich, and InterfaceKinetics::updateROP().

void getDestructionRates ( doublereal *  ddot)
virtual

Species destruction rates [kmol/m^3/s or kmol/m^2/s].

Return the species destruction rates in array ddot, which must be dimensioned at least as large as the total number of species.

See Also
nTotalSpecies.
Parameters
ddotOutput vector of destruction rates. Length: m_kk.

Reimplemented from Kinetics.

Definition at line 387 of file InterfaceKinetics.cpp.

References ReactionStoichMgr::getDestructionRates(), Kinetics::m_kk, InterfaceKinetics::m_rxnstoich, and InterfaceKinetics::updateROP().

void getNetProductionRates ( doublereal *  wdot)
virtual

Species net production rates [kmol/m^3/s or kmol/m^2/s].

Return the species net production rates (creation - destruction) in array wdot, which must be dimensioned at least as large as the total number of species.

See Also
nTotalSpecies.
Parameters
wdotOutput vector of net production rates. Length: m_kk.

Reimplemented from Kinetics.

Definition at line 393 of file InterfaceKinetics.cpp.

References ReactionStoichMgr::getNetProductionRates(), Kinetics::m_kk, InterfaceKinetics::m_rxnstoich, and InterfaceKinetics::updateROP().

Referenced by solveSP::calc_t(), ReactingSurf1D::eval(), solveSP::fun_eval(), solveProb::print_header(), and solveProb::printIteration().

virtual doublereal reactantStoichCoeff ( size_t  k,
size_t  i 
) const
inlinevirtual

Stoichiometric coefficient of species k as a reactant in reaction i.

Parameters
kkinetic species index
ireaction index

Reimplemented from Kinetics.

Definition at line 130 of file InterfaceKinetics.h.

References InterfaceKinetics::m_rrxn.

virtual doublereal productStoichCoeff ( size_t  k,
size_t  i 
) const
inlinevirtual

Stoichiometric coefficient of species k as a product in reaction i.

Parameters
kkinetic species index
ireaction index

Reimplemented from Kinetics.

Definition at line 134 of file InterfaceKinetics.h.

References InterfaceKinetics::m_prxn.

virtual int reactionType ( size_t  i) const
inlinevirtual

Flag specifying the type of reaction.

The legal values and their meaning are specific to the particular kinetics manager.

Parameters
ireaction index

Reimplemented from Kinetics.

Definition at line 138 of file InterfaceKinetics.h.

References InterfaceKinetics::m_index.

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

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

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.

Definition at line 1026 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_ctrxn.

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

References InterfaceKinetics::m_revindex.

virtual std::string reactionString ( size_t  i) const
inlinevirtual

Return a string representing the reaction.

Parameters
ireaction index

Reimplemented from Kinetics.

Definition at line 169 of file InterfaceKinetics.h.

References InterfaceKinetics::m_rxneqn.

void getFwdRateConstants ( doublereal *  kfwd)
virtual

Return the forward rate constants.

length is the number of reactions. units depends on many issues.

Parameters
kfwdOutput vector containing the forward reaction rate constants. Length: m_ii.

Reimplemented from Kinetics.

Definition at line 472 of file InterfaceKinetics.cpp.

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

Referenced by InterfaceKinetics::getRevRateConstants().

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

Return the reverse rate constants.

length is the number of reactions. units depends on many issues. 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 485 of file InterfaceKinetics.cpp.

References InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::getFwdRateConstants(), Kinetics::m_ii, Cantera::multiply_each(), and Kinetics::nReactions().

void getActivationEnergies ( doublereal *  E)
virtual

Return the activation energies in Kelvin.

length is the number of reactions

Parameters
EOuptut vector of activation energies. Length: m_ii.
Deprecated:
To be removed in Cantera 2.2.

Reimplemented from Kinetics.

Definition at line 498 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_E, and Cantera::warn_deprecated().

void addPhase ( thermo_t 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 969 of file InterfaceKinetics.cpp.

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

Referenced by InterfaceKinetics::InterfaceKinetics().

void init ( )
virtual

Prepare the class for the addition of reactions.

This method is called by importKinetics() after all phases have been added but before any reactions have been. The base class method does nothing, but derived classes may use this to perform any initialization (allocating arrays, etc.) that requires knowing the phases and species, but before any reactions are added.

Reimplemented from Kinetics.

Definition at line 976 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_conc, InterfaceKinetics::m_grt, Kinetics::m_kk, InterfaceKinetics::m_mu0, InterfaceKinetics::m_phi, InterfaceKinetics::m_pot, InterfaceKinetics::m_prxn, InterfaceKinetics::m_rrxn, Kinetics::nPhases(), Phase::nSpecies(), and Kinetics::thermo().

void addReaction ( ReactionData r)
virtual

Add a single reaction to the mechanism.

This routine must be called after init() and before finalize().

Parameters
rReference to the ReactionData object for the reaction to be added.

Reimplemented from Kinetics.

Definition at line 739 of file InterfaceKinetics.cpp.

References ReactionData::equation, Kinetics::incrementRxnCount(), InterfaceKinetics::m_rxneqn, InterfaceKinetics::m_rxnPhaseIsProduct, InterfaceKinetics::m_rxnPhaseIsReactant, Kinetics::nPhases(), Kinetics::products(), Kinetics::reactants(), and Kinetics::speciesPhaseIndex().

void finalize ( )
virtual

Finish adding reactions and prepare for use.

This method is called by importKinetics() after all reactions have been entered into the mechanism and before the mechanism is used to calculate reaction rates. The base class method does nothing, but derived classes may use this to perform any initialization (allocating arrays, etc.) that must be done after the reactions are entered.

Reimplemented from Kinetics.

Reimplemented in EdgeKinetics.

Definition at line 991 of file InterfaceKinetics.cpp.

References Kinetics::finalize(), Cantera::int2str(), InterfaceKinetics::m_finalized, Kinetics::m_kk, Kinetics::m_perturb, InterfaceKinetics::m_phaseExists, InterfaceKinetics::m_rwork, InterfaceKinetics::m_surf, Kinetics::m_thermo, Phase::nDim(), Cantera::npos, Kinetics::nReactions(), Kinetics::reactionPhaseIndex(), and Kinetics::thermo().

bool ready ( ) const
virtual

Returns true if the kinetics manager has been properly initialized and finalized.

Reimplemented from Kinetics.

Reimplemented in Interface.

Definition at line 1036 of file InterfaceKinetics.cpp.

References InterfaceKinetics::m_finalized.

void updateROP ( )
void _update_rates_T ( )
void _update_rates_phi ( )

Update properties that depend on the electric potential.

Definition at line 200 of file InterfaceKinetics.cpp.

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

Referenced by InterfaceKinetics::_update_rates_T().

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

References DATA_PTR, ThermoPhase::getActivityConcentrations(), InterfaceKinetics::m_conc, Kinetics::m_start, Kinetics::nPhases(), and Kinetics::thermo().

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

void advanceCoverages ( doublereal  tstep)

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

Definition at line 1042 of file InterfaceKinetics.cpp.

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

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

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

void updateKc ( )

Update the equilibrium constants in molar units for all reversible reactions.

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

References Phase::charge(), DATA_PTR, Cantera::GasConstant, ReactionStoichMgr::getRevReactionDelta(), ThermoPhase::getStandardChemPotentials(), Cantera::int2str(), ThermoPhase::logStandardConc(), Kinetics::m_ii, InterfaceKinetics::m_irrev, InterfaceKinetics::m_mu0, InterfaceKinetics::m_nirrev, InterfaceKinetics::m_nrev, InterfaceKinetics::m_phi, InterfaceKinetics::m_revindex, InterfaceKinetics::m_rxnstoich, Kinetics::m_start, Kinetics::nPhases(), Cantera::npos, Kinetics::nReactions(), Phase::nSpecies(), Phase::temperature(), and Kinetics::thermo().

Referenced by InterfaceKinetics::_update_rates_T().

void registerReaction ( size_t  rxnNumber,
int  type,
size_t  loc 
)
inline

Write values into m_index.

Parameters
rxnNumberreaction number
typereaction type
loclocation ??

Definition at line 293 of file InterfaceKinetics.h.

References InterfaceKinetics::m_index, and InterfaceKinetics::type().

void applyButlerVolmerCorrection ( doublereal *const  kf)

Apply corrections 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
kfVector of forward reaction rate constants on which to have the correction applied

Definition at line 399 of file InterfaceKinetics.cpp.

References Phase::charge(), DATA_PTR, Cantera::fp2str(), Cantera::GasConstant, ReactionStoichMgr::getReactionDelta(), Cantera::int2str(), InterfaceKinetics::m_ctrxn, InterfaceKinetics::m_E, InterfaceKinetics::m_phi, InterfaceKinetics::m_pot, InterfaceKinetics::m_rwork, InterfaceKinetics::m_rxnstoich, Kinetics::nPhases(), Phase::nSpecies(), Phase::temperature(), Kinetics::thermo(), and Cantera::writelog().

Referenced by InterfaceKinetics::_update_rates_T().

void applyExchangeCurrentDensityFormulation ( doublereal *const  kfwd)

When an electrode reaction rate is optionally specified in terms of its exchange current density, extra vectors need to be precalculated.

Definition at line 455 of file InterfaceKinetics.cpp.

References Cantera::GasConstant, InterfaceKinetics::m_ctrxn, InterfaceKinetics::m_ctrxn_ecdf, Phase::temperature(), and Kinetics::thermo().

Referenced by InterfaceKinetics::_update_rates_T().

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

References InterfaceKinetics::m_phaseExists, InterfaceKinetics::m_phaseExistsCheck, InterfaceKinetics::m_phaseIsStable, and Kinetics::m_thermo.

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

References InterfaceKinetics::m_phaseIsStable, and Kinetics::m_thermo.

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

References InterfaceKinetics::m_phaseExists, and Kinetics::m_thermo.

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

References InterfaceKinetics::m_phaseIsStable, and Kinetics::m_thermo.

Member Data Documentation

vector_fp m_grt
protected
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 373 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::isReversible(), InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

Rate1<SurfaceArrhenius> m_rates
protected

Templated class containing the vector of reactions for this interface.

The templated class is described in RateCoeffMgr.h The class SurfaceArrhenius is described in RxnRates.h

Definition at line 380 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

std::map<size_t, std::pair<int, size_t> > m_index
mutableprotected

Vector of information about reactions in the mechanism.

The key is the reaction index (0 < i < m_ii). The first pair is the reactionType of the reaction. The second pair is ...

Definition at line 390 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::operator=(), InterfaceKinetics::reactionType(), and InterfaceKinetics::registerReaction().

std::vector<size_t> m_irrev
protected

Vector of irreversible reaction numbers.

vector containing the reaction numbers of irreversible reactions.

Definition at line 396 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

ReactionStoichMgr m_rxnstoich
protected
size_t m_nirrev
protected

Number of irreversible reactions in the mechanism.

Definition at line 407 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

size_t m_nrev
protected

Number of reversible reactions in the mechanism.

Definition at line 410 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

std::vector<std::map<size_t, doublereal> > m_rrxn
mutableprotected

m_rrxn is a vector of maps, containing the reactant stoichiometric coefficient information

m_rrxn has a length equal to the total number of species in the kinetics object. For each species, there exists a map, with the reaction number being the key, and the reactant stoichiometric coefficient for the species being the value.

HKM -> mutable because search sometimes creates extra entries. To be fixed in future...

Definition at line 423 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::init(), InterfaceKinetics::operator=(), and InterfaceKinetics::reactantStoichCoeff().

std::vector<std::map<size_t, doublereal> > m_prxn
mutableprotected

m_prxn is a vector of maps, containing the reactant stoichiometric coefficient information

m_prxn is a vector of maps. m_prxn has a length equal to the total number of species in the kinetics object. For each species, there exists a map, with the reaction number being the key, and the product stoichiometric coefficient for the species being the value.

Definition at line 433 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::init(), InterfaceKinetics::operator=(), and InterfaceKinetics::productStoichCoeff().

std::vector<std::string> m_rxneqn
protected

String expression for each rxn.

Vector of strings of length m_ii, the number of reactions, containing the string expressions for each reaction (e.g., reactants <=> product1 + product2)

Definition at line 441 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::operator=(), and InterfaceKinetics::reactionString().

vector_fp m_conc
protected

an array of generalized concentrations for each species

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

Referenced by InterfaceKinetics::_update_rates_C(), InterfaceKinetics::_update_rates_T(), InterfaceKinetics::getActivityConcentrations(), InterfaceKinetics::init(), InterfaceKinetics::operator=(), and InterfaceKinetics::updateROP().

vector_fp m_mu0
protected

Vector of standard state chemical potentials.

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

Length = m_k. Units = J/kmol.

Definition at line 464 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::init(), InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

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

Referenced by InterfaceKinetics::_update_rates_phi(), InterfaceKinetics::applyButlerVolmerCorrection(), InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::init(), InterfaceKinetics::operator=(), and InterfaceKinetics::updateKc().

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

Referenced by InterfaceKinetics::applyButlerVolmerCorrection(), InterfaceKinetics::init(), and InterfaceKinetics::operator=().

vector_fp m_rwork
protected

Vector temporary.

Length is number of reactions. It's used to store the voltage contribution to the activation energy.

Definition at line 487 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::applyButlerVolmerCorrection(), EdgeKinetics::finalize(), InterfaceKinetics::finalize(), and InterfaceKinetics::operator=().

vector_fp m_E
protected

Vector of raw activation energies for the reactions.

units are in Kelvin Length is number of reactions.

Definition at line 494 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::applyButlerVolmerCorrection(), InterfaceKinetics::getActivationEnergies(), and InterfaceKinetics::operator=().

SurfPhase* m_surf
protected

Pointer to the single surface phase.

Definition at line 497 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), EdgeKinetics::finalize(), InterfaceKinetics::finalize(), and InterfaceKinetics::operator=().

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 objects's surface problem uncoupled from other surface phases.

Definition at line 505 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::advanceCoverages(), InterfaceKinetics::operator=(), InterfaceKinetics::solvePseudoSteadyStateProblem(), and InterfaceKinetics::~InterfaceKinetics().

std::vector<size_t> m_ctrxn
protected

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

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

  irxn = m_ctrxn[i]

Definition at line 517 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::applyButlerVolmerCorrection(), InterfaceKinetics::applyExchangeCurrentDensityFormulation(), InterfaceKinetics::electrochem_beta(), and InterfaceKinetics::operator=().

vector_int m_ctrxn_ecdf
protected

Vector of booleans indicating whether the charge transfer reaction may be described by an exchange current density expression.

Definition at line 521 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::applyExchangeCurrentDensityFormulation(), and InterfaceKinetics::operator=().

doublereal m_temp
protected

Current temperature of the data.

Definition at line 536 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

doublereal m_logtemp
protected

Current log of the temperature.

Definition at line 538 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

bool m_finalized
protected

boolean indicating whether mechanism has been finalized

Definition at line 543 of file InterfaceKinetics.h.

Referenced by EdgeKinetics::finalize(), InterfaceKinetics::finalize(), InterfaceKinetics::operator=(), and InterfaceKinetics::ready().

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

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

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))

Definition at line 561 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

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.

Definition at line 570 of file InterfaceKinetics.h.

Referenced by InterfaceKinetics::_update_rates_T(), and InterfaceKinetics::operator=().

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

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

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

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

Definition at line 588 of file InterfaceKinetics.h.

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

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

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

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

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

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

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


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