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

Kinetics manager for elementary gas-phase chemistry. More...

#include <GasKinetics.h>

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

Public Member Functions

virtual void resizeReactions ()
 Finalize Kinetics object and associated objects. More...
 
void updateROP ()
 
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 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 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 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 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 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 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 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 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_ddC (double *drop)
 Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant temperature, pressure 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 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 > fwdRatesOfProgress_ddX ()
 Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constant temperature, pressure and molar concentration. 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 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...
 
virtual void update_rates_T ()
 Update temperature-dependent portions of reaction rates and falloff functions. More...
 
virtual void update_rates_C ()
 Update properties that depend on concentrations. More...
 
Constructors and General Information
 GasKinetics (ThermoPhase *thermo=0)
 Constructor. More...
 
virtual std::string kineticsType () const
 Identifies the Kinetics manager type. 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...
 
Reaction Rates Of Progress
virtual void getEquilibriumConstants (doublereal *kc)
 Return a vector of Equilibrium constants. More...
 
virtual void getFwdRateConstants (double *kfwd)
 Return the forward rate constants. More...
 
Reaction Mechanism Setup Routines
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...
 
virtual void invalidateCache ()
 
- Public Member Functions inherited from BulkKinetics
 BulkKinetics (ThermoPhase *thermo=0)
 
virtual bool isReversible (size_t i)
 True if reaction i has been declared to be reversible. More...
 
virtual void getDeltaGibbs (doublereal *deltaG)
 Return the vector of values for the reaction Gibbs 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 getRevRateConstants (double *krev, bool doIrreversible=false)
 Return the reverse rate constants. More...
 
virtual void resizeSpecies ()
 Resize arrays with sizes that depend on the total number of species. More...
 
virtual void setMultiplier (size_t i, double f)
 Set the multiplier for reaction i to f. More...
 
void addThirdBody (shared_ptr< Reaction > r)
 
- Public Member Functions inherited from Kinetics
virtual std::pair< size_t, size_t > checkDuplicates (bool throw_err=true) const
 Check for unmarked duplicate reactions and unmatched marked duplicates. More...
 
void selectPhase (const double *data, const ThermoPhase *phase, double *phase_data)
 
virtual void setRoot (std::shared_ptr< Solution > root)
 Set root Solution holding all phase information. More...
 
virtual double reactionEnthalpy (const Composition &reactants, const Composition &products)
 Calculate the reaction enthalpy of a reaction which has not necessarily been added into the Kinetics object. More...
 
 Kinetics ()
 Default constructor. More...
 
virtual ~Kinetics ()
 
 Kinetics (const Kinetics &)=delete
 Kinetics objects are not copyable or assignable. More...
 
Kineticsoperator= (const Kinetics &)=delete
 
size_t nReactions () const
 Number of reactions in the reaction mechanism. More...
 
void checkReactionIndex (size_t m) const
 Check that the specified reaction index is in range Throws an exception if i is greater than nReactions() More...
 
void checkReactionArraySize (size_t ii) const
 Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions(). More...
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range Throws an exception if k is greater than nSpecies()-1. More...
 
void checkSpeciesArraySize (size_t mm) const
 Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies(). More...
 
size_t nPhases () const
 The number of phases participating in the reaction mechanism. More...
 
void checkPhaseIndex (size_t m) const
 Check that the specified phase index is in range Throws an exception if m is greater than nPhases() More...
 
void checkPhaseArraySize (size_t mm) const
 Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases(). More...
 
size_t phaseIndex (const std::string &ph) const
 Return the phase index of a phase in the list of phases defined within the object. More...
 
size_t surfacePhaseIndex () const
 This returns the integer index of the phase which has ThermoPhase type cSurf. More...
 
size_t reactionPhaseIndex () const
 Phase where the reactions occur. More...
 
ThermoPhasethermo (size_t n=0)
 This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism. More...
 
const ThermoPhasethermo (size_t n=0) const
 
size_t nTotalSpecies () const
 The total number of species in all phases participating in the kinetics mechanism. More...
 
size_t kineticsSpeciesIndex (size_t k, size_t n) const
 The location of species k of phase n in species arrays. More...
 
std::string kineticsSpeciesName (size_t k) const
 Return the name of the kth species in the kinetics manager. More...
 
size_t kineticsSpeciesIndex (const std::string &nm) const
 This routine will look up a species number based on the input std::string nm. More...
 
