Cantera  2.5.1
Reaction.cpp
Go to the documentation of this file.
1 /**
2  * @file Reaction.cpp
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 
11 #include "cantera/base/ctml.h"
12 #include "cantera/base/Array.h"
13 #include "cantera/base/AnyMap.h"
14 #include <sstream>
15 #include <set>
16 
17 #include <boost/algorithm/string.hpp>
18 
19 namespace ba = boost::algorithm;
20 
21 namespace Cantera
22 {
23 
24 Reaction::Reaction(int type)
25  : reaction_type(type)
26  , reversible(true)
27  , duplicate(false)
28  , allow_nonreactant_orders(false)
29  , allow_negative_orders(false)
30 {
31 }
32 
33 Reaction::Reaction(int type, const Composition& reactants_,
34  const Composition& products_)
35  : reaction_type(type)
36  , reactants(reactants_)
37  , products(products_)
38  , reversible(true)
39  , duplicate(false)
40  , allow_nonreactant_orders(false)
41  , allow_negative_orders(false)
42 {
43 }
44 
45 void Reaction::validate()
46 {
47  if (!allow_nonreactant_orders) {
48  for (const auto& order : orders) {
49  if (reactants.find(order.first) == reactants.end()) {
50  throw CanteraError("Reaction::validate", "Reaction order "
51  "specified for non-reactant species '" + order.first + "'");
52  }
53  }
54  }
55 
56  if (!allow_negative_orders) {
57  for (const auto& order : orders) {
58  if (order.second < 0.0) {
59  throw CanteraError("Reaction::validate", "Negative reaction "
60  "order specified for species '" + order.first + "'");
61  }
62  }
63  }
64 }
65 
66 std::string Reaction::reactantString() const
67 {
68  std::ostringstream result;
69  for (auto iter = reactants.begin(); iter != reactants.end(); ++iter) {
70  if (iter != reactants.begin()) {
71  result << " + ";
72  }
73  if (iter->second != 1.0) {
74  result << iter->second << " ";
75  }
76  result << iter->first;
77  }
78  return result.str();
79 }
80 
81 std::string Reaction::productString() const
82 {
83  std::ostringstream result;
84  for (auto iter = products.begin(); iter != products.end(); ++iter) {
85  if (iter != products.begin()) {
86  result << " + ";
87  }
88  if (iter->second != 1.0) {
89  result << iter->second << " ";
90  }
91  result << iter->first;
92  }
93  return result.str();
94 }
95 
96 std::string Reaction::equation() const
97 {
98  if (reversible) {
99  return reactantString() + " <=> " + productString();
100  } else {
101  return reactantString() + " => " + productString();
102  }
103 }
104 
105 ElementaryReaction::ElementaryReaction(const Composition& reactants_,
106  const Composition products_,
107  const Arrhenius& rate_)
108  : Reaction(ELEMENTARY_RXN, reactants_, products_)
109  , rate(rate_)
110  , allow_negative_pre_exponential_factor(false)
111 {
112 }
113 
114 ElementaryReaction::ElementaryReaction()
115  : Reaction(ELEMENTARY_RXN)
116  , allow_negative_pre_exponential_factor(false)
117 {
118 }
119 
121 {
123  if (!allow_negative_pre_exponential_factor &&
124  rate.preExponentialFactor() < 0) {
125  throw CanteraError("ElementaryReaction::validate",
126  "Undeclared negative pre-exponential factor found in reaction '"
127  + equation() + "'");
128  }
129 }
130 
131 ThirdBody::ThirdBody(double default_eff)
132  : default_efficiency(default_eff)
133 {
134 }
135 
136 ThreeBodyReaction::ThreeBodyReaction()
137 {
139 }
140 
141 ThreeBodyReaction::ThreeBodyReaction(const Composition& reactants_,
142  const Composition& products_,
143  const Arrhenius& rate_,
144  const ThirdBody& tbody)
145  : ElementaryReaction(reactants_, products_, rate_)
146  , third_body(tbody)
147 {
149 }
150 
152  return ElementaryReaction::reactantString() + " + M";
153 }
154 
155 std::string ThreeBodyReaction::productString() const {
156  return ElementaryReaction::productString() + " + M";
157 }
158 
159 FalloffReaction::FalloffReaction()
161  , falloff(new Falloff())
162  , allow_negative_pre_exponential_factor(false)
163 {
164 }
165 
166 FalloffReaction::FalloffReaction(
167  const Composition& reactants_, const Composition& products_,
168  const Arrhenius& low_rate_, const Arrhenius& high_rate_,
169  const ThirdBody& tbody)
170  : Reaction(FALLOFF_RXN, reactants_, products_)
171  , low_rate(low_rate_)
172  , high_rate(high_rate_)
173  , third_body(tbody)
174  , falloff(new Falloff())
175 {
176 }
177 
178 std::string FalloffReaction::reactantString() const {
179  if (third_body.default_efficiency == 0 &&
180  third_body.efficiencies.size() == 1) {
181  return Reaction::reactantString() + " (+" +
182  third_body.efficiencies.begin()->first + ")";
183  } else {
184  return Reaction::reactantString() + " (+M)";
185  }
186 }
187 
188 std::string FalloffReaction::productString() const {
189  if (third_body.default_efficiency == 0 &&
190  third_body.efficiencies.size() == 1) {
191  return Reaction::productString() + " (+" +
192  third_body.efficiencies.begin()->first + ")";
193  } else {
194  return Reaction::productString() + " (+M)";
195  }
196 }
197 
200  if (!allow_negative_pre_exponential_factor &&
203  throw CanteraError("FalloffReaction::validate", "Negative "
204  "pre-exponential factor found for reaction '" + equation() + "'");
205  }
207  throw CanteraError("FalloffReaction::validate", "High and "
208  "low rate pre-exponential factors must have the same sign."
209  "Reaction: '{}'", equation());
210  }
211 }
212 
213 ChemicallyActivatedReaction::ChemicallyActivatedReaction()
214 {
216 }
217 
218 ChemicallyActivatedReaction::ChemicallyActivatedReaction(
219  const Composition& reactants_, const Composition& products_,
220  const Arrhenius& low_rate_, const Arrhenius& high_rate_,
221  const ThirdBody& tbody)
222  : FalloffReaction(reactants_, products_, low_rate_, high_rate_, tbody)
223 {
225 }
226 
227 PlogReaction::PlogReaction()
228  : Reaction(PLOG_RXN)
229 {
230 }
231 
232 PlogReaction::PlogReaction(const Composition& reactants_,
233  const Composition& products_, const Plog& rate_)
234  : Reaction(PLOG_RXN, reactants_, products_)
235  , rate(rate_)
236 {
237 }
238 
239 ChebyshevReaction::ChebyshevReaction()
240  : Reaction(CHEBYSHEV_RXN)
241 {
242 }
243 
244 ChebyshevReaction::ChebyshevReaction(const Composition& reactants_,
245  const Composition& products_,
246  const ChebyshevRate& rate_)
247  : Reaction(CHEBYSHEV_RXN, reactants_, products_)
248  , rate(rate_)
249 {
250 }
251 
252 InterfaceReaction::InterfaceReaction()
253  : is_sticking_coefficient(false)
254  , use_motz_wise_correction(false)
255 {
257 }
258 
259 InterfaceReaction::InterfaceReaction(const Composition& reactants_,
260  const Composition& products_,
261  const Arrhenius& rate_,
262  bool isStick)
263  : ElementaryReaction(reactants_, products_, rate_)
264  , is_sticking_coefficient(isStick)
265  , use_motz_wise_correction(false)
266 {
268 }
269 
270 ElectrochemicalReaction::ElectrochemicalReaction()
271  : film_resistivity(0.0)
272  , beta(0.5)
273  , exchange_current_density_formulation(false)
274 {
275 }
276 
277 ElectrochemicalReaction::ElectrochemicalReaction(const Composition& reactants_,
278  const Composition& products_,
279  const Arrhenius& rate_)
280  : InterfaceReaction(reactants_, products_, rate_)
281  , film_resistivity(0.0)
282  , beta(0.5)
283  , exchange_current_density_formulation(false)
284 {
285 }
286 
287 Arrhenius readArrhenius(const XML_Node& arrhenius_node)
288 {
289  return Arrhenius(getFloat(arrhenius_node, "A", "toSI"),
290  getFloat(arrhenius_node, "b"),
291  getFloat(arrhenius_node, "E", "actEnergy") / GasConstant);
292 }
293 
294 Units rateCoeffUnits(const Reaction& R, const Kinetics& kin,
295  int pressure_dependence=0)
296 {
297  if (R.reaction_type == INVALID_RXN) {
298  // If a reaction is invalid because of missing species in the Kinetics
299  // object, determining the units of the rate coefficient is impossible.
300  return Units();
301  } else if (R.reaction_type == INTERFACE_RXN
302  && dynamic_cast<const InterfaceReaction&>(R).is_sticking_coefficient) {
303  // Sticking coefficients are dimensionless
304  return Units();
305  }
306 
307  // Determine the units of the rate coefficient
308  Units rxn_phase_units = kin.thermo(kin.reactionPhaseIndex()).standardConcentrationUnits();
309  Units rcUnits = rxn_phase_units;
310  rcUnits *= Units(1.0, 0, 0, -1);
311  for (const auto& order : R.orders) {
312  const auto& phase = kin.speciesPhase(order.first);
313  rcUnits *= phase.standardConcentrationUnits().pow(-order.second);
314  }
315  for (const auto& stoich : R.reactants) {
316  // Order for each reactant is the reactant stoichiometric coefficient,
317  // unless already overridden by user-specified orders
318  if (stoich.first == "M") {
319  rcUnits *= rxn_phase_units.pow(-1);
320  } else if (R.orders.find(stoich.first) == R.orders.end()) {
321  const auto& phase = kin.speciesPhase(stoich.first);
322  rcUnits *= phase.standardConcentrationUnits().pow(-stoich.second);
323  }
324  }
325 
326  // Incorporate pressure dependence for low-pressure falloff and high-
327  // pressure chemically-activated reaction limits
328  rcUnits *= rxn_phase_units.pow(-pressure_dependence);
329  return rcUnits;
330 }
331 
332 Arrhenius readArrhenius(const Reaction& R, const AnyValue& rate,
333  const Kinetics& kin, const UnitSystem& units,
334  int pressure_dependence=0)
335 {
336  double A, b, Ta;
337  Units rc_units = rateCoeffUnits(R, kin, pressure_dependence);
338  if (rate.is<AnyMap>()) {
339  auto& rate_map = rate.as<AnyMap>();
340  A = units.convert(rate_map["A"], rc_units);
341  b = rate_map["b"].asDouble();
342  Ta = units.convertActivationEnergy(rate_map["Ea"], "K");
343  } else {
344  auto& rate_vec = rate.asVector<AnyValue>(3);
345  A = units.convert(rate_vec[0], rc_units);
346  b = rate_vec[1].asDouble();
347  Ta = units.convertActivationEnergy(rate_vec[2], "K");
348  }
349  return Arrhenius(A, b, Ta);
350 }
351 
352 //! Parse falloff parameters, given a rateCoeff node
353 /*!
354  * @verbatim
355  <falloff type="Troe"> 0.5 73.2 5000. 9999. </falloff>
356  @endverbatim
357 */
358 void readFalloff(FalloffReaction& R, const XML_Node& rc_node)
359 {
360  XML_Node& falloff = rc_node.child("falloff");
361  std::vector<std::string> p;
362  vector_fp falloff_parameters;
363  getStringArray(falloff, p);
364  size_t np = p.size();
365  for (size_t n = 0; n < np; n++) {
366  falloff_parameters.push_back(fpValueCheck(p[n]));
367  }
368 
369  if (caseInsensitiveEquals(falloff["type"], "lindemann")) {
370  if (np != 0) {
371  throw CanteraError("readFalloff", "Lindemann parameterization "
372  "takes no parameters, but {} were given", np);
373  }
374  R.falloff = newFalloff("Lindemann", falloff_parameters);
375  } else if (caseInsensitiveEquals(falloff["type"], "troe")) {
376  if (np != 3 && np != 4) {
377  throw CanteraError("readFalloff", "Troe parameterization takes "
378  "3 or 4 parameters, but {} were given", np);
379  }
380  R.falloff = newFalloff("Troe", falloff_parameters);
381  } else if (caseInsensitiveEquals(falloff["type"], "sri")) {
382  if (np != 3 && np != 5) {
383  throw CanteraError("readFalloff", "SRI parameterization takes "
384  "3 or 5 parameters, but {} were given", np);
385  }
386  R.falloff = newFalloff("SRI", falloff_parameters);
387  } else {
388  throw CanteraError("readFalloff", "Unrecognized falloff type: '{}'",
389  falloff["type"]);
390  }
391 }
392 
393 void readFalloff(FalloffReaction& R, const AnyMap& node)
394 {
395  if (node.hasKey("Troe")) {
396  auto& f = node["Troe"].as<AnyMap>();
397  vector_fp params{
398  f["A"].asDouble(),
399  f["T3"].asDouble(),
400  f["T1"].asDouble()
401  };
402  if (f.hasKey("T2")) {
403  params.push_back(f["T2"].asDouble());
404  }
405  R.falloff = newFalloff("Troe", params);
406  } else if (node.hasKey("SRI")) {
407  auto& f = node["SRI"].as<AnyMap>();
408  vector_fp params{
409  f["A"].asDouble(),
410  f["B"].asDouble(),
411  f["C"].asDouble()
412  };
413  if (f.hasKey("D")) {
414  params.push_back(f["D"].asDouble());
415  }
416  if (f.hasKey("E")) {
417  params.push_back(f["E"].asDouble());
418  }
419  R.falloff = newFalloff("SRI", params);
420  } else {
421  R.falloff = newFalloff("Lindemann", {});
422  }
423 }
424 
425 void readEfficiencies(ThirdBody& tbody, const XML_Node& rc_node)
426 {
427  if (!rc_node.hasChild("efficiencies")) {
428  tbody.default_efficiency = 1.0;
429  return;
430  }
431  const XML_Node& eff_node = rc_node.child("efficiencies");
432  tbody.default_efficiency = fpValue(eff_node["default"]);
433  tbody.efficiencies = parseCompString(eff_node.value());
434 }
435 
436 void readEfficiencies(ThirdBody& tbody, const AnyMap& node)
437 {
438  tbody.default_efficiency = node.getDouble("default-efficiency", 1.0);
439  if (node.hasKey("efficiencies")) {
440  tbody.efficiencies = node["efficiencies"].asMap<double>();
441  }
442 }
443 
444 void setupReaction(Reaction& R, const XML_Node& rxn_node)
445 {
446  // Reactant and product stoichiometries
447  R.reactants = parseCompString(rxn_node.child("reactants").value());
448  R.products = parseCompString(rxn_node.child("products").value());
449 
450  // Non-stoichiometric reaction orders
451  std::vector<XML_Node*> orders = rxn_node.getChildren("order");
452  for (size_t i = 0; i < orders.size(); i++) {
453  R.orders[orders[i]->attrib("species")] = orders[i]->fp_value();
454  }
455 
456  // Flags
457  R.id = rxn_node.attrib("id");
458  R.duplicate = rxn_node.hasAttrib("duplicate");
459  const std::string& rev = rxn_node["reversible"];
460  R.reversible = (rev == "true" || rev == "yes");
461 }
462 
463 void parseReactionEquation(Reaction& R, const AnyValue& equation,
464  const Kinetics& kin) {
465  // Parse the reaction equation to determine participating species and
466  // stoichiometric coefficients
467  std::vector<std::string> tokens;
468  tokenizeString(equation.asString(), tokens);
469  tokens.push_back("+"); // makes parsing last species not a special case
470 
471  size_t last_used = npos; // index of last-used token
472  bool reactants = true;
473  for (size_t i = 1; i < tokens.size(); i++) {
474  if (tokens[i] == "+" || ba::starts_with(tokens[i], "(+") ||
475  tokens[i] == "<=>" || tokens[i] == "=" || tokens[i] == "=>") {
476  std::string species = tokens[i-1];
477 
478  double stoich;
479  if (last_used != npos && tokens[last_used] == "(+") {
480  // Falloff third body with space, e.g. "(+ M)"
481  species = "(+" + species;
482  stoich = -1;
483  } else if (last_used == i-1 && ba::starts_with(species, "(+")
484  && ba::ends_with(species, ")")) {
485  // Falloff 3rd body written without space, e.g. "(+M)"
486  stoich = -1;
487  } else if (last_used == i-2) { // Species with no stoich. coefficient
488  stoich = 1.0;
489  } else if (last_used == i-3) { // Stoich. coefficient and species
490  try {
491  stoich = fpValueCheck(tokens[i-2]);
492  } catch (CanteraError& err) {
493  throw InputFileError("parseReactionEquation", equation,
494  err.getMessage());
495  }
496  } else {
497  throw InputFileError("parseReactionEquation", equation,
498  "Error parsing reaction string '{}'.\n"
499  "Current token: '{}'\nlast_used: '{}'",
500  equation.asString(), tokens[i],
501  (last_used == npos) ? "n/a" : tokens[last_used]);
502  }
503  if (kin.kineticsSpeciesIndex(species) == npos
504  && stoich != -1 && species != "M") {
505  R.reaction_type = INVALID_RXN;
506  }
507 
508  if (reactants) {
509  R.reactants[species] += stoich;
510  } else {
511  R.products[species] += stoich;
512  }
513 
514  last_used = i;
515  }
516 
517  // Tokens after this point are part of the products string
518  if (tokens[i] == "<=>" || tokens[i] == "=") {
519  R.reversible = true;
520  reactants = false;
521  } else if (tokens[i] == "=>") {
522  R.reversible = false;
523  reactants = false;
524  }
525  }
526 }
527 
528 void setupReaction(Reaction& R, const AnyMap& node, const Kinetics& kin)
529 {
530  parseReactionEquation(R, node["equation"], kin);
531  // Non-stoichiometric reaction orders
532  std::map<std::string, double> orders;
533  if (node.hasKey("orders")) {
534  for (const auto& order : node["orders"].asMap<double>()) {
535  R.orders[order.first] = order.second;
536  if (kin.kineticsSpeciesIndex(order.first) == npos) {
537  R.reaction_type = INVALID_RXN;
538  }
539  }
540  }
541 
542  //Flags
543  R.id = node.getString("id", "");
544  R.duplicate = node.getBool("duplicate", false);
545  R.allow_negative_orders = node.getBool("negative-orders", false);
546  R.allow_nonreactant_orders = node.getBool("nonreactant-orders", false);
547 
548  R.input = node;
549 }
550 
551 void setupElementaryReaction(ElementaryReaction& R, const XML_Node& rxn_node)
552 {
553  const XML_Node& rc_node = rxn_node.child("rateCoeff");
554  if (rc_node.hasChild("Arrhenius")) {
555  R.rate = readArrhenius(rc_node.child("Arrhenius"));
556  } else if (rc_node.hasChild("Arrhenius_ExchangeCurrentDensity")) {
557  R.rate = readArrhenius(rc_node.child("Arrhenius_ExchangeCurrentDensity"));
558  } else {
559  throw CanteraError("setupElementaryReaction", "Couldn't find Arrhenius node");
560  }
561  if (rxn_node["negative_A"] == "yes") {
562  R.allow_negative_pre_exponential_factor = true;
563  }
564  if (rxn_node["negative_orders"] == "yes") {
565  R.allow_negative_orders = true;
566  }
567  if (rxn_node["nonreactant_orders"] == "yes") {
568  R.allow_nonreactant_orders = true;
569  }
570  setupReaction(R, rxn_node);
571 }
572 
573 void setupElementaryReaction(ElementaryReaction& R, const AnyMap& node,
574  const Kinetics& kin)
575 {
576  setupReaction(R, node, kin);
577  R.allow_negative_pre_exponential_factor = node.getBool("negative-A", false);
578  R.rate = readArrhenius(R, node["rate-constant"], kin, node.units());
579 }
580 
581 void setupThreeBodyReaction(ThreeBodyReaction& R, const XML_Node& rxn_node)
582 {
583  readEfficiencies(R.third_body, rxn_node.child("rateCoeff"));
584  setupElementaryReaction(R, rxn_node);
585 }
586 
587 void setupThreeBodyReaction(ThreeBodyReaction& R, const AnyMap& node,
588  const Kinetics& kin)
589 {
590  setupElementaryReaction(R, node, kin);
591  if (R.reactants.count("M") != 1 || R.products.count("M") != 1) {
592  throw InputFileError("setupThreeBodyReaction", node["equation"],
593  "Reaction equation '{}' does not contain third body 'M'",
594  node["equation"].asString());
595  }
596  R.reactants.erase("M");
597  R.products.erase("M");
598  readEfficiencies(R.third_body, node);
599 }
600 
601 void setupFalloffReaction(FalloffReaction& R, const XML_Node& rxn_node)
602 {
603  XML_Node& rc_node = rxn_node.child("rateCoeff");
604  std::vector<XML_Node*> rates = rc_node.getChildren("Arrhenius");
605  int nLow = 0;
606  int nHigh = 0;
607  for (size_t i = 0; i < rates.size(); i++) {
608  XML_Node& node = *rates[i];
609  if (node["name"] == "") {
610  R.high_rate = readArrhenius(node);
611  nHigh++;
612  } else if (node["name"] == "k0") {
613  R.low_rate = readArrhenius(node);
614  nLow++;
615  } else {
616  throw CanteraError("setupFalloffReaction", "Found an Arrhenius XML "
617  "node with an unexpected type '" + node["name"] + "'");
618  }
619  }
620  if (nLow != 1 || nHigh != 1) {
621  throw CanteraError("setupFalloffReaction", "Did not find the correct "
622  "number of Arrhenius rate expressions");
623  }
624  if (rxn_node["negative_A"] == "yes") {
625  R.allow_negative_pre_exponential_factor = true;
626  }
627  readFalloff(R, rc_node);
628  readEfficiencies(R.third_body, rc_node);
629  setupReaction(R, rxn_node);
630 }
631 
632 void setupFalloffReaction(FalloffReaction& R, const AnyMap& node,
633  const Kinetics& kin)
634 {
635  setupReaction(R, node, kin);
636  // setupReaction sets the stoichiometric coefficient for the falloff third
637  // body to -1.
638  std::string third_body;
639  for (auto& reactant : R.reactants) {
640  if (reactant.second == -1 && ba::starts_with(reactant.first, "(+")) {
641  third_body = reactant.first;
642  break;
643  }
644  }
645 
646  // Equation must contain a third body, and it must appear on both sides
647  if (third_body == "") {
648  throw InputFileError("setupFalloffReaction", node["equation"],
649  "Reactants for reaction '{}' do not contain a pressure-dependent "
650  "third body", node["equation"].asString());
651  } else if (R.products.count(third_body) == 0) {
652  throw InputFileError("setupFalloffReaction", node["equation"],
653  "Unable to match third body '{}' in reactants and products of "
654  "reaction '{}'", third_body, node["equation"].asString());
655  }
656 
657  // Remove the dummy species
658  R.reactants.erase(third_body);
659  R.products.erase(third_body);
660 
661  R.allow_negative_pre_exponential_factor = node.getBool("negative-A", false);
662  if (third_body == "(+M)") {
663  readEfficiencies(R.third_body, node);
664  } else {
665  // Specific species is listed as the third body
666  R.third_body.default_efficiency = 0;
667  R.third_body.efficiencies[third_body.substr(2, third_body.size() - 3)] = 1.0;
668  }
669 
670  if (node["type"] == "falloff") {
671  R.low_rate = readArrhenius(R, node["low-P-rate-constant"], kin,
672  node.units(), 1);
673  R.high_rate = readArrhenius(R, node["high-P-rate-constant"], kin,
674  node.units());
675  } else { // type == "chemically-activated"
676  R.low_rate = readArrhenius(R, node["low-P-rate-constant"], kin,
677  node.units());
678  R.high_rate = readArrhenius(R, node["high-P-rate-constant"], kin,
679  node.units(), -1);
680  }
681 
682  readFalloff(R, node);
683 }
684 
685 void setupChemicallyActivatedReaction(ChemicallyActivatedReaction& R,
686  const XML_Node& rxn_node)
687 {
688  XML_Node& rc_node = rxn_node.child("rateCoeff");
689  std::vector<XML_Node*> rates = rc_node.getChildren("Arrhenius");
690  int nLow = 0;
691  int nHigh = 0;
692  for (size_t i = 0; i < rates.size(); i++) {
693  XML_Node& node = *rates[i];
694  if (node["name"] == "kHigh") {
695  R.high_rate = readArrhenius(node);
696  nHigh++;
697  } else if (node["name"] == "") {
698  R.low_rate = readArrhenius(node);
699  nLow++;
700  } else {
701  throw CanteraError("setupChemicallyActivatedReaction", "Found an "
702  "Arrhenius XML node with an unexpected type '" + node["name"] + "'");
703  }
704  }
705  if (nLow != 1 || nHigh != 1) {
706  throw CanteraError("setupChemicallyActivatedReaction", "Did not find "
707  "the correct number of Arrhenius rate expressions");
708  }
709  readFalloff(R, rc_node);
710  readEfficiencies(R.third_body, rc_node);
711  setupReaction(R, rxn_node);
712 }
713 
714 void setupPlogReaction(PlogReaction& R, const XML_Node& rxn_node)
715 {
716  XML_Node& rc = rxn_node.child("rateCoeff");
717  std::multimap<double, Arrhenius> rates;
718  for (size_t m = 0; m < rc.nChildren(); m++) {
719  const XML_Node& node = rc.child(m);
720  rates.insert({getFloat(node, "P", "toSI"), readArrhenius(node)});
721  }
722  R.rate = Plog(rates);
723  setupReaction(R, rxn_node);
724 }
725 
726 void setupPlogReaction(PlogReaction& R, const AnyMap& node, const Kinetics& kin)
727 {
728  setupReaction(R, node, kin);
729  std::multimap<double, Arrhenius> rates;
730  for (const auto& rate : node.at("rate-constants").asVector<AnyMap>()) {
731  rates.insert({rate.convert("P", "Pa"),
732  readArrhenius(R, AnyValue(rate), kin, node.units())});
733  }
734  R.rate = Plog(rates);
735 }
736 
738 {
740  rate.validate(equation());
741 }
742 
743 void setupChebyshevReaction(ChebyshevReaction& R, const XML_Node& rxn_node)
744 {
745  XML_Node& rc = rxn_node.child("rateCoeff");
746  const XML_Node& coeff_node = rc.child("floatArray");
747  size_t nP = atoi(coeff_node["degreeP"].c_str());
748  size_t nT = atoi(coeff_node["degreeT"].c_str());
749 
750  vector_fp coeffs_flat;
751  getFloatArray(rc, coeffs_flat, false);
752  Array2D coeffs(nT, nP);
753  for (size_t t = 0; t < nT; t++) {
754  for (size_t p = 0; p < nP; p++) {
755  coeffs(t,p) = coeffs_flat[nP*t + p];
756  }
757  }
758  R.rate = ChebyshevRate(getFloat(rc, "Tmin", "toSI"),
759  getFloat(rc, "Tmax", "toSI"),
760  getFloat(rc, "Pmin", "toSI"),
761  getFloat(rc, "Pmax", "toSI"),
762  coeffs);
763  setupReaction(R, rxn_node);
764 }
765 
766 void setupChebyshevReaction(ChebyshevReaction&R, const AnyMap& node,
767  const Kinetics& kin)
768 {
769  setupReaction(R, node, kin);
770  R.reactants.erase("(+M)"); // remove optional third body notation
771  R.products.erase("(+M)");
772  const auto& T_range = node["temperature-range"].asVector<AnyValue>(2);
773  const auto& P_range = node["pressure-range"].asVector<AnyValue>(2);
774  auto& vcoeffs = node["data"].asVector<vector_fp>();
775  Array2D coeffs(vcoeffs.size(), vcoeffs[0].size());
776  for (size_t i = 0; i < coeffs.nRows(); i++) {
777  if (vcoeffs[i].size() != vcoeffs[0].size()) {
778  throw InputFileError("setupChebyshevReaction", node["data"],
779  "Inconsistent number of coefficients in row {} of matrix", i + 1);
780  }
781  for (size_t j = 0; j < coeffs.nColumns(); j++) {
782  coeffs(i, j) = vcoeffs[i][j];
783  }
784  }
785  const UnitSystem& units = node.units();
786  Units rcUnits = rateCoeffUnits(R, kin);
787  coeffs(0, 0) += std::log10(units.convert(1.0, rcUnits));
788  R.rate = ChebyshevRate(units.convert(T_range[0], "K"),
789  units.convert(T_range[1], "K"),
790  units.convert(P_range[0], "Pa"),
791  units.convert(P_range[1], "Pa"),
792  coeffs);
793 }
794 
795 void setupInterfaceReaction(InterfaceReaction& R, const XML_Node& rxn_node)
796 {
797  if (caseInsensitiveEquals(rxn_node["type"], "global")) {
798  R.reaction_type = GLOBAL_RXN;
799  }
800  XML_Node& arr = rxn_node.child("rateCoeff").child("Arrhenius");
801  if (caseInsensitiveEquals(arr["type"], "stick")) {
802  R.is_sticking_coefficient = true;
803  R.sticking_species = arr["species"];
804 
805  if (caseInsensitiveEquals(arr["motz_wise"], "true")) {
806  R.use_motz_wise_correction = true;
807  } else if (caseInsensitiveEquals(arr["motz_wise"], "false")) {
808  R.use_motz_wise_correction = false;
809  } else {
810  // Default value for all reactions
811  XML_Node* parent = rxn_node.parent();
812  if (parent && parent->name() == "reactionData"
813  && caseInsensitiveEquals((*parent)["motz_wise"], "true")) {
814  R.use_motz_wise_correction = true;
815  }
816  }
817  }
818  std::vector<XML_Node*> cov = arr.getChildren("coverage");
819  for (const auto& node : cov) {
820  CoverageDependency& cdep = R.coverage_deps[node->attrib("species")];
821  cdep.a = getFloat(*node, "a", "toSI");
822  cdep.m = getFloat(*node, "m");
823  cdep.E = getFloat(*node, "e", "actEnergy") / GasConstant;
824  }
825  setupElementaryReaction(R, rxn_node);
826 }
827 
828 void setupInterfaceReaction(InterfaceReaction& R, const AnyMap& node,
829  const Kinetics& kin)
830 {
831  setupReaction(R, node, kin);
832  R.allow_negative_pre_exponential_factor = node.getBool("negative-A", false);
833 
834  if (node.hasKey("rate-constant")) {
835  R.rate = readArrhenius(R, node["rate-constant"], kin, node.units());
836  } else if (node.hasKey("sticking-coefficient")) {
837  R.is_sticking_coefficient = true;
838  R.rate = readArrhenius(R, node["sticking-coefficient"], kin, node.units());
839  R.use_motz_wise_correction = node.getBool("Motz-Wise",
840  kin.thermo().input().getBool("Motz-Wise", false));
841  R.sticking_species = node.getString("sticking-species", "");
842  } else {
843  throw InputFileError("setupInterfaceReaction", node,
844  "Reaction must include either a 'rate-constant' or"
845  " 'sticking-coefficient' node.");
846  }
847 
848  if (node.hasKey("coverage-dependencies")) {
849  for (const auto& item : node["coverage-dependencies"].as<AnyMap>()) {
850  double a, E, m;
851  if (item.second.is<AnyMap>()) {
852  auto& cov_map = item.second.as<AnyMap>();
853  a = cov_map["a"].asDouble();
854  m = cov_map["m"].asDouble();
855  E = node.units().convertActivationEnergy(cov_map["E"], "K");
856  } else {
857  auto& cov_vec = item.second.asVector<AnyValue>(3);
858  a = cov_vec[0].asDouble();
859  m = cov_vec[1].asDouble();
860  E = node.units().convertActivationEnergy(cov_vec[2], "K");
861  }
862  R.coverage_deps[item.first] = CoverageDependency(a, E, m);
863  }
864  }
865 }
866 
867 void setupElectrochemicalReaction(ElectrochemicalReaction& R,
868  const XML_Node& rxn_node)
869 {
870  // Fix reaction_type for some specialized reaction types
871  std::string type = toLowerCopy(rxn_node["type"]);
872  if (type == "butlervolmer") {
873  R.reaction_type = BUTLERVOLMER_RXN;
874  warn_deprecated("reaction type 'ButlerVolmer'",
875  "To be removed after Cantera 2.5.");
876  } else if (type == "butlervolmer_noactivitycoeffs") {
877  R.reaction_type = BUTLERVOLMER_NOACTIVITYCOEFFS_RXN;
878  warn_deprecated("reaction type 'butlervolmer_noactivitycoeffs'",
879  "To be removed after Cantera 2.5.");
880  } else if (type == "surfaceaffinity") {
881  R.reaction_type = SURFACEAFFINITY_RXN;
882  warn_deprecated("reaction type 'surfaceaffinity'",
883  "To be removed after Cantera 2.5.");
884  } else if (type == "global") {
885  R.reaction_type = GLOBAL_RXN;
886  warn_deprecated("reaction type 'global'",
887  "To be removed after Cantera 2.5.");
888  }
889 
890  XML_Node& rc = rxn_node.child("rateCoeff");
891  std::string rc_type = toLowerCopy(rc["type"]);
892  if (rc_type == "exchangecurrentdensity") {
893  R.exchange_current_density_formulation = true;
894  } else if (rc_type != "" && rc_type != "arrhenius") {
895  throw CanteraError("setupElectrochemicalReaction",
896  "Unknown rate coefficient type: '" + rc_type + "'");
897  }
898  if (rc.hasChild("Arrhenius_ExchangeCurrentDensity")) {
899  R.exchange_current_density_formulation = true;
900  }
901 
902  if (rc.hasChild("electrochem") && rc.child("electrochem").hasAttrib("beta")) {
903  R.beta = fpValueCheck(rc.child("electrochem")["beta"]);
904  }
905 
906  if (rxn_node.hasChild("filmResistivity")) {
907  warn_deprecated("reaction filmResistivity",
908  "Not implemented. To be removed after Cantera 2.5.");
909  }
910  getOptionalFloat(rxn_node, "filmResistivity", R.film_resistivity);
911  setupInterfaceReaction(R, rxn_node);
912 
913  // For Butler Volmer reactions, install the orders for the exchange current
914  if (R.reaction_type == BUTLERVOLMER_NOACTIVITYCOEFFS_RXN ||
915  R.reaction_type == BUTLERVOLMER_RXN) {
916  if (!R.reversible) {
917  throw CanteraError("setupElectrochemicalReaction",
918  "A Butler-Volmer reaction must be reversible");
919  }
920 
921  R.orders.clear();
922  // Reaction orders based on species stoichiometric coefficients
923  R.allow_nonreactant_orders = true;
924  for (const auto& sp : R.reactants) {
925  R.orders[sp.first] += sp.second * (1.0 - R.beta);
926  }
927  for (const auto& sp : R.products) {
928  R.orders[sp.first] += sp.second * R.beta;
929  }
930  }
931 
932  // For affinity reactions, fill in the global reaction formulation terms
933  if (rxn_node.hasChild("reactionOrderFormulation")) {
934  warn_deprecated("reactionOrderFormulation",
935  "To be removed after Cantera 2.5.");
936  Composition initial_orders = R.orders;
937  R.orders.clear();
938  R.allow_nonreactant_orders = true;
939  const XML_Node& rof_node = rxn_node.child("reactionOrderFormulation");
940  if (caseInsensitiveEquals(rof_node["model"], "reactantorders")) {
941  R.orders = initial_orders;
942  } else if (caseInsensitiveEquals(rof_node["model"], "zeroorders")) {
943  for (const auto& sp : R.reactants) {
944  R.orders[sp.first] = 0.0;
945  }
946  } else if (caseInsensitiveEquals(rof_node["model"], "butlervolmerorders")) {
947  // Reaction orders based on provided reaction orders
948  for (const auto& sp : R.reactants) {
949  double c = getValue(initial_orders, sp.first, sp.second);
950  R.orders[sp.first] += c * (1.0 - R.beta);
951  }
952  for (const auto& sp : R.products) {
953  double c = getValue(initial_orders, sp.first, sp.second);
954  R.orders[sp.first] += c * R.beta;
955  }
956  } else {
957  throw CanteraError("setupElectrochemicalReaction", "unknown model "
958  "for reactionOrderFormulation XML_Node: '" +
959  rof_node["model"] + "'");
960  }
961  }
962 
963  // Override orders based on the <orders> node
964  if (rxn_node.hasChild("orders")) {
965  Composition orders = parseCompString(rxn_node.child("orders").value());
966  for (const auto& order : orders) {
967  R.orders[order.first] = order.second;
968  }
969  }
970 }
971 
972 void setupElectrochemicalReaction(ElectrochemicalReaction& R,
973  const AnyMap& node, const Kinetics& kin)
974 {
975  setupInterfaceReaction(R, node, kin);
976  R.beta = node.getDouble("beta", 0.5);
977  R.exchange_current_density_formulation = node.getBool(
978  "exchange-current-density-formulation", false);
979 }
980 
981 bool isElectrochemicalReaction(Reaction& R, const Kinetics& kin)
982 {
983  vector_fp e_counter(kin.nPhases(), 0.0);
984 
985  // Find the number of electrons in the products for each phase
986  for (const auto& sp : R.products) {
987  size_t kkin = kin.kineticsSpeciesIndex(sp.first);
988  size_t i = kin.speciesPhaseIndex(kkin);
989  size_t kphase = kin.thermo(i).speciesIndex(sp.first);
990  e_counter[i] += sp.second * kin.thermo(i).charge(kphase);
991  }
992 
993  // Subtract the number of electrons in the reactants for each phase
994  for (const auto& sp : R.reactants) {
995  size_t kkin = kin.kineticsSpeciesIndex(sp.first);
996  size_t i = kin.speciesPhaseIndex(kkin);
997  size_t kphase = kin.thermo(i).speciesIndex(sp.first);
998  e_counter[i] -= sp.second * kin.thermo(i).charge(kphase);
999  }
1000 
1001  // If the electrons change phases then the reaction is electrochemical
1002  for (double delta_e : e_counter) {
1003  if (std::abs(delta_e) > 1e-4) {
1004  return true;
1005  }
1006  }
1007  return false;
1008 }
1009 
1010 shared_ptr<Reaction> newReaction(const XML_Node& rxn_node)
1011 {
1012  std::string type = toLowerCopy(rxn_node["type"]);
1013 
1014  // Modify the reaction type for interface reactions which contain
1015  // electrochemical reaction data
1016  if (rxn_node.child("rateCoeff").hasChild("electrochem")
1017  && (type == "edge" || type == "surface")) {
1018  type = "electrochemical";
1019  }
1020 
1021  // Create a new Reaction object of the appropriate type
1022  if (type == "elementary" || type == "arrhenius" || type == "") {
1023  auto R = make_shared<ElementaryReaction>();
1024  setupElementaryReaction(*R, rxn_node);
1025  return R;
1026  } else if (type == "threebody" || type == "three_body") {
1027  auto R = make_shared<ThreeBodyReaction>();
1028  setupThreeBodyReaction(*R, rxn_node);
1029  return R;
1030  } else if (type == "falloff") {
1031  auto R = make_shared<FalloffReaction>();
1032  setupFalloffReaction(*R, rxn_node);
1033  return R;
1034  } else if (type == "chemact" || type == "chemically_activated") {
1035  auto R = make_shared<ChemicallyActivatedReaction>();
1036  setupChemicallyActivatedReaction(*R, rxn_node);
1037  return R;
1038  } else if (type == "plog" || type == "pdep_arrhenius") {
1039  auto R = make_shared<PlogReaction>();
1040  setupPlogReaction(*R, rxn_node);
1041  return R;
1042  } else if (type == "chebyshev") {
1043  auto R = make_shared<ChebyshevReaction>();
1044  setupChebyshevReaction(*R, rxn_node);
1045  return R;
1046  } else if (type == "interface" || type == "surface" || type == "edge" ||
1047  type == "global") {
1048  auto R = make_shared<InterfaceReaction>();
1049  setupInterfaceReaction(*R, rxn_node);
1050  return R;
1051  } else if (type == "electrochemical" ||
1052  type == "butlervolmer_noactivitycoeffs" ||
1053  type == "butlervolmer" ||
1054  type == "surfaceaffinity") {
1055  auto R = make_shared<ElectrochemicalReaction>();
1056  setupElectrochemicalReaction(*R, rxn_node);
1057  return R;
1058  } else {
1059  throw CanteraError("newReaction",
1060  "Unknown reaction type '" + rxn_node["type"] + "'");
1061  }
1062 }
1063 
1064 unique_ptr<Reaction> newReaction(const AnyMap& node, const Kinetics& kin)
1065 {
1066  std::string type = "elementary";
1067  if (node.hasKey("type")) {
1068  type = node["type"].asString();
1069  }
1070 
1071  if (kin.thermo(kin.reactionPhaseIndex()).nDim() < 3) {
1072  // See if this is an electrochemical reaction
1073  Reaction testReaction(0);
1074  parseReactionEquation(testReaction, node["equation"], kin);
1075  if (isElectrochemicalReaction(testReaction, kin)) {
1076  unique_ptr<ElectrochemicalReaction> R(new ElectrochemicalReaction());
1077  setupElectrochemicalReaction(*R, node, kin);
1078  return unique_ptr<Reaction>(move(R));
1079  } else {
1080  unique_ptr<InterfaceReaction> R(new InterfaceReaction());
1081  setupInterfaceReaction(*R, node, kin);
1082  return unique_ptr<Reaction>(move(R));
1083  }
1084  }
1085 
1086  if (type == "elementary") {
1087  unique_ptr<ElementaryReaction> R(new ElementaryReaction());
1088  setupElementaryReaction(*R, node, kin);
1089  return unique_ptr<Reaction>(move(R));
1090  } else if (type == "three-body") {
1091  unique_ptr<ThreeBodyReaction> R(new ThreeBodyReaction());
1092  setupThreeBodyReaction(*R, node, kin);
1093  return unique_ptr<Reaction>(move(R));
1094  } else if (type == "falloff") {
1095  unique_ptr<FalloffReaction> R(new FalloffReaction());
1096  setupFalloffReaction(*R, node, kin);
1097  return unique_ptr<Reaction>(move(R));
1098  } else if (type == "chemically-activated") {
1099  unique_ptr<ChemicallyActivatedReaction> R(new ChemicallyActivatedReaction());
1100  setupFalloffReaction(*R, node, kin);
1101  return unique_ptr<Reaction>(move(R));
1102  } else if (type == "pressure-dependent-Arrhenius") {
1103  unique_ptr<PlogReaction> R(new PlogReaction());
1104  setupPlogReaction(*R, node, kin);
1105  return unique_ptr<Reaction>(move(R));
1106  } else if (type == "Chebyshev") {
1107  unique_ptr<ChebyshevReaction> R(new ChebyshevReaction());
1108  setupChebyshevReaction(*R, node, kin);
1109  return unique_ptr<Reaction>(move(R));
1110  } else {
1111  throw InputFileError("newReaction", node["type"],
1112  "Unknown reaction type '{}'", type);
1113  }
1114 }
1115 
1116 std::vector<shared_ptr<Reaction> > getReactions(const XML_Node& node)
1117 {
1118  std::vector<shared_ptr<Reaction> > all_reactions;
1119  for (const auto& rxnnode : node.child("reactionData").getChildren("reaction")) {
1120  all_reactions.push_back(newReaction(*rxnnode));
1121  }
1122  return all_reactions;
1123 }
1124 
1125 std::vector<shared_ptr<Reaction>> getReactions(const AnyValue& items,
1126  Kinetics& kinetics)
1127 {
1128  std::vector<shared_ptr<Reaction>> all_reactions;
1129  for (const auto& node : items.asVector<AnyMap>()) {
1130  shared_ptr<Reaction> R(newReaction(node, kinetics));
1131  if (R->reaction_type != INVALID_RXN) {
1132  all_reactions.emplace_back(R);
1133  } else if (!kinetics.skipUndeclaredSpecies()) {
1134  throw InputFileError("getReactions", node,
1135  "Reaction '{}' contains undeclared species.", R->equation());
1136  }
1137  };
1138  return all_reactions;
1139 }
1140 
1141 }
Header file for class Cantera::Array2D.
Parameterizations for reaction falloff functions.
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
A map of string keys to values whose type can vary at runtime.
Definition: AnyMap.h:360
bool hasKey(const std::string &key) const
Returns true if the map contains an item named key.
Definition: AnyMap.cpp:984
A wrapper for a variable whose type is determined at runtime.
Definition: AnyMap.h:77
const std::vector< T > & asVector(size_t nMin=npos, size_t nMax=npos) const
Return the held value, if it is a vector of type T.
Definition: AnyMap.inl.h:73
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:32
Arrhenius reaction rate type depends only on temperature.
Definition: RxnRates.h:32
double preExponentialFactor() const
Return the pre-exponential factor A (in m, kmol, s to powers depending on the reaction order)
Definition: RxnRates.h:82
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 where the rate decreases as pressure increases due to collisional stabilization of a react...
Definition: Reaction.h:164
An interface reaction which involves charged species.
Definition: Reaction.h:242
A reaction which follows mass-action kinetics with a modified Arrhenius reaction rate.
Definition: Reaction.h:84
virtual void validate()
Ensure that the rate constant and other parameters for this reaction are valid.
Definition: Reaction.cpp:120
A reaction that is first-order in [M] at low pressure, like a third-body reaction,...
Definition: Reaction.h:133
virtual std::string reactantString() const
The reactant side of the chemical equation for this reaction.
Definition: Reaction.cpp:178
virtual void validate()
Ensure that the rate constant and other parameters for this reaction are valid.
Definition: Reaction.cpp:198
virtual std::string productString() const
The product side of the chemical equation for this reaction.
Definition: Reaction.cpp:188
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
Base class for falloff function calculators.
Definition: Falloff.h:30
Error thrown for problems processing information contained in an AnyMap or AnyValue.
Definition: AnyMap.h:539
A reaction occurring on an interface (i.e. a SurfPhase or an EdgePhase)
Definition: Reaction.h:213
Public interface for kinetics managers.
Definition: Kinetics.h:111
size_t reactionPhaseIndex() const
Phase where the reactions occur.
Definition: Kinetics.h:214
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
void skipUndeclaredSpecies(bool skip)
Determine behavior when adding a new reaction that contains species not defined in any of the phases ...
Definition: Kinetics.h:740
size_t nDim() const
Returns the number of spatial dimensions (1, 2, or 3)
Definition: Phase.h:651
A pressure-dependent reaction parameterized by logarithmically interpolating between Arrhenius rate e...
Definition: Reaction.h:175
virtual void validate()
Ensure that the rate constant and other parameters for this reaction are valid.
Definition: Reaction.cpp:737
Intermediate class which stores data about a reaction and its rate parameterization so that it can be...
Definition: Reaction.h:24
virtual std::string reactantString() const
The reactant side of the chemical equation for this reaction.
Definition: Reaction.cpp:66
virtual void validate()
Ensure that the rate constant and other parameters for this reaction are valid.
Definition: Reaction.cpp:45
virtual std::string productString() const
The product side of the chemical equation for this reaction.
Definition: Reaction.cpp:81
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
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
virtual std::string reactantString() const
The reactant side of the chemical equation for this reaction.
Definition: Reaction.cpp:151
virtual std::string productString() const
The product side of the chemical equation for this reaction.
Definition: Reaction.cpp:155
Class XML_Node is a tree-based representation of the contents of an XML file.
Definition: xml.h:104
bool hasChild(const std::string &ch) const
Tests whether the current node has a child node with a particular name.
Definition: xml.cpp:528
std::vector< XML_Node * > getChildren(const std::string &name) const
Get a vector of pointers to XML_Node containing all of the children of the current node which match t...
Definition: xml.cpp:711
XML_Node & child(const size_t n) const
Return a changeable reference to the n'th child of the current node.
Definition: xml.cpp:546
CTML ("Cantera Markup Language") is the variant of XML that Cantera uses to store data.
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:54
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:188
std::map< std::string, double > Composition
Map from string names to doubles.
Definition: ct_defs.h:176
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
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
const U & getValue(const std::map< T, U > &m, const T &key, const U &default_val)
Const accessor for a value in a std::map.
Definition: utilities.h:528
const int BUTLERVOLMER_NOACTIVITYCOEFFS_RXN
This is a surface reaction that is formulated using the Butler-Volmer formulation and using concentra...
Definition: reaction_defs.h:83
void getStringArray(const XML_Node &node, std::vector< std::string > &v)
This function interprets the value portion of an XML element as a string.
Definition: ctml.cpp:427
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
std::vector< shared_ptr< Reaction > > getReactions(const XML_Node &node)
Create Reaction objects for all <reaction> nodes in an XML document.
Definition: Reaction.cpp:1116
doublereal getFloat(const XML_Node &parent, const std::string &name, const std::string &type)
Get a floating-point value from a child element.
Definition: ctml.cpp:164
bool getOptionalFloat(const XML_Node &parent, const std::string &name, doublereal &fltRtn, const std::string &type)
Get an optional floating-point value from a child element.
Definition: ctml.cpp:212
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 BUTLERVOLMER_RXN
This is a surface reaction that is formulated using the Butler-Volmer formulation.
Definition: reaction_defs.h:90
shared_ptr< Reaction > newReaction(const XML_Node &rxn_node)
Create a new Reaction object for the reaction defined in rxn_node
Definition: Reaction.cpp:1010
doublereal fpValue(const std::string &val)
Translate a string into one doublereal value.
bool caseInsensitiveEquals(const std::string &input, const std::string &test)
Case insensitive equality predicate.
const int GLOBAL_RXN
A global reaction.
const int THREE_BODY_RXN
A gas-phase reaction that requires a third-body collision partner.
Definition: reaction_defs.h:38
std::string toLowerCopy(const std::string &input)
Convert to lower case.
void readFalloff(FalloffReaction &R, const XML_Node &rc_node)
Parse falloff parameters, given a rateCoeff node.
Definition: Reaction.cpp:358
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 INTERFACE_RXN
A reaction occurring on an interface, e.g a surface or edge.
Definition: reaction_defs.h:77
size_t getFloatArray(const XML_Node &node, vector_fp &v, const bool convert, const std::string &unitsString, const std::string &nodeName)
This function reads the current node or a child node of the current node with the default name,...
Definition: ctml.cpp:256
const int SURFACEAFFINITY_RXN
This is a surface reaction that is formulated using the affinity representation, common in the geoche...
Definition: reaction_defs.h:97
const int CHEMACT_RXN
A chemical activation reaction.
Definition: reaction_defs.h:66
doublereal fpValueCheck(const std::string &val)
Translate a string into one doublereal value, with error checking.
shared_ptr< Falloff > newFalloff(int type, const vector_fp &c)
Return a pointer to a new falloff function calculator.
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
void tokenizeString(const std::string &in_val, std::vector< std::string > &v)
This function separates a string up into tokens according to the location of white space.
compositionMap parseCompString(const std::string &ss, const std::vector< std::string > &names)
Parse a composition string into a map consisting of individual key:composition pairs.
Definition: stringUtils.cpp:60