Public interface for kinetics managers. More...
#include <Kinetics.h>
Public interface for kinetics managers.
This class serves as a base class to derive 'kinetics managers', which are classes that manage homogeneous chemistry within one phase, or heterogeneous chemistry at one interface. The virtual methods of this class are meant to be overloaded in subclasses. The non-virtual methods perform generic functions and are implemented in Kinetics. They should not be overloaded. Only those methods required by a subclass need to be overloaded; the rest will throw exceptions if called.
When the nomenclature "kinetics species index" is used below, this means that the species index ranges over all species in all phases handled by the kinetics manager.
Definition at line 124 of file Kinetics.h.
Public Member Functions | |
virtual pair< size_t, size_t > | checkDuplicates (bool throw_err=true) const |
Check for unmarked duplicate reactions and unmatched marked duplicates. | |
virtual void | setRoot (shared_ptr< Solution > root) |
Set root Solution holding all phase information. | |
shared_ptr< Solution > | root () const |
Get the Solution object containing this Kinetics object and associated ThermoPhase objects. | |
Constructors and General Information about Mechanism | |
Kinetics ()=default | |
Default constructor. | |
Kinetics (const Kinetics &)=delete | |
Kinetics objects are not copyable or assignable. | |
Kinetics & | operator= (const Kinetics &)=delete |
virtual string | kineticsType () const |
Identifies the Kinetics manager type. | |
virtual void | resizeReactions () |
Finalize Kinetics object and associated objects. | |
size_t | nReactions () const |
Number of reactions in the reaction mechanism. | |
void | checkReactionIndex (size_t m) const |
Check that the specified reaction index is in range Throws an exception if i is greater than nReactions() | |
void | checkReactionArraySize (size_t ii) const |
Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions(). | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range Throws an exception if k is greater than nSpecies()-1. | |
void | checkSpeciesArraySize (size_t mm) const |
Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). | |
Information/Lookup Functions about Phases and Species | |
size_t | nPhases () const |
The number of phases participating in the reaction mechanism. | |
void | checkPhaseIndex (size_t m) const |
Check that the specified phase index is in range Throws an exception if m is greater than nPhases() | |
void | checkPhaseArraySize (size_t mm) const |
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). | |
size_t | phaseIndex (const string &ph) const |
Return the phase index of a phase in the list of phases defined within the object. | |
size_t | reactionPhaseIndex () const |
Phase where the reactions occur. | |
shared_ptr< ThermoPhase > | reactionPhase () const |
Return pointer to phase where the reactions occur. | |
ThermoPhase & | thermo (size_t n=0) |
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism. | |
const ThermoPhase & | thermo (size_t n=0) const |
size_t | nTotalSpecies () const |
The total number of species in all phases participating in the kinetics mechanism. | |
size_t | kineticsSpeciesIndex (size_t k, size_t n) const |
The location of species k of phase n in species arrays. | |
string | kineticsSpeciesName (size_t k) const |
Return the name of the kth species in the kinetics manager. | |
size_t | kineticsSpeciesIndex (const string &nm) const |
This routine will look up a species number based on the input string nm. | |
ThermoPhase & | speciesPhase (const string &nm) |
This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides. | |
const ThermoPhase & | speciesPhase (const string &nm) const |
ThermoPhase & | speciesPhase (size_t k) |
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the species' owning ThermoPhase object. | |
size_t | speciesPhaseIndex (size_t k) const |
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the index of the phase owning the species. | |
Reaction Rates Of Progress | |
virtual void | getFwdRatesOfProgress (double *fwdROP) |
Return the forward rates of progress of the reactions. | |
virtual void | getRevRatesOfProgress (double *revROP) |
Return the Reverse rates of progress of the reactions. | |
virtual void | getNetRatesOfProgress (double *netROP) |
Net rates of progress. | |
virtual void | getEquilibriumConstants (double *kc) |
Return a vector of Equilibrium constants. | |
virtual void | getReactionDelta (const double *property, double *deltaProperty) const |
Change in species properties. | |
virtual void | getRevReactionDelta (const double *g, double *dg) const |
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the reversible reactions. | |
virtual void | getDeltaGibbs (double *deltaG) |
Return the vector of values for the reaction Gibbs free energy change. | |
virtual void | getDeltaElectrochemPotentials (double *deltaM) |
Return the vector of values for the reaction electrochemical free energy change. | |
virtual void | getDeltaEnthalpy (double *deltaH) |
Return the vector of values for the reactions change in enthalpy. | |
virtual void | getDeltaEntropy (double *deltaS) |
Return the vector of values for the reactions change in entropy. | |
virtual void | getDeltaSSGibbs (double *deltaG) |
Return the vector of values for the reaction standard state Gibbs free energy change. | |
virtual void | getDeltaSSEnthalpy (double *deltaH) |
Return the vector of values for the change in the standard state enthalpies of reaction. | |
virtual void | getDeltaSSEntropy (double *deltaS) |
Return the vector of values for the change in the standard state entropies for each reaction. | |
virtual void | getThirdBodyConcentrations (double *concm) |
Return a vector of values of effective concentrations of third-body collision partners of any reaction. | |
virtual const vector< double > & | thirdBodyConcentrations () const |
Provide direct access to current third-body concentration values. | |
Species Production Rates | |
virtual void | getCreationRates (double *cdot) |
Species creation rates [kmol/m^3/s or kmol/m^2/s]. | |
virtual void | getDestructionRates (double *ddot) |
Species destruction rates [kmol/m^3/s or kmol/m^2/s]. | |
virtual void | getNetProductionRates (double *wdot) |
Species net production rates [kmol/m^3/s or kmol/m^2/s]. | |
Kinetics derivatives are calculated with respect to temperature, pressure, molar concentrations and species mole fractions for forward/reverse/net rates of progress as well as creation/destruction and net production of species. The following suffixes are used to indicate derivatives:
Source term derivatives are based on a generic rate-of-progress expression for the \( i \)-th reaction \( R_i \), which is a function of temperature \( T \), pressure \( P \) and molar concentrations \( C_j \): \[ R_i = k_{f,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^\prime} - k_{r,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^{\prime\prime}} \] Forward/reverse rate expressions \( k_{f,i} \) and \( k_{r,i} \) are implemented by ReactionRate specializations; forward/reverse stoichiometric coefficients are \( \nu_{ji}^\prime \) and \( \nu_{ji}^{\prime\prime} \). Unless the reaction involves third-body colliders, \( \nu_{M,i} = 0 \). For three-body reactions, effective ThirdBody collider concentrations \( C_M \) are considered with \( \nu_{M,i} = 1 \). For more detailed information on relevant theory, see, for example, Perini, et al. [31] or Niemeyer, et al. [29], although specifics of Cantera's implementation may differ. Partial derivatives are obtained from the product rule, where resulting terms consider reaction rate derivatives, derivatives of the concentration product term, and, if applicable, third-body term derivatives. ReactionRate specializations may implement exact derivatives (example: ArrheniusRate::ddTScaledFromStruct) or approximate them numerically (examples: ReactionData::perturbTemperature, PlogData::perturbPressure, FalloffData::perturbThirdBodies). Derivatives of concentration and third-body terms are based on analytic expressions. Species creation and destruction rates are obtained by multiplying rate-of-progress vectors by stoichiometric coefficient matrices. As this is a linear operation, it is possible to calculate derivatives the same way. All derivatives are calculated for source terms while holding other properties constant, independent of whether equation of state or \( \sum X_k = 1 \) constraints are satisfied. Thus, derivatives deviate from Jacobians and numerical derivatives that implicitly enforce these constraints. Depending on application and equation of state, derivatives can nevertheless be used to obtain Jacobians, for example:
While some applications require exact derivatives, others can tolerate approximate derivatives that neglect terms to increase computational speed and/or improve Jacobian sparsity (example: AdaptivePreconditioner). Derivative evaluations settings are accessible by keyword/value pairs using the methods getDerivativeSettings() and setDerivativeSettings(). For BulkKinetics, the following keyword/value pairs are supported:
For InterfaceKinetics, the following keyword/value pairs are supported:
| |
virtual void | getDerivativeSettings (AnyMap &settings) const |
Retrieve derivative settings. | |
virtual void | setDerivativeSettings (const AnyMap &settings) |
Set/modify derivative settings. | |
virtual void | getFwdRateConstants_ddT (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions. | |
virtual void | getFwdRateConstants_ddP (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions. | |
virtual void | getFwdRateConstants_ddC (double *dkfwd) |
Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions. | |
virtual void | getFwdRatesOfProgress_ddT (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. | |
virtual void | getFwdRatesOfProgress_ddP (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. | |
virtual void | getFwdRatesOfProgress_ddC (double *drop) |
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. | |
virtual Eigen::SparseMatrix< double > | fwdRatesOfProgress_ddX () |
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
virtual Eigen::SparseMatrix< double > | fwdRatesOfProgress_ddCi () |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. | |
virtual void | getRevRatesOfProgress_ddT (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. | |
virtual void | getRevRatesOfProgress_ddP (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. | |
virtual void | getRevRatesOfProgress_ddC (double *drop) |
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. | |
virtual Eigen::SparseMatrix< double > | revRatesOfProgress_ddX () |
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
virtual Eigen::SparseMatrix< double > | revRatesOfProgress_ddCi () |
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. | |
virtual void | getNetRatesOfProgress_ddT (double *drop) |
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. | |
virtual void | getNetRatesOfProgress_ddP (double *drop) |
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. | |
virtual void | getNetRatesOfProgress_ddC (double *drop) |
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. | |
virtual Eigen::SparseMatrix< double > | netRatesOfProgress_ddX () |
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
virtual Eigen::SparseMatrix< double > | netRatesOfProgress_ddCi () |
Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations. | |
void | getCreationRates_ddT (double *dwdot) |
Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions. | |
void | getCreationRates_ddP (double *dwdot) |
Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions. | |
void | getCreationRates_ddC (double *dwdot) |
Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions. | |
Eigen::SparseMatrix< double > | creationRates_ddX () |
Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
Eigen::SparseMatrix< double > | creationRates_ddCi () |
Calculate derivatives for species creation rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. | |
void | getDestructionRates_ddT (double *dwdot) |
Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions. | |
void | getDestructionRates_ddP (double *dwdot) |
Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions. | |
void | getDestructionRates_ddC (double *dwdot) |
Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions. | |
Eigen::SparseMatrix< double > | destructionRates_ddX () |
Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
Eigen::SparseMatrix< double > | destructionRates_ddCi () |
Calculate derivatives for species destruction rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. | |
void | getNetProductionRates_ddT (double *dwdot) |
Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions. | |
void | getNetProductionRates_ddP (double *dwdot) |
Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions. | |
void | getNetProductionRates_ddC (double *dwdot) |
Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions. | |
Eigen::SparseMatrix< double > | netProductionRates_ddX () |
Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration. | |
Eigen::SparseMatrix< double > | netProductionRates_ddCi () |
Calculate derivatives for species net production rates with respect to species concentration at constant temperature, pressure, and concentration of all other species. | |
Reaction Mechanism Informational Query Routines | |
virtual double | reactantStoichCoeff (size_t k, size_t i) const |
Stoichiometric coefficient of species k as a reactant in reaction i. | |
Eigen::SparseMatrix< double > | reactantStoichCoeffs () const |
Stoichiometric coefficient matrix for reactants. | |
virtual double | productStoichCoeff (size_t k, size_t i) const |
Stoichiometric coefficient of species k as a product in reaction i. | |
Eigen::SparseMatrix< double > | productStoichCoeffs () const |
Stoichiometric coefficient matrix for products. | |
Eigen::SparseMatrix< double > | revProductStoichCoeffs () const |
Stoichiometric coefficient matrix for products of reversible reactions. | |
virtual double | reactantOrder (size_t k, size_t i) const |
Reactant order of species k in reaction i. | |
virtual double | productOrder (int k, int i) const |
product Order of species k in reaction i. | |
virtual void | getActivityConcentrations (double *const conc) |
Get the vector of activity concentrations used in the kinetics object. | |
virtual bool | isReversible (size_t i) |
True if reaction i has been declared to be reversible. | |
virtual void | getFwdRateConstants (double *kfwd) |
Return the forward rate constants. | |
virtual void | getRevRateConstants (double *krev, bool doIrreversible=false) |
Return the reverse rate constants. | |
Reaction Mechanism Construction | |
virtual void | addThermo (shared_ptr< ThermoPhase > thermo) |
Add a phase to the kinetics manager object. | |
virtual void | init () |
Prepare the class for the addition of reactions, after all phases have been added. | |
AnyMap | parameters () |
Return the parameters for a phase definition which are needed to reconstruct an identical object using the newKinetics function. | |
virtual void | resizeSpecies () |
Resize arrays with sizes that depend on the total number of species. | |
virtual bool | addReaction (shared_ptr< Reaction > r, bool resize=true) |
Add a single reaction to the mechanism. | |
virtual void | modifyReaction (size_t i, shared_ptr< Reaction > rNew) |
Modify the rate expression associated with a reaction. | |
shared_ptr< Reaction > | reaction (size_t i) |
Return the Reaction object for reaction i. | |
shared_ptr< const Reaction > | reaction (size_t i) const |
void | skipUndeclaredSpecies (bool skip) |
Determine behavior when adding a new reaction that contains species not defined in any of the phases associated with this kinetics manager. | |
bool | skipUndeclaredSpecies () const |
void | skipUndeclaredThirdBodies (bool skip) |
Determine behavior when adding a new reaction that contains third-body efficiencies for species not defined in any of the phases associated with this kinetics manager. | |
bool | skipUndeclaredThirdBodies () const |
Altering Reaction Rates | |
These methods alter reaction rates. They are designed primarily for carrying out sensitivity analysis, but may be used for any purpose requiring dynamic alteration of rate constants. For each reaction, a real-valued multiplier may be defined that multiplies the reaction rate coefficient. The multiplier may be set to zero to completely remove a reaction from the mechanism. | |
double | multiplier (size_t i) const |
The current value of the multiplier for reaction i. | |
virtual void | setMultiplier (size_t i, double f) |
Set the multiplier for reaction i to f. | |
virtual void | invalidateCache () |
Protected Member Functions | |
virtual void | updateROP () |
double | checkDuplicateStoich (map< int, double > &r1, map< int, double > &r2) const |
Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reactions are duplicates of one another, and 0.0 otherwise. | |
Protected Attributes | |
ValueCache | m_cache |
Cache for saved calculations within each Kinetics object. | |
bool | m_ready = false |
Boolean indicating whether Kinetics object is fully configured. | |
size_t | m_kk = 0 |
The number of species in all of the phases that participate in this kinetics mechanism. | |
vector< double > | m_perturb |
Vector of perturbation factors for each reaction's rate of progress vector. | |
vector< shared_ptr< Reaction > > | m_reactions |
Vector of Reaction objects represented by this Kinetics manager. | |
vector< shared_ptr< ThermoPhase > > | m_thermo |
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator | |
vector< size_t > | m_start |
m_start is a vector of integers specifying the beginning position for the species vector for the n'th phase in the kinetics class. | |
map< string, size_t > | m_phaseindex |
Mapping of the phase name to the position of the phase within the kinetics object. | |
size_t | m_mindim = 4 |
number of spatial dimensions of lowest-dimensional phase. | |
vector< double > | m_rfn |
Forward rate constant for each reaction. | |
vector< double > | m_delta_gibbs0 |
Delta G^0 for all reactions. | |
vector< double > | m_rkcn |
Reciprocal of the equilibrium constant in concentration units. | |
vector< double > | m_ropf |
Forward rate-of-progress for each reaction. | |
vector< double > | m_ropr |
Reverse rate-of-progress for each reaction. | |
vector< double > | m_ropnet |
Net rate-of-progress for each reaction. | |
vector< double > | m_dH |
The enthalpy change for each reaction to calculate Blowers-Masel rates. | |
vector< double > | m_rbuf |
Buffer used for storage of intermediate reaction-specific results. | |
bool | m_skipUndeclaredSpecies = false |
See skipUndeclaredSpecies() | |
bool | m_skipUndeclaredThirdBodies = false |
See skipUndeclaredThirdBodies() | |
bool | m_hasUndeclaredThirdBodies = false |
Flag indicating whether reactions include undeclared third bodies. | |
std::weak_ptr< Solution > | m_root |
reference to Solution | |
Stoichiometry management | |
These objects and functions handle turning reaction extents into species production rates and also handle turning thermo properties into reaction thermo properties. | |
StoichManagerN | m_reactantStoich |
Stoichiometry manager for the reactants for each reaction. | |
StoichManagerN | m_productStoich |
Stoichiometry manager for the products for each reaction. | |
StoichManagerN | m_revProductStoich |
Stoichiometry manager for the products of reversible reactions. | |
Eigen::SparseMatrix< double > | m_stoichMatrix |
Net stoichiometry (products - reactants) | |
|
default |
Default constructor.
|
inlinevirtual |
Identifies the Kinetics manager type.
Each class derived from Kinetics should override this method to return a meaningful identifier.
Reimplemented in BulkKinetics, EdgeKinetics, and InterfaceKinetics.
Definition at line 144 of file Kinetics.h.
|
virtual |
Finalize Kinetics object and associated objects.
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 35 of file Kinetics.cpp.
|
inline |
Number of reactions in the reaction mechanism.
Definition at line 152 of file Kinetics.h.
void checkReactionIndex | ( | size_t | m | ) | const |
Check that the specified reaction index is in range Throws an exception if i is greater than nReactions()
Definition at line 27 of file Kinetics.cpp.
void checkReactionArraySize | ( | size_t | ii | ) | const |
Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions().
Used before calls which take an array pointer.
Definition at line 54 of file Kinetics.cpp.
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.
Definition at line 88 of file Kinetics.cpp.
void checkSpeciesArraySize | ( | size_t | mm | ) | const |
Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies().
Used before calls which take an array pointer.
Definition at line 95 of file Kinetics.cpp.
|
inline |
The number of phases participating in the reaction mechanism.
For a homogeneous reaction mechanism, this will always return 1, but for a heterogeneous mechanism it will return the total number of phases in the mechanism.
Definition at line 184 of file Kinetics.h.
void checkPhaseIndex | ( | size_t | m | ) | const |
Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
Definition at line 62 of file Kinetics.cpp.
void checkPhaseArraySize | ( | size_t | mm | ) | const |
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases().
Used before calls which take an array pointer.
Definition at line 69 of file Kinetics.cpp.
|
inline |
Return the phase index of a phase in the list of phases defined within the object.
ph | string name of the phase |
If a -1 is returned, then the phase is not defined in the Kinetics object.
Definition at line 206 of file Kinetics.h.
size_t reactionPhaseIndex | ( | ) | const |
Phase where the reactions occur.
For heterogeneous mechanisms, one of the phases in the list of phases represents the 2D interface or 1D edge at which the reactions take place. This method returns the index of the phase with the smallest spatial dimension (1, 2, or 3) among the list of phases. If there is more than one, the index of the first one is returned. For homogeneous mechanisms, the value 0 is returned.
Definition at line 76 of file Kinetics.cpp.
shared_ptr< ThermoPhase > reactionPhase | ( | ) | const |
Return pointer to phase where the reactions occur.
Definition at line 83 of file Kinetics.cpp.
|
inline |
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
It is typically used so that member functions of the ThermoPhase object may be called. For homogeneous mechanisms, there is only one object, and this method can be called without an argument to access it.
n | Index of the ThermoPhase being sought. |
Definition at line 242 of file Kinetics.h.
|
inline |
Definition at line 245 of file Kinetics.h.
|
inline |
The total number of species in all phases participating in the kinetics mechanism.
This is useful to dimension arrays for use in calls to methods that return the species production rates, for example.
Definition at line 254 of file Kinetics.h.
|
inline |
The location of species k of phase n in species arrays.
Kinetics manager classes return species production rates in flat arrays, with the species of each phases following one another, in the order the phases were added. This method is useful to find the value for a particular species of a particular phase in arrays returned from methods like getCreationRates that return an array of species-specific quantities.
Example: suppose a heterogeneous mechanism involves three phases. The first contains 12 species, the second 26, and the third 3. Then species arrays must have size at least 41, and positions 0 - 11 are the values for the species in the first phase, positions 12 - 37 are the values for the species in the second phase, etc. Then kineticsSpeciesIndex(7, 0) = 7, kineticsSpeciesIndex(4, 1) = 16, and kineticsSpeciesIndex(2, 2) = 40.
k | species index |
n | phase index for the species |
Definition at line 276 of file Kinetics.h.
string kineticsSpeciesName | ( | size_t | k | ) | const |
Return the name of the kth species in the kinetics manager.
k is an integer from 0 to ktot - 1, where ktot is the number of species in the kinetics manager, which is the sum of the number of species in all phases participating in the kinetics manager. If k is out of bounds, the string "<unknown>" is returned.
k | species index |
Definition at line 238 of file Kinetics.cpp.
size_t kineticsSpeciesIndex | ( | const string & | nm | ) | const |
This routine will look up a species number based on the input string nm.
The lookup of species will occur for all phases listed in the kinetics object.
return
nm | Input string name of the species |
Definition at line 248 of file Kinetics.cpp.
ThermoPhase & speciesPhase | ( | const string & | nm | ) |
This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides.
Will throw an error if the species doesn't match.
nm | String containing the name of the species. |
Definition at line 260 of file Kinetics.cpp.
const ThermoPhase & speciesPhase | ( | const string & | nm | ) | const |
Definition at line 271 of file Kinetics.cpp.
|
inline |
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.
k | Species index |
Definition at line 321 of file Kinetics.h.
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.
k | Species index |
Definition at line 281 of file Kinetics.cpp.
|
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.
fwdROP | Output vector containing forward rates of progress of the reactions. Length: nReactions(). |
Definition at line 302 of file Kinetics.cpp.
|
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.
revROP | Output vector containing reverse rates of progress of the reactions. Length: nReactions(). |
Definition at line 308 of file Kinetics.cpp.
|
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.
netROP | Output vector of the net ROP. Length: nReactions(). |
Definition at line 314 of file Kinetics.cpp.
|
inlinevirtual |
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.
\[ Kc_i = \exp [ \Delta G_{ss,i} ] \prod(Cs_k) \exp(\sum_k \nu_{k,i} F \phi_n) \]
kc | Output vector containing the equilibrium constants. Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 381 of file Kinetics.h.
|
virtual |
Change in species properties.
Given an array of molar species property values \( z_k, k = 1, \dots, K \), return the array of reaction values
\[ \Delta Z_i = \sum_k \nu_{k,i} z_k, i = 1, \dots, I. \]
For example, if this method is called with the array of standard-state molar Gibbs free energies for the species, then the values returned in array deltaProperty
would be the standard-state Gibbs free energies of reaction for each reaction.
property | Input vector of property value. Length: m_kk. |
deltaProperty | Output vector of deltaRxn. Length: nReactions(). |
Definition at line 320 of file Kinetics.cpp.
|
virtual |
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the reversible reactions.
Array 'g' must have a length at least as great as the number of species, and array 'dg' must have a length as great as the total number of reactions. This method only computes 'dg' for the reversible reactions, and the entries of 'dg' for the irreversible reactions are unaltered. This is primarily designed for use in calculating reverse rate coefficients from thermochemistry for reversible reactions.
Definition at line 329 of file Kinetics.cpp.
|
inlinevirtual |
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
deltaG | Output vector of deltaG's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 423 of file Kinetics.h.
|
inlinevirtual |
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
deltaM | Output vector of deltaM's for reactions Length: nReactions(). |
Reimplemented in InterfaceKinetics.
Definition at line 438 of file Kinetics.h.
|
inlinevirtual |
Return the vector of values for the reactions change in enthalpy.
These values depend upon the concentration of the solution.
units = J kmol-1
deltaH | Output vector of deltaH's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 451 of file Kinetics.h.
|
inlinevirtual |
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
deltaS | Output vector of deltaS's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 464 of file Kinetics.h.
|
inlinevirtual |
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
deltaG | Output vector of ss deltaG's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 478 of file Kinetics.h.
|
inlinevirtual |
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
deltaH | Output vector of ss deltaH's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 492 of file Kinetics.h.
|
inlinevirtual |
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
deltaS | Output vector of ss deltaS's for reactions Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 506 of file Kinetics.h.
|
inlinevirtual |
Return a vector of values of effective concentrations of third-body collision partners of any reaction.
Entries for reactions not involving third-body collision partners are not defined and set to not-a-number.
concm | Output vector of effective third-body concentrations. Length: nReactions(). |
Reimplemented in BulkKinetics.
Definition at line 518 of file Kinetics.h.
|
inlinevirtual |
Provide direct access to current third-body concentration values.
Reimplemented in BulkKinetics.
Definition at line 528 of file Kinetics.h.
|
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.
cdot | Output vector of creation rates. Length: m_kk. |
Definition at line 338 of file Kinetics.cpp.
|
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.
ddot | Output vector of destruction rates. Length: m_kk. |
Definition at line 352 of file Kinetics.cpp.
|
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.
wdot | Output vector of net production rates. Length: m_kk. |
Definition at line 363 of file Kinetics.cpp.
|
virtual |
Stoichiometric coefficient of species k as a reactant in reaction i.
k | kinetic species index |
i | reaction index |
Definition at line 292 of file Kinetics.cpp.
|
inline |
Stoichiometric coefficient matrix for reactants.
Definition at line 1140 of file Kinetics.h.
|
virtual |
Stoichiometric coefficient of species k as a product in reaction i.
k | kinetic species index |
i | reaction index |
Definition at line 297 of file Kinetics.cpp.
|
inline |
Stoichiometric coefficient matrix for products.
Definition at line 1155 of file Kinetics.h.
|
inline |
Stoichiometric coefficient matrix for products of reversible reactions.
Definition at line 1162 of file Kinetics.h.
|
inlinevirtual |
Reactant order of species k in reaction i.
This is the nominal order of the activity concentration in determining the forward rate of progress of the reaction
k | kinetic species index |
i | reaction index |
Definition at line 1174 of file Kinetics.h.
|
inlinevirtual |
product Order of species k in reaction i.
This is the nominal order of the activity concentration of species k in determining the reverse rate of progress of the reaction i
For irreversible reactions, this will all be zero.
k | kinetic species index |
i | reaction index |
Definition at line 1188 of file Kinetics.h.
|
inlinevirtual |
Get the vector of activity concentrations used in the kinetics object.
[out] | conc | Vector of activity concentrations. Length is equal to the number of species in the kinetics object |
Reimplemented in InterfaceKinetics.
Definition at line 1197 of file Kinetics.h.
|
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.
i | reaction index |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 1208 of file Kinetics.h.
|
inlinevirtual |
Return the forward rate constants.
The computed values include all temperature-dependent and pressure-dependent contributions. By default, third-body concentrations are only considered if they are part of the reaction rate definition; for a legacy implementation that includes third-body concentrations see Cantera::use_legacy_rate_constants(). Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.
kfwd | Output vector containing the forward reaction rate constants. Length: nReactions(). |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 1225 of file Kinetics.h.
|
inlinevirtual |
Return the reverse rate constants.
The computed values include all temperature-dependent and pressure-dependent contributions. By default, third-body concentrations are only considered if they are part of the reaction rate definition; for a legacy implementation that includes third-body concentrations see Cantera::use_legacy_rate_constants(). Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction. Note, this routine will return rate constants for irreversible reactions if the default for doIrreversible
is overridden.
krev | Output vector of reverse rate constants |
doIrreversible | boolean indicating whether irreversible reactions should be included. |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 1245 of file Kinetics.h.
|
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. The following fields are updated:
thermo | Reference to the ThermoPhase to be added. |
Reimplemented in InterfaceKinetics.
Definition at line 520 of file Kinetics.cpp.
|
inlinevirtual |
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 in InterfaceKinetics.
Definition at line 1280 of file Kinetics.h.
AnyMap parameters | ( | ) |
Return the parameters for a phase definition which are needed to reconstruct an identical object using the newKinetics function.
This excludes the reaction definitions, which are handled separately.
Definition at line 537 of file Kinetics.cpp.
|
virtual |
Resize arrays with sizes that depend on the total number of species.
Automatically called before adding each Reaction and Phase.
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 553 of file Kinetics.cpp.
|
virtual |
Add a single reaction to the mechanism.
Derived classes should call the base class method in addition to handling their own specialized behavior.
r | Pointer to the Reaction object to be added. |
resize | If true , resizeReactions is called after reaction is added. |
true
if the reaction is added or false
if it was skipped Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 565 of file Kinetics.cpp.
|
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.
i | Index of the reaction to be modified |
rNew | Reaction with the new rate expressions |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 653 of file Kinetics.cpp.
shared_ptr< Reaction > reaction | ( | size_t | i | ) |
Return the Reaction object for reaction i.
Changes to this object do not affect the Kinetics object until the modifyReaction function is called.
Definition at line 684 of file Kinetics.cpp.
shared_ptr< const Reaction > reaction | ( | size_t | i | ) | const |
Definition at line 690 of file Kinetics.cpp.
|
inline |
Determine behavior when adding a new reaction that contains species not defined in any of the phases associated with this kinetics manager.
If set to true, the reaction will silently be ignored. If false, (the default) an exception will be raised.
Definition at line 1326 of file Kinetics.h.
|
inline |
Definition at line 1329 of file Kinetics.h.
|
inline |
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.
If set to true, the given third-body efficiency will be ignored. If false, (the default) an exception will be raised.
Definition at line 1338 of file Kinetics.h.
|
inline |
Definition at line 1341 of file Kinetics.h.
|
inline |
The current value of the multiplier for reaction i.
i | index of the reaction |
Definition at line 1360 of file Kinetics.h.
|
inlinevirtual |
Set the multiplier for reaction i to f.
i | index of the reaction |
f | value of the multiplier. |
Reimplemented in BulkKinetics, and InterfaceKinetics.
Definition at line 1369 of file Kinetics.h.
|
inlinevirtual |
Definition at line 1373 of file Kinetics.h.
|
virtual |
Check for unmarked duplicate reactions and unmatched marked duplicates.
If throw_err
is true, then an exception will be thrown if either any unmarked duplicate reactions are found, or if any marked duplicate reactions do not have a matching duplicate reaction. If throw_err
is false, the indices of the first pair of duplicate reactions found will be returned, or the index of the unmatched duplicate will be returned as both elements of the pair. If no unmarked duplicates or unmatched marked duplicate reactions are found, returns (npos, npos)
.
Definition at line 102 of file Kinetics.cpp.
|
inlinevirtual |
Set root Solution holding all phase information.
Definition at line 1391 of file Kinetics.h.
|
inline |
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition at line 1397 of file Kinetics.h.
|
inlineprotectedvirtual |
Reimplemented in InterfaceKinetics.
Definition at line 1406 of file Kinetics.h.
|
protected |
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.
r1
and r2
are maps of species key to stoichiometric coefficient, one for each reaction, where the species key is 1+k
for reactants and -1-k
for products and k
is the species index.
Definition at line 195 of file Kinetics.cpp.
|
protected |
Cache for saved calculations within each Kinetics object.
Definition at line 1403 of file Kinetics.h.
|
protected |
Stoichiometry manager for the reactants for each reaction.
Definition at line 1431 of file Kinetics.h.
|
protected |
Stoichiometry manager for the products for each reaction.
Definition at line 1434 of file Kinetics.h.
|
protected |
Stoichiometry manager for the products of reversible reactions.
Definition at line 1437 of file Kinetics.h.
|
protected |
Net stoichiometry (products - reactants)
Definition at line 1440 of file Kinetics.h.
|
protected |
Boolean indicating whether Kinetics object is fully configured.
Definition at line 1444 of file Kinetics.h.
|
protected |
The number of species in all of the phases that participate in this kinetics mechanism.
Definition at line 1448 of file Kinetics.h.
|
protected |
Vector of perturbation factors for each reaction's rate of progress vector.
It is initialized to one.
Definition at line 1452 of file Kinetics.h.
|
protected |
Vector of Reaction objects represented by this Kinetics manager.
Definition at line 1455 of file Kinetics.h.
|
protected |
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
For homogeneous kinetics applications, this vector will only have one entry. For interfacial reactions, this vector will consist of multiple entries; some of them will be surface phases, and the other ones will be bulk phases. The order that the objects are listed determines the order in which the species comprising each phase are listed in the source term vector, originating from the reaction mechanism.
Definition at line 1467 of file Kinetics.h.
|
protected |
m_start is a vector of integers specifying the beginning position for the species vector for the n'th phase in the kinetics class.
Definition at line 1473 of file Kinetics.h.
|
protected |
Mapping of the phase name to the position of the phase within the kinetics object.
Positions start with the value of 1. The member function, phaseIndex() decrements by one before returning the index value, so that missing phases return -1.
Definition at line 1481 of file Kinetics.h.
|
protected |
number of spatial dimensions of lowest-dimensional phase.
Definition at line 1484 of file Kinetics.h.
|
protected |
Forward rate constant for each reaction.
Definition at line 1487 of file Kinetics.h.
|
protected |
Delta G^0 for all reactions.
Definition at line 1490 of file Kinetics.h.
|
protected |
Reciprocal of the equilibrium constant in concentration units.
Definition at line 1493 of file Kinetics.h.
|
protected |
Forward rate-of-progress for each reaction.
Definition at line 1496 of file Kinetics.h.
|
protected |
Reverse rate-of-progress for each reaction.
Definition at line 1499 of file Kinetics.h.
|
protected |
Net rate-of-progress for each reaction.
Definition at line 1502 of file Kinetics.h.
|
protected |
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition at line 1505 of file Kinetics.h.
|
protected |
Buffer used for storage of intermediate reaction-specific results.
Definition at line 1508 of file Kinetics.h.
|
protected |
Definition at line 1511 of file Kinetics.h.
|
protected |
See skipUndeclaredThirdBodies()
Definition at line 1514 of file Kinetics.h.
|
protected |
Flag indicating whether reactions include undeclared third bodies.
Definition at line 1517 of file Kinetics.h.
|
protected |
reference to Solution
Definition at line 1520 of file Kinetics.h.