48static scalar
minimaxs(
const scalar P)
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;
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;
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))));
69 return P < 0.5 ? -
s :
s;
73static scalar
Sn(
const scalar a,
const scalar
x)
80 for (
int i=1; i<100; ++i)
102 <<
"The parameter (i.e. a) cannot be negative or zero"
106 return std::numeric_limits<scalar>::infinity();
112 <<
"The domain of the parameter (i.e. P) should be limited to [0,1]"
116 return std::numeric_limits<scalar>::infinity();
120 const scalar Q = 1 - P;
128 const scalar Ga = tgamma(a);
129 const scalar
B = Q*Ga;
131 if (
B > 0.6 || (
B >= 0.45 && a >= 0.3))
135 (
B*Q > 1
e-8) ?
pow(P*Ga*a, 1/a) :
exp((-Q/a) -
Eu);
137 return u/(1 - (u/(a + 1)));
139 else if (a < 0.3 && B >= 0.35)
142 const scalar t =
exp(-
Eu -
B);
143 const scalar u = t*
exp(t);
147 else if (
B > 0.15 || a >= 0.3)
150 const scalar
y = -
log(
B);
151 const scalar u =
y - (1 - a)*
log(
y);
153 return y - (1 - a)*
log(u) -
log(1 + (1 - a)/(1 + u));
158 const scalar
y = -
log(
B);
159 const scalar u =
y - (1 - a)*
log(
y);
165 (
sqr(u) + 2*(3 - a)*u + (2 - a)*(3 - a))
166 /(
sqr(u) + (5 - a)*u + 2)
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);
187 + (11*a2 - 46*a + 47)/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;
199 return y + c1 + (c2/
y) + (c3/y2) + (c4/y3) + (c5/y4);
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);
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);
219 if (a >= 500 &&
mag(1 - w/a) < 1
e-6)
231 const scalar
D =
max(scalar(2), scalar(a*(a - 1)));
232 const scalar lnGa =
lgamma(a);
233 const scalar lnB =
log(Q) + lnGa;
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;
246 const scalar c2 = (a - 1)*(1 + c1);
260 + (11*a2 - 46*a + 47)/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
272 const scalar y2 =
y*
y;
273 const scalar y3 = y2*
y;
274 const scalar y4 = y2*y2;
276 return y + c1 + (c2/
y) + (c3/y2) + (c4/y3) + (c5/y4);
282 -lnB + (a - 1)*
log(w) -
log(1 + (1 - a)/(1 + w));
284 return -lnB + (a - 1)*
log(u) -
log(1 + (1 - a)/(1 + u));
291 const scalar ap1 = a + 1;
296 const scalar ap2 = a + 2;
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);
307 if (z <= 0.01*ap1 || z > 0.7*ap1)
314 const scalar lnSn =
log(
Sn(a, z));
316 z =
exp((v + z - lnSn)/a);
318 return z*(1 - (a*
log(z) - z - v + lnSn)/(a - z));
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
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))
#define WarningInFunction
Report a warning using Foam::Warning.
scalar invIncGamma(const scalar a, const scalar P)
Inverse of regularised lower incomplete gamma function.
constexpr scalar Eu(0.57721566490153286060651209)
Euler's constant.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static scalar minimaxs(const scalar P)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static scalar Sn(const scalar a, const scalar x)
const dimensionedScalar & D