Cantera  3.1.0a1
BulkKinetics Class Reference

Specialization of Kinetics for chemistry in a single bulk phase. More...

#include <BulkKinetics.h>

Inheritance diagram for BulkKinetics:
[legend]

Detailed Description

Specialization of Kinetics for chemistry in a single bulk phase.

Definition at line 20 of file BulkKinetics.h.

Public Member Functions

Constructors and General Information
 BulkKinetics ()
 
string kineticsType () const override
 Identifies the Kinetics manager type. More...
 
bool isReversible (size_t i) override
 True if reaction i has been declared to be reversible. More...
 
Reaction Mechanism Setup Routines
bool addReaction (shared_ptr< Reaction > r, bool resize=true) override
 Add a single reaction to the mechanism. More...
 
void addThirdBody (shared_ptr< Reaction > r)
 
void modifyReaction (size_t i, shared_ptr< Reaction > rNew) override
 Modify the rate expression associated with a reaction. More...
 
void resizeSpecies () override
 Resize arrays with sizes that depend on the total number of species. More...
 
void resizeReactions () override
 Finalize Kinetics object and associated objects. More...
 
void setMultiplier (size_t i, double f) override
 Set the multiplier for reaction i to f. More...
 
void invalidateCache () override
 
Reaction rate constants, rates of progress, and thermodynamic properties
void getFwdRateConstants (double *kfwd) override
 Return the forward rate constants. More...
 
void getEquilibriumConstants (double *kc) override
 Return a vector of Equilibrium constants. More...
 
void getRevRateConstants (double *krev, bool doIrreversible=false) override
 Return the reverse rate constants. More...
 
void getDeltaGibbs (double *deltaG) override
 Return the vector of values for the reaction Gibbs free energy change. More...
 
void getDeltaEnthalpy (double *deltaH) override
 Return the vector of values for the reactions change in enthalpy. More...
 
void getDeltaEntropy (double *deltaS) override
 Return the vector of values for the reactions change in entropy. More...
 
void getDeltaSSGibbs (double *deltaG) override
 Return the vector of values for the reaction standard state Gibbs free energy change. More...
 
void getDeltaSSEnthalpy (double *deltaH) override
 Return the vector of values for the change in the standard state enthalpies of reaction. More...
 
void getDeltaSSEntropy (double *deltaS) override
 Return the vector of values for the change in the standard state entropies for each reaction. More...
 
Derivatives of rate constants and rates of progress
void getDerivativeSettings (AnyMap &settings) const override
 Retrieve derivative settings. More...
 
void setDerivativeSettings (const AnyMap &settings) override
 Set/modify derivative settings. More...
 
void getFwdRateConstants_ddT (double *dkfwd) override
 Calculate derivatives for forward rate constants with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getFwdRatesOfProgress_ddT (double *drop) override
 Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getRevRatesOfProgress_ddT (double *drop) override
 Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getNetRatesOfProgress_ddT (double *drop) override
 Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getFwdRateConstants_ddP (double *dkfwd) override
 Calculate derivatives for forward rate constants with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getFwdRatesOfProgress_ddP (double *drop) override
 Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getRevRatesOfProgress_ddP (double *drop) override
 Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getNetRatesOfProgress_ddP (double *drop) override
 Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getFwdRateConstants_ddC (double *dkfwd) override
 Calculate derivatives for forward rate constants with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
void getFwdRatesOfProgress_ddC (double *drop) override
 Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
void getRevRatesOfProgress_ddC (double *drop) override
 Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
void getNetRatesOfProgress_ddC (double *drop) override
 Calculate derivatives for net rates-of-progress with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX () override
 Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
Eigen::SparseMatrix< double > revRatesOfProgress_ddX () override
 Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
Eigen::SparseMatrix< double > netRatesOfProgress_ddX () override
 Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi () override
 Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. More...
 
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi () override
 Calculate derivatives for forward rates-of-progress with respect to species concentration at constant temperature, pressure and remaining species concentrations. More...
 
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi () override
 Calculate derivatives for net rates-of-progress with respect to species concentration at constant temperature, pressure, and remaining species concentrations. More...
 
Rate calculation intermediate methods
void updateROP () override
 
void getThirdBodyConcentrations (double *concm) override
 Return a vector of values of effective concentrations of third-body collision partners of any reaction. More...
 