size_t kineticsSpeciesIndex (const std::string &nm, const std::string &ph) const
 This routine will look up a species number based on the input std::string nm. More...
 
ThermoPhasespeciesPhase (const std::string &nm)
 This function looks up the name of a species and returns a reference to the ThermoPhase object of the phase where the species resides. More...
 
const ThermoPhasespeciesPhase (const std::string &nm) const
 
ThermoPhasespeciesPhase (size_t k)
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the species' owning ThermoPhase object. More...
 
size_t speciesPhaseIndex (size_t k) const
 This function takes as an argument the kineticsSpecies index (that is, the list index in the list of species in the kinetics manager) and returns the index of the phase owning the species. More...
 
virtual void getFwdRatesOfProgress (doublereal *fwdROP)
 Return the forward rates of progress of the reactions. More...
 
virtual void getRevRatesOfProgress (doublereal *revROP)
 Return the Reverse rates of progress of the reactions. More...
 
virtual void getNetRatesOfProgress (doublereal *netROP)
 Net rates of progress. More...
 
virtual void getReactionDelta (const double *property, double *deltaProperty) const
 Change in species properties. More...
 
virtual void getRevReactionDelta (const double *g, double *dg) const
 Given an array of species properties 'g', return in array 'dg' the change in this quantity in the reversible reactions. More...
 
virtual void getDeltaElectrochemPotentials (doublereal *deltaM)
 Return the vector of values for the reaction electrochemical free energy change. More...
 
virtual void getCreationRates (doublereal *cdot)
 Species creation rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getDestructionRates (doublereal *ddot)
 Species destruction rates [kmol/m^3/s or kmol/m^2/s]. More...
 
virtual void getNetProductionRates (doublereal *wdot)
 Species net production rates [kmol/m^3/s or kmol/m^2/s]. More...
 
void getCreationRates_ddT (double *dwdot)
 Calculate derivatives for species creation rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getCreationRates_ddP (double *dwdot)
 Calculate derivatives for species creation rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getCreationRates_ddC (double *dwdot)
 Calculate derivatives for species creation rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > creationRates_ddX ()
 Calculate derivatives for species creation rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getDestructionRates_ddT (double *dwdot)
 Calculate derivatives for species destruction rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getDestructionRates_ddP (double *dwdot)
 Calculate derivatives for species destruction rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getDestructionRates_ddC (double *dwdot)
 Calculate derivatives for species destruction rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > destructionRates_ddX ()
 Calculate derivatives for species destruction rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
void getNetProductionRates_ddT (double *dwdot)
 Calculate derivatives for species net production rates with respect to temperature at constant pressure, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddP (double *dwdot)
 Calculate derivatives for species net production rates with respect to pressure at constant temperature, molar concentration and mole fractions. More...
 
void getNetProductionRates_ddC (double *dwdot)
 Calculate derivatives for species net production rates with respect to molar concentration at constant temperature, pressure and mole fractions. More...
 
Eigen::SparseMatrix< double > netProductionRates_ddX ()
 Calculate derivatives for species net production rates with respect to species mole fractions at constant temperature, pressure and molar concentration. More...
 
virtual double reactantStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a reactant in reaction i. More...
 
Eigen::SparseMatrix< double > reactantStoichCoeffs () const
 Stoichiometric coefficient matrix for reactants. More...
 
virtual double productStoichCoeff (size_t k, size_t i) const
 Stoichiometric coefficient of species k as a product in reaction i. More...
 
Eigen::SparseMatrix< double > productStoichCoeffs () const
 Stoichiometric coefficient matrix for products. More...
 
Eigen::SparseMatrix< double > revProductStoichCoeffs () const
 Stoichiometric coefficient matrix for products of reversible reactions. More...
 
virtual doublereal reactantOrder (size_t k, size_t i) const
 Reactant order of species k in reaction i. More...
 
virtual doublereal productOrder (int k, int i) const
 product Order of species k in reaction i. More...
 
virtual 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...
 
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 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...
 
shared_ptr< Reactionreaction (size_t i)
 Return the Reaction object for reaction i. More...
 
shared_ptr< const Reactionreaction (size_t i) const
 
void skipUndeclaredSpecies (bool skip)
 Determine behavior when adding a new reaction that contains species not defined in any of the phases associated with this kinetics manager. More...
 
