/// <summary> /// Returns the sum of a hypergeometric series 2F1. /// <para>Sum2F1 = addend + Σ(( (a1)_n * (a2)_n)/(b1)_n * (z^n/n!)) n={0,∞}</para> /// </summary> /// <param name="a1">First numerator</param> /// <param name="a2">Second numerator</param> /// <param name="b1">First denominator. Requires b1 > 0</param> /// <param name="z">Argument</param> /// <param name="addend">The addend to the series</param> /// <returns> /// The sum of the series, or NaN if the series does not converge /// within the policy maximum number of iterations /// </returns> public static double Sum2F1(double a1, double a2, double b1, double z, double addend = 0) { if ((double.IsNaN(a1) || double.IsInfinity(a1)) || (double.IsNaN(a2) || double.IsInfinity(a2)) || (double.IsNaN(z) || double.IsInfinity(z)) || (double.IsNaN(b1) || double.IsInfinity(b1)) || double.IsNaN(addend)) { Policies.ReportDomainError("Sum2F1(a1: {0}, a2: {1}, b1: {2}, z: {3}, addend: {4}): Requires finite arguments", a1, a2, b1, z, addend); return double.NaN; } if (b1 <= 0 && Math2.IsInteger(b1)) { Policies.ReportDomainError("Sum2F1(a1: {0}, a2: {1}, b1: {2}, z: {3}, addend: {4}): Requires b1 is not zero or a negative integer", a1, a2, b1, z, addend); return double.NaN; } double sum = 1.0 + addend; if (z == 0) return sum; if (a1 == b1) return Sum1F0(a2, z, addend); if (a2 == b1) return Sum1F0(a1, z, addend); if (a1 == 1) return Sum2F1_A1(a2, b1, z, addend); if (a2 == 1) return Sum2F1_A1(a1, b1, z, addend); double term = 1.0; int n = 0; for (; n < Policies.MaxSeriesIterations; n++) { double prevSum = sum; term *= (z / (n + 1)) * ((a1 + n) / (b1 + n)) * (a2 + n); sum += term; #if CMP_FLOAT if (Math.Abs(term) <= Math.Abs(prevSum) * Policies.SeriesTolerance) return sum; #else if (sum == prevSum) return sum; #endif } Policies.ReportConvergenceError("Sum2F1(a1: {0}, a2: {1}, b1: {2}, z: {3}, addend: {4}): No convergence after {5} iterations", a1, a2, b1, z, addend, n); return double.NaN; }
/// <summary> /// Compute Γ(a,x) = Γ(a) - γ(a,x). Use this function when x is small /// </summary> /// <param name="a"></param> /// <param name="x"></param> /// <returns></returns> static double TgammaMinusLowerSeries(double a, double x) { double lowerSeriesSum = LowerSeries(a, x); if (a <= DoubleLimits.MaxGamma) { return(Math2.Tgamma(a) - Prefix(a, x) / a * lowerSeriesSum); } // Since we could overflow too soon, // or get (inf - inf) = NaN, instead of inf // use logs: Log(Γ(a)*(1-P(a,x))) return(Math.Exp(Math2.Lgamma(a) + Math2.Log1p(-PrefixRegularized(a, x) / a * lowerSeriesSum))); }
/// <summary> /// Compute γ(a,x) = Γ(a) - Γ(a,x). Use this function when x is large /// </summary> /// <param name="a"></param> /// <param name="x"></param> /// <returns></returns> static double TgammaMinusUpperFraction(double a, double x) { double upperFraction = UpperFraction(a, x); if (a <= DoubleLimits.MaxGamma) { return(Math2.Tgamma(a) - Prefix(a, x) * upperFraction); } // Since we could overflow too soon, // or get (inf - inf) = NaN, instead of inf // use logs: Log(Γ(a)*(1-Q(a,x))) return(Math.Exp(Math2.Lgamma(a) + Math2.Log1p(-PrefixRegularized(a, x) * upperFraction))); }
// unreferenced, but kept for future reference static double DidonatoFN(double p, double a, double x, uint N, double tolerance) { // // Computation of the Incomplete Gamma Function Ratios and their Inverse // ARMIDO R. DIDONATO and ALFRED H. MORRIS, JR. // ACM Transactions on Mathematical Software, Vol. 12, No. 4, // December 1986, Pages 377-393. // // See equation 34. // double u = Math.Log(p) + Math2.Lgamma(a + 1); return(Math.Exp((u + x - Math.Log(didonato_SN(a, x, N, tolerance))) / a)); }
/// <summary> /// Returns the value of the Associated Laguerre polynomial of degree <paramref name="n"/> and order <paramref name="m"/> at point <paramref name="x"/> /// </summary> /// <param name="n">Polynomial degree. Requires <paramref name="n"/> ≥ 0. Limited to <paramref name="n"/> ≤ 127</param> /// <param name="m">Polynomial order. Requires <paramref name="m"/> ≥ 0</param> /// <param name="x">Polynomial argument</param> public static double LaguerreL(int n, int m, double x) { if ((n < 0) || (m < 0)) { Policies.ReportDomainError("LaguerreL(n: {0}, m: {1}, x: {2}): Requires n,m >= 0", n, m, x); return(double.NaN); } if (n >= 128) { Policies.ReportNotImplementedError("LaguerreL(n: {0}, m: {1}, x: {2}): Not implemented for n >= 128", n, m, x); return(double.NaN); } // Special cases: if (m == 0) { return(Math2.LaguerreL(n, x)); } if (n == 0) { return(1); } if (double.IsNaN(x)) { Policies.ReportDomainError("LaguerreL(n: {0}, m: {1}, x: {2}): NaN not allowed", n, m, x); return(double.NaN); } if (double.IsInfinity(x)) { return((x > 0 && IsOdd(n)) ? double.NegativeInfinity : double.PositiveInfinity); } // use recurrence // p0 - current // p1 - next double p0 = 1; double p1 = m + 1 - x; for (uint k = 1; k < (uint)n; k++) { double next = ((2 * k + (uint)m + 1 - x) * p1 - (k + (uint)m) * p0) / (k + 1); p0 = p1; p1 = next; } return(p1); }
/// <summary> /// Calculates the specified Gamma function when a is an integer /// <para>Q(a,x) = e^-x * (Σ(x^k/k!) k={0,a-1})</para> /// </summary> /// <param name="r"></param> /// <param name="a"></param> /// <param name="x"></param> /// <returns></returns> static double IntegerA(Routine r, double a, double x) { if (r == Routine.Q || r == Routine.Upper) { double q = Math.Exp(-x) * ExponentialSeries((int)(a - 1), x); return((r == Routine.Q) ? q : q *Math2.Tgamma(a)); } else { // avoid this routine for small x, because there's too much cancellation // unsuccessfuly tried: //p = -Math2.Expm1(-x) - Math.Exp(-x)*ExponentialSeries(a - 1, x, -1); double p = 1.0 - Math.Exp(-x) * ExponentialSeries((int)(a - 1), x); return((r == Routine.P) ? p : p *Math2.Tgamma(a)); } }
/// <summary> /// Returns the binomial coefficient /// <para>BinomialCoefficient(n,k) = n! / (k!(n-k)!)</para> /// </summary> /// <param name="n">Requires n ≥ 0</param> /// <param name="k">Requires k in [0,n]</param> public static double BinomialCoefficient(int n, int k) { if ((n < 0) || (k < 0 || k > n)) { Policies.ReportDomainError("BinomialCoefficient(n: {0},k: {1}): Requires n >= 0; k in [0,n]", n, k); return(double.NaN); } if (k == 0 || k == n) { return(1); } if (k == 1 || k == n - 1) { return(n); } double result; if (n < Math2.FactorialTable.Length) { // Use fast table lookup: result = (FactorialTable[n] / FactorialTable[n - k]) / FactorialTable[k]; } else { // Use the beta function: if (k < n - k) { result = k * Math2.Beta(k, n - k + 1); } else { result = (n - k) * Math2.Beta(k + 1, n - k); } if (result == 0) { return(double.PositiveInfinity); } result = 1 / result; } // convert to nearest integer: return(Math.Ceiling(result - 0.5)); }
/// <summary> /// Returns the Inverse Hyperbolic Sine /// <para>Asinh(x) = ln(x + Sqrt(x^2 + 1))</para> /// </summary> /// <param name="x">Asinh argument</param> public static double Asinh(double x) { if (double.IsNaN(x)) { Policies.ReportDomainError("Asinh(x: {0}): NaNs not allowed", x); return(double.NaN); } if (double.IsInfinity(x)) { return(x); } // odd function if (x < 0) { return(-Math2.Asinh(-x)); } // For |x| < eps^(1/4), use a taylor series // Asinh(x) = z - z^3/6 + 3 * z^5/40 ... // http://functions.wolfram.com/ElementaryFunctions/ArcSinh/06/01/03/01/ // // otherwise use varying forms of // Ln(x + Sqrt(x * x + 1)) // http://functions.wolfram.com/ElementaryFunctions/ArcSinh/02/ if (x < 0.75) { if (x < DoubleLimits.RootMachineEpsilon._4) { return(x * (1.0 - x * x / 6)); } // rearrangment of log(x + sqrt(x^2 + 1)) to preserve digits: return(Math2.Log1p(x + Math2.Sqrt1pm1(x * x))); } // for large x, result = Math.Log(2*x), but watch for overflows if (x > 1 / DoubleLimits.RootMachineEpsilon._2) { return(Constants.Ln2 + Math.Log(x)); } return(Math.Log(x + Math.Sqrt(x * x + 1))); }
/// <summary> /// Returns Log(Beta(a,b)) using Stirling approximations /// </summary> /// <param name="a"></param> /// <param name="b"></param> /// <returns></returns> public static double LogBeta(double a, double b) { Debug.Assert(a >= LowerLimit && b >= LowerLimit); double c = a + b; // Beta(a, b) = // (a/c)^a * (b/c)^b * Sqrt(2πc/(a*b)) * (GammaSeries(a)*GammaSeries(b)/GammaSeries(c)) double logValue = 0; double factor = GammaSeries(a) * GammaSeries(b) / GammaSeries(c); factor *= Constants.Sqrt2PI * Math.Sqrt(1 / a + 1 / b); // calculate (a/c)^a or the equivalent(1-b/c)^a // choose the latter method if agh/cgh is sufficiently close to 1 double t1 = b / c; if (t1 < 0.5) { logValue += a * Math2.Log1p(-t1); } else { logValue += a * Math.Log(a / c); } // calculate (b/c)^b or the equivalent (1 - a/c)^b // choose the latter method if b/c is sufficiently close to 1 double t2 = a / c; if (t2 < 0.5) { logValue += b * Math2.Log1p(-t2); } else { logValue += b * Math.Log(b / c); } logValue += Math.Log(factor); return(logValue); }
/// <summary> /// Returns Γ(z)/Γ(z+delta) using Stirling approximations /// </summary> public static double TgammaDeltaRatio(double x, double delta) { Debug.Assert(x >= LowerLimit && delta + x >= LowerLimit); // Compute Γ(x)/Γ(x+Δ), which reduces to: // // (x/(x+Δ))^(x-0.5) * (e/(x+Δ))^Δ * (GammaSeries(x)/GammaSeries(x+Δ)) // OR // (1+Δ/x)^-(x+Δ-0.5) * (e/x)^Δ * (GammaSeries(x)/GammaSeries(x+Δ)) double h = delta / x; double factor = (Math.Abs(delta) < 1) ? Math.Exp(LgammaSeries(x) - LgammaSeries(x + delta)) : GammaSeries(x) / GammaSeries(x + delta); factor *= Math.Sqrt(1 + h); double d = delta / (x + delta); if (Math.Abs(d) <= 0.5) { double l = (x + delta) * Math2.Log1pmx(-d); double deltalogx = -delta *Math.Log(x); if (deltalogx > DoubleLimits.MinLogValue && deltalogx < DoubleLimits.MaxLogValue) { return(factor * Math.Pow(x, -delta) * Math.Exp(l)); } return(factor * Math.Exp(l + deltalogx)); } double result; if (Math.Abs(h) <= 0.5) { result = Math.Exp(-(x + delta) * Math2.Log1p(h)); } else { result = Math.Pow(x / (x + delta), (x + delta)); } result *= factor; result *= Math.Pow(Math.E / x, delta); return(result); }
/// <summary> /// Computes the specified Gamma function when a is small. Uses: /// <para>LowerSmallASeriesSum = Σ((-x)^k/((a+k)*k!)) k={1, inf}</para> /// <para>γ(a,x) = z^a * LowerSmallASeriesSum(a, x, 1/a)</para> /// </summary> /// <param name="r"></param> /// <param name="a"></param> /// <param name="x"></param> /// <returns></returns> static double SmallA(Routine r, double a, double x) { if (r == Routine.P || r == Routine.Lower) { double l = Math.Pow(x, a) * LowerSmallASeries(a, x, 1.0 / a); return((r == Routine.Lower) ? l : l / Math2.Tgamma(a)); } else { double p = Math2.Powm1(x, a); double g = Math2.Tgamma1pm1(a); double initValue = (g - p) / a; double u = initValue - (p + 1) * LowerSmallASeries(a, x, 0.0); return((r == Routine.Upper) ? u : u *a / (g + 1)); } }
/// <summary> /// Returns the Inverse Hyperbolic Cosine /// <para>Acosh(x) = ln(x + Sqrt(x^2 - 1))</para> /// </summary> /// <param name="x">Acosh argument. Requires x ≥ 1</param> public static double Acosh(double x) { if (!(x >= 1)) { Policies.ReportDomainError("Acosh(x:{0}): Requires x >= 1", x); return(double.NaN); } if (double.IsInfinity(x)) { return(double.PositiveInfinity); } // for 1 < x < 1 + sqrt(eps), use a taylor series // Acosh(x) = Sqrt(2 * y) * (1 - y / 12 + 3 * y^2 / 160 ...) where y = x-1 // http://functions.wolfram.com/ElementaryFunctions/ArcCosh/06/01/04/01/ // // otherwise use variations of the standard form // Acosh(x) = ln(x + sqrt(x^2 - 1)) // http://functions.wolfram.com/ElementaryFunctions/ArcCosh/02/ if (x < 1.5) { double y = x - 1; if (y < DoubleLimits.RootMachineEpsilon._2) { const double c0 = 1.0; const double c1 = -1 / 12.0; const double c2 = 3 / 160.0; return(Math.Sqrt(2 * y) * (c0 + y * (c1 + y * c2))); } // This is just a rearrangement of the standard form // devised to minimize loss of precision when x ~ 1: return(Math2.Log1p(y + Math.Sqrt(y * (y + 2)))); } // for large x, result = Math.Log(2*x), but watch for overflows if (x > 1 / DoubleLimits.RootMachineEpsilon._2) { return(Constants.Ln2 + Math.Log(x)); } return(Math.Log(x + Math.Sqrt(x * x - 1))); }
/// <summary> /// Returns Y{n}(z) for positive integer order and z <= sqrt(eps) /// </summary> /// <param name="n"></param> /// <param name="z"></param> /// <returns></returns> static DoubleX YN_SmallArg(int n, double z) { Debug.Assert(n >= 0); Debug.Assert(z <= DoubleLimits.RootMachineEpsilon._2); // // See http://functions.wolfram.com/Bessel-TypeFunctions/BesselY/06/01/04/01/02/ // // Note that when called we assume that x < epsilon and n is a positive integer. // if (n == 0) { return(new DoubleX((2 / Math.PI) * (Math.Log(z / 2) + Constants.EulerMascheroni))); } if (n == 1) { //return (z / Math.PI) * Math.Log(z / 2) - 2 / (Math.PI * z) - (z / (2 * Math.PI)) * (1 - 2 * Constants.EulerMascheroni); return(new DoubleX((z / 2 * (2 * Math.Log(z / 2) - (1 - 2 * Constants.EulerMascheroni)) - 2 / z) / Math.PI)); } if (n == 2) { double z2 = z * z; return(new DoubleX((z2 / 4 * (Math.Log(z / 2) - (3 - 4 * Constants.EulerMascheroni)) - (4 / z2)) / Math.PI)); } // To double check the maximum x value on Mathematica // BesselYDiff[n_, x_] := Block[{ a = BesselY[n, x], b = -(Factorial[n - 1]/Pi)*(x/2)^(-n)}, (b - a)/a] // Table[Block[{eps = 2^-52}, N[BesselYDiff[n, Sqrt[eps]], 20]], {n, 3, 100}] double num = -((Math2.Factorial(n - 1) / Math.PI)); double den = Math.Pow(z / 2, n); if (double.IsInfinity(num) || den == 0) { return(DoubleX.NegativeInfinity); } DoubleX result = new DoubleX(num) / new DoubleX(den); return(result); }
/// <summary> /// Returns the value of the Legendre polynomial that is the second solution to the Legendre differential equation /// </summary> /// <param name="n">Polynomial degree. Requires n ≥ 0. Limited to n ≤ 127</param> /// <param name="x">Requires |x| ≤ 1</param> public static double LegendreQ(int n, double x) { if ((n < 0) || (!(x >= -1 && x <= 1))) { Policies.ReportDomainError("LegendreQ(n: {0}, x: {1}): Requires n >= 0; |x| <= 1", n, x); return(double.NaN); } if (x == 1) { return(double.PositiveInfinity); } if (x == -1) { return(IsOdd(n) ? double.PositiveInfinity : double.NegativeInfinity); } if (n >= 128) { Policies.ReportNotImplementedError("LegendreQ(n: {0}, x: {1}): Not implemented for n >= 128", n, x); return(double.NaN); } // Implement Legendre Q polynomials via recurrance // p0 - current // p1 - next double p0 = Math2.Atanh(x); double p1 = x * p0 - 1; if (n == 0) { return(p0); } for (uint k = 1; k < (uint)n; k++) { double next = ((2 * k + 1) * x * p1 - k * p0) / (k + 1); p0 = p1; p1 = next; } return(p1); }
/// <summary> /// Compute (z^a)(e^-z)/Γ(a) used in the regularized incomplete gammas /// </summary> /// <param name="a"></param> /// <param name="z"></param> /// <returns></returns> /// <remarks>Note: most of the error occurs in this function</remarks> public static double PrefixRegularized(double a, double z) { Debug.Assert(a > 0); Debug.Assert(z >= 0); // if we can't compute Tgamma directly, use Stirling if (a > DoubleLimits.MaxGamma || (a >= StirlingGamma.LowerLimit && z > -DoubleLimits.MinLogValue)) { return(StirlingGamma.PowExpDivGamma(a, z)); } // avoid overflows from dividing small Gamma(a) // for z > 1, e^(eps*ln(z)-z) ~= e^-z, so underflow is ok if (a < DoubleLimits.MachineEpsilon / Constants.EulerMascheroni) { return(a * Math.Pow(z, a) * Math.Exp(-z)); } double alz = a * Math.Log(z); double g = Math2.Tgamma(a); double prefix; if ((alz > DoubleLimits.MinLogValue && alz < DoubleLimits.MaxLogValue) && (z < -DoubleLimits.MinLogValue)) { prefix = Math.Pow(z, a) * Math.Exp(-z) / g; } else if ((alz > 2 * DoubleLimits.MinLogValue && alz < 2 * DoubleLimits.MaxLogValue) && (z < -2 * DoubleLimits.MinLogValue)) { double rp = Math.Pow(z, a / 2) * Math.Exp(-z / 2); prefix = (rp / g) * rp; } else if ((alz > 4 * DoubleLimits.MinLogValue && alz < 4 * DoubleLimits.MaxLogValue) && (z < -4 * DoubleLimits.MinLogValue)) { double rp = Math.Pow(z, a / 4) * Math.Exp(-z / 4); prefix = (rp / g) * rp * rp * rp; } else { prefix = Math.Exp(alz - z - Math.Log(g)); } return(prefix); }
/// <summary> /// Returns Sqrt(x+1)-1 with improved accuracy for |x| ≤ 0.75 /// </summary> /// <param name="x">The function argument</param> public static double Sqrt1pm1(double x) { if (!(x >= -1)) { Policies.ReportDomainError("Sqrt1pm1(x: {0}): Requires x >= -1", x); return(double.NaN); } if (double.IsInfinity(x)) { return(x); } if (Math.Abs(x) <= 0.75) { return(Math2.Expm1(Math2.Log1p(x) / 2)); } return(Math.Sqrt(1 + x) - 1); }
/// <summary> /// Computes the Jacobi elliptic functions using the arithmetic-geometric mean /// </summary> /// <param name="k">The modulus</param> /// <param name="u">The argument</param> /// <param name="anm1"></param> /// <param name="bnm1"></param> /// <param name="N"></param> /// <param name="pTn"></param> /// <returns></returns> /// <see href="http://dlmf.nist.gov/22.20#ii"/> private static double Jacobi_Recurse(double k, double u, double anm1, double bnm1, int N, out double pTn) { ++N; double Tn; double cn = (anm1 - bnm1) / 2; double an = (anm1 + bnm1) / 2; if (cn < DoubleLimits.MachineEpsilon) { Tn = Math2.Ldexp(1, N) * u * an; } else { Tn = Jacobi_Recurse(k, u, an, Math.Sqrt(anm1 * bnm1), N, out double unused); } pTn = Tn; return((Tn + Math.Asin((cn / an) * Math.Sin(Tn))) / 2); }
public static (double Result, int Iterations) RefineQ_SmallA_LargeZ(double z, double p, double q, double aGuess) { Debug.Assert(z >= 0, "Requires z >= 0"); Debug.Assert(aGuess >= 0 && aGuess <= double.MaxValue, "Requires finite aGuess >= 0"); Debug.Assert(q <= 0.5, "Requires q <= 0.5"); Debug.Assert(aGuess < z + 1); // f = Log(z^a/Gamma(a))-Log(q); // f/f' = (Log(z^a/Gamma(a))-Log(q))/ (Log(z) - Polygamma(0, a)) // // http://functions.wolfram.com/GammaBetaErf/GammaRegularized/06/02/02/0001/ // For z->Infinity, Q(a,z) ~= z^(a-1)*e^-z/Gamma(a) * (1 - (1-a)/z ...) // Find a quick estimate of a // where Log[q] = Log[z^(a-1)*e^-z/Gamma(a)] to get into the ballpark double lz = Math.Log(z); double lq = Math.Log(q); int i = 0; for (; i < 1; i++) { const double eps = 1.0 / 2; double f = (aGuess - 1) * lz - z - Math2.Lgamma(aGuess) - lq; if (f == 0) { break; } var delta = f / (lz - Math2.Digamma(aGuess)); double aLast = aGuess; aGuess -= delta; if (Math.Abs(delta) <= eps * Math.Abs(aLast)) { break; } } return(aGuess, i); }
/// <summary> /// Returns the next representable value after x in the direction of y. /// </summary> /// <param name="x">the value to get the next of</param> /// <param name="y">direction</param> /// <returns> /// The next floating point value in the direction of y /// <para>If x == NaN || x == Infinity, returns NaN</para> /// <para>If y == NaN, returns NaN</para> /// <para>If y > x, returns FloatNext(x)</para> /// <para>If y = x, returns x</para> /// <para>If y < x, returns FloatPrior(x)</para> /// </returns> public static double Nextafter(double x, double y) { if ((double.IsNaN(x) || double.IsInfinity(x)) || (double.IsNaN(y))) { Policies.ReportDomainError("Nextafter(x: {0}, y: {1}): Requires finite x, y not NaN", x, y); return(double.NaN); } if (y > x) { return(Math2.FloatNext(x)); } if (y == x) { return(x); } return(Math2.FloatPrior(x)); }
/// <summary> /// Returns Γ(z)/Γ(z+delta) using the Lanczos approximation /// </summary> public static double TgammaDeltaRatio(double x, double delta) { Debug.Assert(x > 0 && delta + x > 0); // Compute Γ(x)/Γ(x+Δ), which reduces to: // // ((x+g-0.5)/(x+Δ+g-0.5))^(x-0.5) * (x+Δ+g-0.5)^-Δ * e^Δ * (Sum(x)/Sum(x+Δ)) // = (1+Δ/(x+g-0.5))^-(x-0.5) * (x+Δ+g-0.5)^-Δ * e^Δ * (Sum(x)/Sum(x+Δ)) double xgh = x + (Lanczos.G - 0.5); double xdgh = x + delta + (Lanczos.G - 0.5); double mu = delta / xgh; double result = Lanczos.Series(x) / Lanczos.Series(x + delta); result *= (Math.Abs(mu) < 0.5) ? Math.Exp((0.5 - x) * Math2.Log1p(mu)) : Math.Pow(xgh / xdgh, x - 0.5); result *= Math.Pow(Math.E / xdgh, delta); return(result); }
/// <summary> /// Returns Γ(x + 1) - 1 with accuracy improvements in [-0.5, 2] /// </summary> /// <param name="x">Tgamma1pm1 function argument</param> public static double Tgamma1pm1(double x) { if (double.IsNaN(x)) { Policies.ReportDomainError("Tgamma1pm1(x: {0}): NaN not allowed", x); return(double.NaN); } if (double.IsInfinity(x)) { if (x < 0) { Policies.ReportDomainError("Tgamma1pm1(x: {0}): Requires x != -Infinity", x); return(double.NaN); } return(double.PositiveInfinity); } double result; if (x < -0.5 || x >= 2) { return(Math2.Tgamma(1 + x) - 1); // Best method is simply to subtract 1 from tgamma: } if (Math.Abs(x) < GammaSmall.UpperLimit) { return(GammaSmall.Tgamma1pm1(x)); } if (x < 0) { // compute exp( log( Γ(x+2)/(1+x) ) ) - 1 result = Math2.Expm1(-Math2.Log1p(x) + _Lgamma.Rational_1_3(x + 2, x + 1, x)); } else { // compute exp( log(Γ(x+1)) ) - 1 result = Math2.Expm1(_Lgamma.Rational_1_3(x + 1, x, x - 1)); } return(result); }
/// <summary> /// Series evaluation for Jv(x), for small x /// </summary> /// <param name="v"></param> /// <param name="x"></param> /// <returns></returns> public static double J_SmallArg(double v, double x) { Debug.Assert(v >= 0); // See http://functions.wolfram.com/Bessel-TypeFunctions/BesselJ/06/01/04/01/01/0003/ // Converges rapidly for all x << v. // Most of the error occurs calculating the prefix double prefix; if (v <= DoubleLimits.MaxGamma - 1) prefix = Math.Pow(x / 2, v) / Math2.Tgamma(v + 1); else prefix = StirlingGamma.PowDivGamma(x / 2, v)/v; if (prefix == 0) return prefix; double series = HypergeometricSeries.Sum0F1(v + 1, -x * x / 4); return prefix * series; }
/// <summary> /// Returns the double factorial = n!! /// <para>If n is odd, returns: n * (n-2) * ... * 5 * 3 * 1</para> /// <para>If n is even, returns: n * (n-2) * ... * 6 * 4 * 2</para> /// <para>If n = 0, -1, returns 1</para> /// </summary> /// <param name="n">Requires n ≥ -1</param> public static double Factorial2(int n) { if (n < -1) { Policies.ReportDomainError("Factorial2(n: {0}): Requires n >= -1", n); return(double.NaN); } // common values if (n == -1 || n == 0) { return(1); } double result; if (IsOdd(n)) { // odd n: if (n < FactorialTable.Length) { int k = (n - 1) / 2; return(Math.Ceiling(Math2.Ldexp(FactorialTable[n] / FactorialTable[k], -k) - 0.5)); } // // Fallthrough: n is too large to use table lookup, try the // gamma function instead. // result = Constants.RecipSqrtPI * Math2.Tgamma(0.5 * n + 1.0); result = Math.Ceiling(Math2.Ldexp(result, (n + 1) / 2) - 0.5); } else { // even n: int k = n / 2; result = Math2.Ldexp(Factorial(k), k); } return(result); }
/// <summary> /// Returns the Inverse Hyperbolic Tangent /// <para>Atanh(x) = 1/2 * ln((1+x)/(1-x))</para> /// </summary> /// <param name="x">Atanh argument. Requires |x| <= 1</param> public static double Atanh(double x) { if (!(x >= -1 && x <= 1)) { Policies.ReportDomainError("Atanh(x: {0}): Requires |x| <= 1", x); return(double.NaN); } if (x == -1) { return(double.NegativeInfinity); } if (x == 1) { return(double.PositiveInfinity); } // For |x| < eps^(1/4), use a taylor series // Atanh(x) = z + z^3/3 + z^5/5 ... // http://functions.wolfram.com/ElementaryFunctions/ArcTanh/06/01/03/01/ // // for eps^(1/4) <= |x| < 0.5 // Atanh(x) = 1/2 * ( ln(1+x) - ln(1-x) ) // http://functions.wolfram.com/ElementaryFunctions/ArcTanh/02/ // // otherwise use // Atanh(x) = 1/2 * ln((1+x)/(1-x)) double absX = Math.Abs(x); if (absX < 0.5) { if (absX < DoubleLimits.RootMachineEpsilon._4) { return(x * (1.0 + x * x / 3)); } return((Math2.Log1p(x) - Math2.Log1p(-x)) / 2); } return(Math.Log((1 + x) / (1 - x)) / 2); }
public static (double Result, int Iterations) SolveA(double z, double p, double q, double aGuess, double step, double min, double max) { Debug.Assert(aGuess >= min && aGuess <= max, "Guess out of range"); Func <double, double> f; if (p < q) { f = a => Math2.GammaP(a, z) - p; } else { f = a => Math2.GammaQ(a, z) - q; } // // Use our generic derivative-free root finding procedure. // We could use Newton steps here, taking the PDF of the // Poisson distribution as our derivative, but that's // even worse performance-wise than the generic method :-( // // dQ/da is increasing, dP/da is decreasing FunctionShape fShape = (p < q) ? FunctionShape.Decreasing : FunctionShape.Increasing; RootResults rr = RootFinder.Toms748Bracket(f, aGuess, step, fShape, min, max); if (rr == null) { Policies.ReportRootNotFoundError($"Invalid Parameter: PQInvA(z: {z}, p: {p}, q: {q}) [{min}, {max}] g: {aGuess}"); return(double.NaN, 0); } if (!rr.Success) { Policies.ReportRootNotFoundError("Root not found after {0} iterations", rr.Iterations); return(double.NaN, rr.Iterations); } return(rr.SolutionX, rr.Iterations); }
/// <summary> /// Returns the Log(Beta(a,b)) using the Lanczos approximation /// </summary> /// <param name="a"></param> /// <param name="b"></param> /// <returns></returns> public static double LogBeta(double a, double b) { // B(a,b) == B(b, a) if (a < b) { Utility.Swap(ref a, ref b); } // from this point a >= b // Lanczos calculation. Compute: // ln((agh/cgh)^a * (bgh/cgh)^b * Sqrt((cgh * e)/(agh * bgh))) // = ln((1-b/cgh)^a * (1 - a/cgh)^b * Sqrt((cgh * e)/(agh * bgh))) double c = a + b; double agh = a + (Lanczos.G - 0.5); double bgh = b + (Lanczos.G - 0.5); double cgh = c + (Lanczos.G - 0.5); double logValue = 0; double factor = (SeriesExpGScaled(a) / SeriesExpGScaled(c)) * SeriesExpGScaled(b); factor *= Constants.SqrtE * Math.Sqrt((cgh / agh) / bgh); // calculate a*ln(agh/cgh) or the equivalent a*ln(1 - b/cgh) // calculate b*ln(bgh/cgh) or the equivalent b*ln(1 - a/cgh) double t1 = b / cgh; logValue += (t1 <= 0.5) ? a * Math2.Log1p(-t1) : a *Math.Log(agh / cgh); double t2 = a / cgh; logValue += (t2 <= 0.5) ? b * Math2.Log1p(-t2) : b *Math.Log(bgh / cgh); logValue += Math.Log(factor); return(logValue); }
/// <summary> /// Returns the Incomplete Beta Complement /// <para>Betac(a, b, x) = B(a,b) - B<sub>x</sub>(a,b)</para> /// </summary> /// <param name="a">Requires a ≥ 0 and not (a == b == 0)</param> /// <param name="b">Requires b ≥ 0 and not (a == b == 0)</param> /// <param name="x">0 ≤ x ≤ 1</param> public static double Betac(double a, double b, double x) { if ((!(a > 0) || double.IsInfinity(a)) || (!(b > 0) || double.IsInfinity(b)) || !(x >= 0 && x <= 1)) { Policies.ReportDomainError("Betac(a: {0}, b: {1}, x: {2}): Requires finite a,b > 0; x in [0,1]", a, b, x); return(double.NaN); } // special values if (x == 0) { return(Math2.Beta(a, b)); } if (x == 1) { return(0); } return(_Ibeta.BetaImp(a, b, x, true)); }
static double InverseNegativeBinomialCornishFisher(double n, double sf, double sfc, double p, double q) { // mean: double m = n * (sfc) / sf; double t = Math.Sqrt(n * (sfc)); // standard deviation: double sigma = t / sf; // skewness double sk = (1 + sfc) / t; // kurtosis: double k = (6 - sf * (5 + sfc)) / (n * (sfc)); // Get the inverse of a std normal distribution: double x = Math2.ErfcInv(p > q ? 2 * q : 2 * p) * Constants.Sqrt2; // Set the sign: if (p < 0.5) { x = -x; } double x2 = x * x; // w is correction term due to skewness double w = x + sk * (x2 - 1) / 6; // // Add on correction due to kurtosis. // if (n >= 10) { w += k * x * (x2 - 3) / 24 + sk * sk * x * (2 * x2 - 5) / -36; } w = m + sigma * w; if (w < DoubleLimits.MinNormalValue) { return(DoubleLimits.MinNormalValue); } return(w); }
/// <summary> /// Returns the complete elliptic integral of the second kind E(k) = E(π/2,k) /// <para>E(k) = ∫ sqrt(1-k^2*sin^2(θ)) dθ, θ={0,π/2}</para> /// </summary> /// <param name="k">The modulus. Requires |k| ≤ 1</param> public static double EllintE(double k) { if (!(k >= -1 && k <= 1)) { Policies.ReportDomainError("EllintE(k: {0}): Requires |k| <= 1", k); return(double.NaN); } // Small series at k == 0 // Note that z = k^2 in the following link // http://functions.wolfram.com/EllipticIntegrals/EllipticE/06/01/03/01/0004/ double k2 = k * k; if (k2 < DoubleLimits.RootMachineEpsilon._13) { if (k2 < 4 * DoubleLimits.MachineEpsilon) { return(Math.PI / 2); // E(0) = Pi/2 } return((Math.PI / 2) * HypergeometricSeries.Sum2F1(-0.5, 0.5, 1, k2)); } if (Math.Abs(k) == 1) { return(1); } double x = 0; double y = 1 - k2; double z = 1; double value = 2 * Math2.EllintRG(x, y, z); return(value); }
/// <summary> /// Returns Ln(|Γ(x)|) for |x| >= <c>LowerLimit</c> using the Stirling asymptotic series /// </summary> /// <param name="x">The argument</param> /// <param name="sign">Sets sign = Sign(Γ(x))</param> /// <returns></returns> public static double Lgamma(double x, out int sign) { Debug.Assert(Math.Abs(x) >= LowerLimit); // (Ln(2*Pi)-1)/2 = Ln(Sqrt(2*Pi/E)) const double Half_Ln2PIm1 = 0.418938533204672741780329736405617639861397473637783412817151; sign = 0; if (x > 0) { //return (x - 0.5) * Math.Log(x) - x + Constants.Ln2PI / 2 + LgammaSeries(x); sign = 1; return((x - 0.5) * (Math.Log(x) - 1) + (Half_Ln2PIm1 + LgammaSeries(x))); } // Use the reflection formula: Γ(-x) = -π/(x * Γ(x) * Sin(π * x)) x = -x; double sin = Math2.SinPI(x); if (sin == 0) { return(double.NaN); } if (sin < 0) { sin = -sin; sign = 1; } else { sign = -1; } sin *= (x / Math.PI); return(-Math.Log(sin) - Lgamma(x)); }