const vector< double > & thirdBodyConcentrations () const override
 Provide direct access to current third-body concentration values. More...
 
- Public Member Functions inherited from Kinetics
virtual pair< size_t, size_t > checkDuplicates (bool throw_err=true) const
 Check for unmarked duplicate reactions and unmatched marked duplicates. More...
 
virtual void setRoot (shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
shared_ptr< Solutionroot () const
 Get the Solution object containing this Kinetics object and associated ThermoPhase objects. More...
 
 Kinetics ()=default
 Default constructor. More...
 
virtual ~Kinetics ()=default
 
 Kinetics (const Kinetics &)=delete
 Kinetics objects are not copyable or assignable. More...
 
Kineticsoperator= (const Kinetics &)=delete
 
size_t nReactions () const
 Number of reactions in the reaction mechanism. More...
 
void checkReactionIndex (size_t m) const
 Check that the specified reaction index is in range Throws an exception if i is greater than nReactions() More...
 
void checkReactionArraySize (size_t ii) const
 Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions(). More...
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range Throws an exception if k is greater than nSpecies()-1. More...
 
void checkSpeciesArraySize (size_t mm) const
 Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). More...
 
size_t nPhases () const
 The number of phases participating in the reaction mechanism. More...
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases() More...
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). More...
 
size_t phaseIndex (const string &ph) const
 Return the phase index of a phase in the list of phases defined within the object. More...
 
size_t reactionPhaseIndex () const
 Phase where the reactions occur. More...
 
shared_ptr< ThermoPhasereactionPhase () const
 Return pointer to 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...
 
string kineticsSpeciesName (size_t k) const
 Return the name of the kth species in the kinetics manager. More...
 
size_t kineticsSpeciesIndex (const string &nm) const
 This routine will look up a species number based on the input string nm. More...
 
ThermoPhasespeciesPhase (const string &nm)
 This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides. More...
 
const ThermoPhasespeciesPhase (const string &nm) const
 
ThermoPhasespeciesPhase (size_t k)
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the species' owning ThermoPhase object. More...
 
size_t speciesPhaseIndex (size_t k) const
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the index of the phase owning the species. More...
 
virtual void getFwdRatesOfProgress (double *fwdROP)
 Return the forward rates of progress of the reactions. More...
 
virtual void getRevRatesOfProgress (double *revROP)
 Return the Reverse rates of progress of the reactions. More...
 
virtual void getNetRatesOfProgress (double *netROP)
 Net rates of progress. More...
 
virtual void getReactionDelta (const double *property, double *deltaProperty) const
 Change in species properties. More...
 
virtual void getRevReactionDelta (const double *g, double *dg) const
 Given an array of species properties 'g', return in array 'dg' the change in this quantity in the reversible reactions. More...
 
virtual void getDeltaElectrochemPotentials (double *deltaM)
 Return the vector of values for the reaction electrochemical free energy change. More...
 
virtual void getCreationRates (double *cdot)
 Species creation rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getDestructionRates (double *ddot)
 Species destruction rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getNetProductionRates (double *wdot)
 Species net production rates [kmol/m^3/s or kmol/m^2/s]. 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...
 
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. 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...
 
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. 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...
 
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. More...
 
virtual double reactantStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a reactant in reaction i. More...
 
Eigen::SparseMatrix< double > reactantStoichCoeffs () const
 Stoichiometric coefficient matrix for reactants. More...
 
virtual double productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i. More...
 
Eigen::SparseMatrix< double > productStoichCoeffs () const
 Stoichiometric coefficient matrix for products. More...
 
Eigen::SparseMatrix< double > revProductStoichCoeffs () const
 Stoichiometric coefficient matrix for products of reversible reactions. More...
 
virtual double reactantOrder (size_t k, size_t i) const
 Reactant order of species k in reaction i. More...
 
virtual double productOrder (int k, int i) const
 product Order of species k in reaction i. More...
 
virtual void getActivityConcentrations (double *const conc)
 Get the vector of activity concentrations used in the kinetics object. More...
 
