45 alpha_(dict.
lookup<scalar>(
"alpha")),
46 Tsss_(dict.
lookup<scalar>(
"Tsss")),
47 Ts_(dict.
lookup<scalar>(
"Ts")),
48 Tss_(dict.
lookup<scalar>(
"Tss"))
54 inline Foam::scalar Foam::TroeFallOffFunction::operator()
60 scalar logFcent =
log10 64 (1 - alpha_)*
exp(-T/Tsss_) + alpha_*
exp(-T/Ts_) +
exp(-Tss_/T),
69 scalar
c = -0.4 - 0.67*logFcent;
70 static const scalar d = 0.14;
71 scalar
n = 0.75 - 1.27*logFcent;
73 scalar logPr =
log10(
max(Pr, small));
74 return pow(10, logFcent/(1 +
sqr((logPr + c)/(n - d*(logPr + c)))));
86 scalar logPr =
log10(
max(Pr, small));
87 scalar logTen =
log(10.0);
92 (1 - alpha_)*
exp(-T/Tsss_) + alpha_*
exp(-T/Ts_) +
exp(-Tss_/T),
96 scalar logFcent =
log10(Fcent);
100 (alpha_ - 1)*
exp(-T/Tsss_)/Tsss_
101 - alpha_*
exp(-T/Ts_)/Ts_
102 + Tss_*
exp(-Tss_/T)/
sqr(T)
106 scalar dlogFcentdT = dFcentdT/Fcent/logTen;
107 scalar
c = -0.4 - 0.67*logFcent;
108 scalar dcdT = -0.67*dlogFcentdT;
109 scalar
n = 0.75 - 1.27*logFcent;
110 scalar dndT = -1.27*dlogFcentdT;
112 scalar dlogPrdT = dPrdT/Pr/logTen;
115 2.0*(logPr +
c)/
sqr(n - d*(logPr + c))
118 - (logPr +
c)*(dndT - d*(dlogPrdT + dcdT))/(n - d*(logPr + c))
125 dlogFcentdT/(1.0 +
sqr((logPr + c)/(n - d*(logPr + c))))
126 - logFcent*dParentdT/
sqr(1.0 +
sqr((logPr + c)/(n - d*(logPr + c))))
140 scalar logPr =
log10(
max(Pr, small));
141 scalar logTen =
log(10.0);
142 scalar logFcent =
log10 146 (1 - alpha_)*
exp(-T/Tsss_) + alpha_*
exp(-T/Ts_) +
exp(-Tss_/T),
151 scalar dlogPrdc = dPrdc/Pr/logTen;
153 scalar
c = -0.4 - 0.67*logFcent;
154 scalar
n = 0.75 - 1.27*logFcent;
157 2.0*(logPr +
c)/
sqr(n - d*(logPr + c))
160 - (logPr +
c)*(-d*(dlogPrdc))/(n - d*(logPr + c))
167 - logFcent*dParentdc/
sqr(1.0 +
sqr((logPr + c)/(n - d*(logPr + c))))
virtual Ostream & write(const char)=0
Write character.
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
The Troe fall-off function.
void write(Ostream &os) const
Write to stream.
const dimensionedScalar c
Speed of light in a vacuum.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
TroeFallOffFunction(const scalar alpha, const scalar Tsss, const scalar Ts, const scalar Tss)
Construct from components.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
scalar ddc(const scalar Pr, const scalar F, const scalar dPrdc, const scalar T) const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar ddT(const scalar Pr, const scalar F, const scalar dPrdT, const scalar T) const
dimensionedScalar log10(const dimensionedScalar &ds)