incGamma.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Global
28 Foam::Math::incGamma
29
30Description
31 Implementation of the incomplete gamma functions.
32
33\*---------------------------------------------------------------------------*/
34
35#include "MathFunctions.H"
37#include "error.H"
38#include <cmath>
39#include <limits>
40
41// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
42
43namespace Foam
45
46// (DM:Eq. 13)
47static scalar calcQE11(const scalar a, const scalar x, const int e = 30)
48{
49 scalar a_2n = 0;
50 scalar b_2n = 1;
51
52 scalar a_2np1 = 1;
53 scalar b_2np1 = x;
54
55 int n = 1;
56 for (n = 1; (2*n) <= e; n++)
57 {
58 const scalar a_2nm1 = a_2np1;
59 const scalar b_2nm1 = b_2np1;
60
61 a_2n = a_2nm1 + (n - a)*a_2n;
62 b_2n = b_2nm1 + (n - a)*b_2n;
63
64 a_2np1 = x*a_2n + n*a_2nm1;
65 b_2np1 = x*b_2n + n*b_2nm1;
66 }
67
68 if (2*(n - 1) < e)
69 {
70 return a_2np1/b_2np1;
71 }
72 else
73 {
74 return a_2n/b_2n;
75 }
76}
78
79// (DM:Eq. 15)
80static scalar calcPE15(const scalar a, const scalar x, const int nmax = 20)
81{
82 scalar prod = 1;
83 scalar sum = 0;
84
85 for (int n = 1; n <= nmax; n++)
86 {
87 prod *= (a + n);
88 sum += pow(x, n)/prod;
89 }
90
91 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
92
93 return R/a*(1 + sum);
94}
96
97// (DM:Eq. 16)
98static scalar calcQE16(const scalar a, const scalar x, const int N = 20)
99{
100 scalar an = 1;
101 scalar sum = 0;
102
103 for (int n = 1; n <= (N - 1); n++)
104 {
105 an *= (a - n);
106 sum += an/pow(x, n);
107 }
108
109 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
110
111 return R/x*(1 + sum);
112}
114
115// (DM:Eq. 18)
116static scalar calcTE18
117(
118 const scalar a,
119 const scalar e0,
120 const scalar x,
121 const scalar lambda,
122 const scalar sigma,
123 const scalar phi
124)
125{
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;
133 // unused: constexpr scalar D0_7 = -0.218544851067999E-05;
134 // unused: constexpr scalar D0_8 = -0.185406221071516E-05;
135 // unused: constexpr scalar D0_9 = 0.829671134095309E-06;
136 // unused: constexpr scalar D0_10 = -0.176659527368261E-06;
137 // unused: constexpr scalar D0_11 = 0.670785354340150E-08;
138 // unused: constexpr scalar D0_12 = 0.102618097842403E-07;
139 // unused: constexpr scalar D0_13 = -0.438203601845335E-08;
140
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;
146 // unused: constexpr scalar D1_5 = -0.401877572016461E-06;
147 // unused: constexpr scalar D1_6 = -0.180985503344900E-04;
148 // unused: constexpr scalar D1_7 = 0.764916091608111E-05;
149 // unused: constexpr scalar D1_8 = -0.161209008945634E-05;
150 // unused: constexpr scalar D1_9 = 0.464712780280743E-08;
151 // unused: constexpr scalar D1_10 = 0.137863344691572E-06;
152 // unused: constexpr scalar D1_11 = -0.575254560351770E-07;
153 // unused: constexpr scalar D1_12 = 0.119516285997781E-07;
154
155 constexpr scalar D2_0 = 0.413359788359788E-02;
156 constexpr scalar D2_1 = -0.268132716049383E-02;
157 // unused: constexpr scalar D2_2 = 0.771604938271605E-03;
158 // unused: constexpr scalar D2_3 = 0.200938786008230E-05;
159 // unused: constexpr scalar D2_4 = -0.107366532263652E-03;
160 // unused: constexpr scalar D2_5 = 0.529234488291201E-04;
161 // unused: constexpr scalar D2_6 = -0.127606351886187E-04;
162 // unused: constexpr scalar D2_7 = 0.342357873409614E-07;
163 // unused: constexpr scalar D2_8 = 0.137219573090629E-05;
164 // unused: constexpr scalar D2_9 = -0.629899213838006E-06;
165 // unused: constexpr scalar D2_10 = 0.142806142060642E-06;
166
167 const scalar u = 1/a;
168 scalar z = sqrt(2*phi);
169
170 if (lambda < 1)
171 {
172 z = -z;
173 }
174
175 if (sigma > (e0/sqrt(a)))
176 {
177 const scalar C0 =
178 D0_6*pow6(z) + D0_5*pow5(z) + D0_4*pow4(z)
179 + D0_3*pow3(z) + D0_2*sqr(z) + D0_1*z + D0_0;
180
181 const scalar C1 =
182 D1_4*pow4(z) + D1_3*pow3(z) + D1_2*sqr(z) + D1_1*z + D1_0;
183
184 const scalar C2 = D2_1*z + D2_0;
185
186 return C2*sqr(u) + C1*u + C0;
187 }
188 else
189 {
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;
193
194 return C2*sqr(u) + C1*u + C0;
195 }
196}
197
198} // End namespace Foam
199
201// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202
203Foam::scalar Foam::Math::incGammaRatio_Q(const scalar a, const scalar x)
204{
205 using namespace Foam::constant::mathematical;
206
207 #ifdef FULLDEBUG
208 if (a <= 0)
209 {
211 << "The parameter (i.e. a) cannot be negative or zero"
212 << " a = " << a
213 << endl;
214
215 return std::numeric_limits<scalar>::infinity();
216 }
217
218 if (x < 0)
219 {
221 << "The parameter (i.e. x) cannot be negative"
222 << " x = " << x
223 << endl;
224
225 return std::numeric_limits<scalar>::infinity();
226 }
227 #endif
228
229 constexpr scalar BIG = 14;
230 constexpr scalar x0 = 17;
231 constexpr scalar e0 = 0.025;
232
233 if (a < 1)
234 {
235 if (a == 0.5)
236 {
237 // (DM:Eq. 8)
238 if (x < 0.25)
239 {
240 return 1 - erf(sqrt(x));
241 }
242 else
243 {
244 return erfc(sqrt(x));
245 }
246 }
247 else if ( x < 1.1)
248 {
249 // (DM:Eq. 12)
250 scalar alpha = x/2.59;
251
252 if (x < 0.5)
253 {
254 alpha = log(sqrt(0.765))/log(x);
255 }
256
257 scalar sum = 0;
258
259 for (label n = 1; n <= 10; ++n)
260 {
261 sum += pow((-x), n)/((a + n)*factorial(n));
262 }
263
264 const scalar J = -a*sum;
265
266 if (a > alpha || a == alpha)
267 {
268 // (DM:Eq. 9)
269 return 1 - (pow(x, a)*(1 - J))/tgamma(a + 1);
270 }
271 else
272 {
273 // (DM:Eq. 10)
274 const scalar L = exp(a*log(x)) - 1;
275 const scalar H = 1/(tgamma(a + 1)) - 1;
276
277 return (pow(x, a)*J - L)/tgamma(a + 1) - H;
278 }
279 }
280 else
281 {
282 // (DM:Eq. 11)
283 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
284
285 return R*calcQE11(a, x);
286 }
287 }
288 else if (a >= BIG)
289 {
290 const scalar sigma = fabs(1 - x/a);
291
292 if (sigma <= e0/sqrt(a))
293 {
294 // (DM:Eq. 19)
295 const scalar lambda = x/a;
296 const scalar phi = lambda - 1 - log(lambda);
297 const scalar y = a*phi;
298
299 const scalar E = 0.5 - (1 - y/3)*sqrt(y/pi);
300
301 if (lambda <= 1)
302 {
303 return
304 1
305 - (
306 E
307 - (1 - y)/sqrt(2*pi*a)
308 *calcTE18(a, e0, x, lambda, sigma, phi)
309 );
310 }
311 else
312 {
313 return
314 E
315 + (1 - y)/sqrt(2*pi*a)
316 *calcTE18(a, e0, x, lambda, sigma, phi);
317 }
318 }
319 else
320 {
321 if (sigma <= 0.4)
322 {
323 // (DM:Eq. 17)
324 const scalar lambda = x/a;
325 const scalar phi = lambda - 1 - log(lambda);
326 const scalar y = a*phi;
327
328 if (lambda <= 1)
329 {
330 return
331 1
332 - (0.5*erfc(sqrt(y))
333 - exp(-y)/sqrt(2*pi*a)
334 *calcTE18(a, e0, x, lambda, sigma, phi));
335 }
336 else
337 {
338 return
339 0.5*erfc(sqrt(y))
340 + exp(-y)/sqrt(2*pi*a)
341 *calcTE18(a, e0, x, lambda, sigma, phi);
342 }
343 }
344 else
345 {
346 if (x <= max(a, log(10.0)))
347 {
348 // (DM:Eq. 15)
349 return 1 - calcPE15(a, x);
350 }
351 else if (x < x0)
352 {
353 // (DM:Eq. 11)
354 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
355
356 return R*calcQE11(a, x);
357 }
358 else
359 {
360 // (DM:Eq. 16)
361 return calcQE16(a, x);
362 }
363 }
364 }
365 }
366 else
367 {
368 if (a > x || x >= x0)
369 {
370 if (x <= max(a, log(10.0)))
371 {
372 // (DM:Eq. 15)
373 return 1 - calcPE15(a, x);
374 }
375 else if ( x < x0)
376 {
377 // (DM:Eq. 11)
378 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
379
380 return R*calcQE11(a, x);
381 }
382 else
383 {
384 // (DM:Eq. 16)
385 return calcQE16(a, x);
386 }
387 }
388 else
389 {
390 if (floor(2*a) == 2*a)
391 {
392 // (DM:Eq. 14)
393 if (floor(a) == a)
394 {
395 scalar sum = 0;
396
397 for (label n = 0; n <= (a - 1); ++n)
398 {
399 sum += pow(x, n)/factorial(n);
400 }
401
402 return exp(-x)*sum;
403 }
404 else
405 {
406 int i = a - 0.5;
407 scalar prod = 1;
408 scalar sum = 0;
409
410 for (int n = 1; n <= i; n++)
411 {
412 prod *= (n - 0.5);
413 sum += pow(x, n)/prod;
414 }
415
416 return erfc(sqrt(x)) + exp(-x)/sqrt(pi*x)*sum;
417 }
418 }
419 else if (x <= max(a, log(10.0)))
420 {
421 // (DM:Eq. 15)
422 return 1 - calcPE15(a, x);
423 }
424 else if (x < x0)
425 {
426 // (DM:Eq. 11)
427 const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
428
429 return R*calcQE11(a, x);
430 }
431 else
432 {
433 // (DM:Eq. 16)
434 return calcQE16(a, x);
435 }
436 }
437 }
439
440
441Foam::scalar Foam::Math::incGammaRatio_P(const scalar a, const scalar x)
442{
443 return 1 - incGammaRatio_Q(a, x);
445
446
447Foam::scalar Foam::Math::incGamma_Q(const scalar a, const scalar x)
448{
449 return incGammaRatio_Q(a, x)*tgamma(a);
451
452
453Foam::scalar Foam::Math::incGamma_P(const scalar a, const scalar x)
454{
455 return incGammaRatio_P(a, x)*tgamma(a);
456}
457
458
459// ************************************************************************* //
scalar y
#define R(A, B, C, D, E, F, K, M)
label n
Y[inertIndex] max(0.0)
surfaceScalarField & phi
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.
Definition: incGamma.C:450
scalar incGamma_Q(const scalar a, const scalar x)
Upper incomplete gamma function.
Definition: incGamma.C:444
scalar incGammaRatio_P(const scalar a, const scalar x)
Regularised lower incomplete gamma function.
Definition: incGamma.C:438
scalar incGammaRatio_Q(const scalar a, const scalar x)
Regularised upper incomplete gamma function.
Definition: incGamma.C:200
Mathematical constants.
Namespace for OpenFOAM.
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)
Definition: incGamma.C:44
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)
Definition: incGamma.C:95
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static scalar calcTE18(const scalar a, const scalar e0, const scalar x, const scalar lambda, const scalar sigma, const scalar phi)
Definition: incGamma.C:114
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
label factorial(label n)
Evaluate n! : 0 < n <= 12.
Definition: label.C:145
static scalar calcPE15(const scalar a, const scalar x, const int nmax=20)
Definition: incGamma.C:77
volScalarField & alpha
volScalarField & e
Definition: createFields.H:11
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
const vector L(dict.get< vector >("L"))
const Vector< label > N(dict.get< Vector< label > >("N"))