virtual void addThermo (shared_ptr< 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...
 
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
 
double multiplier (size_t i) const
 The current value of the multiplier for reaction i. More...
 

Protected Member Functions

Internal service methods
Note
These methods are for internal use, and seek to avoid code duplication while evaluating terms used for rate constants, rates of progress, and their derivatives.
void processThirdBodies (double *rop)
 Multiply rate with third-body collider concentrations. More...
 
void applyEquilibriumConstants (double *rop)
 Multiply rate with inverse equilibrium constant. More...
 
void applyEquilibriumConstants_ddT (double *drkcn)
 Multiply rate with scaled temperature derivatives of the inverse equilibrium constant. More...
 
void process_ddT (const vector< double > &in, double *drop)
 Process temperature derivative. More...
 
void process_ddP (const vector< double > &in, double *drop)
 Process pressure derivative. More...
 
void process_ddC (StoichManagerN &stoich, const vector< double > &in, double *drop, bool mass_action=true)
 Process concentration (molar density) derivative. More...
 
Eigen::SparseMatrix< double > calculateCompositionDerivatives (StoichManagerN &stoich, const vector< double > &in, bool ddX=true)
 Process derivatives. More...
 
void assertDerivativesValid (const string &name)
 Helper function ensuring that all rate derivatives can be calculated. More...
 
- Protected Member Functions inherited from Kinetics
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. More...
 

Protected Attributes

vector< unique_ptr< MultiRateBase > > m_bulk_rates
 Vector of rate handlers. More...
 
map< string, size_t > m_bulk_types
 Mapping of rate handlers. More...
 
vector< size_t > m_revindex
 Indices of reversible reactions. More...
 
vector< size_t > m_irrev
 Indices of irreversible reactions. More...
 
vector< double > m_dn
 Difference between the global reactants order and the global products order. More...
 
ThirdBodyCalc m_multi_concm
 used with MultiRate evaluator More...
 
vector< double > m_concm
 Third body concentrations. More...
 
vector< double > m_act_conc
 Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations. More...
 
vector< double > m_phys_conc
 Physical concentrations, as calculated by ThermoPhase::getConcentrations. More...
 
bool m_jac_skip_third_bodies
 Derivative settings. More...
 
bool m_jac_skip_falloff
 
double m_jac_rtol_delta
 
bool m_ROP_ok = false
 
vector< double > m_rbuf0
 Buffers for partial rop results with length nReactions() More...
 
vector< double > m_rbuf1
 
vector< double > m_rbuf2
 
vector< double > m_kf0
 Forward rate constants without perturbation. More...
 
vector< double > m_sbuf0
 
vector< double > m_state
 
vector< double > m_grt
 Standard chemical potentials for each species. More...
 
- Protected Attributes inherited from Kinetics
ValueCache m_cache
 Cache for saved calculations within each Kinetics object. More...
 
bool m_ready = false
 Boolean indicating whether Kinetics object is fully configured. More...
 
size_t m_kk = 0
 The number of species in all of the phases that participate in this kinetics mechanism. More...
 
vector< double > m_perturb
 Vector of perturbation factors for each reaction's rate of progress vector. More...
 
vector< shared_ptr< Reaction > > m_reactions
 Vector of Reaction objects represented by this Kinetics manager. More...
 
vector< shared_ptr< ThermoPhase > > m_thermo
 m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator More...
 
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...
 
map< string, size_t > m_phaseindex
 Mapping of the phase name to the position of the phase within the kinetics object. More...
 
size_t m_mindim = 4
 number of spatial dimensions of lowest-dimensional phase. More...
 
vector< double > m_rfn
 Forward rate constant for each reaction. More...
 
vector< double > m_delta_gibbs0
 Delta G^0 for all reactions. More...
 
vector< double > m_rkcn
 Reciprocal of the equilibrium constant in concentration units. More...
 
vector< double > m_ropf
 Forward rate-of-progress for each reaction. More...
 
vector< double > m_ropr
 Reverse rate-of-progress for each reaction. More...
 
vector< double > m_ropnet
 Net rate-of-progress for each reaction. More...
 
vector< double > m_dH
 The enthalpy change for each reaction to calculate Blowers-Masel rates. More...
 
vector< double > m_rbuf
 Buffer used for storage of intermediate reaction-specific results. More...
 
bool m_skipUndeclaredSpecies = false
 See skipUndeclaredSpecies() More...
 
bool m_skipUndeclaredThirdBodies = false
 See skipUndeclaredThirdBodies() More...
 
bool m_hasUndeclaredThirdBodies = false
 Flag indicating whether reactions include undeclared third bodies. More...
 
std::weak_ptr< Solutionm_root
 reference to Solution More...
 
StoichManagerN m_reactantStoich
 Stoichiometry manager for the reactants for each reaction. More...
 
StoichManagerN m_productStoich
 Stoichiometry manager for the products for each reaction. More...
 
StoichManagerN m_revProductStoich
 Stoichiometry manager for the products of reversible reactions. More...
 
Eigen::SparseMatrix< double > m_stoichMatrix
 Net stoichiometry (products - reactants) More...
 

Member Function Documentation

◆ kineticsType()

string kineticsType ( ) const
inlineoverridevirtual

Identifies the Kinetics manager type.

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

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

Reimplemented from Kinetics.

Definition at line 27 of file BulkKinetics.h.

◆ isReversible()

bool isReversible ( size_t  i)
overridevirtual

True if reaction i has been declared to be reversible.

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

Parameters
ireaction index

Reimplemented from Kinetics.

Definition at line 15 of file BulkKinetics.cpp.

◆ addReaction()

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

Add a single reaction to the mechanism.

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

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

Reimplemented from Kinetics.

Definition at line 19 of file BulkKinetics.cpp.

◆ modifyReaction()

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

Modify the rate expression associated with a reaction.

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

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

Reimplemented from Kinetics.

Definition at line 93 of file BulkKinetics.cpp.

◆ resizeSpecies()

void resizeSpecies ( )
overridevirtual

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

Automatically called before adding each Reaction and Phase.

Reimplemented from Kinetics.

Definition at line 118 of file BulkKinetics.cpp.

◆ resizeReactions()

void resizeReactions ( )
overridevirtual

Finalize Kinetics object and associated objects.

Reimplemented from Kinetics.

Definition at line 129 of file BulkKinetics.cpp.

◆ setMultiplier()

void setMultiplier ( size_t  i,
double  f 
)
overridevirtual

Set the multiplier for reaction i to f.

Parameters
iindex of the reaction
fvalue of the multiplier.

Reimplemented from Kinetics.

Definition at line 147 of file BulkKinetics.cpp.

◆ getFwdRateConstants()

void getFwdRateConstants ( double *  kfwd)
overridevirtual

Return the forward rate constants.

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

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

Reimplemented from Kinetics.

Definition at line 159 of file BulkKinetics.cpp.

◆ getEquilibriumConstants()

void getEquilibriumConstants ( double *  kc)
overridevirtual

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 from Kinetics.

Definition at line 168 of file BulkKinetics.cpp.

◆ getRevRateConstants()

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

Return the reverse rate constants.

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

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

Reimplemented from Kinetics.

Definition at line 185 of file BulkKinetics.cpp.

◆ getDeltaGibbs()

void getDeltaGibbs ( double *  deltaG)
overridevirtual

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

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

units = J kmol-1

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

Reimplemented from Kinetics.

Definition at line 204 of file BulkKinetics.cpp.

◆ getDeltaEnthalpy()

void getDeltaEnthalpy ( double *  deltaH)
overridevirtual

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

These values depend upon the concentration of the solution.

units = J kmol-1

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

Reimplemented from Kinetics.

Definition at line 212 of file BulkKinetics.cpp.

◆ getDeltaEntropy()

void getDeltaEntropy ( double *  deltaS)
overridevirtual

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

These values depend upon the concentration of the solution.

units = J kmol-1 Kelvin-1

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

Reimplemented from Kinetics.

Definition at line 220 of file BulkKinetics.cpp.

◆ getDeltaSSGibbs()

void getDeltaSSGibbs ( double *  deltaG)
overridevirtual

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

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

units = J kmol-1

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

Reimplemented from Kinetics.

Definition at line 228 of file BulkKinetics.cpp.

◆ getDeltaSSEnthalpy()

void getDeltaSSEnthalpy ( double *  deltaH)
overridevirtual

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

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

units = J kmol-1

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

Reimplemented from Kinetics.

Definition at line 239 of file BulkKinetics.cpp.

◆ getDeltaSSEntropy()

void getDeltaSSEntropy ( double *  deltaS)
overridevirtual

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

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

units = J kmol-1 Kelvin-1

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

Reimplemented from Kinetics.

Definition at line 250 of file BulkKinetics.cpp.

◆ getDerivativeSettings()

void getDerivativeSettings ( AnyMap settings) const
overridevirtual

Retrieve derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented from Kinetics.

Definition at line 263 of file BulkKinetics.cpp.

◆ setDerivativeSettings()

void setDerivativeSettings ( const AnyMap settings)
overridevirtual

Set/modify derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented from Kinetics.

Definition at line 270 of file BulkKinetics.cpp.

◆ getFwdRateConstants_ddT()

void getFwdRateConstants_ddT ( double *  dkfwd)
overridevirtual

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 from Kinetics.

Definition at line 284 of file BulkKinetics.cpp.

◆ getFwdRatesOfProgress_ddT()

void getFwdRatesOfProgress_ddT ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 291 of file BulkKinetics.cpp.

◆ getRevRatesOfProgress_ddT()

void getRevRatesOfProgress_ddT ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 298 of file BulkKinetics.cpp.

◆ getNetRatesOfProgress_ddT()

void getNetRatesOfProgress_ddT ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 312 of file BulkKinetics.cpp.

◆ getFwdRateConstants_ddP()

void getFwdRateConstants_ddP ( double *  dkfwd)
overridevirtual

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 from Kinetics.

Definition at line 326 of file BulkKinetics.cpp.

◆ getFwdRatesOfProgress_ddP()

void getFwdRatesOfProgress_ddP ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 333 of file BulkKinetics.cpp.

◆ getRevRatesOfProgress_ddP()

void getRevRatesOfProgress_ddP ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 340 of file BulkKinetics.cpp.

◆ getNetRatesOfProgress_ddP()

void getNetRatesOfProgress_ddP ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 347 of file BulkKinetics.cpp.

◆ getFwdRateConstants_ddC()

void getFwdRateConstants_ddC ( double *  dkfwd)
overridevirtual

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 from Kinetics.

Definition at line 354 of file BulkKinetics.cpp.

◆ getFwdRatesOfProgress_ddC()

void getFwdRatesOfProgress_ddC ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 361 of file BulkKinetics.cpp.

◆ getRevRatesOfProgress_ddC()

void getRevRatesOfProgress_ddC ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 368 of file BulkKinetics.cpp.

◆ getNetRatesOfProgress_ddC()

void getNetRatesOfProgress_ddC ( double *  drop)
overridevirtual

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 from Kinetics.

Definition at line 375 of file BulkKinetics.cpp.

◆ fwdRatesOfProgress_ddX()

Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX ( )
overridevirtual

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 from Kinetics.

Definition at line 387 of file BulkKinetics.cpp.

◆ revRatesOfProgress_ddX()

Eigen::SparseMatrix< double > revRatesOfProgress_ddX ( )
overridevirtual

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 from Kinetics.

Definition at line 397 of file BulkKinetics.cpp.

◆ netRatesOfProgress_ddX()

Eigen::SparseMatrix< double > netRatesOfProgress_ddX ( )
overridevirtual

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 from Kinetics.

Definition at line 408 of file BulkKinetics.cpp.

◆ fwdRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi ( )
overridevirtual

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

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

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

Reimplemented from Kinetics.

Definition at line 422 of file BulkKinetics.cpp.

◆ revRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > revRatesOfProgress_ddCi ( )
overridevirtual

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

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

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

Reimplemented from Kinetics.

Definition at line 432 of file BulkKinetics.cpp.

◆ netRatesOfProgress_ddCi()

Eigen::SparseMatrix< double > netRatesOfProgress_ddCi ( )
overridevirtual

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

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

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

Reimplemented from Kinetics.

Definition at line 443 of file BulkKinetics.cpp.

◆ getThirdBodyConcentrations()

void getThirdBodyConcentrations ( double *  concm)
overridevirtual

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.

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

Reimplemented from Kinetics.

Definition at line 541 of file BulkKinetics.cpp.

◆ thirdBodyConcentrations()

const vector<double>& thirdBodyConcentrations ( ) const
inlineoverridevirtual

Provide direct access to current third-body concentration values.

See also
getThirdBodyConcentrations.

Reimplemented from Kinetics.

Definition at line 90 of file BulkKinetics.h.

◆ processThirdBodies()

void processThirdBodies ( double *  rop)
protected

Multiply rate with third-body collider concentrations.

Definition at line 547 of file BulkKinetics.cpp.

◆ applyEquilibriumConstants()

void applyEquilibriumConstants ( double *  rop)
protected

Multiply rate with inverse equilibrium constant.

Definition at line 555 of file BulkKinetics.cpp.

◆ applyEquilibriumConstants_ddT()

void applyEquilibriumConstants_ddT ( double *  drkcn)
protected

Multiply rate with scaled temperature derivatives of the inverse equilibrium constant.

This (scaled) derivative is handled by a finite difference.

Definition at line 564 of file BulkKinetics.cpp.

◆ process_ddT()

void process_ddT ( const vector< double > &  in,
double *  drop 
)
protected

Process temperature derivative.

Parameters
inrate expression used for the derivative calculation
droppointer to output buffer

Definition at line 599 of file BulkKinetics.cpp.

◆ process_ddP()

void process_ddP ( const vector< double > &  in,
double *  drop 
)
protected

Process pressure derivative.

Parameters
inrate expression used for the derivative calculation
droppointer to output buffer

Definition at line 608 of file BulkKinetics.cpp.

◆ process_ddC()

void process_ddC ( StoichManagerN stoich,
const vector< double > &  in,
double *  drop,
bool  mass_action = true 
)
protected

Process concentration (molar density) derivative.

Parameters
stoichstoichiometry manager
inrate expression used for the derivative calculation
droppointer to output buffer
mass_actionboolean indicating whether law of mass action applies

Definition at line 617 of file BulkKinetics.cpp.

◆ calculateCompositionDerivatives()

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

Process derivatives.

Parameters
stoichstoichiometry manager
inrate expression used for the derivative calculation
ddXtrue: w.r.t mole fractions false: w.r.t species concentrations
Returns
a sparse matrix of derivative contributions for each reaction of dimensions nTotalReactions by nTotalSpecies

Definition at line 652 of file BulkKinetics.cpp.

◆ assertDerivativesValid()

void assertDerivativesValid ( const string &  name)
protected

Helper function ensuring that all rate derivatives can be calculated.

Parameters
namemethod name used for error output
Exceptions
CanteraErrorif ideal gas assumption does not hold

Definition at line 694 of file BulkKinetics.cpp.

Member Data Documentation

◆ m_bulk_rates

vector<unique_ptr<MultiRateBase> > m_bulk_rates
protected

Vector of rate handlers.

Definition at line 152 of file BulkKinetics.h.

◆ m_bulk_types

map<string, size_t> m_bulk_types
protected

Mapping of rate handlers.

Definition at line 153 of file BulkKinetics.h.

◆ m_revindex

vector<size_t> m_revindex
protected

Indices of reversible reactions.

Definition at line 155 of file BulkKinetics.h.

◆ m_irrev

vector<size_t> m_irrev
protected

Indices of irreversible reactions.

Definition at line 156 of file BulkKinetics.h.

◆ m_dn

vector<double> m_dn
protected

Difference between the global reactants order and the global products order.

Of type "double" to account for the fact that we can have real- valued stoichiometries.

Definition at line 161 of file BulkKinetics.h.

◆ m_multi_concm

ThirdBodyCalc m_multi_concm
protected

used with MultiRate evaluator

Definition at line 163 of file BulkKinetics.h.

◆ m_concm

vector<double> m_concm
protected

Third body concentrations.

Definition at line 166 of file BulkKinetics.h.

◆ m_act_conc

vector<double> m_act_conc
protected

Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations.

Definition at line 169 of file BulkKinetics.h.

◆ m_phys_conc

vector<double> m_phys_conc
protected

Physical concentrations, as calculated by ThermoPhase::getConcentrations.

Definition at line 172 of file BulkKinetics.h.

◆ m_jac_skip_third_bodies

bool m_jac_skip_third_bodies
protected

Derivative settings.

Definition at line 175 of file BulkKinetics.h.

◆ m_rbuf0

vector<double> m_rbuf0
protected

Buffers for partial rop results with length nReactions()

Definition at line 182 of file BulkKinetics.h.

◆ m_kf0

vector<double> m_kf0
protected

Forward rate constants without perturbation.

Definition at line 185 of file BulkKinetics.h.

◆ m_grt

vector<double> m_grt
protected

Standard chemical potentials for each species.

Definition at line 188 of file BulkKinetics.h.


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