28 template<
class ReactionRate,
class ChemicallyActivationFunction>
32 ChemicallyActivationFunction
33 >::ChemicallyActivatedReactionRate
35 const ReactionRate& k0,
36 const ReactionRate& kInf,
37 const ChemicallyActivationFunction& F,
44 thirdBodyEfficiencies_(tbes)
53 template<
class ReactionRate,
class ChemicallyActivationFunction>
57 ChemicallyActivationFunction
67 thirdBodyEfficiencies_(species, dict)
69 forAll(thirdBodyEfficiencies_, i)
76 thirdBodyEfficiencies_[i]
85 template<
class ReactionRate,
class ChemicallyActivationFunction>
89 ChemicallyActivationFunction
90 >::preEvaluate()
const 97 template<
class ReactionRate,
class ChemicallyActivationFunction>
101 ChemicallyActivationFunction
102 >::postEvaluate()
const 105 kInf_.postEvaluate();
109 template<
class ReactionRate,
class ChemicallyActivationFunction>
113 ChemicallyActivationFunction
122 const scalar k0 = k0_(p, T, c, li);
123 const scalar kInf = kInf_(p, T, c, li);
124 const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
126 return k0*(1/(1 + Pr))*F_(T, Pr);
130 template<
class ReactionRate,
class ChemicallyActivationFunction>
134 ChemicallyActivationFunction
143 const scalar k0 = k0_(p, T, c, li);
144 const scalar kInf = kInf_(p, T, c, li);
145 const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
147 return (1/(1 + Pr))*F_(T, Pr)*k0_.ddT(p, T, c, li);
151 template<
class ReactionRate,
class ChemicallyActivationFunction>
155 ChemicallyActivationFunction
165 const scalar
M = thirdBodyEfficiencies_.M(c);
169 const scalar k0 = k0_(p, T, c, li);
170 const scalar kInf = kInf_(p, T, c, li);
172 const scalar Pr = k0*M/kInf;
173 const scalar F = F_(T, Pr);
177 const scalar dPrdci = -beta_[i].second()*k0/kInf;
178 const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
179 dcidc[i] = (-dPrdci/(1 + Pr) + dFdci/F);
192 template<
class ReactionRate,
class ChemicallyActivationFunction>
196 ChemicallyActivationFunction
205 const scalar
M = thirdBodyEfficiencies_.M(c);
209 const scalar k0 = k0_(p, T, c, li);
210 const scalar kInf = kInf_(p, T, c, li);
212 const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
213 const scalar F = F_(T, Pr);
215 Pr*(k0_.ddT(p, T, c, li)/k0 - kInf_.ddT(p, T, c, li)/kInf - 1/
T);
216 const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
218 return (-dPrdT/(1 + Pr) + dFdT/F);
227 template<
class ReactionRate,
class ChemicallyActivationFunction>
231 ChemicallyActivationFunction
237 thirdBodyEfficiencies_.write(os);
241 template<
class ReactionRate,
class ChemicallyActivationFunction>
246 <ReactionRate, ChemicallyActivationFunction>& carr
#define forAll(list, i)
Loop across all elements in list.
virtual Ostream & write(const char)=0
Write character.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A list of keyword definitions, which are a keyword followed by any number of values (e...
A 2-tuple for storing two objects of different types.
const dimensionedScalar c
Speed of light in a vacuum.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
void preEvaluate() const
Pre-evaluation hook.
General class for handling chemically-activated bimolecular reactions.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
A wordList with hashed indices for faster lookup by name.
void postEvaluate() const
Post-evaluation hook.