47 alpha_(
dict.lookup<scalar>(
"alpha")),
48 Tsss_(
dict.lookup<scalar>(
"Tsss")),
49 Ts_(
dict.lookup<scalar>(
"Ts")),
50 Tss_(
dict.lookup<scalar>(
"Tss"))
56 inline Foam::scalar Foam::TroeFallOffFunction::operator()
62 const scalar logPr =
log10(
max(Pr, small));
65 (1 - alpha_)*
exp(-
T/Tsss_)
69 const scalar logFcent =
log10(
max(Fcent, small));
71 const scalar
c = -0.4 - 0.67*logFcent;
72 static const scalar d = 0.14;
73 const scalar
n = 0.75 - 1.27*logFcent;
75 const scalar x1 =
n - d*(logPr +
c);
76 const scalar x2 = (logPr +
c)/x1;
77 const scalar x3 = 1 +
sqr(x2);
78 const scalar x4 = logFcent/x3;
91 static const scalar logTen =
log(scalar(10));
93 const scalar logPr =
log10(
max(Pr, small));
96 (1 - alpha_)*
exp(-
T/Tsss_)
99 const scalar dFcentdT =
100 - (1 - alpha_)/Tsss_*
exp(-
T/Tsss_)
101 - alpha_/Ts_*
exp(-
T/Ts_)
104 const scalar logFcent =
log10(
max(Fcent, small));
105 const scalar dlogFcentdT = Fcent >= small ? dFcentdT/Fcent/logTen : 0;
107 const scalar
c = -0.4 - 0.67*logFcent;
108 const scalar dcdT = -0.67*dlogFcentdT;
109 static const scalar d = 0.14;
110 const scalar
n = 0.75 - 1.27*logFcent;
111 const scalar dndT = - 1.27*dlogFcentdT;
113 const scalar x1 =
n - d*(logPr +
c);
114 const scalar dx1dT = dndT - d*dcdT;
115 const scalar x2 = (logPr +
c)/x1;
116 const scalar dx2dT = (dcdT - x2*dx1dT)/x1;
117 const scalar x3 = 1 +
sqr(x2);
118 const scalar dx3dT = 2*x2*dx2dT;
119 const scalar x4 = logFcent/x3;
120 const scalar dx4dT = (dlogFcentdT - x4*dx3dT)/x3;
122 return logTen*
F*dx4dT;
133 static const scalar logTen =
log(scalar(10));
135 const scalar logPr =
log10(
max(Pr, small));
136 const scalar dlogPrdPr = Pr >= small ? 1/(logTen*Pr) : 0;
139 (1 - alpha_)*
exp(-
T/Tsss_)
143 const scalar logFcent =
log10(
max(Fcent, small));
145 const scalar
c = -0.4 - 0.67*logFcent;
146 static const scalar d = 0.14;
147 const scalar
n = 0.75 - 1.27*logFcent;
149 const scalar x1 =
n - d*(logPr +
c);
150 const scalar dx1dPr = -d*dlogPrdPr;
151 const scalar x2 = (logPr +
c)/x1;
152 const scalar dx2dPr = (dlogPrdPr - x2*dx1dPr)/x1;
153 const scalar x3 = 1 +
sqr(x2);
154 const scalar dx3dPr = 2*x2*dx2dPr;
155 const scalar x4 = logFcent/x3;
156 const scalar dx4dPr = -x4*dx3dPr/x3;
158 return logTen*
F*dx4dPr;
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual Ostream & write(const char)=0
Write character.
The Troe fall-off function.
scalar ddT(const scalar T, const scalar Pr, const scalar F) const
TroeFallOffFunction(const scalar alpha, const scalar Tsss, const scalar Ts, const scalar Tss)
Construct from components.
void write(Ostream &os) const
Write to stream.
scalar ddPr(const scalar T, const scalar Pr, const scalar F) const
A list of keyword definitions, which are a keyword followed by any number of values (e....
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)