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 const scalar logPr =
log10(
max(Pr, small));
63 (1 - alpha_)*
exp(-T/Tsss_)
67 const scalar logFcent =
log10(
max(Fcent, small));
69 const scalar
c = -0.4 - 0.67*logFcent;
70 static const scalar d = 0.14;
71 const scalar
n = 0.75 - 1.27*logFcent;
73 const scalar x1 = n - d*(logPr +
c);
74 const scalar x2 = (logPr +
c)/x1;
75 const scalar x3 = 1 +
sqr(x2);
76 const scalar x4 = logFcent/x3;
89 static const scalar logTen =
log(scalar(10));
91 const scalar logPr =
log10(
max(Pr, small));
94 (1 - alpha_)*
exp(-T/Tsss_)
97 const scalar dFcentdT =
98 - (1 - alpha_)/Tsss_*
exp(-T/Tsss_)
99 - alpha_/Ts_*
exp(-T/Ts_)
100 + Tss_/
sqr(T)*
exp(-Tss_/T);
102 const scalar logFcent =
log10(
max(Fcent, small));
103 const scalar dlogFcentdT = Fcent >= small ? dFcentdT/Fcent/logTen : 0;
105 const scalar
c = -0.4 - 0.67*logFcent;
106 const scalar dcdT = -0.67*dlogFcentdT;
107 static const scalar d = 0.14;
108 const scalar
n = 0.75 - 1.27*logFcent;
109 const scalar dndT = - 1.27*dlogFcentdT;
111 const scalar x1 = n - d*(logPr +
c);
112 const scalar dx1dT = dndT - d*dcdT;
113 const scalar x2 = (logPr +
c)/x1;
114 const scalar dx2dT = (dcdT - x2*dx1dT)/x1;
115 const scalar x3 = 1 +
sqr(x2);
116 const scalar dx3dT = 2*x2*dx2dT;
117 const scalar x4 = logFcent/x3;
118 const scalar dx4dT = (dlogFcentdT - x4*dx3dT)/x3;
120 return logTen*F*dx4dT;
131 static const scalar logTen =
log(scalar(10));
133 const scalar logPr =
log10(
max(Pr, small));
134 const scalar dlogPrdPr = Pr >= small ? 1/(logTen*Pr) : 0;
137 (1 - alpha_)*
exp(-T/Tsss_)
141 const scalar logFcent =
log10(
max(Fcent, small));
143 const scalar
c = -0.4 - 0.67*logFcent;
144 static const scalar d = 0.14;
145 const scalar
n = 0.75 - 1.27*logFcent;
147 const scalar x1 = n - d*(logPr +
c);
148 const scalar dx1dPr = -d*dlogPrdPr;
149 const scalar x2 = (logPr +
c)/x1;
150 const scalar dx2dPr = (dlogPrdPr - x2*dx1dPr)/x1;
151 const scalar x3 = 1 +
sqr(x2);
152 const scalar dx3dPr = 2*x2*dx2dPr;
153 const scalar x4 = logFcent/x3;
154 const scalar dx4dPr = -x4*dx3dPr/x3;
156 return logTen*F*dx4dPr;
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
scalar ddPr(const scalar T, const scalar Pr, const scalar F) const
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
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 ddT(const scalar T, const scalar Pr, const scalar F) const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar log10(const dimensionedScalar &ds)