Cantera  2.5.1
GasKinetics.cpp
Go to the documentation of this file.
1 /**
2  * @file GasKinetics.cpp Homogeneous kinetics in ideal gases
3  */
4 
5 // This file is part of Cantera. See License.txt in the top-level directory or
6 // at https://cantera.org/license.txt for license and copyright information.
7 
9 
10 using namespace std;
11 
12 namespace Cantera
13 {
14 GasKinetics::GasKinetics(thermo_t* thermo) :
15  BulkKinetics(thermo),
16  m_logp_ref(0.0),
17  m_logc_ref(0.0),
18  m_logStandConc(0.0),
19  m_pres(0.0)
20 {
21 }
22 
24 {
25  doublereal T = thermo().temperature();
26  doublereal P = thermo().pressure();
27  m_logStandConc = log(thermo().standardConcentration());
28  doublereal logT = log(T);
29 
30  if (T != m_temp) {
31  if (!m_rfn.empty()) {
32  m_rates.update(T, logT, m_rfn.data());
33  }
34 
35  if (!m_rfn_low.empty()) {
36  m_falloff_low_rates.update(T, logT, m_rfn_low.data());
37  m_falloff_high_rates.update(T, logT, m_rfn_high.data());
38  }
39  if (!falloff_work.empty()) {
40  m_falloffn.updateTemp(T, falloff_work.data());
41  }
42  updateKc();
43  m_ROP_ok = false;
44  }
45 
46  if (T != m_temp || P != m_pres) {
47  if (m_plog_rates.nReactions()) {
48  m_plog_rates.update(T, logT, m_rfn.data());
49  m_ROP_ok = false;
50  }
51 
52  if (m_cheb_rates.nReactions()) {
53  m_cheb_rates.update(T, logT, m_rfn.data());
54  m_ROP_ok = false;
55  }
56  }
57  m_pres = P;
58  m_temp = T;
59 }
60 
62 {
63  thermo().getActivityConcentrations(m_conc.data());
64  doublereal ctot = thermo().molarDensity();
65 
66  // 3-body reactions
67  if (!concm_3b_values.empty()) {
68  m_3b_concm.update(m_conc, ctot, concm_3b_values.data());
69  }
70 
71  // Falloff reactions
72  if (!concm_falloff_values.empty()) {
73  m_falloff_concm.update(m_conc, ctot, concm_falloff_values.data());
74  }
75 
76  // P-log reactions
77  if (m_plog_rates.nReactions()) {
78  double logP = log(thermo().pressure());
79  m_plog_rates.update_C(&logP);
80  }
81 
82  // Chebyshev reactions
83  if (m_cheb_rates.nReactions()) {
84  double log10P = log10(thermo().pressure());
85  m_cheb_rates.update_C(&log10P);
86  }
87 
88  m_ROP_ok = false;
89 }
90 
92 {
93  thermo().getStandardChemPotentials(m_grt.data());
94  fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
95 
96  // compute Delta G^0 for all reversible reactions
97  getRevReactionDelta(m_grt.data(), m_rkcn.data());
98 
99  doublereal rrt = 1.0 / thermo().RT();
100  for (size_t i = 0; i < m_revindex.size(); i++) {
101  size_t irxn = m_revindex[i];
102  m_rkcn[irxn] = std::min(exp(m_rkcn[irxn]*rrt - m_dn[irxn]*m_logStandConc),
103  BigNumber);
104  }
105 
106  for (size_t i = 0; i != m_irrev.size(); ++i) {
107  m_rkcn[ m_irrev[i] ] = 0.0;
108  }
109 }
110 
112 {
113  update_rates_T();
114  thermo().getStandardChemPotentials(m_grt.data());
115  fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
116 
117  // compute Delta G^0 for all reactions
118  getReactionDelta(m_grt.data(), m_rkcn.data());
119 
120  doublereal rrt = 1.0 / thermo().RT();
121  for (size_t i = 0; i < nReactions(); i++) {
122  kc[i] = exp(-m_rkcn[i]*rrt + m_dn[i]*m_logStandConc);
123  }
124 
125  // force an update of T-dependent properties, so that m_rkcn will
126  // be updated before it is used next.
127  m_temp = 0.0;
128 }
129 
130 void GasKinetics::processFalloffReactions()
131 {
132  // use m_ropr for temporary storage of reduced pressure
133  vector_fp& pr = m_ropr;
134 
135  for (size_t i = 0; i < m_falloff_low_rates.nReactions(); i++) {
136  pr[i] = concm_falloff_values[i] * m_rfn_low[i] / (m_rfn_high[i] + SmallNumber);
137  AssertFinite(pr[i], "GasKinetics::processFalloffReactions",
138  "pr[{}] is not finite.", i);
139  }
140 
141  m_falloffn.pr_to_falloff(pr.data(), falloff_work.data());
142 
143  for (size_t i = 0; i < m_falloff_low_rates.nReactions(); i++) {
144  if (reactionType(m_fallindx[i]) == FALLOFF_RXN) {
145  pr[i] *= m_rfn_high[i];
146  } else { // CHEMACT_RXN
147  pr[i] *= m_rfn_low[i];
148  }
149  m_ropf[m_fallindx[i]] = pr[i];
150  }
151 }
152 
153 void GasKinetics::updateROP()
154 {
155  update_rates_C();
156  update_rates_T();
157  if (m_ROP_ok) {
158  return;
159  }
160 
161  // copy rate coefficients into ropf
162  m_ropf = m_rfn;
163 
164  // multiply ropf by enhanced 3b conc for all 3b rxns
165  if (!concm_3b_values.empty()) {
166  m_3b_concm.multiply(m_ropf.data(), concm_3b_values.data());
167  }
168 
169  if (m_falloff_high_rates.nReactions()) {
170  processFalloffReactions();
171  }
172 
173  for (size_t i = 0; i < nReactions(); i++) {
174  // Scale the forward rate coefficient by the perturbation factor
175  m_ropf[i] *= m_perturb[i];
176  // For reverse rates computed from thermochemistry, multiply the forward
177  // rate coefficients by the reciprocals of the equilibrium constants
178  m_ropr[i] = m_ropf[i] * m_rkcn[i];
179  }
180 
181  // multiply ropf by concentration products
182  m_reactantStoich.multiply(m_conc.data(), m_ropf.data());
183 
184  // for reversible reactions, multiply ropr by concentration products
185  m_revProductStoich.multiply(m_conc.data(), m_ropr.data());
186 
187  for (size_t j = 0; j != nReactions(); ++j) {
188  m_ropnet[j] = m_ropf[j] - m_ropr[j];
189  }
190 
191  for (size_t i = 0; i < m_rfn.size(); i++) {
192  AssertFinite(m_rfn[i], "GasKinetics::updateROP",
193  "m_rfn[{}] is not finite.", i);
194  AssertFinite(m_ropf[i], "GasKinetics::updateROP",
195  "m_ropf[{}] is not finite.", i);
196  AssertFinite(m_ropr[i], "GasKinetics::updateROP",
197  "m_ropr[{}] is not finite.", i);
198  }
199  m_ROP_ok = true;
200 }
201 
202 void GasKinetics::getFwdRateConstants(doublereal* kfwd)
203 {
204  update_rates_C();
205  update_rates_T();
206 
207  // copy rate coefficients into ropf
208  m_ropf = m_rfn;
209 
210  // multiply ropf by enhanced 3b conc for all 3b rxns
211  if (!concm_3b_values.empty()) {
212  m_3b_concm.multiply(m_ropf.data(), concm_3b_values.data());
213  }
214 
215  if (m_falloff_high_rates.nReactions()) {
216  processFalloffReactions();
217  }
218 
219  for (size_t i = 0; i < nReactions(); i++) {
220  // multiply by perturbation factor
221  kfwd[i] = m_ropf[i] * m_perturb[i];
222  }
223 }
224 
225 bool GasKinetics::addReaction(shared_ptr<Reaction> r)
226 {
227  // operations common to all reaction types
228  bool added = BulkKinetics::addReaction(r);
229  if (!added) {
230  return false;
231  }
232 
233  switch (r->reaction_type) {
234  case ELEMENTARY_RXN:
235  addElementaryReaction(dynamic_cast<ElementaryReaction&>(*r));
236  break;
237  case THREE_BODY_RXN:
238  addThreeBodyReaction(dynamic_cast<ThreeBodyReaction&>(*r));
239  break;
240  case FALLOFF_RXN:
241  case CHEMACT_RXN:
242  addFalloffReaction(dynamic_cast<FalloffReaction&>(*r));
243  break;
244  case PLOG_RXN:
245  addPlogReaction(dynamic_cast<PlogReaction&>(*r));
246  break;
247  case CHEBYSHEV_RXN:
248  addChebyshevReaction(dynamic_cast<ChebyshevReaction&>(*r));
249  break;
250  default:
251  throw CanteraError("GasKinetics::addReaction",
252  "Unknown reaction type specified: {}", r->reaction_type);
253  }
254  return true;
255 }
256 
257 void GasKinetics::addFalloffReaction(FalloffReaction& r)
258 {
259  // install high and low rate coeff calculators and extend the high and low
260  // rate coeff value vectors
261  size_t nfall = m_falloff_high_rates.nReactions();
262  m_falloff_high_rates.install(nfall, r.high_rate);
263  m_rfn_high.push_back(0.0);
264  m_falloff_low_rates.install(nfall, r.low_rate);
265  m_rfn_low.push_back(0.0);
266 
267  // add this reaction number to the list of falloff reactions
268  m_fallindx.push_back(nReactions()-1);
269  m_rfallindx[nReactions()-1] = nfall;
270 
271  // install the enhanced third-body concentration calculator
272  map<size_t, double> efficiencies;
273  for (const auto& eff : r.third_body.efficiencies) {
274  size_t k = kineticsSpeciesIndex(eff.first);
275  if (k != npos) {
276  efficiencies[k] = eff.second;
277  } else if (!m_skipUndeclaredThirdBodies) {
278  throw CanteraError("GasKinetics::addFalloffReaction", "Found "
279  "third-body efficiency for undefined species '" + eff.first +
280  "' while adding reaction '" + r.equation() + "'");
281  }
282  }
283  m_falloff_concm.install(nfall, efficiencies,
285  concm_falloff_values.resize(m_falloff_concm.workSize());
286 
287  // install the falloff function calculator for this reaction
288  m_falloffn.install(nfall, r.reaction_type, r.falloff);
289  falloff_work.resize(m_falloffn.workSize());
290 }
291 
292 void GasKinetics::addThreeBodyReaction(ThreeBodyReaction& r)
293 {
294  m_rates.install(nReactions()-1, r.rate);
295  map<size_t, double> efficiencies;
296  for (const auto& eff : r.third_body.efficiencies) {
297  size_t k = kineticsSpeciesIndex(eff.first);
298  if (k != npos) {
299  efficiencies[k] = eff.second;
300  } else if (!m_skipUndeclaredThirdBodies) {
301  throw CanteraError("GasKinetics::addThreeBodyReaction", "Found "
302  "third-body efficiency for undefined species '" + eff.first +
303  "' while adding reaction '" + r.equation() + "'");
304  }
305  }
306  m_3b_concm.install(nReactions()-1, efficiencies,
307  r.third_body.default_efficiency);
308  concm_3b_values.resize(m_3b_concm.workSize());
309 }
310 
311 void GasKinetics::addPlogReaction(PlogReaction& r)
312 {
313  m_plog_rates.install(nReactions()-1, r.rate);
314 }
315 
316 void GasKinetics::addChebyshevReaction(ChebyshevReaction& r)
317 {
318  m_cheb_rates.install(nReactions()-1, r.rate);
319 }
320 
321 void GasKinetics::modifyReaction(size_t i, shared_ptr<Reaction> rNew)
322 {
323  // operations common to all reaction types
325 
326  switch (rNew->reaction_type) {
327  case ELEMENTARY_RXN:
328  modifyElementaryReaction(i, dynamic_cast<ElementaryReaction&>(*rNew));
329  break;
330  case THREE_BODY_RXN:
331  modifyThreeBodyReaction(i, dynamic_cast<ThreeBodyReaction&>(*rNew));
332  break;
333  case FALLOFF_RXN:
334  case CHEMACT_RXN:
335  modifyFalloffReaction(i, dynamic_cast<FalloffReaction&>(*rNew));
336  break;
337  case PLOG_RXN:
338  modifyPlogReaction(i, dynamic_cast<PlogReaction&>(*rNew));
339  break;
340  case CHEBYSHEV_RXN:
341  modifyChebyshevReaction(i, dynamic_cast<ChebyshevReaction&>(*rNew));
342  break;
343  default:
344  throw CanteraError("GasKinetics::modifyReaction",
345  "Unknown reaction type specified: {}", rNew->reaction_type);
346  }
347 
348  // invalidate all cached data
349  m_ROP_ok = false;
350  m_temp += 0.1234;
351  m_pres += 0.1234;
352 }
353 
354 void GasKinetics::modifyThreeBodyReaction(size_t i, ThreeBodyReaction& r)
355 {
356  m_rates.replace(i, r.rate);
357 }
358 
359 void GasKinetics::modifyFalloffReaction(size_t i, FalloffReaction& r)
360 {
361  size_t iFall = m_rfallindx[i];
362  m_falloff_high_rates.replace(iFall, r.high_rate);
363  m_falloff_low_rates.replace(iFall, r.low_rate);
364  m_falloffn.replace(iFall, r.falloff);
365 }
366 
367 void GasKinetics::modifyPlogReaction(size_t i, PlogReaction& r)
368 {
369  m_plog_rates.replace(i, r.rate);
370 }
371 
372 void GasKinetics::modifyChebyshevReaction(size_t i, ChebyshevReaction& r)
373 {
374  m_cheb_rates.replace(i, r.rate);
375 }
376 
378 {
380  m_logp_ref = log(thermo().refPressure()) - log(GasConstant);
381 }
382 
383 void GasKinetics::invalidateCache()
384 {
385  BulkKinetics::invalidateCache();
386  m_pres += 0.13579;
387 }
388 
389 }
Partial specialization of Kinetics for chemistry in a single bulk phase.
Definition: BulkKinetics.h:22
std::vector< size_t > m_revindex
Indices of reversible reactions.
Definition: BulkKinetics.h:50
vector_fp m_dn
Difference between the global reactants order and the global products order.
Definition: BulkKinetics.h:56
virtual bool addReaction(shared_ptr< Reaction > r)
Add a single reaction to the mechanism.
std::vector< size_t > m_irrev
Indices of irreversible reactions.
Definition: BulkKinetics.h:51
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
A pressure-dependent reaction parameterized by a bi-variate Chebyshev polynomial in temperature and p...
Definition: Reaction.h:187
A reaction which follows mass-action kinetics with a modified Arrhenius reaction rate.
Definition: Reaction.h:84
size_t workSize()
Size of the work array required to store intermediate results.
Definition: FalloffMgr.h:59
void pr_to_falloff(doublereal *values, const doublereal *work)
Given a vector of reduced pressures for each falloff reaction, replace each entry by the value of the...
Definition: FalloffMgr.h:79
void install(size_t rxn, int reactionType, shared_ptr< Falloff > f)
Install a new falloff function calculator.
Definition: FalloffMgr.h:39
void replace(size_t rxn, shared_ptr< Falloff > f)
Definition: FalloffMgr.h:54
void updateTemp(doublereal t, doublereal *work)
Update the cached temperature-dependent intermediate results for all installed falloff functions.
Definition: FalloffMgr.h:69
A reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Definition: Reaction.h:133
ThirdBody third_body
Relative efficiencies of third-body species in enhancing the reaction rate.
Definition: Reaction.h:150
Arrhenius low_rate
The rate constant in the low-pressure limit.
Definition: Reaction.h:144
Arrhenius high_rate
The rate constant in the high-pressure limit.
Definition: Reaction.h:147
shared_ptr< Falloff > falloff
Falloff function which determines how low_rate and high_rate are combined to determine the rate const...
Definition: Reaction.h:154
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
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:76
doublereal m_pres
Last pressure at which rates were evaluated.
Definition: GasKinetics.h:100
Rate1< Arrhenius > m_falloff_high_rates
Rate expressions for falloff reactions at the high-pressure limit.
Definition: GasKinetics.h:82
virtual void getFwdRateConstants(doublereal *kfwd)
Return the forward rate constants.
virtual bool addReaction(shared_ptr< Reaction > r)
Add a single reaction to the mechanism.
void updateKc()
Update the equilibrium constants in molar units.
Definition: GasKinetics.cpp:91
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
virtual void update_rates_C()
Update properties that depend on concentrations.
Definition: GasKinetics.cpp:61
virtual void update_rates_T()
Update temperature-dependent portions of reaction rates and falloff functions.
Definition: GasKinetics.cpp:23
std::vector< size_t > m_fallindx
Reaction index of each falloff reaction.
Definition: GasKinetics.h:72
Rate1< Arrhenius > m_falloff_low_rates
Rate expressions for falloff reactions at the low-pressure limit.
Definition: GasKinetics.h:79
virtual void getEquilibriumConstants(doublereal *kc)
Return a vector of Equilibrium constants.
vector_fp m_ropnet
Net rate-of-progress for each reaction.
Definition: Kinetics.h:936
thermo_t & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition: Kinetics.h:227
virtual int reactionType(size_t i) const
Flag specifying the type of reaction.
Definition: Kinetics.h:603
vector_fp m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition: Kinetics.h:876
virtual void getReactionDelta(const doublereal *property, doublereal *deltaProperty)
Change in species properties.
Definition: Kinetics.cpp:388
vector_fp m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition: Kinetics.h:927
vector_fp m_ropf
Forward rate-of-progress for each reaction.
Definition: Kinetics.h:930
virtual void getRevReactionDelta(const doublereal *g, doublereal *dg)
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition: Kinetics.cpp:398
vector_fp m_ropr
Reverse rate-of-progress for each reaction.
Definition: Kinetics.h:933
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition: Kinetics.h:700
bool m_skipUndeclaredThirdBodies
Definition: Kinetics.h:942
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition: Kinetics.cpp:588
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition: Kinetics.h:864
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition: Kinetics.h:135
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:261
vector_fp m_rfn
Forward rate constant for each reaction.
Definition: Kinetics.h:924
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition: Kinetics.h:861
double molarDensity() const
Molar density (kmol/m^3).
Definition: Phase.cpp:700
doublereal temperature() const
Temperature (K).
Definition: Phase.h:667
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Definition: Phase.h:679
A pressure-dependent reaction parameterized by logarithmically interpolating between Arrhenius rate e...
Definition: Reaction.h:175
int reaction_type
Type of the reaction.
Definition: Reaction.h:46
std::string equation() const
The chemical equation for this reaction.
Definition: Reaction.cpp:96
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
doublereal RT() const
Return the Gas Constant multiplied by the current temperature.
Definition: ThermoPhase.h:776
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
Definition: ThermoPhase.h:398
virtual void getStandardChemPotentials(doublereal *mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
Definition: ThermoPhase.h:574
double default_efficiency
The default third body efficiency for species not listed in efficiencies.
Definition: Reaction.h:111
Composition efficiencies
Map of species to third body efficiency.
Definition: Reaction.h:107
A reaction with a non-reacting third body "M" that acts to add or remove energy from the reacting spe...
Definition: Reaction.h:117
#define AssertFinite(expr, procedure,...)
Throw an exception if the specified exception is not a finite number.
Definition: ctexceptions.h:272
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:188
const double SmallNumber
smallest number to compare to zero.
Definition: ct_defs.h:149
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:180
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition: ct_defs.h:109
const double BigNumber
largest number to compare to inf.
Definition: ct_defs.h:151
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
const int CHEBYSHEV_RXN
A general gas-phase pressure-dependent reaction where k(T,P) is defined in terms of a bivariate Cheby...
Definition: reaction_defs.h:58
const int PLOG_RXN
A pressure-dependent rate expression consisting of several Arrhenius rate expressions evaluated at di...
Definition: reaction_defs.h:52
const int THREE_BODY_RXN
A gas-phase reaction that requires a third-body collision partner.
Definition: reaction_defs.h:38
const int ELEMENTARY_RXN
A reaction with a rate coefficient that depends only on temperature and voltage that also obeys mass-...
Definition: reaction_defs.h:32
const int CHEMACT_RXN
A chemical activation reaction.
Definition: reaction_defs.h:66
const int FALLOFF_RXN
The general form for a gas-phase association or dissociation reaction, with a pressure-dependent rate...
Definition: reaction_defs.h:44