Cantera  2.1.2
GasKinetics.cpp
Go to the documentation of this file.
1 /**
2  * @file GasKinetics.cpp
3  *
4  * Homogeneous kinetics in ideal gases
5  */
6 
7 // Copyright 2001 California Institute of Technology
8 
10 
15 
16 using namespace std;
17 
18 namespace Cantera
19 {
20 GasKinetics::
21 GasKinetics(thermo_t* thermo) :
22  Kinetics(),
23  m_nfall(0),
24  m_nirrev(0),
25  m_nrev(0),
26  m_logp_ref(0.0),
27  m_logc_ref(0.0),
28  m_logStandConc(0.0),
29  m_ROP_ok(false),
30  m_temp(0.0),
31  m_pres(0.0),
32  m_finalized(false)
33 {
34  if (thermo != 0) {
35  addPhase(*thermo);
36  }
37  m_temp = 0.0;
38 }
39 
41  Kinetics(),
42  m_nfall(0),
43  m_nirrev(0),
44  m_nrev(0),
45  m_logp_ref(0.0),
46  m_logc_ref(0.0),
47  m_logStandConc(0.0),
48  m_ROP_ok(false),
49  m_temp(0.0),
50  m_pres(0.0),
51  m_finalized(false)
52 {
53  m_temp = 0.0;
54  *this = right;
55 }
56 
58 {
59  if (this == &right) {
60  return *this;
61  }
62 
63  Kinetics::operator=(right);
64 
65  m_nfall = right.m_nfall;
66  m_fallindx = right.m_fallindx;
67  m_falloff_low_rates = right.m_falloff_low_rates;
68  m_falloff_high_rates = right.m_falloff_high_rates;
69  m_rates = right.m_rates;
70  m_index = right.m_index;
71  m_falloffn = right.m_falloffn;
72  m_3b_concm = right.m_3b_concm;
73  m_falloff_concm = right.m_falloff_concm;
74  m_irrev = right.m_irrev;
75  m_plog_rates = right.m_plog_rates;
76  m_cheb_rates = right.m_cheb_rates;
77 
78  m_rxnstoich = right.m_rxnstoich;
79 
80  m_fwdOrder = right.m_fwdOrder;
81  m_nirrev = right.m_nirrev;
82  m_nrev = right.m_nrev;
83  m_rgroups = right.m_rgroups;
84  m_pgroups = right.m_pgroups;
85  m_rxntype = right.m_rxntype;
86  m_rrxn = right.m_rrxn;
87  m_prxn = right.m_prxn;
88  m_dn = right.m_dn;
89  m_revindex = right.m_revindex;
90  m_rxneqn = right.m_rxneqn;
91 
92  m_logp_ref = right.m_logp_ref;
93  m_logc_ref = right.m_logc_ref;
94  m_logStandConc = right.m_logStandConc;
95  m_ropf = right.m_ropf;
96  m_ropr = right.m_ropr;
97  m_ropnet = right.m_ropnet;
98  m_rfn_low = right.m_rfn_low;
99  m_rfn_high = right.m_rfn_high;
100  m_ROP_ok = right.m_ROP_ok;
101  m_temp = right.m_temp;
102  m_rfn = right.m_rfn;
103  falloff_work = right.falloff_work;
104  concm_3b_values = right.concm_3b_values;
105  concm_falloff_values = right.concm_falloff_values;
106  m_rkcn = right.m_rkcn;
107 
108  m_conc = right.m_conc;
109  m_grt = right.m_grt;
110  m_finalized = right.m_finalized;
111 
112  throw CanteraError("GasKinetics::operator=()",
113  "Unfinished implementation");
114 
115  return *this;
116 }
117 
118 Kinetics* GasKinetics::duplMyselfAsKinetics(const std::vector<thermo_t*> & tpVector) const
119 {
120  GasKinetics* gK = new GasKinetics(*this);
121  gK->assignShallowPointers(tpVector);
122  return gK;
123 }
124 
126 {
127  doublereal T = thermo().temperature();
128  doublereal P = thermo().pressure();
129  m_logStandConc = log(thermo().standardConcentration());
130  doublereal logT = log(T);
131 
132  if (T != m_temp) {
133  if (!m_rfn.empty()) {
134  m_rates.update(T, logT, &m_rfn[0]);
135  }
136 
137  if (!m_rfn_low.empty()) {
138  m_falloff_low_rates.update(T, logT, &m_rfn_low[0]);
139  m_falloff_high_rates.update(T, logT, &m_rfn_high[0]);
140  }
141  if (!falloff_work.empty()) {
142  m_falloffn.updateTemp(T, &falloff_work[0]);
143  }
144  updateKc();
145  m_ROP_ok = false;
146  }
147 
148  if (T != m_temp || P != m_pres) {
149  if (m_plog_rates.nReactions()) {
150  m_plog_rates.update(T, logT, &m_rfn[0]);
151  m_ROP_ok = false;
152  }
153 
154  if (m_cheb_rates.nReactions()) {
155  m_cheb_rates.update(T, logT, &m_rfn[0]);
156  m_ROP_ok = false;
157  }
158  }
159  m_pres = P;
160  m_temp = T;
161 }
162 
164 {
165  thermo().getActivityConcentrations(&m_conc[0]);
166  doublereal ctot = thermo().molarDensity();
167 
168  // 3-body reactions
169  if (!concm_3b_values.empty()) {
170  m_3b_concm.update(m_conc, ctot, &concm_3b_values[0]);
171  }
172 
173  // Falloff reactions
174  if (!concm_falloff_values.empty()) {
175  m_falloff_concm.update(m_conc, ctot, &concm_falloff_values[0]);
176  }
177 
178  // P-log reactions
179  if (m_plog_rates.nReactions()) {
180  double logP = log(thermo().pressure());
181  m_plog_rates.update_C(&logP);
182  }
183 
184  // Chebyshev reactions
185  if (m_cheb_rates.nReactions()) {
186  double log10P = log10(thermo().pressure());
187  m_cheb_rates.update_C(&log10P);
188  }
189 
190  m_ROP_ok = false;
191 }
192 
194 {
195  thermo().getStandardChemPotentials(&m_grt[0]);
196  fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
197 
198  // compute Delta G^0 for all reversible reactions
199  m_rxnstoich.getRevReactionDelta(m_ii, &m_grt[0], &m_rkcn[0]);
200 
201  doublereal rrt = 1.0/(GasConstant * thermo().temperature());
202  for (size_t i = 0; i < m_nrev; i++) {
203  size_t irxn = m_revindex[i];
204  m_rkcn[irxn] = std::min(exp(m_rkcn[irxn]*rrt - m_dn[irxn]*m_logStandConc),
205  BigNumber);
206  }
207 
208  for (size_t i = 0; i != m_nirrev; ++i) {
209  m_rkcn[ m_irrev[i] ] = 0.0;
210  }
211 }
212 
214 {
215  update_rates_T();
216  thermo().getStandardChemPotentials(&m_grt[0]);
217  fill(m_rkcn.begin(), m_rkcn.end(), 0.0);
218 
219  // compute Delta G^0 for all reactions
220  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], &m_rkcn[0]);
221 
222  doublereal rrt = 1.0/(GasConstant * thermo().temperature());
223  for (size_t i = 0; i < m_ii; i++) {
224  kc[i] = exp(-m_rkcn[i]*rrt + m_dn[i]*m_logStandConc);
225  }
226 
227  // force an update of T-dependent properties, so that m_rkcn will
228  // be updated before it is used next.
229  m_temp = 0.0;
230 }
231 
232 void GasKinetics::getDeltaGibbs(doublereal* deltaG)
233 {
234  /*
235  * Get the chemical potentials of the species in the
236  * ideal gas solution.
237  */
238  thermo().getChemPotentials(&m_grt[0]);
239  /*
240  * Use the stoichiometric manager to find deltaG for each
241  * reaction.
242  */
243  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaG);
244 }
245 
246 void GasKinetics::getDeltaEnthalpy(doublereal* deltaH)
247 {
248  /*
249  * Get the partial molar enthalpy of all species in the
250  * ideal gas.
251  */
252  thermo().getPartialMolarEnthalpies(&m_grt[0]);
253  /*
254  * Use the stoichiometric manager to find deltaG for each
255  * reaction.
256  */
257  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaH);
258 }
259 
260 void GasKinetics::getDeltaEntropy(doublereal* deltaS)
261 {
262  /*
263  * Get the partial molar entropy of all species in the
264  * solid solution.
265  */
266  thermo().getPartialMolarEntropies(&m_grt[0]);
267  /*
268  * Use the stoichiometric manager to find deltaS for each
269  * reaction.
270  */
271  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaS);
272 }
273 
274 void GasKinetics::getDeltaSSGibbs(doublereal* deltaG)
275 {
276  /*
277  * Get the standard state chemical potentials of the species.
278  * This is the array of chemical potentials at unit activity
279  * We define these here as the chemical potentials of the pure
280  * species at the temperature and pressure of the solution.
281  */
282  thermo().getStandardChemPotentials(&m_grt[0]);
283  /*
284  * Use the stoichiometric manager to find deltaG for each
285  * reaction.
286  */
287  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaG);
288 }
289 
290 void GasKinetics::getDeltaSSEnthalpy(doublereal* deltaH)
291 {
292  /*
293  * Get the standard state enthalpies of the species.
294  * This is the array of chemical potentials at unit activity
295  * We define these here as the enthalpies of the pure
296  * species at the temperature and pressure of the solution.
297  */
298  thermo().getEnthalpy_RT(&m_grt[0]);
299  doublereal RT = thermo().temperature() * GasConstant;
300  for (size_t k = 0; k < m_kk; k++) {
301  m_grt[k] *= RT;
302  }
303  /*
304  * Use the stoichiometric manager to find deltaG for each
305  * reaction.
306  */
307  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaH);
308 }
309 
310 void GasKinetics::getDeltaSSEntropy(doublereal* deltaS)
311 {
312  /*
313  * Get the standard state entropy of the species.
314  * We define these here as the entropies of the pure
315  * species at the temperature and pressure of the solution.
316  */
317  thermo().getEntropy_R(&m_grt[0]);
318  doublereal R = GasConstant;
319  for (size_t k = 0; k < m_kk; k++) {
320  m_grt[k] *= R;
321  }
322  /*
323  * Use the stoichiometric manager to find deltaS for each
324  * reaction.
325  */
326  m_rxnstoich.getReactionDelta(m_ii, &m_grt[0], deltaS);
327 }
328 
330 {
331  updateROP();
332  m_rxnstoich.getNetProductionRates(m_kk, &m_ropnet[0], net);
333 }
334 
335 void GasKinetics::getCreationRates(doublereal* cdot)
336 {
337  updateROP();
338  m_rxnstoich.getCreationRates(m_kk, &m_ropf[0], &m_ropr[0], cdot);
339 }
340 
341 void GasKinetics::getDestructionRates(doublereal* ddot)
342 {
343  updateROP();
344  m_rxnstoich.getDestructionRates(m_kk, &m_ropf[0], &m_ropr[0], ddot);
345 }
346 
347 void GasKinetics::processFalloffReactions()
348 {
349  // use m_ropr for temporary storage of reduced pressure
350  vector_fp& pr = m_ropr;
351 
352  for (size_t i = 0; i < m_nfall; i++) {
353  pr[i] = concm_falloff_values[i] * m_rfn_low[i] / (m_rfn_high[i] + SmallNumber);
354  AssertFinite(pr[i], "GasKinetics::processFalloffReactions",
355  "pr[" + int2str(i) + "] is not finite.");
356  }
357 
358  double* work = (falloff_work.empty()) ? 0 : &falloff_work[0];
359  m_falloffn.pr_to_falloff(&pr[0], work);
360 
361  for (size_t i = 0; i < m_nfall; i++) {
362  if (m_rxntype[m_fallindx[i]] == FALLOFF_RXN) {
363  pr[i] *= m_rfn_high[i];
364  } else { // CHEMACT_RXN
365  pr[i] *= m_rfn_low[i];
366  }
367  }
368 
369  scatter_copy(pr.begin(), pr.begin() + m_nfall,
370  m_ropf.begin(), m_fallindx.begin());
371 }
372 
373 void GasKinetics::updateROP()
374 {
375  update_rates_C();
376  update_rates_T();
377 
378  if (m_ROP_ok) {
379  return;
380  }
381 
382  // copy rate coefficients into ropf
383  copy(m_rfn.begin(), m_rfn.end(), m_ropf.begin());
384 
385  // multiply ropf by enhanced 3b conc for all 3b rxns
386  if (!concm_3b_values.empty()) {
387  m_3b_concm.multiply(&m_ropf[0], &concm_3b_values[0]);
388  }
389 
390  if (m_nfall) {
391  processFalloffReactions();
392  }
393 
394  // multiply by perturbation factor
395  multiply_each(m_ropf.begin(), m_ropf.end(), m_perturb.begin());
396 
397  // copy the forward rates to the reverse rates
398  copy(m_ropf.begin(), m_ropf.end(), m_ropr.begin());
399 
400  // for reverse rates computed from thermochemistry, multiply
401  // the forward rates copied into m_ropr by the reciprocals of
402  // the equilibrium constants
403  multiply_each(m_ropr.begin(), m_ropr.end(), m_rkcn.begin());
404 
405  // multiply ropf by concentration products
406  m_rxnstoich.multiplyReactants(&m_conc[0], &m_ropf[0]);
407  //m_reactantStoich.multiply(m_conc.begin(), ropf.begin());
408 
409  // for reversible reactions, multiply ropr by concentration
410  // products
411  m_rxnstoich.multiplyRevProducts(&m_conc[0], &m_ropr[0]);
412  //m_revProductStoich.multiply(m_conc.begin(), ropr.begin());
413 
414  for (size_t j = 0; j != m_ii; ++j) {
415  m_ropnet[j] = m_ropf[j] - m_ropr[j];
416  }
417 
418  for (size_t i = 0; i < m_rfn.size(); i++) {
419  AssertFinite(m_rfn[i], "GasKinetics::updateROP",
420  "m_rfn[" + int2str(i) + "] is not finite.");
421  AssertFinite(m_ropf[i], "GasKinetics::updateROP",
422  "m_ropf[" + int2str(i) + "] is not finite.");
423  AssertFinite(m_ropr[i], "GasKinetics::updateROP",
424  "m_ropr[" + int2str(i) + "] is not finite.");
425  }
426 
427  m_ROP_ok = true;
428 }
429 
430 void GasKinetics::
431 getFwdRateConstants(doublereal* kfwd)
432 {
433  update_rates_C();
434  update_rates_T();
435 
436  // copy rate coefficients into ropf
437  copy(m_rfn.begin(), m_rfn.end(), m_ropf.begin());
438 
439  // multiply ropf by enhanced 3b conc for all 3b rxns
440  if (!concm_3b_values.empty()) {
441  m_3b_concm.multiply(&m_ropf[0], &concm_3b_values[0]);
442  }
443 
444  /*
445  * This routine is hardcoded to replace some of the values
446  * of the ropf vector.
447  */
448  if (m_nfall) {
449  processFalloffReactions();
450  }
451 
452  // multiply by perturbation factor
453  multiply_each(m_ropf.begin(), m_ropf.end(), m_perturb.begin());
454 
455  for (size_t i = 0; i < m_ii; i++) {
456  kfwd[i] = m_ropf[i];
457  }
458 }
459 
460 void GasKinetics::
461 getRevRateConstants(doublereal* krev, bool doIrreversible)
462 {
463  /*
464  * go get the forward rate constants. -> note, we don't
465  * really care about speed or redundancy in these
466  * informational routines.
467  */
468  getFwdRateConstants(krev);
469 
470  if (doIrreversible) {
471  getEquilibriumConstants(&m_ropnet[0]);
472  for (size_t i = 0; i < m_ii; i++) {
473  krev[i] /= m_ropnet[i];
474  }
475  } else {
476  // m_rkcn[] is zero for irreversible reactions
477  for (size_t i = 0; i < m_ii; i++) {
478  krev[i] *= m_rkcn[i];
479  }
480  }
481 }
482 
483 void GasKinetics::
485 {
486  switch (r.reactionType) {
487  case ELEMENTARY_RXN:
488  addElementaryReaction(r);
489  break;
490  case THREE_BODY_RXN:
491  addThreeBodyReaction(r);
492  break;
493  case FALLOFF_RXN:
494  case CHEMACT_RXN:
495  addFalloffReaction(r);
496  break;
497  case PLOG_RXN:
498  addPlogReaction(r);
499  break;
500  case CHEBYSHEV_RXN:
501  addChebyshevReaction(r);
502  break;
503  default:
504  throw CanteraError("GasKinetics::addReaction", "Invalid reaction type specified");
505  }
506 
507  // operations common to all reaction types
508  installReagents(r);
509  installGroups(reactionNumber(), r.rgroups, r.pgroups);
511  m_rxneqn.push_back(r.equation);
512  m_rxntype.push_back(r.reactionType);
513 }
514 
515 void GasKinetics::
516 addFalloffReaction(ReactionData& r)
517 {
518  // install high and low rate coeff calculators
519  // and add constant terms to high and low rate coeff value vectors
520  size_t iloc = m_falloff_high_rates.install(m_nfall, r);
521  m_rfn_high.push_back(r.rateCoeffParameters[0]);
523  m_falloff_low_rates.install(m_nfall, r);
524  m_rfn_low.push_back(r.rateCoeffParameters[0]);
525 
526  // add a dummy entry in m_rf, where computed falloff
527  // rate coeff will be put
528  m_rfn.push_back(0.0);
529 
530  // add this reaction number to the list of
531  // falloff reactions
532  m_fallindx.push_back(reactionNumber());
533 
534  // install the enhanced third-body concentration
535  // calculator for this reaction
536  m_falloff_concm.install(m_nfall, r.thirdBodyEfficiencies,
537  r.default_3b_eff);
538 
539  // install the falloff function calculator for
540  // this reaction
541  m_falloffn.install(m_nfall, r.falloffType, r.reactionType,
543 
544  // forward rxn order equals number of reactants, since rate
545  // coeff is defined in terms of the high-pressure limit
546  m_fwdOrder.push_back(r.reactants.size());
547 
548  // increment the falloff reaction counter
549  ++m_nfall;
550  registerReaction(reactionNumber(), r.reactionType, iloc);
551 }
552 
553 void GasKinetics::
554 addElementaryReaction(ReactionData& r)
555 {
556  // install rate coeff calculator
557  size_t iloc = m_rates.install(reactionNumber(), r);
558 
559  // add constant term to rate coeff value vector
560  m_rfn.push_back(r.rateCoeffParameters[0]);
561 
562  // forward rxn order equals number of reactants
563  m_fwdOrder.push_back(r.reactants.size());
564  registerReaction(reactionNumber(), ELEMENTARY_RXN, iloc);
565 }
566 
567 void GasKinetics::
568 addThreeBodyReaction(ReactionData& r)
569 {
570  // install rate coeff calculator
571  size_t iloc = m_rates.install(reactionNumber(), r);
572 
573  // add constant term to rate coeff value vector
574  m_rfn.push_back(r.rateCoeffParameters[0]);
575 
576  // forward rxn order equals number of reactants + 1
577  m_fwdOrder.push_back(r.reactants.size() + 1);
578 
579  m_3b_concm.install(reactionNumber(), r.thirdBodyEfficiencies,
580  r.default_3b_eff);
581  registerReaction(reactionNumber(), THREE_BODY_RXN, iloc);
582 }
583 
584 void GasKinetics::addPlogReaction(ReactionData& r)
585 {
586  // install rate coefficient calculator
587  size_t iloc = m_plog_rates.install(reactionNumber(), r);
588 
589  // add a dummy entry in m_rfn, where computed rate coeff will be put
590  m_rfn.push_back(0.0);
591 
592  m_fwdOrder.push_back(r.reactants.size());
593  registerReaction(reactionNumber(), PLOG_RXN, iloc);
594 }
595 
596 void GasKinetics::addChebyshevReaction(ReactionData& r)
597 {
598  // install rate coefficient calculator
599  size_t iloc = m_cheb_rates.install(reactionNumber(), r);
600 
601  // add a dummy entry in m_rfn, where computed rate coeff will be put
602  m_rfn.push_back(0.0);
603 
604  m_fwdOrder.push_back(r.reactants.size());
605  registerReaction(reactionNumber(), CHEBYSHEV_RXN, iloc);
606 }
607 
608 void GasKinetics::installReagents(const ReactionData& r)
609 {
610  m_ropf.push_back(0.0); // extend by one for new rxn
611  m_ropr.push_back(0.0);
612  m_ropnet.push_back(0.0);
613  size_t n, ns, m;
614  doublereal nsFlt;
615  doublereal reactantGlobalOrder = 0.0;
616  doublereal productGlobalOrder = 0.0;
617  size_t rnum = reactionNumber();
618 
619  std::vector<size_t> rk;
620  size_t nr = r.reactants.size();
621  for (n = 0; n < nr; n++) {
622  nsFlt = r.rstoich[n];
623  reactantGlobalOrder += nsFlt;
624  ns = (size_t) nsFlt;
625  if ((doublereal) ns != nsFlt) {
626  if (ns < 1) {
627  ns = 1;
628  }
629  }
630  if (r.rstoich[n] != 0.0) {
631  m_rrxn[r.reactants[n]][rnum] += r.rstoich[n];
632  }
633  for (m = 0; m < ns; m++) {
634  rk.push_back(r.reactants[n]);
635  }
636  }
637  m_reactants.push_back(rk);
638 
639  std::vector<size_t> pk;
640  size_t np = r.products.size();
641  for (n = 0; n < np; n++) {
642  nsFlt = r.pstoich[n];
643  productGlobalOrder += nsFlt;
644  ns = (size_t) nsFlt;
645  if ((double) ns != nsFlt) {
646  if (ns < 1) {
647  ns = 1;
648  }
649  }
650  if (r.pstoich[n] != 0.0) {
651  m_prxn[r.products[n]][rnum] += r.pstoich[n];
652  }
653  for (m = 0; m < ns; m++) {
654  pk.push_back(r.products[n]);
655  }
656  }
657  m_products.push_back(pk);
658  m_rkcn.push_back(0.0);
659  m_rxnstoich.add(reactionNumber(), r);
660 
661  if (r.reversible) {
662  m_dn.push_back(productGlobalOrder - reactantGlobalOrder);
663  m_revindex.push_back(reactionNumber());
664  m_nrev++;
665  } else {
666  m_dn.push_back(productGlobalOrder - reactantGlobalOrder);
667  m_irrev.push_back(reactionNumber());
668  m_nirrev++;
669  }
670 }
671 
672 void GasKinetics::installGroups(size_t irxn,
673  const vector<grouplist_t>& r, const vector<grouplist_t>& p)
674 {
675  if (!r.empty()) {
676  writelog("installing groups for reaction "+int2str(reactionNumber()));
677  m_rgroups[reactionNumber()] = r;
678  m_pgroups[reactionNumber()] = p;
679  }
680 }
681 
683 {
684  m_kk = thermo().nSpecies();
685  m_rrxn.resize(m_kk);
686  m_prxn.resize(m_kk);
687  m_conc.resize(m_kk);
688  m_grt.resize(m_kk);
689  m_logp_ref = log(thermo().refPressure()) - log(GasConstant);
690 }
691 
693 {
694  if (!m_finalized) {
695  falloff_work.resize(m_falloffn.workSize());
696  concm_3b_values.resize(m_3b_concm.workSize());
697  concm_falloff_values.resize(m_falloff_concm.workSize());
698  m_finalized = true;
699 
700  // Guarantee that these arrays can be converted to double* even in the
701  // special case where there are no reactions defined.
702  if (!m_ii) {
703  m_perturb.resize(1, 1.0);
704  m_ropf.resize(1, 0.0);
705  m_ropr.resize(1, 0.0);
706  m_ropnet.resize(1, 0.0);
707  m_rkcn.resize(1, 0.0);
708  }
709  }
710 }
711 
712 bool GasKinetics::ready() const
713 {
714  return m_finalized;
715 }
716 
717 }
std::vector< grouplist_t > rgroups
Optional data used in reaction path diagrams.
Definition: ReactionData.h:64
virtual void getNetProductionRates(size_t nsp, const doublereal *ropnet, doublereal *w)
Species net production rates.
const int PLOG_RXN
A pressure-dependent rate expression consisting of several Arrhenius rate expressions evaluated at di...
Definition: reaction_defs.h:46
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:82
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
void incrementRxnCount()
Increment the number of reactions in the mechanism by one.
Definition: Kinetics.h:859
std::vector< std::vector< size_t > > m_reactants
This is a vector of vectors containing the reactants for each reaction.
Definition: Kinetics.h:908
vector_fp falloffParameters
Values used in the falloff parameterization.
Definition: ReactionData.h:104
virtual void getDestructionRates(size_t nSpecies, const doublereal *fwdRatesOfProgress, const doublereal *revRatesOfProgress, doublereal *destructionRates)
Species destruction rates.
void install(size_t rxn, int falloffType, int reactionType, const vector_fp &c)
Install a new falloff function calculator.
Definition: FalloffMgr.h:51
int reactionType
Type of the reaction.
Definition: ReactionData.h:39
virtual void assignShallowPointers(const std::vector< thermo_t * > &tpVector)
Reassign the pointers within the Kinetics object.
Definition: Kinetics.cpp:144
virtual void getEnthalpy_RT(doublereal *hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
Definition: ThermoPhase.h:721
vector_fp m_dn
Difference between the input global reactants order and the input global products order...
Definition: GasKinetics.h:199
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:288
vector_fp auxRateCoeffParameters
Vector of auxiliary rate coefficient parameters.
Definition: ReactionData.h:96
virtual void getCreationRates(size_t nSpecies, const doublereal *fwdRatesOfProgress, const doublereal *revRatesOfProgress, doublereal *creationRates)
Species creation rates.
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition: Kinetics.h:891
virtual void getRevReactionDelta(size_t nr, const doublereal *g, doublereal *dg)
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
virtual void getFwdRateConstants(doublereal *kfwd)
Return the forward rate constants.
int falloffType
Type of falloff parameterization to use.
Definition: ReactionData.h:100
const int CHEBYSHEV_RXN
A general pressure-dependent reaction where k(T,P) is defined in terms of a bivariate Chebyshev polyn...
Definition: reaction_defs.h:52
virtual void getDeltaSSEntropy(doublereal *deltaS)
Return the vector of values for the change in the standard state entropies for each reaction...
virtual void getReactionDelta(size_t nReactions, const doublereal *g, doublereal *dg)
Calculates the change of a molar species property in a reaction.
virtual void init()
Prepare the class for the addition of reactions.
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:711
doublereal molarDensity() const
Molar density (kmol/m^3).
Definition: Phase.cpp:597
const int CHEMACT_RXN
A chemical activation reaction.
Definition: reaction_defs.h:60
virtual void getNetProductionRates(doublereal *net)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
void updateTemp(doublereal t, doublereal *work)
Update the cached temperature-dependent intermediate results for all installed falloff functions...
Definition: FalloffMgr.h:72
size_t workSize()
Size of the work array required to store intermediate results.
Definition: FalloffMgr.h:62
Kinetics manager for elementary gas-phase chemistry.
Definition: GasKinetics.h:35
#define AssertFinite(expr, procedure, message)
Throw an exception if the specified exception is not a finite number.
Definition: ctexceptions.h:254
void multiply_each(OutputIter x_begin, OutputIter x_end, InputIter y_begin)
Multiply each entry in x by the corresponding entry in y.
Definition: utilities.h:227
vector_fp rateCoeffParameters
Vector of rate coefficient parameters.
Definition: ReactionData.h:91
std::map< size_t, doublereal > thirdBodyEfficiencies
Map of species index to third body efficiency.
Definition: ReactionData.h:68
virtual void update_rates_T()
Update temperature-dependent portions of reaction rates and falloff functions.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
const int FALLOFF_RXN
The general form for an association or dissociation reaction, with a pressure-dependent rate...
Definition: reaction_defs.h:38
virtual void getEquilibriumConstants(doublereal *kc)
Return a vector of Equilibrium constants.
std::string equation
The reaction equation. Used only for display purposes.
Definition: ReactionData.h:109
virtual void multiplyReactants(const doublereal *C, doublereal *R)
Given an array of concentrations C, multiply the entries in array R by the concentration products for...
std::vector< grouplist_t > pgroups
Optional data used in reaction path diagrams.
Definition: ReactionData.h:65
virtual bool ready() const
Returns true if the kinetics manager has been properly initialized and finalized. ...
virtual void getDestructionRates(doublereal *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
virtual void getDeltaEntropy(doublereal *deltaS)
Return the vector of values for the reactions change in entropy.
void update_C(const doublereal *c)
Update the concentration-dependent parts of the rate coefficient, if any.
Definition: RateCoeffMgr.h:63
virtual void getDeltaSSEnthalpy(doublereal *deltaH)
Return the vector of values for the change in the standard state enthalpies of reaction.
virtual void getCreationRates(doublereal *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
const doublereal BigNumber
largest number to compare to inf.
Definition: ct_defs.h:141
virtual void getChemPotentials(doublereal *mu) const
Get the species chemical potentials. Units: J/kmol.
Definition: ThermoPhase.h:595
void updateKc()
Update the equilibrium constants in molar units.
vector_fp m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition: Kinetics.h:895
Public interface for kinetics managers.
Definition: Kinetics.h:131
doublereal default_3b_eff
The default third body efficiency for species not listed in thirdBodyEfficiencies.
Definition: ReactionData.h:113
Intermediate class which stores data about a reaction and its rate parameterization before adding the...
Definition: ReactionData.h:16
virtual void getPartialMolarEntropies(doublereal *sbar) const
Returns an array of partial molar entropies of the species in the solution.
Definition: ThermoPhase.h:638
GasKinetics(thermo_t *thermo=0)
Constructor.
Definition: GasKinetics.cpp:21
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
size_t install(size_t rxnNumber, const ReactionData &rdata)
Install a rate coefficient calculator.
Definition: RateCoeffMgr.h:38
virtual Kinetics * duplMyselfAsKinetics(const std::vector< thermo_t * > &tpVector) const
Duplication routine for objects which inherit from Kinetics.
void update(doublereal T, doublereal logT, doublereal *values)
Write the rate coefficients into array values.
Definition: RateCoeffMgr.h:77
Kinetics & operator=(const Kinetics &right)
Assignment operator.
Definition: Kinetics.cpp:62
const int THREE_BODY_RXN
A reaction that requires a third-body collision partner.
Definition: reaction_defs.h:32
std::vector< size_t > reactants
Indices of reactant species.
Definition: ReactionData.h:44
virtual void getRevRateConstants(doublereal *krev, bool doIrreversible=false)
Return the reverse rate constants.
virtual void update_rates_C()
Update properties that depend on concentrations.
GasKinetics & operator=(const GasKinetics &right)
Assignment operator.
Definition: GasKinetics.cpp:57
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:252
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
Definition: ThermoPhase.h:314
virtual void getEntropy_R(doublereal *sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
Definition: ThermoPhase.h:731
doublereal temperature() const
Temperature (K).
Definition: Phase.h:528
virtual void getActivityConcentrations(doublereal *c) const
This method returns an array of generalized concentrations.
Definition: ThermoPhase.h:465
virtual void getDeltaSSGibbs(doublereal *deltaG)
Return the vector of values for the reaction standard state gibbs free energy change.
virtual void getDeltaGibbs(doublereal *deltaG)
Return the vector of values for the reaction gibbs free energy change.
const doublereal SmallNumber
smallest number to compare to zero.
Definition: ct_defs.h:139
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:165
std::vector< std::vector< size_t > > m_products
This is a vector of vectors containing the products for each reaction.
Definition: Kinetics.h:921
virtual void addReaction(ReactionData &r)
Add a single reaction to the mechanism.
doublereal m_pres
Last pressure at which rates were evaluated.
Definition: GasKinetics.h:217
void scatter_copy(InputIter begin, InputIter end, OutputIter result, IndexIter index)
Copies a contiguous range in a sequence to indexed positions in another sequence. ...
Definition: utilities.h:430
const doublereal GasConstant
Universal Gas Constant. [J/kmol/K].
Definition: ct_defs.h:66
size_t m_ii
Number of reactions in the mechanism.
Definition: Kinetics.h:887
const int ELEMENTARY_RXN
A reaction with a rate coefficient that depends only on temperature.
Definition: reaction_defs.h:26
virtual void multiplyRevProducts(const doublereal *c, doublereal *r)
Given an array of concentrations C, multiply the entries in array R by the concentration products for...
virtual void getDeltaEnthalpy(doublereal *deltaH)
Return the vector of values for the reactions change in enthalpy.
virtual void getPartialMolarEnthalpies(doublereal *hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
Definition: ThermoPhase.h:628
void writelog(const std::string &msg)
Write a message to the screen.
Definition: global.cpp:43
virtual void addPhase(thermo_t &thermo)
Add a phase to the kinetics manager object.
Definition: Kinetics.cpp:255
virtual void finalize()
Finish adding reactions and prepare for use.
virtual void add(size_t rxn, const std::vector< size_t > &reactants, const std::vector< size_t > &products, bool reversible)
Add a reaction with mass-action kinetics.