bool skipUndeclaredSpecies () const
 
void skipUndeclaredThirdBodies (bool skip)
 Determine behavior when adding a new reaction that contains third-body efficiencies for species not defined in any of the phases associated with this kinetics manager. More...
 
bool skipUndeclaredThirdBodies () const
 
doublereal multiplier (size_t i) const
 The current value of the multiplier for reaction i. More...
 

Protected Member Functions

void processFalloffReactions (double *ropf)
 
void addThreeBodyReaction (ThreeBodyReaction2 &r)
 
void addFalloffReaction (FalloffReaction2 &r)
 
void addPlogReaction (PlogReaction2 &r)
 
void addChebyshevReaction (ChebyshevReaction2 &r)
 
void modifyThreeBodyReaction (size_t i, ThreeBodyReaction2 &r)
 
void modifyFalloffReaction (size_t i, FalloffReaction2 &r)
 
void modifyPlogReaction (size_t i, PlogReaction2 &r)
 
void modifyChebyshevReaction (size_t i, ChebyshevReaction2 &r)
 
void updateKc ()
 Update the equilibrium constants in molar units. More...
 
Internal service methods
void processFwdRateCoefficients (double *ropf)
 Calculate rate coefficients. More...
 
void processThirdBodies (double *rop)
 Multiply rate with third-body collider concentrations. More...
 
void processEquilibriumConstants (double *rop)
 Multiply rate with inverse equilibrium constant. More...
 
void processEquilibriumConstants_ddT (double *drkcn)
 Multiply rate with scaled temperature derivatives of the inverse equilibrium constant. More...
 
void process_ddT (const vector_fp &in, double *drop)
 Process temperature derivative. More...
 
void process_ddP (const vector_fp &in, double *drop)
 Process pressure derivative. More...
 
void process_ddC (StoichManagerN &stoich, const vector_fp &in, double *drop, bool mass_action=true)
 Process concentration (molar density) derivative. More...
 
Eigen::SparseMatrix< double > process_ddX (StoichManagerN &stoich, const vector_fp &in)
 Process mole fraction derivative. More...
 
void assertDerivativesValid (const std::string &name)
 Helper function ensuring that all rate derivatives can be calculated. More...
 
- Protected Member Functions inherited from BulkKinetics
virtual void addElementaryReaction (ElementaryReaction2 &r)
 
virtual void modifyElementaryReaction (size_t i, ElementaryReaction2 &rNew)
 
- Protected Member Functions inherited from Kinetics
double checkDuplicateStoich (std::map< int, double > &r1, std::map< int, double > &r2) const
 Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reactions are duplicates of one another, and 0.0 otherwise. More...
 
void checkReactionBalance (const Reaction &R)
 Check that the specified reaction is balanced (same number of atoms for each element in the reactants and products). More...
 

Protected Attributes

std::vector< size_t > m_fallindx
 Reaction index of each falloff reaction. More...
 
std::vector< size_t > m_legacy
 Reaction index of each legacy reaction (old framework) More...
 
std::map< size_t, size_t > m_rfallindx
 Map of reaction index to falloff reaction index (i.e indices in m_falloff_low_rates and m_falloff_high_rates) More...
 
Rate1< Arrhenius2m_falloff_low_rates
 Rate expressions for falloff reactions at the low-pressure limit. More...
 
Rate1< Arrhenius2m_falloff_high_rates
 Rate expressions for falloff reactions at the high-pressure limit. More...
 
FalloffMgr m_falloffn
 
ThirdBodyCalc m_3b_concm
 
ThirdBodyCalc m_falloff_concm
 
Rate1< Plogm_plog_rates
 
Rate1< ChebyshevRatem_cheb_rates
 
vector_fp m_rbuf0
 Buffers for partial rop results with length nReactions() More...
 
vector_fp m_rbuf1
 
vector_fp m_rbuf2
 
vector_fp m_sbuf0
 
vector_fp m_state
 
bool m_jac_skip_third_bodies
 Derivative settings. More...
 
bool m_jac_skip_falloff
 
double m_jac_rtol_delta
 
Reaction rate data
doublereal m_logStandConc
 
vector_fp m_rfn_low
 
vector_fp m_rfn_high
 
doublereal m_pres
 Last pressure at which rates were evaluated. More...
 
vector_fp falloff_work
 
