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

Public interface for kinetics managers. More...

#include <Kinetics.h>

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

Public Member Functions

virtual std::pair< size_t, size_t > checkDuplicates (bool throw_err=true) const
 Check for unmarked duplicate reactions and unmatched marked duplicates. More...
 
void selectPhase (const double *data, const ThermoPhase *phase, double *phase_data)
 
virtual void setRoot (std::shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
virtual double reactionEnthalpy (const Composition &reactants, const Composition &products)
 Calculate the reaction enthalpy of a reaction which has not necessarily been added into the Kinetics object. More...
 
Constructors and General Information about Mechanism
 Kinetics ()
 Default constructor. More...
 
virtual ~Kinetics ()
 
 Kinetics (const Kinetics &)=delete
 Kinetics objects are not copyable or assignable. More...
 
Kineticsoperator= (const Kinetics &)=delete
 
virtual std::string kineticsType () const
 Identifies the Kinetics manager type. More...
 
virtual void resizeReactions ()
 Finalize Kinetics object and associated objects. 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...
 
Information/Lookup Functions about Phases and Species
size_t nPhases () const
 The number of phases participating in the reaction mechanism. More...
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases() More...
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). More...
 
size_t phaseIndex (const std::string &ph) const
 Return the phase index of a phase in the list of phases defined within the object. More...
 
size_t surfacePhaseIndex () const
 This returns the integer index of the phase which has ThermoPhase type cSurf. More...
 
size_t reactionPhaseIndex () const
 Phase where the reactions occur. More...
 
ThermoPhasethermo (size_t n=0)
 This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism. More...
 
const ThermoPhasethermo (size_t n=0) const
 
size_t nTotalSpecies () const
 The total number of species in all phases participating in the kinetics mechanism. More...
 
size_t kineticsSpeciesIndex (size_t k, size_t n) const
 The location of species k of phase n in species arrays. More...
 
std::string kineticsSpeciesName (size_t k) const
 Return the name of the kth species in the kinetics manager. More...
 
size_t kineticsSpeciesIndex (const std::string &nm) const
 This routine will look up a species number based on the input std::string nm. More...
 
size_t kineticsSpeciesIndex (const std::string &nm, const std::string &ph) const
 This routine will look up a species number based on the input std::string nm. More...
 
ThermoPhasespeciesPhase (const std::string &nm)
 This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides. More...
 
const ThermoPhasespeciesPhase (const std::string &nm) const
 
ThermoPhasespeciesPhase (size_t k)
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the species' owning ThermoPhase object. More...
 
size_t speciesPhaseIndex (size_t k) const
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the index of the phase owning the species. More...
 
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...
 
virtual void getReactionDelta (const double *property, double *deltaProperty) const
 Change in species properties. More...
 
virtual void getRevReactionDelta (const double *g, double *dg) const
 Given an array of species properties 'g', return in array 'dg' the change in this quantity in the reversible reactions. More...
 
virtual void 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...
 
virtual void getThirdBodyConcentrations (double *concm)
 Return a vector of values of effective concentrations of third-body collision partners of any reaction. More...
 
virtual const vector_fpthirdBodyConcentrations () const
 Provide direct access to current third-body concentration values. More...
 
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 *wdot)
 Species net production rates [kmol/m^3/s or kmol/m^2/s]. More...
 
Routines to Calculate Derivatives (Jacobians)
virtual void getDerivativeSettings (AnyMap &settings) const
 Retrieve derivative settings. More...
 
virtual void setDerivativeSettings (const AnyMap &settings)
 Set/modify derivative settings. More...
 
