9#ifndef CT_GASKINETICS_H
10#define CT_GASKINETICS_H
56 virtual bool addReaction(shared_ptr<Reaction> r,
bool resize=
true);
58 virtual void invalidateCache();
133 double* drop,
bool mass_action=
true);
174 doublereal m_logStandConc;
194 bool m_jac_skip_falloff;
195 double m_jac_rtol_delta;
A map of string keys to values whose type can vary at runtime.
Partial specialization of Kinetics for chemistry in a single bulk phase.
vector_fp m_concm
Third body concentrations.
A pressure-dependent reaction parameterized by a bi-variate Chebyshev polynomial in temperature and p...
A falloff manager that implements any set of falloff functions.
A reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Kinetics manager for elementary gas-phase chemistry.
void processFalloffReactions(double *ropf)
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Eigen::SparseMatrix< double > process_ddX(StoichManagerN &stoich, const vector_fp &in)
Process mole fraction derivative.
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
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_hig...
virtual void getFwdRateConstants(double *kfwd)
Return the forward rate constants.
virtual void getDerivativeSettings(AnyMap &settings) const
Retrieve derivative settings.
virtual void getFwdRateConstants_ddC(double *dkfwd)
Calculate derivatives for forward rate constants with respect to molar concentration at constant temp...
Rate1< Arrhenius2 > m_falloff_low_rates
Rate expressions for falloff reactions at the low-pressure limit.
Rate1< ChebyshevRate > m_cheb_rates
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
doublereal m_pres
Last pressure at which rates were evaluated.
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX()
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
void addFalloffReaction(FalloffReaction2 &r)
void processFwdRateCoefficients(double *ropf)
Calculate rate coefficients.
Rate1< Plog > m_plog_rates
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
void addPlogReaction(PlogReaction2 &r)
void modifyFalloffReaction(size_t i, FalloffReaction2 &r)
Rate1< Arrhenius2 > m_falloff_high_rates
Rate expressions for falloff reactions at the high-pressure limit.
void modifyPlogReaction(size_t i, PlogReaction2 &r)
GasKinetics(ThermoPhase *thermo=0)
Constructor.
bool m_jac_skip_third_bodies
Derivative settings.
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
void process_ddP(const vector_fp &in, double *drop)
Process pressure derivative.
void assertDerivativesValid(const std::string &name)
Helper function ensuring that all rate derivatives can be calculated.
void processEquilibriumConstants(double *rop)
Multiply rate with inverse equilibrium constant.
void updateKc()
Update the equilibrium constants in molar units.
void processThirdBodies(double *rop)
Multiply rate with third-body collider concentrations.
void process_ddC(StoichManagerN &stoich, const vector_fp &in, double *drop, bool mass_action=true)
Process concentration (molar density) derivative.
vector_fp concm_falloff_values
virtual void getFwdRateConstants_ddT(double *dkfwd)
Calculate derivatives for forward rate constants with respect to temperature at constant pressure,...
std::vector< size_t > m_legacy
Reaction index of each legacy reaction (old framework)
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
virtual void getThirdBodyConcentrations(double *concm)
Return a vector of values of effective concentrations of third-body collision partners of any reactio...
virtual void update_rates_C()
Update properties that depend on concentrations.
virtual void update_rates_T()
Update temperature-dependent portions of reaction rates and falloff functions.
virtual void getFwdRateConstants_ddP(double *dkfwd)
Calculate derivatives for forward rate constants with respect to pressure at constant temperature,...
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
std::vector< size_t > m_fallindx
Reaction index of each falloff reaction.
virtual void getEquilibriumConstants(doublereal *kc)
Return a vector of Equilibrium constants.
vector_fp m_rbuf0
Buffers for partial rop results with length nReactions()
void process_ddT(const vector_fp &in, double *drop)
Process temperature derivative.
virtual const vector_fp & thirdBodyConcentrations() const
Provide direct access to current third-body concentration values.
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
virtual void setDerivativeSettings(const AnyMap &settings)
Set/modify derivative settings.
void modifyThreeBodyReaction(size_t i, ThreeBodyReaction2 &r)
virtual std::string kineticsType() const
Identifies the Kinetics manager type.
void modifyChebyshevReaction(size_t i, ChebyshevReaction2 &r)
void addChebyshevReaction(ChebyshevReaction2 &r)
void addThreeBodyReaction(ThreeBodyReaction2 &r)
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
vector_fp concm_3b_values
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
void processEquilibriumConstants_ddT(double *drkcn)
Multiply rate with scaled temperature derivatives of the inverse equilibrium constant.
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
ThirdBodyCalc m_falloff_concm
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
A pressure-dependent reaction parameterized by logarithmically interpolating between Arrhenius rate e...
This rate coefficient manager supports one parameterization of the rate constant of any type.
Base class for a phase with thermodynamic properties.
Calculate and apply third-body effects on reaction rates, including non- unity third-body efficiencie...
A reaction with a non-reacting third body "M" that acts to add or remove energy from the reacting spe...
Namespace for the Cantera kernel.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.