vector_fp concm_3b_values
 
vector_fp concm_falloff_values
 
- Protected Attributes inherited from BulkKinetics
std::vector< unique_ptr< MultiRateBase > > m_bulk_rates
 Vector of rate handlers. More...
 
std::map< std::string, size_t > m_bulk_types
 Mapping of rate handlers. More...
 
Rate1< Arrhenius2m_rates
 
std::vector< size_t > m_revindex
 Indices of reversible reactions. More...
 
std::vector< size_t > m_irrev
 Indices of irreversible reactions. More...
 
vector_fp m_dn
 Difference between the global reactants order and the global products order. More...
 
ThirdBodyCalc3 m_multi_concm
 used with MultiRate evaluator More...
 
vector_fp m_concm
 Third body concentrations. More...
 
vector_fp m_act_conc
 Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations. More...
 
vector_fp m_phys_conc
 Physical concentrations, as calculated by ThermoPhase::getConcentrations. More...
 
vector_fp m_grt
 
bool m_ROP_ok
 
doublereal m_temp
 
- Protected Attributes inherited from Kinetics
ValueCache m_cache
 Cache for saved calculations within each Kinetics object. More...
 
bool m_ready
 Boolean indicating whether Kinetics object is fully configured. More...
 
size_t m_kk
 The number of species in all of the phases that participate in this kinetics mechanism. More...
 
vector_fp m_perturb
 Vector of perturbation factors for each reaction's rate of progress vector. More...
 
std::vector< shared_ptr< Reaction > > m_reactions
 Vector of Reaction objects represented by this Kinetics manager. More...
 
std::vector< ThermoPhase * > m_thermo
 m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator More...
 
std::vector< size_t > m_start
 m_start is a vector of integers specifying the beginning position for the species vector for the n'th phase in the kinetics class. More...
 
std::map< std::string, size_t > m_phaseindex
 Mapping of the phase name to the position of the phase within the kinetics object. More...
 
size_t m_surfphase
 Index in the list of phases of the one surface phase. More...
 
size_t m_rxnphase
 Phase Index where reactions are assumed to be taking place. More...
 
size_t m_mindim
 number of spatial dimensions of lowest-dimensional phase. More...
 
vector_fp m_rfn
 Forward rate constant for each reaction. More...
 
vector_fp m_delta_gibbs0
 Delta G^0 for all reactions. More...
 
vector_fp m_rkcn
 Reciprocal of the equilibrium constant in concentration units. More...
 
vector_fp m_ropf
 Forward rate-of-progress for each reaction. More...
 
vector_fp m_ropr
 Reverse rate-of-progress for each reaction. More...
 
vector_fp m_ropnet
 Net rate-of-progress for each reaction. More...
 
vector_fp m_dH
 The enthalpy change for each reaction to calculate Blowers-Masel rates. More...
 
vector_fp m_rbuf
 Buffer used for storage of intermediate reaction-specific results. More...
 
bool m_skipUndeclaredSpecies
 
bool m_skipUndeclaredThirdBodies
 
std::weak_ptr< Solutionm_root
 reference to Solution More...
 
StoichManagerN m_reactantStoich
 Stoichiometry manager for the reactants for each reaction. More...
 
StoichManagerN m_productStoich
 Stoichiometry manager for the products for each reaction. More...
 
StoichManagerN m_revProductStoich
 Stoichiometry manager for the products of reversible reactions. More...
 
Eigen::SparseMatrix< double > m_stoichMatrix
 Net stoichiometry (products - reactants) More...
 

Detailed Description

Kinetics manager for elementary gas-phase chemistry.

This kinetics manager implements standard mass-action reaction rate expressions for low-density gases.

Definition at line 25 of file GasKinetics.h.

Constructor & Destructor Documentation

◆ GasKinetics()

GasKinetics ( ThermoPhase thermo = 0)

Constructor.

Parameters
thermoPointer to the gas ThermoPhase (optional)

Definition at line 15 of file GasKinetics.cpp.

References GasKinetics::setDerivativeSettings().

Member Function Documentation

◆ kineticsType()

virtual std::string kineticsType ( ) const
inlinevirtual

Identifies the Kinetics manager type.

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

Reimplemented from Kinetics.

Definition at line 37 of file GasKinetics.h.

◆ getThirdBodyConcentrations()

