Cantera 2.6.0
GasKinetics.h
Go to the documentation of this file.
1/**
2 * @file GasKinetics.h
3 * @ingroup chemkinetics
4 */
5
6// This file is part of Cantera. See License.txt in the top-level directory or
7// at https://cantera.org/license.txt for license and copyright information.
8
9#ifndef CT_GASKINETICS_H
10#define CT_GASKINETICS_H
11
12#include "BulkKinetics.h"
13#include "FalloffMgr.h"
14#include "Reaction.h"
15
16namespace Cantera
17{
18
19/**
20 * Kinetics manager for elementary gas-phase chemistry. This kinetics manager
21 * implements standard mass-action reaction rate expressions for low-density
22 * gases.
23 * @ingroup kinetics
24 */
26{
27public:
28 //! @name Constructors and General Information
29 //! @{
30
31 //! Constructor.
32 /*!
33 * @param thermo Pointer to the gas ThermoPhase (optional)
34 */
36
37 virtual std::string kineticsType() const {
38 return "Gas";
39 }
40
41 virtual void getThirdBodyConcentrations(double* concm);
42 virtual const vector_fp& thirdBodyConcentrations() const {
43 return m_concm;
44 }
45
46 //! @}
47 //! @name Reaction Rates Of Progress
48 //! @{
49
50 virtual void getEquilibriumConstants(doublereal* kc);
51 virtual void getFwdRateConstants(double* kfwd);
52
53 //! @}
54 //! @name Reaction Mechanism Setup Routines
55 //! @{
56 virtual bool addReaction(shared_ptr<Reaction> r, bool resize=true);
57 virtual void modifyReaction(size_t i, shared_ptr<Reaction> rNew);
58 virtual void invalidateCache();
59 //! @}
60
61 virtual void resizeReactions();
62 void updateROP();
63
64 virtual void getDerivativeSettings(AnyMap& settings) const;
65 virtual void setDerivativeSettings(const AnyMap& settings);
66 virtual void getFwdRateConstants_ddT(double* dkfwd);
67 virtual void getFwdRatesOfProgress_ddT(double* drop);
68 virtual void getRevRatesOfProgress_ddT(double* drop);
69 virtual void getNetRatesOfProgress_ddT(double* drop);
70 virtual void getFwdRateConstants_ddP(double* dkfwd);
71 virtual void getFwdRatesOfProgress_ddP(double* drop);
72 virtual void getRevRatesOfProgress_ddP(double* drop);
73 virtual void getNetRatesOfProgress_ddP(double* drop);
74 virtual void getFwdRateConstants_ddC(double* dkfwd);
75 virtual void getFwdRatesOfProgress_ddC(double* drop);
76 virtual void getRevRatesOfProgress_ddC(double* drop);
77 virtual void getNetRatesOfProgress_ddC(double* drop);
78 virtual Eigen::SparseMatrix<double> fwdRatesOfProgress_ddX();
79 virtual Eigen::SparseMatrix<double> revRatesOfProgress_ddX();
80 virtual Eigen::SparseMatrix<double> netRatesOfProgress_ddX();
81
82 //! Update temperature-dependent portions of reaction rates and falloff
83 //! functions.
84 virtual void update_rates_T();
85
86 //! Update properties that depend on concentrations.
87 //! Currently the enhanced collision partner concentrations are updated
88 //! here, as well as the pressure-dependent portion of P-log and Chebyshev
89 //! reactions.
90 virtual void update_rates_C();
91
92protected:
93 //! @name Internal service methods
94 /*!
95 * These methods are for @internal use, and seek to avoid code duplication
96 * while evaluating terms used for rate constants, rates of progress, and
97 * their derivatives.
98 */
99 //! @{
100
101 //! Calculate rate coefficients
102 void processFwdRateCoefficients(double* ropf);
103
104 //! Multiply rate with third-body collider concentrations
105 void processThirdBodies(double* rop);
106
107 //! Multiply rate with inverse equilibrium constant
108 void processEquilibriumConstants(double* rop);
109
110 //! Multiply rate with scaled temperature derivatives of the inverse
111 //! equilibrium constant
112 /*!
113 * This (scaled) derivative is handled by a finite difference.
114 */
115 void processEquilibriumConstants_ddT(double* drkcn);
116
117 //! Process temperature derivative
118 //! @param in rate expression used for the derivative calculation
119 //! @param drop pointer to output buffer
120 void process_ddT(const vector_fp& in, double* drop);
121
122 //! Process pressure derivative
123 //! @param in rate expression used for the derivative calculation
124 //! @param drop pointer to output buffer
125 void process_ddP(const vector_fp& in, double* drop);
126
127 //! Process concentration (molar density) derivative
128 //! @param stoich stoichiometry manager
129 //! @param in rate expression used for the derivative calculation
130 //! @param drop pointer to output buffer
131 //! @param mass_action boolean indicating whether law of mass action applies
132 void process_ddC(StoichManagerN& stoich, const vector_fp& in,
133 double* drop, bool mass_action=true);
134
135 //! Process mole fraction derivative
136 //! @param stoich stoichiometry manager
137 //! @param in rate expression used for the derivative calculation
138 Eigen::SparseMatrix<double> process_ddX(StoichManagerN& stoich,
139 const vector_fp& in);
140
141 //! Helper function ensuring that all rate derivatives can be calculated
142 //! @param name method name used for error output
143 //! @throw CanteraError if legacy rates are present
144 void assertDerivativesValid(const std::string& name);
145
146 //! @}
147
148 //! Reaction index of each falloff reaction
149 std::vector<size_t> m_fallindx; //!< @deprecated (legacy only)
150
151 //! Reaction index of each legacy reaction (old framework)
152 std::vector<size_t> m_legacy;
153
154 //! Map of reaction index to falloff reaction index (i.e indices in
155 //! #m_falloff_low_rates and #m_falloff_high_rates)
156 std::map<size_t, size_t> m_rfallindx; //!< @deprecated (legacy only)
157
158 //! Rate expressions for falloff reactions at the low-pressure limit
159 Rate1<Arrhenius2> m_falloff_low_rates; //!< @deprecated (legacy only)
160
161 //! Rate expressions for falloff reactions at the high-pressure limit
162 Rate1<Arrhenius2> m_falloff_high_rates; //!< @deprecated (legacy only)
163
164 FalloffMgr m_falloffn; //!< @deprecated (legacy only)
165
166 ThirdBodyCalc m_3b_concm; //!< @deprecated (legacy only)
167 ThirdBodyCalc m_falloff_concm; //!< @deprecated (legacy only)
168
169 Rate1<Plog> m_plog_rates; //!< @deprecated (legacy only)
170 Rate1<ChebyshevRate> m_cheb_rates; //!< @deprecated (legacy only)
171
172 //! @name Reaction rate data
173 //!@{
174 doublereal m_logStandConc;
175 vector_fp m_rfn_low; //!< @deprecated (legacy only)
176 vector_fp m_rfn_high; //!< @deprecated (legacy only)
177
178 doublereal m_pres; //!< Last pressure at which rates were evaluated
179 vector_fp falloff_work; //!< @deprecated (legacy only)
180 vector_fp concm_3b_values; //!< @deprecated (legacy only)
181 vector_fp concm_falloff_values; //!< @deprecated (legacy only)
182
183 //!@}
184
185 //! Buffers for partial rop results with length nReactions()
187 vector_fp m_rbuf1;
188 vector_fp m_rbuf2;
189 vector_fp m_sbuf0;
190 vector_fp m_state;
191
192 //! Derivative settings
194 bool m_jac_skip_falloff;
195 double m_jac_rtol_delta;
196
197 // functions marked as deprecated below are only used for XML import and
198 // transitional reaction types that are marked as '-legacy'
199
200 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
201 void processFalloffReactions(double* ropf);
202
203 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
205 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
207 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
209 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
211
212 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
214 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
215 void modifyFalloffReaction(size_t i, FalloffReaction2& r);
216 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
217 void modifyPlogReaction(size_t i, PlogReaction2& r);
218 //! @deprecated To be removed after Cantera 2.6 (replaced by MultiRate approach)
220
221 //! Update the equilibrium constants in molar units.
222 void updateKc();
223};
224
225}
226
227#endif
A map of string keys to values whose type can vary at runtime.
Definition: AnyMap.h:399
Partial specialization of Kinetics for chemistry in a single bulk phase.
Definition: BulkKinetics.h:24
vector_fp m_concm
Third body concentrations.
Definition: BulkKinetics.h:73
A pressure-dependent reaction parameterized by a bi-variate Chebyshev polynomial in temperature and p...
Definition: Reaction.h:386
A falloff manager that implements any set of falloff functions.
Definition: FalloffMgr.h:29
A reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Definition: Reaction.h:297
Kinetics manager for elementary gas-phase chemistry.
Definition: GasKinetics.h:26
void processFalloffReactions(double *ropf)
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition: GasKinetics.cpp:23
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...
Definition: GasKinetics.h:156
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.
Definition: GasKinetics.h:159
Rate1< ChebyshevRate > m_cheb_rates
Definition: GasKinetics.h:170
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.
Definition: GasKinetics.h:178
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
Definition: GasKinetics.h:169
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.
Definition: GasKinetics.h:162
void modifyPlogReaction(size_t i, PlogReaction2 &r)
GasKinetics(ThermoPhase *thermo=0)
Constructor.
Definition: GasKinetics.cpp:15
bool m_jac_skip_third_bodies
Derivative settings.
Definition: GasKinetics.h:193
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
Definition: GasKinetics.h:181
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)
Definition: GasKinetics.h:152
ThirdBodyCalc m_3b_concm
Definition: GasKinetics.h:166
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...
Definition: GasKinetics.cpp:34
virtual void update_rates_C()
Update properties that depend on concentrations.
Definition: GasKinetics.cpp:91
virtual void update_rates_T()
Update temperature-dependent portions of reaction rates and falloff functions.
Definition: GasKinetics.cpp:40
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.
Definition: GasKinetics.h:149
virtual void getEquilibriumConstants(doublereal *kc)
Return a vector of Equilibrium constants.
vector_fp m_rbuf0
Buffers for partial rop results with length nReactions()
Definition: GasKinetics.h:186
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.
Definition: GasKinetics.h:42
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.
Definition: GasKinetics.h:37
FalloffMgr m_falloffn
Definition: GasKinetics.h:164
void modifyChebyshevReaction(size_t i, ChebyshevReaction2 &r)
vector_fp falloff_work
Definition: GasKinetics.h:179
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
Definition: GasKinetics.h:180
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
Definition: GasKinetics.h:167
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition: Kinetics.h:231
A pressure-dependent reaction parameterized by logarithmically interpolating between Arrhenius rate e...
Definition: Reaction.h:365
This rate coefficient manager supports one parameterization of the rate constant of any type.
Definition: RateCoeffMgr.h:28
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
Calculate and apply third-body effects on reaction rates, including non- unity third-body efficiencie...
Definition: ThirdBodyCalc.h:20
A reaction with a non-reacting third body "M" that acts to add or remove energy from the reacting spe...
Definition: Reaction.h:269
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184