virtual void getFwdRateConstants_ddT (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
virtual void getFwdRateConstants_ddP (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
virtual void getFwdRateConstants_ddC (double *dkfwd)
 Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
virtual void getFwdRatesOfProgress_ddT (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
virtual void getFwdRatesOfProgress_ddP (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
virtual void getFwdRatesOfProgress_ddC (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX ()
 Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
virtual void getRevRatesOfProgress_ddT (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
virtual void getRevRatesOfProgress_ddP (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
virtual void getRevRatesOfProgress_ddC (double *drop)
 Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX ()
 Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
virtual void getNetRatesOfProgress_ddT (double *drop)
 Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
virtual void getNetRatesOfProgress_ddP (double *drop)
 Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
virtual void getNetRatesOfProgress_ddC (double *drop)
 Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX ()
 Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getCreationRates_ddT (double *dwdot)
 Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getCreationRates_ddP (double *dwdot)
 Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getCreationRates_ddC (double *dwdot)
 Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > creationRates_ddX ()
 Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getDestructionRates_ddT (double *dwdot)
 Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getDestructionRates_ddP (double *dwdot)
 Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getDestructionRates_ddC (double *dwdot)
 Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > destructionRates_ddX ()
 Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getNetProductionRates_ddT (double *dwdot)
 Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddP (double *dwdot)
 Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddC (double *dwdot)
 Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > netProductionRates_ddX ()
 Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
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. More...
 
Eigen::SparseMatrix< double > reactantStoichCoeffs () const
 Stoichiometric coefficient matrix for reactants. More...
 
virtual double productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i. More...
 
Eigen::SparseMatrix< double > productStoichCoeffs () const
 Stoichiometric coefficient matrix for products. More...
 
Eigen::SparseMatrix< double > revProductStoichCoeffs () const
 Stoichiometric coefficient matrix for products of reversible reactions. More...
 
virtual doublereal reactantOrder (size_t k, size_t i) const
 Reactant order of species k in reaction i. More...
 
virtual doublereal productOrder (int k, int i) const
 product Order of species k in reaction i. More...
 
virtual void getActivityConcentrations (doublereal *const conc)
 Get the vector of activity concentrations used in the kinetics object. More...
 
virtual int reactionType (size_t i) const
 Flag specifying the type of reaction. More...
 
virtual std::string reactionTypeStr (size_t i) const
 String specifying the type of reaction. More...
 
virtual bool isReversible (size_t i)
 True if reaction i has been declared to be reversible. More...
 
std::string reactionString (size_t i) const
 Return a string representing the reaction. More...
 
std::string reactantString (size_t i) const
 Returns a string containing the reactants side of the reaction equation. More...
 
std::string productString (size_t i) const
 Returns a string containing the products side of the reaction equation. More...
 
virtual void getFwdRateConstants (double *kfwd)
 Return the forward rate constants. More...
 
virtual void getRevRateConstants (double *krev, bool doIrreversible=false)
 Return the reverse rate constants. More...
 
Reaction Mechanism Construction
virtual void addPhase (ThermoPhase &thermo)
 Add a phase to the kinetics manager object. More...
 
virtual void init ()
 Prepare the class for the addition of reactions, after all phases have been added. More...
 
AnyMap parameters ()
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newKinetics function. More...
 
virtual void resizeSpecies ()
 Resize arrays with sizes that depend on the total number of species. More...
 
virtual bool addReaction (shared_ptr< Reaction > r, bool resize=true)
 Add a single reaction to the mechanism. More...
 
virtual void modifyReaction (size_t i, shared_ptr< Reaction > rNew)
 Modify the rate expression associated with a reaction. More...
 
shared_ptr< Reactionreaction (size_t i)
 Return the Reaction object for reaction i. More...
 
shared_ptr< const Reactionreaction (size_t i) const
 
void skipUndeclaredSpecies (bool skip)
 Determine behavior when adding a new reaction that contains species not defined in any of the phases associated with this kinetics manager. More...
 
bool skipUndeclaredSpecies () const
 
void skipUndeclaredThirdBodies (bool skip)
 Determine behavior when adding a new reaction that contains third-body efficiencies for species not defined in any of the phases associated with this kinetics manager. More...
 
bool skipUndeclaredThirdBodies () const
 
Altering Reaction Rates
doublereal multiplier (size_t i) const
 The current value of the multiplier for reaction i. More...
 
virtual void setMultiplier (size_t i, doublereal f)
 Set the multiplier for reaction i to f. More...
 
virtual void invalidateCache ()
 

Protected Member Functions

virtual void updateROP ()
 
double checkDuplicateStoich (std::map< int, double > &r1, std::map< int, double > &r2) const
 Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reactions are duplicates of one another, and 0.0 otherwise. More...
 
void checkReactionBalance (const Reaction &R)
 Check that the specified reaction is balanced (same number of atoms for each element in the reactants and products). More...
 

Protected Attributes

ValueCache m_cache
 Cache for saved calculations within each Kinetics object. More...
 
bool m_ready
 Boolean indicating whether Kinetics object is fully configured. More...
 
size_t m_kk
 The number of species in all of the phases that participate in this kinetics mechanism. More...
 
vector_fp m_perturb
 Vector of perturbation factors for each reaction's rate of progress vector. More...
 
std::vector< shared_ptr< Reaction > > m_reactions
 Vector of Reaction objects represented by this Kinetics manager. More...
 
std::vector< ThermoPhase * > m_thermo
 m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator More...
 
std::vector< size_t > m_start
 m_start is a vector of integers specifying the beginning position for the species vector for the n'th phase in the kinetics class. More...
 
std::map< std::string, size_t > m_phaseindex
 Mapping of the phase name to the position of the phase within the kinetics object. More...
 
size_t m_surfphase
 Index in the list of phases of the one surface phase. More...
 
size_t m_rxnphase
 Phase Index where reactions are assumed to be taking place. More...
 
size_t m_mindim
 number of spatial dimensions of lowest-dimensional phase. More...
 
vector_fp m_rfn
 Forward rate constant for each reaction. More...
 
vector_fp m_delta_gibbs0
 Delta G^0 for all reactions. More...
 
vector_fp m_rkcn
 Reciprocal of the equilibrium constant in concentration units. More...
 
vector_fp m_ropf
 Forward rate-of-progress for each reaction. More...
 
vector_fp m_ropr
 Reverse rate-of-progress for each reaction. More...
 
vector_fp m_ropnet
 Net rate-of-progress for each reaction. More...
 
vector_fp m_dH
 The enthalpy change for each reaction to calculate Blowers-Masel rates. More...
 
vector_fp m_rbuf
 Buffer used for storage of intermediate reaction-specific results. More...
 
bool m_skipUndeclaredSpecies
 
bool m_skipUndeclaredThirdBodies
 
std::weak_ptr< Solutionm_root
 reference to Solution More...
 
Stoichiometry management
StoichManagerN m_reactantStoich
 Stoichiometry manager for the reactants for each reaction. More...
 
StoichManagerN m_productStoich
 Stoichiometry manager for the products for each reaction. More...
 
StoichManagerN m_revProductStoich
 Stoichiometry manager for the products of reversible reactions. More...
 
Eigen::SparseMatrix< double > m_stoichMatrix
 Net stoichiometry (products - reactants) More...
 

Detailed Description

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 113 of file Kinetics.h.

Constructor & Destructor Documentation

◆ Kinetics() [1/2]

Kinetics ( )

Default constructor.

Definition at line 26 of file Kinetics.cpp.

◆ ~Kinetics()

~Kinetics ( )
virtual

Definition at line 38 of file Kinetics.cpp.

◆ Kinetics() [2/2]

Kinetics ( const Kinetics )
delete

Kinetics objects are not copyable or assignable.

Member Function Documentation

◆ kineticsType()

virtual std::string kineticsType ( ) const
inlinevirtual

◆ resizeReactions()

void resizeReactions ( )
virtual

◆ nReactions()

size_t nReactions ( ) const
inline

Number of reactions in the reaction mechanism.

Definition at line 139 of file Kinetics.h.

References Kinetics::m_reactions.

Referenced by GasKinetics::addChebyshevReaction(), GasKinetics::addFalloffReaction(), GasKinetics::addPlogReaction(), GasKinetics::addReaction(), InterfaceKinetics::addReaction(), Kinetics::addReaction(), GasKinetics::addThreeBodyReaction(), Kinetics::checkReactionArraySize(), Kinetics::checkReactionIndex(), Kinetics::getCreationRates_ddC(), Kinetics::getCreationRates_ddP(), Kinetics::getCreationRates_ddT(), InterfaceKinetics::getDeltaGibbs(), Kinetics::getDestructionRates_ddC(), Kinetics::getDestructionRates_ddP(), Kinetics::getDestructionRates_ddT(), GasKinetics::getEquilibriumConstants(), InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::getFwdRateConstants(), Kinetics::getNetProductionRates_ddC(), Kinetics::getNetProductionRates_ddP(), Kinetics::getNetProductionRates_ddT(), GasKinetics::getNetRatesOfProgress_ddC(), GasKinetics::getNetRatesOfProgress_ddT(), Kinetics::getReactionDelta(), InterfaceKinetics::getRevRateConstants(), GasKinetics::getRevRatesOfProgress_ddT(), Kinetics::getRevReactionDelta(), Kinetics::parameters(), GasKinetics::process_ddC(), GasKinetics::process_ddX(), GasKinetics::processEquilibriumConstants(), GasKinetics::processFwdRateCoefficients(), GasKinetics::resizeReactions(), InterfaceKinetics::resizeReactions(), Kinetics::resizeReactions(), InterfaceKinetics::resizeSpecies(), InterfaceKinetics::updateExchangeCurrentQuantities(), InterfaceKinetics::updateKc(), and InterfaceKinetics::updateROP().

◆ checkReactionIndex()

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 40 of file Kinetics.cpp.

References Kinetics::nReactions().

Referenced by Kinetics::modifyReaction(), and Kinetics::reaction().

◆ checkReactionArraySize()

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 67 of file Kinetics.cpp.

References Kinetics::nReactions().

◆ checkSpeciesIndex()

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 89 of file Kinetics.cpp.

References Kinetics::m_kk.

◆ checkSpeciesArraySize()

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 96 of file Kinetics.cpp.

References Kinetics::m_kk.

◆ nPhases()

size_t nPhases ( ) const
inline

◆ checkPhaseIndex()

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 75 of file Kinetics.cpp.

References Kinetics::nPhases().

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

◆ checkPhaseArraySize()

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 82 of file Kinetics.cpp.

References Kinetics::nPhases().

◆ phaseIndex()

size_t phaseIndex ( const std::string &  ph) const
inline

Return the phase index of a phase in the list of phases defined within the object.

Parameters
phstd::string name of the phase

If a -1 is returned, then the phase is not defined in the Kinetics object.

Definition at line 193 of file Kinetics.h.

References Kinetics::m_phaseindex, and Cantera::npos.

Referenced by Cantera::importKinetics().

◆ surfacePhaseIndex()

size_t surfacePhaseIndex ( ) const
inline

This returns the integer index of the phase which has ThermoPhase type cSurf.

For heterogeneous mechanisms, this identifies the one surface phase. For homogeneous mechanisms, this returns -1.

Definition at line 206 of file Kinetics.h.

References Kinetics::m_surfphase.

Referenced by solveSP::calc_t(), Reaction::checkBalance(), solveSP::fun_eval(), StickingCoverage::setContext(), ReactorSurface::setKinetics(), and solveSP::solveSP().

◆ reactionPhaseIndex()

size_t reactionPhaseIndex ( ) const
inline

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 218 of file Kinetics.h.

References Kinetics::m_rxnphase.

Referenced by InterfaceKinetics::applyVoltageKfwdCorrection(), InterfaceKinetics::buildSurfaceArrhenius(), Reaction::calculateRateCoeffUnits(), ThreeBodyReaction2::calculateRateCoeffUnits(), FalloffReaction2::calculateRateCoeffUnits(), ChemicallyActivatedReaction2::calculateRateCoeffUnits(), Reaction::calculateRateCoeffUnits3(), Reaction::checkBalance(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), InterfaceKinetics::getDeltaSSEnthalpy(), InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::init(), Reaction::Reaction(), StickingCoverage::setContext(), ImplicitSurfChem::solvePseudoSteadyStateProblem(), InterfaceKinetics::updateKc(), and InterfaceKinetics::updateMu0().

◆ thermo() [1/2]

ThermoPhase & thermo ( size_t  n = 0)
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.

Parameters
nIndex of the ThermoPhase being sought.

Definition at line 231 of file Kinetics.h.

References Kinetics::m_thermo.

Referenced by InterfaceKinetics::_update_rates_phi(), InterfaceKinetics::addPhase(), Kinetics::addPhase(), InterfaceKinetics::applyVoltageKfwdCorrection(), GasKinetics::assertDerivativesValid(), InterfaceKinetics::buildSurfaceArrhenius(), solveSP::calc_t(), Reaction::calculateRateCoeffUnits(), ThreeBodyReaction2::calculateRateCoeffUnits(), FalloffReaction2::calculateRateCoeffUnits(), ChemicallyActivatedReaction2::calculateRateCoeffUnits(), Reaction::calculateRateCoeffUnits3(), Reaction::checkBalance(), InterfaceKinetics::convertExchangeCurrentDensityFormulation(), ReactingSurf1D::eval(), InterfaceKinetics::getDeltaElectrochemPotentials(), InterfaceKinetics::getDeltaEnthalpy(), InterfaceKinetics::getDeltaEntropy(), InterfaceKinetics::getDeltaSSEnthalpy(), InterfaceKinetics::getDeltaSSEntropy(), InterfaceKinetics::getDeltaSSGibbs(), GasKinetics::getEquilibriumConstants(), InterfaceKinetics::getEquilibriumConstants(), InterfaceKinetics::init(), InterfaceKinetics::InterfaceKinetics(), Kinetics::kineticsSpeciesIndex(), Kinetics::kineticsSpeciesName(), GasKinetics::process_ddC(), GasKinetics::process_ddX(), Reaction::Reaction(), Kinetics::reactionEnthalpy(), GasKinetics::resizeReactions(), InterfaceRateBase::setContext(), StickingCoverage::setContext(), ReactorSurface::setKinetics(), ImplicitSurfChem::solvePseudoSteadyStateProblem(), solveSP::solveSP(), Kinetics::speciesPhase(), InterfaceData::update(), GasKinetics::update_rates_C(), InterfaceKinetics::updateExchangeCurrentQuantities(), GasKinetics::updateKc(), InterfaceKinetics::updateKc(), solveSP::updateMFKinSpecies(), InterfaceKinetics::updateMu0(), and Reaction::usesElectrochemistry().

◆ thermo() [2/2]

const ThermoPhase & thermo ( size_t  n = 0) const
inline

Definition at line 234 of file Kinetics.h.

◆ nTotalSpecies()

size_t nTotalSpecies ( ) const
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 243 of file Kinetics.h.

References Kinetics::m_kk.

Referenced by ReactingSurf1D::init(), Kinetics::reactionEnthalpy(), GasKinetics::resizeReactions(), and InterfaceKinetics::resizeReactions().

◆ kineticsSpeciesIndex() [1/3]

size_t kineticsSpeciesIndex ( size_t  k,
size_t  n 
) const
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.

Parameters
kspecies index
nphase index for the species

Definition at line 265 of file Kinetics.h.

References Kinetics::m_start.

Referenced by GasKinetics::addFalloffReaction(), InterfaceKinetics::addReaction(), Kinetics::addReaction(), GasKinetics::addThreeBodyReaction(), InterfaceKinetics::buildSurfaceArrhenius(), solveSP::calc_t(), Cantera::checkElectrochemReaction(), ReactingSurf1D::eval(), solveSP::fun_eval(), Kinetics::kineticsSpeciesIndex(), Kinetics::reactionEnthalpy(), InterfaceRateBase::setContext(), StickingCoverage::setContext(), BlowersMaselRate::setContext(), Reaction::setParameters(), solveSP::solveSP(), solveSP::updateMFKinSpecies(), and Reaction::usesElectrochemistry().

◆ kineticsSpeciesName()

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.

Parameters
kspecies index

Definition at line 310 of file Kinetics.cpp.

References Kinetics::m_start, Cantera::npos, Phase::speciesName(), and Kinetics::thermo().

◆ kineticsSpeciesIndex() [2/3]

size_t kineticsSpeciesIndex ( const std::string &  nm) const

This routine will look up a species number based on the input std::string nm.

The lookup of species will occur for all phases listed in the kinetics object.

return

  • If a match is found, the position in the species list is returned.
  • If no match is found, the value -1 is returned.
Parameters
nmInput string name of the species

Definition at line 320 of file Kinetics.cpp.

References Kinetics::m_start, Kinetics::m_thermo, Cantera::npos, Phase::speciesIndex(), and Kinetics::thermo().

◆ kineticsSpeciesIndex() [3/3]

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.

The lookup of species will occur in the specified phase of the object, or all phases if ph is "<any>".

return

  • If a match is found, the position in the species list is returned.
  • If no match is found, the value npos (-1) is returned.
Parameters
nmInput string name of the species
phInput string name of the phase.

Definition at line 332 of file Kinetics.cpp.

References Kinetics::kineticsSpeciesIndex(), Kinetics::m_start, Kinetics::m_thermo, Phase::name(), Cantera::npos, Phase::speciesIndex(), and Kinetics::thermo().

◆ speciesPhase() [1/3]

ThermoPhase & speciesPhase ( 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.

Will throw an error if the species doesn't match.

Parameters
nmString containing the name of the species.

Definition at line 352 of file Kinetics.cpp.

References Kinetics::m_thermo, Cantera::npos, Phase::speciesIndex(), and Kinetics::thermo().

Referenced by Reaction::calculateRateCoeffUnits(), Reaction::calculateRateCoeffUnits3(), Reaction::checkBalance(), and Cantera::checkElectrochemReaction().

◆ speciesPhase() [2/3]

const ThermoPhase & speciesPhase ( const std::string &  nm) const

Definition at line 363 of file Kinetics.cpp.

◆ speciesPhase() [3/3]

ThermoPhase & speciesPhase ( size_t  k)
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.

Parameters
kSpecies index

Definition at line 325 of file Kinetics.h.

References Kinetics::speciesPhaseIndex(), and Kinetics::thermo().

◆ speciesPhaseIndex()

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.

Parameters
kSpecies index

Definition at line 373 of file Kinetics.cpp.

References Kinetics::m_start, and Cantera::npos.

Referenced by InterfaceKinetics::addReaction(), InterfaceKinetics::buildSurfaceArrhenius(), InterfaceRateBase::setContext(), StickingCoverage::setContext(), Kinetics::speciesPhase(), and Reaction::usesElectrochemistry().

◆ getFwdRatesOfProgress()

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: nReactions().

Definition at line 423 of file Kinetics.cpp.

References Kinetics::m_ropf.

◆ getRevRatesOfProgress()

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: nReactions().

Definition at line 429 of file Kinetics.cpp.

References Kinetics::m_ropr.

◆ getNetRatesOfProgress()

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: nReactions().

Definition at line 435 of file Kinetics.cpp.

References Kinetics::m_ropnet.

◆ getEquilibriumConstants()

virtual void getEquilibriumConstants ( doublereal *  kc)
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) ] \]

Parameters
kcOutput vector containing the equilibrium constants. Length: nReactions().

Reimplemented in GasKinetics, and InterfaceKinetics.

Definition at line 385 of file Kinetics.h.

◆ getReactionDelta()

void getReactionDelta ( const double *  property,
double *  deltaProperty 
) const
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.

Parameters
propertyInput vector of property value. Length: m_kk.
deltaPropertyOutput vector of deltaRxn. Length: nReactions().

Definition at line 441 of file Kinetics.cpp.

References Kinetics::m_productStoich, Kinetics::m_reactantStoich, and Kinetics::nReactions().

Referenced by InterfaceKinetics::applyVoltageKfwdCorrection(), InterfaceKinetics::getDeltaElectrochemPotentials(), InterfaceKinetics::getDeltaEnthalpy(), InterfaceKinetics::getDeltaEntropy(), InterfaceKinetics::getDeltaGibbs(), InterfaceKinetics::getDeltaSSEnthalpy(), InterfaceKinetics::getDeltaSSEntropy(), InterfaceKinetics::getDeltaSSGibbs(), GasKinetics::getEquilibriumConstants(), InterfaceKinetics::getEquilibriumConstants(), and InterfaceKinetics::updateExchangeCurrentQuantities().

◆ getRevReactionDelta()

void getRevReactionDelta ( const double *  g,
double *  dg 
) const
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 450 of file Kinetics.cpp.

References Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, and Kinetics::nReactions().

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

◆ getDeltaGibbs()

virtual void getDeltaGibbs ( doublereal *  deltaG)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 427 of file Kinetics.h.

◆ getDeltaElectrochemPotentials()

virtual void getDeltaElectrochemPotentials ( doublereal *  deltaM)
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

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

Reimplemented in InterfaceKinetics.

Definition at line 442 of file Kinetics.h.

◆ getDeltaEnthalpy()

virtual void getDeltaEnthalpy ( doublereal *  deltaH)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 455 of file Kinetics.h.

◆ getDeltaEntropy()

virtual void getDeltaEntropy ( doublereal *  deltaS)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 468 of file Kinetics.h.

◆ getDeltaSSGibbs()

virtual void getDeltaSSGibbs ( doublereal *  deltaG)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 482 of file Kinetics.h.

◆ getDeltaSSEnthalpy()

virtual void getDeltaSSEnthalpy ( doublereal *  deltaH)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 496 of file Kinetics.h.

◆ getDeltaSSEntropy()

virtual void getDeltaSSEntropy ( doublereal *  deltaS)
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

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

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 510 of file Kinetics.h.

◆ getThirdBodyConcentrations()

virtual void getThirdBodyConcentrations ( double *  concm)
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 parters are not defined and set to not-a-number.

Parameters
concmOutput vector of effective third-body concentrations. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 522 of file Kinetics.h.

References Kinetics::kineticsType().

◆ thirdBodyConcentrations()

virtual const vector_fp & thirdBodyConcentrations ( ) const
inlinevirtual

Provide direct access to current third-body concentration values.

See also
getThirdBodyConcentrations.

Reimplemented in GasKinetics.

Definition at line 532 of file Kinetics.h.

References Kinetics::kineticsType().

◆ getCreationRates()

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.

Definition at line 459 of file Kinetics.cpp.

References Kinetics::m_kk, Kinetics::m_productStoich, Kinetics::m_reactantStoich, Kinetics::m_ropf, and Kinetics::m_ropr.

◆ getDestructionRates()

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.

Definition at line 473 of file Kinetics.cpp.

References Kinetics::m_kk, Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::m_ropf, and Kinetics::m_ropr.

◆ getNetProductionRates()

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.

Definition at line 484 of file Kinetics.cpp.

References Kinetics::m_kk, Kinetics::m_productStoich, Kinetics::m_reactantStoich, and Kinetics::m_ropnet.

Referenced by solveSP::calc_t(), ConstPressureReactor::eval(), IdealGasConstPressureReactor::eval(), IdealGasReactor::eval(), ReactingSurf1D::eval(), solveSP::fun_eval(), and StFlow::getWdot().

◆ getDerivativeSettings()

virtual void getDerivativeSettings ( AnyMap settings) const
inlinevirtual

Retrieve derivative settings.

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:

  • _ddT: derivative with respect to temperature (a vector)
  • _ddP: derivative with respect to pressure (a vector)
  • _ddC: derivative with respect to molar concentration (a vector)
  • _ddX: derivative with respect to species mole fractions (a matrix)

Settings for derivative evaluation are set by keyword/value pairs using the methods

See also
getDerivativeSettings and
setDerivativeSettings.

For GasKinetics, the following keyword/value pairs are supported:

  • skip-third-bodies (boolean) ... if false (default), third body concentrations are considered for the evaluation of jacobians
  • skip-falloff (boolean) ... if false (default), third-body effects on rate constants are considered for the evaluation of derivatives.
  • rtol-delta (double) ... relative tolerance used to perturb properties when calculating numerical derivatives. The default value is 1e-8.
Warning
The calculation of derivatives is an experimental part of the Cantera API and may be changed or removed without notice.
Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented in GasKinetics.

Definition at line 604 of file Kinetics.h.

References Kinetics::kineticsType().

◆ setDerivativeSettings()

virtual void setDerivativeSettings ( const AnyMap settings)
inlinevirtual

Set/modify derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented in GasKinetics.

Definition at line 615 of file Kinetics.h.

References Kinetics::kineticsType().

◆ getFwdRateConstants_ddT()

virtual void getFwdRateConstants_ddT ( double *  dkfwd)
inlinevirtual

Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dkfwdOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 626 of file Kinetics.h.

References Kinetics::kineticsType().

◆ getFwdRateConstants_ddP()

virtual void getFwdRateConstants_ddP ( double *  dkfwd)
inlinevirtual

Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dkfwdOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 637 of file Kinetics.h.

References Kinetics::kineticsType().

◆ getFwdRateConstants_ddC()

virtual void getFwdRateConstants_ddC ( double *  dkfwd)
inlinevirtual

Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dkfwdOutput vector of derivatives. Length: nReactions().
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Reimplemented in GasKinetics.

Definition at line 651 of file Kinetics.h.

References Kinetics::kineticsType().

◆ getFwdRatesOfProgress_ddT()

virtual void getFwdRatesOfProgress_ddT ( double *  drop)
inlinevirtual

Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 662 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddT(), and Kinetics::getDestructionRates_ddT().

◆ getFwdRatesOfProgress_ddP()

virtual void getFwdRatesOfProgress_ddP ( double *  drop)
inlinevirtual

Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 673 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddP(), and Kinetics::getDestructionRates_ddP().

◆ getFwdRatesOfProgress_ddC()

virtual void getFwdRatesOfProgress_ddC ( double *  drop)
inlinevirtual

Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Reimplemented in GasKinetics.

Definition at line 687 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddC(), and Kinetics::getDestructionRates_ddC().

◆ fwdRatesOfProgress_ddX()

virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX ( )
inlinevirtual

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

The method returns a matrix with nReactions rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Reimplemented in GasKinetics.

Definition at line 704 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::creationRates_ddX(), and Kinetics::destructionRates_ddX().

◆ getRevRatesOfProgress_ddT()

virtual void getRevRatesOfProgress_ddT ( double *  drop)
inlinevirtual

Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 715 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddT(), and Kinetics::getDestructionRates_ddT().

◆ getRevRatesOfProgress_ddP()

virtual void getRevRatesOfProgress_ddP ( double *  drop)
inlinevirtual

Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 726 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddP(), and Kinetics::getDestructionRates_ddP().

◆ getRevRatesOfProgress_ddC()

virtual void getRevRatesOfProgress_ddC ( double *  drop)
inlinevirtual

Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Reimplemented in GasKinetics.

Definition at line 740 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getCreationRates_ddC(), and Kinetics::getDestructionRates_ddC().

◆ revRatesOfProgress_ddX()

virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX ( )
inlinevirtual

Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration.

The method returns a matrix with nReactions rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Reimplemented in GasKinetics.

Definition at line 757 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::creationRates_ddX(), and Kinetics::destructionRates_ddX().

◆ getNetRatesOfProgress_ddT()

virtual void getNetRatesOfProgress_ddT ( double *  drop)
inlinevirtual

Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 768 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getNetProductionRates_ddT().

◆ getNetRatesOfProgress_ddP()

virtual void getNetRatesOfProgress_ddP ( double *  drop)
inlinevirtual

Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().

Reimplemented in GasKinetics.

Definition at line 779 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getNetProductionRates_ddP().

◆ getNetRatesOfProgress_ddC()

virtual void getNetRatesOfProgress_ddC ( double *  drop)
inlinevirtual

Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dropOutput vector of derivatives. Length: nReactions().
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Reimplemented in GasKinetics.

Definition at line 793 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::getNetProductionRates_ddC().

◆ netRatesOfProgress_ddX()

virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX ( )
inlinevirtual

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

The method returns a matrix with nReactions rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Reimplemented in GasKinetics.

Definition at line 810 of file Kinetics.h.

References Kinetics::kineticsType().

Referenced by Kinetics::netProductionRates_ddX().

◆ getCreationRates_ddT()

void getCreationRates_ddT ( double *  dwdot)

Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 495 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddT(), Kinetics::getRevRatesOfProgress_ddT(), Kinetics::m_kk, Kinetics::m_productStoich, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ getCreationRates_ddP()

void getCreationRates_ddP ( double *  dwdot)

Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 507 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddP(), Kinetics::getRevRatesOfProgress_ddP(), Kinetics::m_kk, Kinetics::m_productStoich, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ getCreationRates_ddC()

void getCreationRates_ddC ( double *  dwdot)

Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Definition at line 519 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddC(), Kinetics::getRevRatesOfProgress_ddC(), Kinetics::m_kk, Kinetics::m_productStoich, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ creationRates_ddX()

Eigen::SparseMatrix< double > creationRates_ddX ( )

Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration.

The method returns a matrix with nTotalSpecies rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Definition at line 531 of file Kinetics.cpp.

References Kinetics::fwdRatesOfProgress_ddX(), Kinetics::m_productStoich, Kinetics::m_reactantStoich, Kinetics::revRatesOfProgress_ddX(), and StoichManagerN::stoichCoeffs().

◆ getDestructionRates_ddT()

void getDestructionRates_ddT ( double *  dwdot)

Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 541 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddT(), Kinetics::getRevRatesOfProgress_ddT(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ getDestructionRates_ddP()

void getDestructionRates_ddP ( double *  dwdot)

Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 553 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddP(), Kinetics::getRevRatesOfProgress_ddP(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ getDestructionRates_ddC()

void getDestructionRates_ddC ( double *  dwdot)

Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Definition at line 565 of file Kinetics.cpp.

References Kinetics::getFwdRatesOfProgress_ddC(), Kinetics::getRevRatesOfProgress_ddC(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::nReactions(), and StoichManagerN::stoichCoeffs().

◆ destructionRates_ddX()

Eigen::SparseMatrix< double > destructionRates_ddX ( )

Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration.

The method returns a matrix with nTotalSpecies rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Definition at line 577 of file Kinetics.cpp.

References Kinetics::fwdRatesOfProgress_ddX(), Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::revRatesOfProgress_ddX(), and StoichManagerN::stoichCoeffs().

◆ getNetProductionRates_ddT()

void getNetProductionRates_ddT ( double *  dwdot)

Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 587 of file Kinetics.cpp.

References Kinetics::getNetRatesOfProgress_ddT(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_stoichMatrix, and Kinetics::nReactions().

◆ getNetProductionRates_ddP()

void getNetProductionRates_ddP ( double *  dwdot)

Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.

Definition at line 595 of file Kinetics.cpp.

References Kinetics::getNetRatesOfProgress_ddP(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_stoichMatrix, and Kinetics::nReactions().

◆ getNetProductionRates_ddC()

void getNetProductionRates_ddC ( double *  dwdot)

Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions.

Parameters
[out]dwdotOutput vector of derivatives. Length: m_kk.
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Definition at line 603 of file Kinetics.cpp.

References Kinetics::getNetRatesOfProgress_ddC(), Kinetics::m_kk, Kinetics::m_rbuf, Kinetics::m_stoichMatrix, and Kinetics::nReactions().

◆ netProductionRates_ddX()

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.

The method returns a matrix with nTotalSpecies rows and nTotalSpecies columns. For a derivative with respect to \(X_i\), all other \(X_j\) are held constant, rather than enforcing \(\sum X_j = 1\).

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

Definition at line 611 of file Kinetics.cpp.

References Kinetics::m_stoichMatrix, and Kinetics::netRatesOfProgress_ddX().

◆ reactantStoichCoeff()

double reactantStoichCoeff ( size_t  k,
size_t  i 
) const
virtual

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

Parameters
kkinetic species index
ireaction index

Definition at line 384 of file Kinetics.cpp.

References Kinetics::m_reactantStoich, and StoichManagerN::stoichCoeffs().

◆ reactantStoichCoeffs()

Eigen::SparseMatrix< double > reactantStoichCoeffs ( ) const
inline

Stoichiometric coefficient matrix for reactants.

Definition at line 942 of file Kinetics.h.

References Kinetics::m_reactantStoich, and StoichManagerN::stoichCoeffs().

◆ productStoichCoeff()

double productStoichCoeff ( size_t  k,
size_t  i 
) const
virtual

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

Parameters
kkinetic species index
ireaction index

Definition at line 389 of file Kinetics.cpp.

References Kinetics::m_productStoich, and StoichManagerN::stoichCoeffs().

◆ productStoichCoeffs()

Eigen::SparseMatrix< double > productStoichCoeffs ( ) const
inline

Stoichiometric coefficient matrix for products.

Definition at line 957 of file Kinetics.h.

References Kinetics::m_productStoich, and StoichManagerN::stoichCoeffs().

◆ revProductStoichCoeffs()

Eigen::SparseMatrix< double > revProductStoichCoeffs ( ) const
inline

Stoichiometric coefficient matrix for products of reversible reactions.

Definition at line 964 of file Kinetics.h.

References Kinetics::m_revProductStoich, and StoichManagerN::stoichCoeffs().

◆ reactantOrder()

virtual doublereal reactantOrder ( size_t  k,
size_t  i 
) const
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

Parameters
kkinetic species index
ireaction index

Definition at line 976 of file Kinetics.h.

◆ productOrder()

virtual doublereal productOrder ( int  k,
int  i 
) const
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.

Parameters
kkinetic species index
ireaction index

Definition at line 990 of file Kinetics.h.

◆ getActivityConcentrations()

virtual void getActivityConcentrations ( doublereal *const  conc)
inlinevirtual

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 in InterfaceKinetics.

Definition at line 999 of file Kinetics.h.

◆ reactionType()

int reactionType ( size_t  i) const
virtual

Flag specifying the type of reaction.

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

Parameters
ireaction index
Deprecated:
To be changed after Cantera 2.6.

Definition at line 394 of file Kinetics.cpp.

References Kinetics::m_reactions, and Cantera::warn_deprecated().

◆ reactionTypeStr()

std::string reactionTypeStr ( size_t  i) const
virtual

String specifying the type of reaction.

Parameters
ireaction index

Definition at line 402 of file Kinetics.cpp.

References Kinetics::m_reactions.

Referenced by GasKinetics::processFalloffReactions().

◆ isReversible()

virtual bool isReversible ( size_t  i)
inlinevirtual

True if reaction i has been declared to be reversible.

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

Parameters
ireaction index

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 1027 of file Kinetics.h.

◆ reactionString()

std::string reactionString ( size_t  i) const

Return a string representing the reaction.

Parameters
ireaction index

Definition at line 406 of file Kinetics.cpp.

References Kinetics::m_reactions.

◆ reactantString()

std::string reactantString ( size_t  i) const

Returns a string containing the reactants side of the reaction equation.

Definition at line 412 of file Kinetics.cpp.

References Kinetics::m_reactions.

◆ productString()

std::string productString ( size_t  i) const

Returns a string containing the products side of the reaction equation.

Definition at line 418 of file Kinetics.cpp.

References Kinetics::m_reactions.

◆ getFwdRateConstants()

virtual void getFwdRateConstants ( double *  kfwd)
inlinevirtual

Return the forward rate constants.

The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Parameters
kfwdOutput vector containing the forward reaction rate constants. Length: nReactions().
Deprecated:
Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin and Glarborg, 'Chemically Reacting Flow', Wiley Interscience, 2003). For new behavior, set 'Cantera::use_legacy_rate_constants(false)'.

Reimplemented in GasKinetics, and InterfaceKinetics.

Definition at line 1064 of file Kinetics.h.

◆ getRevRateConstants()

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

Return the reverse rate constants.

The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Length is the number of reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction. Note, this routine will return rate constants for irreversible reactions if the default for doIrreversible is overridden.

Parameters
krevOutput vector of reverse rate constants
doIrreversibleboolean indicating whether irreversible reactions should be included.
Deprecated:
Behavior to change after Cantera 2.6; for Cantera 2.6, rate constants of three-body reactions are multiplied with third-body concentrations (no change to legacy behavior). After Cantera 2.6, results will no longer include third-body concentrations and be consistent with conventional definitions (see Eq. 9.75 in Kee, Coltrin and Glarborg, 'Chemically Reacting Flow', Wiley Interscience, 2003). For new behavior, set 'Cantera::use_legacy_rate_constants(false)'.

Reimplemented in BulkKinetics, and InterfaceKinetics.

Definition at line 1090 of file Kinetics.h.

◆ addPhase()

void addPhase ( ThermoPhase thermo)
virtual

Add a phase to the kinetics manager object.

This must be done before the function init() is called or before any reactions are input. The following fields are updated:

  • m_start -> vector of integers, containing the starting position of the species for each phase in the kinetics mechanism.
  • m_surfphase -> index of the surface phase.
  • m_thermo -> vector of pointers to ThermoPhase phases that participate in the kinetics mechanism.
  • m_phaseindex -> map containing the std::string id of each ThermoPhase phase as a key and the index of the phase within the kinetics manager object as the value.
Parameters
thermoReference to the ThermoPhase to be added.

Reimplemented in InterfaceKinetics.

Definition at line 616 of file Kinetics.cpp.

References Kinetics::kineticsType(), Kinetics::m_mindim, Kinetics::m_phaseindex, Kinetics::m_rxnphase, Kinetics::m_surfphase, Kinetics::m_thermo, Phase::nDim(), Kinetics::nPhases(), Kinetics::resizeSpecies(), Kinetics::thermo(), and ThermoPhase::type().

Referenced by InterfaceKinetics::addPhase(), and Cantera::importKinetics().

◆ init()

virtual void init ( )
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 1125 of file Kinetics.h.

Referenced by Kinetics::addReaction(), and Cantera::importKinetics().

◆ parameters()

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 635 of file Kinetics.cpp.

References Kinetics::kineticsType(), and Kinetics::nReactions().

◆ resizeSpecies()

void resizeSpecies ( )
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 648 of file Kinetics.cpp.

References Kinetics::m_kk, Kinetics::m_start, Kinetics::m_thermo, and Kinetics::nPhases().

Referenced by Kinetics::addPhase(), Kinetics::addReaction(), and InterfaceKinetics::resizeSpecies().

◆ addReaction()

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

Add a single reaction to the mechanism.

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

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

Reimplemented in BulkKinetics, GasKinetics, and InterfaceKinetics.

Definition at line 660 of file Kinetics.cpp.

References StoichManagerN::add(), Kinetics::init(), Kinetics::kineticsSpeciesIndex(), Kinetics::m_delta_gibbs0, Kinetics::m_dH, Kinetics::m_kk, Kinetics::m_perturb, Kinetics::m_productStoich, Kinetics::m_reactantStoich, Kinetics::m_reactions, Kinetics::m_ready, Kinetics::m_revProductStoich, Kinetics::m_rfn, Kinetics::m_rkcn, Kinetics::m_ropf, Kinetics::m_ropnet, Kinetics::m_ropr, Kinetics::nReactions(), Kinetics::resizeReactions(), and Kinetics::resizeSpecies().

Referenced by InterfaceKinetics::addReaction(), and Cantera::installReactionArrays().

◆ modifyReaction()

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

Modify the rate expression associated with a reaction.

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

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

Reimplemented in BulkKinetics, GasKinetics, and InterfaceKinetics.

Definition at line 748 of file Kinetics.cpp.

References Kinetics::checkReactionIndex(), and Kinetics::m_reactions.

Referenced by InterfaceKinetics::modifyReaction().

◆ reaction() [1/2]

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 781 of file Kinetics.cpp.

References Kinetics::checkReactionIndex(), and Kinetics::m_reactions.

◆ reaction() [2/2]

shared_ptr< const Reaction > reaction ( size_t  i) const

Definition at line 787 of file Kinetics.cpp.

◆ skipUndeclaredSpecies() [1/2]

void skipUndeclaredSpecies ( bool  skip)
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 1171 of file Kinetics.h.

References Kinetics::m_skipUndeclaredSpecies.

Referenced by Cantera::addReactions(), Reaction::checkSpecies(), and Cantera::installReactionArrays().

◆ skipUndeclaredSpecies() [2/2]

bool skipUndeclaredSpecies ( ) const
inline

Definition at line 1174 of file Kinetics.h.

◆ skipUndeclaredThirdBodies() [1/2]

void skipUndeclaredThirdBodies ( bool  skip)
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 1183 of file Kinetics.h.

References Kinetics::m_skipUndeclaredThirdBodies.

Referenced by Cantera::addReactions(), Reaction::checkSpecies(), and Cantera::installReactionArrays().

◆ skipUndeclaredThirdBodies() [2/2]

bool skipUndeclaredThirdBodies ( ) const
inline

Definition at line 1186 of file Kinetics.h.

◆ multiplier()

doublereal multiplier ( size_t  i) const
inline

The current value of the multiplier for reaction i.

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.

Parameters
iindex of the reaction

Definition at line 1206 of file Kinetics.h.

References Kinetics::m_perturb.

Referenced by InterfaceKinetics::buildSurfaceArrhenius().

◆ setMultiplier()

virtual void setMultiplier ( size_t  i,
doublereal  f 
)
inlinevirtual

Set the multiplier for reaction i to f.

Parameters
iindex of the reaction
fvalue of the multiplier.

Reimplemented in BulkKinetics.

Definition at line 1215 of file Kinetics.h.

References Kinetics::m_perturb.

◆ invalidateCache()

virtual void invalidateCache ( )
inlinevirtual

Definition at line 1219 of file Kinetics.h.

◆ checkDuplicates()

std::pair< size_t, size_t > checkDuplicates ( bool  throw_err = true) const
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 103 of file Kinetics.cpp.

References Kinetics::m_reactions.

Referenced by Cantera::installReactionArrays().

◆ selectPhase()

void selectPhase ( const double *  data,
const ThermoPhase phase,
double *  phase_data 
)

Takes as input an array of properties for all species in the mechanism and copies those values belonging to a particular phase to the output array.

Parameters
dataInput data array.
phasePointer to one of the phase objects participating in this reaction mechanism
phase_dataOutput array where the values for the the specified phase are to be written.

Definition at line 296 of file Kinetics.cpp.

References Kinetics::m_start, Kinetics::m_thermo, Kinetics::nPhases(), and Phase::nSpecies().

◆ setRoot()

virtual void setRoot ( std::shared_ptr< Solution root)
inlinevirtual

Set root Solution holding all phase information.

Definition at line 1248 of file Kinetics.h.

References Kinetics::m_root.

◆ reactionEnthalpy()

double reactionEnthalpy ( const Composition reactants,
const Composition products 
)
virtual

Calculate the reaction enthalpy of a reaction which has not necessarily been added into the Kinetics object.

Definition at line 793 of file Kinetics.cpp.

References ThermoPhase::getPartialMolarEnthalpies(), Kinetics::kineticsSpeciesIndex(), Kinetics::m_start, Kinetics::nPhases(), Kinetics::nTotalSpecies(), and Kinetics::thermo().

◆ updateROP()

virtual void updateROP ( )
inlineprotectedvirtual

Reimplemented in InterfaceKinetics.

Definition at line 1262 of file Kinetics.h.

◆ checkReactionBalance()

void checkReactionBalance ( const Reaction R)
protected

Check that the specified reaction is balanced (same number of atoms for each element in the reactants and products).

Raises an exception if the reaction is not balanced.

Deprecated:
To be removed in Cantera 2.6. Replaceable by Reaction::checkBalance.

Definition at line 289 of file Kinetics.cpp.

Member Data Documentation

◆ m_cache

ValueCache m_cache
protected

Cache for saved calculations within each Kinetics object.

Definition at line 1259 of file Kinetics.h.

◆ m_reactantStoich

StoichManagerN m_reactantStoich
protected

◆ m_productStoich

StoichManagerN m_productStoich
protected

◆ m_revProductStoich

StoichManagerN m_revProductStoich
protected

◆ m_stoichMatrix

Eigen::SparseMatrix<double> m_stoichMatrix
protected

◆ m_ready

bool m_ready
protected

Boolean indicating whether Kinetics object is fully configured.

Definition at line 1313 of file Kinetics.h.

Referenced by Kinetics::addReaction(), and Kinetics::resizeReactions().

◆ m_kk

size_t m_kk
protected

◆ m_perturb

vector_fp m_perturb
protected

Vector of perturbation factors for each reaction's rate of progress vector.

It is initialized to one.

Definition at line 1321 of file Kinetics.h.

Referenced by Kinetics::addReaction(), InterfaceKinetics::getFwdRateConstants(), Kinetics::multiplier(), GasKinetics::processFwdRateCoefficients(), Kinetics::setMultiplier(), and InterfaceKinetics::updateROP().

◆ m_reactions

std::vector<shared_ptr<Reaction> > m_reactions
protected

◆ m_thermo

std::vector<ThermoPhase*> m_thermo
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.

Note that this kinetics object doesn't own these ThermoPhase objects and is not responsible for creating or deleting them.

Definition at line 1339 of file Kinetics.h.

Referenced by InterfaceKinetics::_update_rates_C(), Kinetics::addPhase(), InterfaceKinetics::getDeltaGibbs(), Kinetics::kineticsSpeciesIndex(), Kinetics::nPhases(), Kinetics::resizeSpecies(), Kinetics::selectPhase(), Kinetics::speciesPhase(), and Kinetics::thermo().

◆ m_start

std::vector<size_t> m_start
protected

◆ m_phaseindex

std::map<std::string, size_t> m_phaseindex
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 1353 of file Kinetics.h.

Referenced by Kinetics::addPhase(), and Kinetics::phaseIndex().

◆ m_surfphase

size_t m_surfphase
protected

Index in the list of phases of the one surface phase.

Definition at line 1356 of file Kinetics.h.

Referenced by Kinetics::addPhase(), and Kinetics::surfacePhaseIndex().

◆ m_rxnphase

size_t m_rxnphase
protected

Phase Index where reactions are assumed to be taking place.

We calculate this by assuming that the phase with the lowest dimensionality is the phase where reactions are taking place.

Definition at line 1363 of file Kinetics.h.

Referenced by Kinetics::addPhase(), and Kinetics::reactionPhaseIndex().

◆ m_mindim

size_t m_mindim
protected

number of spatial dimensions of lowest-dimensional phase.

Definition at line 1366 of file Kinetics.h.

Referenced by Kinetics::addPhase().

◆ m_rfn

vector_fp m_rfn
protected

◆ m_delta_gibbs0

vector_fp m_delta_gibbs0
protected

Delta G^0 for all reactions.

Definition at line 1372 of file Kinetics.h.

Referenced by Kinetics::addReaction(), and GasKinetics::updateKc().

◆ m_rkcn

vector_fp m_rkcn
protected

◆ m_ropf

vector_fp m_ropf
protected

◆ m_ropr

vector_fp m_ropr
protected

◆ m_ropnet

vector_fp m_ropnet
protected

◆ m_dH

vector_fp m_dH
protected

The enthalpy change for each reaction to calculate Blowers-Masel rates.

Definition at line 1387 of file Kinetics.h.

Referenced by Kinetics::addReaction().

◆ m_rbuf

vector_fp m_rbuf
protected

◆ m_skipUndeclaredSpecies

bool m_skipUndeclaredSpecies
protected
See also
skipUndeclaredSpecies()

Definition at line 1393 of file Kinetics.h.

Referenced by Kinetics::skipUndeclaredSpecies().

◆ m_skipUndeclaredThirdBodies

bool m_skipUndeclaredThirdBodies
protected
See also
skipUndeclaredThirdBodies()

Definition at line 1396 of file Kinetics.h.

Referenced by Kinetics::skipUndeclaredThirdBodies().

◆ m_root

std::weak_ptr<Solution> m_root
protected

reference to Solution

Definition at line 1399 of file Kinetics.h.

Referenced by Kinetics::setRoot().


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