30 #define CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN false
34 template<
class ReactionRate,
class ChemicallyActivationFunction>
38 ChemicallyActivationFunction
41 const ReactionRate& k0,
42 const ReactionRate& kInf,
43 const ChemicallyActivationFunction&
F,
50 thirdBodyEfficiencies_(tbes)
54 template<
class ReactionRate,
class ChemicallyActivationFunction>
58 ChemicallyActivationFunction
66 k0_(species, dims,
dict.subDict(
"k0")),
67 kInf_(species, dims,
dict.subDict(
"kInf")),
68 F_(
dict.subDict(
"F")),
69 thirdBodyEfficiencies_(species,
dict.subDict(
"thirdBodyEfficiencies"))
75 template<
class ReactionRate,
class ChemicallyActivationFunction>
79 ChemicallyActivationFunction
80 >::preEvaluate()
const
87 template<
class ReactionRate,
class ChemicallyActivationFunction>
91 ChemicallyActivationFunction
92 >::postEvaluate()
const
99 template<
class ReactionRate,
class ChemicallyActivationFunction>
103 ChemicallyActivationFunction
112 const scalar k0 = k0_(
p,
T,
c, li);
113 const scalar kInf = kInf_(
p,
T,
c, li);
114 const scalar
M = thirdBodyEfficiencies_.M(
c);
115 const scalar Pr = k0/kInf*
M;
116 const scalar
F = F_(
T, Pr);
118 return k0/(1 + Pr)*
F;
122 template<
class ReactionRate,
class ChemicallyActivationFunction>
126 ChemicallyActivationFunction
135 const scalar k0 = k0_(
p,
T,
c, li);
136 const scalar kInf = kInf_(
p,
T,
c, li);
137 const scalar
M = thirdBodyEfficiencies_.M(
c);
138 const scalar Pr = k0/kInf*
M;
139 const scalar
F = F_(
T, Pr);
141 const scalar dk0dT = k0_.ddT(
p,
T,
c, li);
143 #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
145 const scalar dkInfdT = kInf_.ddT(
p,
T,
c, li);
146 const scalar dPrdT = (
M*dk0dT - Pr*dkInfdT)/kInf;
147 const scalar dFdT = F_.ddT(
T, Pr,
F) + F_.ddPr(
T, Pr,
F)*dPrdT;
151 - k0*dPrdT/
sqr(1 + Pr)*
F
156 return dk0dT/(1 + Pr)*
F;
162 template<
class ReactionRate,
class ChemicallyActivationFunction>
166 ChemicallyActivationFunction
173 template<
class ReactionRate,
class ChemicallyActivationFunction>
177 ChemicallyActivationFunction
187 const scalar k0 = k0_(
p,
T,
c, li);
188 const scalar kInf = kInf_(
p,
T,
c, li);
189 const scalar
M = thirdBodyEfficiencies_.M(
c);
190 const scalar Pr = k0/kInf*
M;
191 const scalar
F = F_(
T, Pr);
193 #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
196 k0_.ddc(
p,
T,
c, li, dk0dc);
198 kInf_.ddc(
p,
T,
c, li, dkInfdc);
201 scalarField dPrdc((
M*dk0dc - Pr*dkInfdc + k0*dMdc)/kInf);
202 const scalar dFdPr = F_.ddPr(
T, Pr,
F);
206 - k0*dPrdc/
sqr(1 + Pr)*
F
207 + k0/(1 + Pr)*dFdPr*dPrdc;
211 k0_.ddc(
p,
T,
c, li, ddc);
219 template<
class ReactionRate,
class ChemicallyActivationFunction>
223 ChemicallyActivationFunction
247 os <<
indent <<
"thirdBodyEfficiencies" <<
nl;
250 thirdBodyEfficiencies_.
write(os);
256 template<
class ReactionRate,
class ChemicallyActivationFunction>
261 <ReactionRate, ChemicallyActivationFunction>& carr
General class for handling chemically-activated bimolecular reactions.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Dimension set for the base types.
A wordList with hashed indices for faster lookup by name.
A class for managing temporary objects.
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
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.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Ostream & indent(Ostream &os)
Indent stream.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)