void getThirdBodyConcentrations ( double *  concm)
virtual

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

Definition at line 34 of file GasKinetics.cpp.

References BulkKinetics::m_concm.

◆ thirdBodyConcentrations()

virtual const vector_fp & thirdBodyConcentrations ( ) const
inlinevirtual

Provide direct access to current third-body concentration values.

See also
getThirdBodyConcentrations.

Reimplemented from Kinetics.

Definition at line 42 of file GasKinetics.h.

References BulkKinetics::m_concm.

◆ getEquilibriumConstants()

void getEquilibriumConstants ( doublereal *  kc)
virtual

Return a vector of Equilibrium constants.

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

\[ 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 188 of file GasKinetics.cpp.

References Kinetics::getReactionDelta(), BulkKinetics::m_dn, GasKinetics::m_rbuf0, Kinetics::nReactions(), ThermoPhase::RT(), Kinetics::thermo(), and GasKinetics::update_rates_T().

◆ getFwdRateConstants()

void getFwdRateConstants ( double *  kfwd)
virtual

Return the forward rate constants.

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

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

Reimplemented from Kinetics.

Definition at line 254 of file GasKinetics.cpp.

References Cantera::legacy_rate_constants_used(), Kinetics::m_ropf, GasKinetics::processFwdRateCoefficients(), GasKinetics::processThirdBodies(), and Cantera::warn_deprecated().

◆ 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 from BulkKinetics.

Definition at line 574 of file GasKinetics.cpp.

References GasKinetics::addChebyshevReaction(), GasKinetics::addFalloffReaction(), GasKinetics::addPlogReaction(), BulkKinetics::addReaction(), GasKinetics::addThreeBodyReaction(), GasKinetics::m_legacy, and Kinetics::nReactions().

◆ modifyReaction()

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

Modify the rate expression associated with a reaction.

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

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

Reimplemented from BulkKinetics.

Definition at line 662 of file GasKinetics.cpp.

References GasKinetics::modifyChebyshevReaction(), GasKinetics::modifyFalloffReaction(), GasKinetics::modifyPlogReaction(), BulkKinetics::modifyReaction(), and GasKinetics::modifyThreeBodyReaction().

◆ invalidateCache()

void invalidateCache ( )
virtual

Reimplemented from BulkKinetics.

Definition at line 716 of file GasKinetics.cpp.

◆ resizeReactions()

void resizeReactions ( )
virtual

Finalize Kinetics object and associated objects.

Reimplemented from BulkKinetics.

Definition at line 23 of file GasKinetics.cpp.

References GasKinetics::m_rbuf0, Kinetics::nReactions(), Kinetics::nTotalSpecies(), BulkKinetics::resizeReactions(), and Kinetics::thermo().

◆ updateROP()

void updateROP ( )
virtual

Reimplemented from Kinetics.

Definition at line 227 of file GasKinetics.cpp.

◆ getDerivativeSettings()

void getDerivativeSettings ( AnyMap settings) const
virtual

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

Definition at line 273 of file GasKinetics.cpp.

References GasKinetics::m_jac_skip_third_bodies.

◆ setDerivativeSettings()

void setDerivativeSettings ( const AnyMap settings)
virtual

Set/modify derivative settings.

Parameters
settingsAnyMap containing settings determining derivative evaluation.

Reimplemented from Kinetics.

Definition at line 280 of file GasKinetics.cpp.

References AnyMap::empty(), AnyMap::getBool(), AnyMap::getDouble(), AnyMap::hasKey(), and GasKinetics::m_jac_skip_third_bodies.

Referenced by GasKinetics::GasKinetics().

◆ getFwdRateConstants_ddT()

void getFwdRateConstants_ddT ( double *  dkfwd)
virtual

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 351 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_rfn, and GasKinetics::process_ddT().

◆ getFwdRatesOfProgress_ddT()

void getFwdRatesOfProgress_ddT ( double *  drop)
virtual

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 358 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropf, and GasKinetics::process_ddT().

◆ getRevRatesOfProgress_ddT()

void getRevRatesOfProgress_ddT ( double *  drop)
virtual

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 365 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropr, Kinetics::nReactions(), GasKinetics::process_ddT(), and GasKinetics::processEquilibriumConstants_ddT().

◆ getNetRatesOfProgress_ddT()

void getNetRatesOfProgress_ddT ( double *  drop)
virtual

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 379 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropnet, Kinetics::m_ropr, Kinetics::nReactions(), GasKinetics::process_ddT(), and GasKinetics::processEquilibriumConstants_ddT().

◆ getFwdRateConstants_ddP()

void getFwdRateConstants_ddP ( double *  dkfwd)
virtual

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 402 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_rfn, and GasKinetics::process_ddP().

◆ getFwdRatesOfProgress_ddP()

void getFwdRatesOfProgress_ddP ( double *  drop)
virtual

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 409 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropf, and GasKinetics::process_ddP().

◆ getRevRatesOfProgress_ddP()

void getRevRatesOfProgress_ddP ( double *  drop)
virtual

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 416 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropr, and GasKinetics::process_ddP().

◆ getNetRatesOfProgress_ddP()

void getNetRatesOfProgress_ddP ( double *  drop)
virtual

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 423 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_ropnet, and GasKinetics::process_ddP().

◆ getFwdRateConstants_ddC()

void getFwdRateConstants_ddC ( double *  dkfwd)
virtual

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 466 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_reactantStoich, Kinetics::m_rfn, and GasKinetics::process_ddC().

◆ getFwdRatesOfProgress_ddC()

void getFwdRatesOfProgress_ddC ( double *  drop)
virtual

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 473 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_reactantStoich, Kinetics::m_ropf, and GasKinetics::process_ddC().

◆ getRevRatesOfProgress_ddC()

void getRevRatesOfProgress_ddC ( double *  drop)
virtual

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 480 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_revProductStoich, Kinetics::m_ropr, and GasKinetics::process_ddC().

◆ getNetRatesOfProgress_ddC()

void getNetRatesOfProgress_ddC ( double *  drop)
virtual

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 487 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, Kinetics::m_ropf, Kinetics::m_ropr, Kinetics::nReactions(), and GasKinetics::process_ddC().

◆ fwdRatesOfProgress_ddX()

Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX ( )
virtual

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 539 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), GasKinetics::m_rbuf0, Kinetics::m_reactantStoich, GasKinetics::process_ddX(), and GasKinetics::processFwdRateCoefficients().

