invIncGamma.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) 2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Global
28  Foam::Math::invIncGamma
29 
30 Description
31  Implementation of the inverse incomplete gamma function.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "MathFunctions.H"
36 #include "mathematicalConstants.H"
37 #include "error.H"
38 #include <cmath>
39 #include <limits>
40 
41 using namespace Foam::constant::mathematical;
42 
43 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 static scalar minimaxs(const scalar P)
49 {
50  // (DM:Eq. 32)
51 
52  constexpr scalar a_0 = 3.31125922108741;
53  constexpr scalar a_1 = 11.6616720288968;
54  constexpr scalar a_2 = 4.28342155967104;
55  constexpr scalar a_3 = 0.213623493715853;
56 
57  constexpr scalar b_0 = 6.61053765625462;
58  constexpr scalar b_1 = 6.40691597760039;
59  constexpr scalar b_2 = 1.27364489782223;
60  constexpr scalar b_3 = 0.03611708101884203;
61 
62  const scalar t = P < 0.5 ? sqrt(-2*log(P)) : sqrt(-2*log(1 - P));
63 
64  const scalar s =
65  t
66  - (a_0 + t*(a_1 + t*(a_2 + t*a_3)))
67  /(1 + t*(b_0 + t*(b_1 + t*(b_2 + t*b_3))));
68 
69  return P < 0.5 ? -s : s;
70 }
71 
72 
73 static scalar Sn(const scalar a, const scalar x)
74 {
75  // (DM:Eq. 34)
76 
77  scalar Sn = 1;
78  scalar Si = 1;
79 
80  for (int i=1; i<100; ++i)
81  {
82  Si *= x/(a + i);
83  Sn += Si;
84 
85  if (Si < 1e-4) break;
86  }
87 
88  return Sn;
89 }
90 
91 } // End namespace Foam
92 
93 
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96 Foam::scalar Foam::Math::invIncGamma(const scalar a, const scalar P)
97 {
98  #ifdef FULLDEBUG
99  if (a <= 0)
100  {
102  << "The parameter (i.e. a) cannot be negative or zero"
103  << " a = " << a
104  << endl;
105 
106  return std::numeric_limits<scalar>::infinity();
107  }
108 
109  if (P < 0 || P > 1)
110  {
112  << "The domain of the parameter (i.e. P) should be limited to [0,1]"
113  << " P = " << P
114  << endl;
115 
116  return std::numeric_limits<scalar>::infinity();
117  }
118  #endif
119 
120  const scalar Q = 1 - P;
121 
122  if (a == 1)
123  {
124  return -log(Q);
125  }
126  else if (a < 1)
127  {
128  const scalar Ga = tgamma(a);
129  const scalar B = Q*Ga;
130 
131  if (B > 0.6 || (B >= 0.45 && a >= 0.3))
132  {
133  // (DM:Eq. 21)
134  const scalar u =
135  (B*Q > 1e-8) ? pow(P*Ga*a, 1/a) : exp((-Q/a) - Eu);
136 
137  return u/(1 - (u/(a + 1)));
138  }
139  else if (a < 0.3 && B >= 0.35)
140  {
141  // (DM:Eq. 22)
142  const scalar t = exp(-Eu - B);
143  const scalar u = t*exp(t);
144 
145  return t*exp(u);
146  }
147  else if (B > 0.15 || a >= 0.3)
148  {
149  // (DM:Eq. 23)
150  const scalar y = -log(B);
151  const scalar u = y - (1 - a)*log(y);
152 
153  return y - (1 - a)*log(u) - log(1 + (1 - a)/(1 + u));
154  }
155  else if (B > 0.1)
156  {
157  // (DM:Eq. 24)
158  const scalar y = -log(B);
159  const scalar u = y - (1 - a)*log(y);
160 
161  return y
162  - (1 - a)*log(u)
163  - log
164  (
165  (sqr(u) + 2*(3 - a)*u + (2 - a)*(3 - a))
166  /(sqr(u) + (5 - a)*u + 2)
167  );
168  }
169  else
170  {
171  // (DM:Eq. 25)
172  const scalar y = -log(B);
173  const scalar c1 = (a - 1)*log(y);
174  const scalar c12 = c1*c1;
175  const scalar c13 = c12*c1;
176  const scalar c14 = c12*c12;
177  const scalar a2 = a*a;
178  const scalar a3 = a2*a;
179  const scalar c2 = (a - 1)*(1 + c1);
180  const scalar c3 = (a - 1)*(-(c12/2) + (a - 2)*c1 + (3*a - 5)/2);
181  const scalar c4 =
182  (a - 1)
183  *(
184  (c13/3)
185  - (3*a - 5)*c12/2
186  + (a2 - 6*a + 7)*c1
187  + (11*a2 - 46*a + 47)/6
188  );
189  const scalar c5 =
190  (a - 1)*(-(c14/4)
191  + (11*a - 17)*c13/6
192  + (-3*a2 + 13*a - 13)*c12
193  + (2*a3 - 25*a2 + 72*a - 61)*c1/2
194  + (25*a3 - 195*a2 + 477*a - 379)/12);
195  const scalar y2 = y*y;
196  const scalar y3 = y2*y;
197  const scalar y4 = y2*y2;
198 
199  return y + c1 + (c2/y) + (c3/y2) + (c4/y3) + (c5/y4);
200  }
201  }
202  else
203  {
204  // (DM:Eq. 31)
205  scalar s = minimaxs(P);
206 
207  const scalar s2 = sqr(s);
208  const scalar s3 = s*s2;
209  const scalar s4 = s2*s2;
210  const scalar s5 = s*s4;
211  const scalar sqrta = sqrt(a);
212 
213  const scalar w =
214  a + s*sqrta + (s2 - 1)/3
215  + (s3 - 7*s)/(36*sqrta)
216  - (3*s4 + 7*s2 - 16)/(810*a)
217  + (9*s5 + 256*s3 - 433*s)/(38880*a*sqrta);
218 
219  if (a >= 500 && mag(1 - w/a) < 1e-6)
220  {
221  return w;
222  }
223  else if (P > 0.5)
224  {
225  if (w < 3*a)
226  {
227  return w;
228  }
229  else
230  {
231  const scalar D = max(scalar(2), scalar(a*(a - 1)));
232  const scalar lnGa = lgamma(a);
233  const scalar lnB = log(Q) + lnGa;
234 
235  if (lnB < -2.3*D)
236  {
237  // (DM:Eq. 25)
238  const scalar y = -lnB;
239  const scalar c1 = (a - 1)*log(y);
240  const scalar c12 = c1*c1;
241  const scalar c13 = c12*c1;
242  const scalar c14 = c12*c12;
243  const scalar a2 = a*a;
244  const scalar a3 = a2*a;
245 
246  const scalar c2 = (a - 1)*(1 + c1);
247  const scalar c3 =
248  (a - 1)
249  *(
250  - (c12/2)
251  + (a - 2)*c1
252  + (3*a - 5)/2
253  );
254  const scalar c4 =
255  (a - 1)
256  *(
257  (c13/3)
258  - (3*a - 5)*c12/2
259  + (a2 - 6*a + 7)*c1
260  + (11*a2 - 46*a + 47)/6
261  );
262  const scalar c5 =
263  (a - 1)
264  *(
265  - (c14/4)
266  + (11*a - 17)*c13/6
267  + (-3*a2 + 13*a - 13)*c12
268  + (2*a3 - 25*a2 + 72*a - 61)*c1/2
269  + (25*a3 - 195*a2 + 477*a - 379)/12
270  );
271 
272  const scalar y2 = y*y;
273  const scalar y3 = y2*y;
274  const scalar y4 = y2*y2;
275 
276  return y + c1 + (c2/y) + (c3/y2) + (c4/y3) + (c5/y4);
277  }
278  else
279  {
280  // (DM:Eq. 33)
281  const scalar u =
282  -lnB + (a - 1)*log(w) - log(1 + (1 - a)/(1 + w));
283 
284  return -lnB + (a - 1)*log(u) - log(1 + (1 - a)/(1 + u));
285  }
286  }
287  }
288  else
289  {
290  scalar z = w;
291  const scalar ap1 = a + 1;
292 
293  if (w < 0.15*ap1)
294  {
295  // (DM:Eq. 35)
296  const scalar ap2 = a + 2;
297  const scalar v = log(P) + lgamma(ap1);
298  z = exp((v + w)/a);
299  s = log1p(z/ap1*(1 + z/ap2));
300  z = exp((v + z - s)/a);
301  s = log1p(z/ap1*(1 + z/ap2));
302  z = exp((v + z - s)/a);
303  s = log1p(z/ap1*(1 + z/ap2*(1 + z/(a + 3))));
304  z = exp((v + z - s)/a);
305  }
306 
307  if (z <= 0.01*ap1 || z > 0.7*ap1)
308  {
309  return z;
310  }
311  else
312  {
313  // (DM:Eq. 36)
314  const scalar lnSn = log(Sn(a, z));
315  const scalar v = log(P) + lgamma(ap1);
316  z = exp((v + z - lnSn)/a);
317 
318  return z*(1 - (a*log(z) - z - v + lnSn)/(a - z));
319  }
320  }
321  }
322 }
323 
324 
325 // ************************************************************************* //
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static scalar Sn(const scalar a, const scalar x)
Definition: invIncGamma.C:64
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static scalar minimaxs(const scalar P)
Definition: invIncGamma.C:39
scalar y
constexpr scalar Eu(0.57721566490153286060651209)
Euler&#39;s constant.
dimensionedScalar exp(const dimensionedScalar &ds)
Mathematical constants.
scalar invIncGamma(const scalar a, const scalar P)
Inverse of regularised lower incomplete gamma function.
Definition: invIncGamma.C:87
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
const dimensionedScalar & D
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
constexpr scalar e(M_E)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Namespace for OpenFOAM.