47static scalar
calcQE11(
const scalar a,
const scalar
x,
const int e = 30)
56 for (
n = 1; (2*
n) <=
e;
n++)
58 const scalar a_2nm1 = a_2np1;
59 const scalar b_2nm1 = b_2np1;
61 a_2n = a_2nm1 + (
n - a)*a_2n;
62 b_2n = b_2nm1 + (
n - a)*b_2n;
64 a_2np1 =
x*a_2n +
n*a_2nm1;
65 b_2np1 =
x*b_2n +
n*b_2nm1;
80static scalar
calcPE15(
const scalar a,
const scalar
x,
const int nmax = 20)
85 for (
int n = 1;
n <= nmax;
n++)
91 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
98static scalar
calcQE16(
const scalar a,
const scalar
x,
const int N = 20)
103 for (
int n = 1;
n <= (
N - 1);
n++)
109 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
111 return R/
x*(1 +
sum);
126 constexpr scalar D0_0 = -0.333333333333333E-00;
127 constexpr scalar D0_1 = 0.833333333333333E-01;
128 constexpr scalar D0_2 = -0.148148148148148E-01;
129 constexpr scalar D0_3 = 0.115740740740741E-02;
130 constexpr scalar D0_4 = 0.352733686067019E-03;
131 constexpr scalar D0_5 = -0.178755144032922E-03;
132 constexpr scalar D0_6 = 0.391926317852244E-04;
141 constexpr scalar D1_0 = -0.185185185185185E-02;
142 constexpr scalar D1_1 = -0.347222222222222E-02;
143 constexpr scalar D1_2 = 0.264550264550265E-02;
144 constexpr scalar D1_3 = -0.990226337448560E-03;
145 constexpr scalar D1_4 = 0.205761316872428E-03;
155 constexpr scalar D2_0 = 0.413359788359788E-02;
156 constexpr scalar D2_1 = -0.268132716049383E-02;
167 const scalar u = 1/a;
175 if (sigma > (e0/
sqrt(a)))
179 + D0_3*
pow3(z) + D0_2*
sqr(z) + D0_1*z + D0_0;
182 D1_4*
pow4(z) + D1_3*
pow3(z) + D1_2*
sqr(z) + D1_1*z + D1_0;
184 const scalar C2 = D2_1*z + D2_0;
186 return C2*
sqr(u) + C1*u + C0;
190 const scalar C0 = D0_2*
sqr(z) + D0_1*z + D0_0;
191 const scalar C1 = D1_1*z + D1_0;
192 const scalar C2 = D2_1*z + D2_0;
194 return C2*
sqr(u) + C1*u + C0;
211 <<
"The parameter (i.e. a) cannot be negative or zero"
215 return std::numeric_limits<scalar>::infinity();
221 <<
"The parameter (i.e. x) cannot be negative"
225 return std::numeric_limits<scalar>::infinity();
229 constexpr scalar BIG = 14;
230 constexpr scalar x0 = 17;
231 constexpr scalar e0 = 0.025;
259 for (label
n = 1;
n <= 10; ++
n)
264 const scalar J = -a*
sum;
269 return 1 - (
pow(
x, a)*(1 - J))/tgamma(a + 1);
275 const scalar
H = 1/(tgamma(a + 1)) - 1;
277 return (
pow(
x, a)*J -
L)/tgamma(a + 1) -
H;
283 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
290 const scalar sigma = fabs(1 -
x/a);
292 if (sigma <= e0/
sqrt(a))
297 const scalar
y = a*
phi;
299 const scalar E = 0.5 - (1 -
y/3)*
sqrt(
y/pi);
307 - (1 -
y)/
sqrt(2*pi*a)
315 + (1 -
y)/
sqrt(2*pi*a)
326 const scalar
y = a*
phi;
354 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
368 if (a >
x ||
x >= x0)
378 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
390 if (floor(2*a) == 2*a)
397 for (label
n = 0;
n <= (a - 1); ++
n)
410 for (
int n = 1;
n <= i;
n++)
419 else if (
x <=
max(a,
log(10.0)))
427 const scalar
R = (
exp(-
x)*
pow(
x, a))/tgamma(a);
443 return 1 - incGammaRatio_Q(a,
x);
449 return incGammaRatio_Q(a,
x)*tgamma(a);
455 return incGammaRatio_P(a,
x)*tgamma(a);
#define R(A, B, C, D, E, F, K, M)
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
#define WarningInFunction
Report a warning using Foam::Warning.
scalar incGamma_P(const scalar a, const scalar x)
Lower incomplete gamma function.
scalar incGamma_Q(const scalar a, const scalar x)
Upper incomplete gamma function.
scalar incGammaRatio_P(const scalar a, const scalar x)
Regularised lower incomplete gamma function.
scalar incGammaRatio_Q(const scalar a, const scalar x)
Regularised upper incomplete gamma function.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
static scalar calcQE11(const scalar a, const scalar x, const int e=30)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar erf(const dimensionedScalar &ds)
static scalar calcQE16(const scalar a, const scalar x, const int N=20)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
static scalar calcTE18(const scalar a, const scalar e0, const scalar x, const scalar lambda, const scalar sigma, const scalar phi)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
label factorial(label n)
Evaluate n! : 0 < n <= 12.
static scalar calcPE15(const scalar a, const scalar x, const int nmax=20)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
const vector L(dict.get< vector >("L"))
const Vector< label > N(dict.get< Vector< label > >("N"))