◆ revRatesOfProgress_ddX()

Eigen::SparseMatrix< double > revRatesOfProgress_ddX ( )
virtual

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 549 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), GasKinetics::m_rbuf0, Kinetics::m_revProductStoich, GasKinetics::process_ddX(), GasKinetics::processEquilibriumConstants(), and GasKinetics::processFwdRateCoefficients().

◆ netRatesOfProgress_ddX()

Eigen::SparseMatrix< double > netRatesOfProgress_ddX ( )
virtual

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 560 of file GasKinetics.cpp.

References GasKinetics::assertDerivativesValid(), GasKinetics::m_rbuf0, Kinetics::m_reactantStoich, Kinetics::m_revProductStoich, GasKinetics::process_ddX(), GasKinetics::processEquilibriumConstants(), and GasKinetics::processFwdRateCoefficients().

◆ update_rates_T()

void update_rates_T ( )
virtual

Update temperature-dependent portions of reaction rates and falloff functions.

Definition at line 40 of file GasKinetics.cpp.

Referenced by GasKinetics::getEquilibriumConstants(), and GasKinetics::processFwdRateCoefficients().

◆ update_rates_C()

void update_rates_C ( )
virtual

◆ processFwdRateCoefficients()

void processFwdRateCoefficients ( double *  ropf)
protected

Calculate rate coefficients.

These methods are for

use, and seek to avoid code duplication while evaluating terms used for rate constants, rates of progress, and their derivatives.

Definition at line 148 of file GasKinetics.cpp.

References GasKinetics::m_falloff_high_rates, Kinetics::m_perturb, Kinetics::m_rfn, Kinetics::nReactions(), GasKinetics::processFalloffReactions(), GasKinetics::update_rates_C(), and GasKinetics::update_rates_T().

Referenced by GasKinetics::fwdRatesOfProgress_ddX(), GasKinetics::getFwdRateConstants(), GasKinetics::netRatesOfProgress_ddX(), and GasKinetics::revRatesOfProgress_ddX().

◆ processThirdBodies()

void processThirdBodies ( double *  rop)
protected

◆ processEquilibriumConstants()

void processEquilibriumConstants ( double *  rop)
protected

