57 static scalar
calcQE11(
const scalar a,
const scalar
x,
const int e = 30)
66 for (
n = 1; (2*
n) <=
e;
n++)
68 const scalar a_2nm1 = a_2np1;
69 const scalar b_2nm1 = b_2np1;
71 a_2n = a_2nm1 + (
n - a)*a_2n;
72 b_2n = b_2nm1 + (
n - a)*b_2n;
74 a_2np1 =
x*a_2n +
n*a_2nm1;
75 b_2np1 =
x*b_2n +
n*b_2nm1;
90 static scalar
calcPE15(
const scalar a,
const scalar
x,
const int nmax = 20)
95 for (
int n = 1;
n <= nmax;
n++)
101 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
103 return R/a*(1 +
sum);
108 static scalar
calcQE16(
const scalar a,
const scalar
x,
const int N = 20)
113 for (
int n = 1;
n <= (N - 1);
n++)
119 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
121 return R/
x*(1 +
sum);
136 static const scalar D0[] =
138 -0.333333333333333E-00,
139 0.833333333333333E-01,
140 -0.148148148148148E-01,
141 0.115740740740741E-02,
142 0.352733686067019E-03,
143 -0.178755144032922E-03,
144 0.391926317852244E-04,
145 -0.218544851067999E-05,
146 -0.185406221071516E-05,
147 0.829671134095309E-06,
148 -0.176659527368261E-06,
149 0.670785354340150E-08,
150 0.102618097842403E-07,
151 -0.438203601845335E-08
154 static const scalar D1[] =
156 -0.185185185185185E-02,
157 -0.347222222222222E-02,
158 0.264550264550265E-02,
159 -0.990226337448560E-03,
160 0.205761316872428E-03,
161 -0.401877572016461E-06,
162 -0.180985503344900E-04,
163 0.764916091608111E-05,
164 -0.161209008945634E-05,
165 0.464712780280743E-08,
166 0.137863344691572E-06,
167 -0.575254560351770E-07,
168 0.119516285997781E-07
171 static const scalar D2[] =
173 0.413359788359788E-02,
174 -0.268132716049383E-02,
175 0.771604938271605E-03,
176 0.200938786008230E-05,
177 -0.107366532263652E-03,
178 0.529234488291201E-04,
179 -0.127606351886187E-04,
180 0.342357873409614E-07,
181 0.137219573090629E-05,
182 -0.629899213838006E-06,
183 0.142806142060642E-06
186 const scalar u = 1/a;
187 scalar z =
sqrt(2*phi);
198 + D0[3]*
pow3(z) + D0[2]*
sqr(z) + D0[1]*z + D0[0];
201 D1[4]*
pow4(z) + D1[3]*
pow3(z) + D1[2]*
sqr(z) + D1[1]*z + D1[0];
203 const scalar C2 = D2[1]*z + D2[0];
205 return C2*
sqr(u) + C1*u + C0;
209 const scalar C0 = D0[2]*
sqr(z) + D0[1]*z + D0[0];
210 const scalar C1 = D1[1]*z + D1[0];
211 const scalar C2 = D2[1]*z + D2[0];
213 return C2*
sqr(u) + C1*u + C0;
226 const scalar BIG = 14;
227 const scalar x0 = 17;
228 const scalar e0 = 0.025;
256 for (
int n = 1;
n <= 10;
n++)
261 const scalar J = -a*
sum;
266 return 1 - (
pow(
x, a)*(1 - J))/tgamma(a + 1);
271 const scalar L =
exp(a*
log(
x)) - 1;
272 const scalar H = 1/(tgamma(a + 1)) - 1;
274 return (
pow(
x, a)*J - L)/tgamma(a + 1) - H;
280 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
287 const scalar
sigma = fabs(1 -
x/a);
294 const scalar
y = a*phi;
296 const scalar E = 0.5 - (1 -
y/3)*
sqrt(
y/
pi);
323 const scalar
y = a*phi;
351 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
365 if (a >
x ||
x >= x0)
375 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
387 if (floor(2*a) == 2*a)
394 for (
int n = 0;
n <= (a - 1);
n++)
407 for (
int n = 1;
n <= i;
n++)
416 else if (
x <=
max(a,
log(10.0)))
424 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
dimensionedScalar lambda(viscosity->lookup("lambda"))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static scalar calcQE11(const scalar a, const scalar x, const int e=30)
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar incGamma_P(const scalar a, const scalar x)
Lower incomplete gamma function.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
static scalar calcPE15(const scalar a, const scalar x, const int nmax=20)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
label factorial(label n)
Return n! : 0 < n <= 12.
scalar incGamma_Q(const scalar a, const scalar x)
Upper incomplete gamma function.
scalar incGammaRatio_P(const scalar a, const scalar x)
Normalised lower incomplete gamma function.
static scalar calcQE16(const scalar a, const scalar x, const int N=20)
static scalar calcTE18(const scalar a, const scalar e0, const scalar x, const scalar lambda, const scalar sigma, const scalar phi)
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalised upper incomplete gamma function.