Multiply rate with inverse equilibrium constant.

Definition at line 179 of file GasKinetics.cpp.

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

Referenced by GasKinetics::netRatesOfProgress_ddX(), and GasKinetics::revRatesOfProgress_ddX().

◆ processEquilibriumConstants_ddT()

void processEquilibriumConstants_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 307 of file GasKinetics.cpp.

Referenced by GasKinetics::getNetRatesOfProgress_ddT(), and GasKinetics::getRevRatesOfProgress_ddT().

◆ process_ddT()

void process_ddT ( const vector_fp in,
double *  drop 
)
protected

Process temperature derivative.

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

Definition at line 342 of file GasKinetics.cpp.

References BulkKinetics::m_bulk_rates, and Kinetics::m_rfn.

Referenced by GasKinetics::getFwdRateConstants_ddT(), GasKinetics::getFwdRatesOfProgress_ddT(), GasKinetics::getNetRatesOfProgress_ddT(), and GasKinetics::getRevRatesOfProgress_ddT().

◆ process_ddP()

void process_ddP ( const vector_fp in,
double *  drop 
)
protected

Process pressure derivative.

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

Definition at line 393 of file GasKinetics.cpp.

References BulkKinetics::m_bulk_rates, and Kinetics::m_rfn.

Referenced by GasKinetics::getFwdRateConstants_ddP(), GasKinetics::getFwdRatesOfProgress_ddP(), GasKinetics::getNetRatesOfProgress_ddP(), and GasKinetics::getRevRatesOfProgress_ddP().

◆ process_ddC()

void process_ddC ( StoichManagerN stoich,
const vector_fp 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 430 of file GasKinetics.cpp.

References ThirdBodyCalc3::empty(), BulkKinetics::m_bulk_rates, BulkKinetics::m_concm, GasKinetics::m_jac_skip_third_bodies, BulkKinetics::m_multi_concm, Kinetics::m_rfn, Phase::molarDensity(), Kinetics::nReactions(), StoichManagerN::scale(), ThirdBodyCalc3::scale(), ThirdBodyCalc3::scaleM(), and Kinetics::thermo().

Referenced by GasKinetics::getFwdRateConstants_ddC(), GasKinetics::getFwdRatesOfProgress_ddC(), GasKinetics::getNetRatesOfProgress_ddC(), and GasKinetics::getRevRatesOfProgress_ddC().

◆ process_ddX()

Eigen::SparseMatrix< double > process_ddX ( StoichManagerN stoich,
const vector_fp in 
)
protected

◆ assertDerivativesValid()

void assertDerivativesValid ( const std::string &  name)
protected

◆ processFalloffReactions()

void processFalloffReactions ( double *  ropf)
protected

◆ addThreeBodyReaction()

void addThreeBodyReaction ( ThreeBodyReaction2 r)
protected

◆ addFalloffReaction()

void addFalloffReaction ( FalloffReaction2 r)
protected

◆ addPlogReaction()

void addPlogReaction ( PlogReaction2 r)
protected
Deprecated:
To be removed after Cantera 2.6 (replaced by MultiRate approach)

Definition at line 652 of file GasKinetics.cpp.

References GasKinetics::m_plog_rates, and Kinetics::nReactions().

Referenced by GasKinetics::addReaction().

◆ addChebyshevReaction()

void addChebyshevReaction ( ChebyshevReaction2 r)
protected
Deprecated:
To be removed after Cantera 2.6 (replaced by MultiRate approach)

Definition at line 657 of file GasKinetics.cpp.

References GasKinetics::m_cheb_rates, and Kinetics::nReactions().

Referenced by GasKinetics::addReaction().

◆ modifyThreeBodyReaction()

void modifyThreeBodyReaction ( size_t  i,
ThreeBodyReaction2 r 
)
protected
Deprecated:
To be removed after Cantera 2.6 (replaced by MultiRate approach)

Definition at line 693 of file GasKinetics.cpp.

References BulkKinetics::m_rates.

Referenced by GasKinetics::modifyReaction().

◆ modifyFalloffReaction()

void modifyFalloffReaction ( size_t  i,
FalloffReaction2 r 
)
protected

◆ modifyPlogReaction()

void modifyPlogReaction ( size_t  i,
PlogReaction2 r 
)
protected
Deprecated:
To be removed after Cantera 2.6 (replaced by MultiRate approach)

Definition at line 706 of file GasKinetics.cpp.

References GasKinetics::m_plog_rates.

Referenced by GasKinetics::modifyReaction().

◆ modifyChebyshevReaction()

void modifyChebyshevReaction ( size_t  i,
ChebyshevReaction2 r 
)
protected
Deprecated:
To be removed after Cantera 2.6 (replaced by MultiRate approach)

Definition at line 711 of file GasKinetics.cpp.

References GasKinetics::m_cheb_rates.

Referenced by GasKinetics::modifyReaction().

◆ updateKc()

void updateKc ( )
protected

Member Data Documentation

◆ m_fallindx

std::vector<size_t> m_fallindx
protected

Reaction index of each falloff reaction.

Deprecated:
(legacy only)

Definition at line 149 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::processFalloffReactions().

◆ m_legacy

std::vector<size_t> m_legacy
protected

Reaction index of each legacy reaction (old framework)

Definition at line 152 of file GasKinetics.h.

Referenced by GasKinetics::addReaction(), and GasKinetics::assertDerivativesValid().

◆ m_rfallindx

std::map<size_t, size_t> m_rfallindx
protected

Map of reaction index to falloff reaction index (i.e indices in m_falloff_low_rates and m_falloff_high_rates)

Deprecated:
(legacy only)

Definition at line 156 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::modifyFalloffReaction().

◆ m_falloff_low_rates

Rate1<Arrhenius2> m_falloff_low_rates
protected

Rate expressions for falloff reactions at the low-pressure limit.

Deprecated:
(legacy only)

Definition at line 159 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), GasKinetics::modifyFalloffReaction(), and GasKinetics::processFalloffReactions().

◆ m_falloff_high_rates

Rate1<Arrhenius2> m_falloff_high_rates
protected

Rate expressions for falloff reactions at the high-pressure limit.

Deprecated:
(legacy only)

Definition at line 162 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), GasKinetics::modifyFalloffReaction(), and GasKinetics::processFwdRateCoefficients().

◆ m_falloffn

FalloffMgr m_falloffn
protected

◆ m_3b_concm

ThirdBodyCalc m_3b_concm
protected

◆ m_falloff_concm

ThirdBodyCalc m_falloff_concm
protected
Deprecated:
(legacy only)

Definition at line 167 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::update_rates_C().

◆ m_plog_rates

Rate1<Plog> m_plog_rates
protected

◆ m_cheb_rates

Rate1<ChebyshevRate> m_cheb_rates
protected

◆ m_logStandConc

doublereal m_logStandConc
protected

Definition at line 174 of file GasKinetics.h.

◆ m_rfn_low

vector_fp m_rfn_low
protected
Deprecated:
(legacy only)

Definition at line 175 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::processFalloffReactions().

◆ m_rfn_high

vector_fp m_rfn_high
protected
Deprecated:
(legacy only)

Definition at line 176 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::processFalloffReactions().

◆ m_pres

doublereal m_pres
protected

Last pressure at which rates were evaluated.

Definition at line 178 of file GasKinetics.h.

◆ falloff_work

vector_fp falloff_work
protected
Deprecated:
(legacy only)

Definition at line 179 of file GasKinetics.h.

Referenced by GasKinetics::addFalloffReaction(), and GasKinetics::processFalloffReactions().

◆ concm_3b_values

vector_fp concm_3b_values
protected

◆ concm_falloff_values

vector_fp concm_falloff_values
protected

◆ m_rbuf0

vector_fp m_rbuf0
protected

◆ m_rbuf1

vector_fp m_rbuf1
protected

Definition at line 187 of file GasKinetics.h.

◆ m_rbuf2

vector_fp m_rbuf2
protected

Definition at line 188 of file GasKinetics.h.

◆ m_sbuf0

vector_fp m_sbuf0
protected

Definition at line 189 of file GasKinetics.h.

◆ m_state

vector_fp m_state
protected

Definition at line 190 of file GasKinetics.h.

◆ m_jac_skip_third_bodies

bool m_jac_skip_third_bodies
protected

◆ m_jac_skip_falloff

bool m_jac_skip_falloff
protected

Definition at line 194 of file GasKinetics.h.

◆ m_jac_rtol_delta

double m_jac_rtol_delta
protected

Definition at line 195 of file GasKinetics.h.


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