include!(concat!(env!("OUT_DIR"), "/simd_lanes.rs"));
use std::{
f64::consts::PI,
simd::{Simd, StdFloat},
};
use crate::kernels::scientific::{
distributions::shared::constants::*,
erf::{erfc, erfc_inv},
};
#[inline(always)]
pub fn ln_gamma(x: f64) -> f64 {
if x.is_nan() {
return f64::NAN;
}
if x.is_infinite() && x.is_sign_positive() {
return f64::INFINITY;
}
if x <= 0.0 && (x.fract().abs() < 1e-14) {
return f64::INFINITY;
}
if x < 0.5 {
return std::f64::consts::PI.ln()
- (std::f64::consts::PI * x).sin().abs().ln()
- ln_gamma(1.0 - x);
}
let z = x - 1.0; let mut a = COF[0];
for (i, &c) in COF.iter().enumerate().skip(1) {
a += c / (z + i as f64);
}
let t = z + 7.5; HALF_LOG_TWO_PI + (z + 0.5) * t.ln() - t + a.ln()
}
#[inline(always)]
pub fn ln_gamma_plus1(k: f64) -> f64 {
ln_gamma(k + 1.0)
}
#[inline(always)]
pub fn ln_gamma_simd<const N: usize>(x: Simd<f64, N>) -> Simd<f64, N>
where
{
let z = x - Simd::splat(1.0); let mut a = Simd::splat(COF[0]); for (i, &c) in COF.iter().enumerate().skip(1) {
a += Simd::splat(c) / (z + Simd::splat(i as f64));
}
let t = z + Simd::splat(7.5); let half_ln_two_pi = Simd::splat(0.9189385332046727_f64); half_ln_two_pi + (z + Simd::splat(0.5)) * t.ln() - t + a.ln()
}
#[inline(always)]
pub fn inv_reg_lower_gamma(a: f64, p: f64) -> f64 {
if !(a.is_finite() && p.is_finite()) || a <= 0.0 {
return f64::NAN;
}
if p <= 0.0 {
return 0.0;
}
if p >= 1.0 {
return f64::INFINITY;
}
let mut x = if p < 1e-8 {
(p * gamma_func(a + 1.0)).powf(1.0 / a).max(1e-300)
} else if a > 1.0 {
let t = (-2.0 * p.ln()).sqrt();
let s = (t - (2.0 * a - 1.0).sqrt()) / (2.0 * (a - 1.0).sqrt());
(a * (1.0 - s + (s * s) / 3.0)).max(1e-300)
} else {
(p * gamma_func(a + 1.0)).powf(1.0 / a).max(1e-300)
};
const MAX_ITERS: usize = 80;
const REL_TOL: f64 = 1e-15;
const ABS_FLOOR: f64 = 5e-16;
for _ in 0..MAX_ITERS {
let f = reg_lower_gamma(a, x) - p; let fp = gamma_pdf_scalar(a, x); if !fp.is_finite() || fp.abs() < 1e-300 {
break;
}
let mut dx = f / fp;
let max_step = x.max(1.0); if dx.abs() > max_step {
dx = max_step * dx.signum();
}
let x_new = (x - dx).max(1e-300);
let tol = (REL_TOL * (1.0 + x_new.abs())).max(ABS_FLOOR);
if (x_new - x).abs() <= tol {
x = x_new;
break;
}
x = x_new;
}
let f = reg_lower_gamma(a, x) - p;
let fp = gamma_pdf_scalar(a, x);
if fp.is_finite() && fp.abs() > 0.0 && x.is_finite() && x > 0.0 {
let fpp = fp * ((a - 1.0) / x - 1.0); let h = f / fp;
let denom = 1.0 - 0.5 * h * (fpp / fp);
if denom.is_finite() && denom != 0.0 {
let xh = (x - h / denom).max(1e-300);
if (reg_lower_gamma(a, xh) - p).abs() <= f.abs() {
x = xh;
}
}
}
x
}
#[inline(always)]
pub fn inv_reg_upper_gamma(a: f64, q: f64) -> f64 {
if !(a.is_finite() && q.is_finite()) || a <= 0.0 {
return f64::NAN;
}
if q <= 0.0 {
return f64::INFINITY;
}
if q >= 1.0 {
return 0.0;
}
let p = 1.0 - q;
let mut x = if q > 1e-8 && a > 1.0 {
let t = (-2.0 * q.ln()).sqrt();
let s = (t - (2.0 * a - 1.0).sqrt()) / (2.0 * (a - 1.0).sqrt());
(a * (1.0 - s + (s * s) / 3.0)).max(1e-300)
} else {
(p * gamma_func(a + 1.0)).powf(1.0 / a).max(1e-300)
};
const MAX_ITERS: usize = 80;
const REL_TOL: f64 = 1e-15;
const ABS_FLOOR: f64 = 5e-16;
for _ in 0..MAX_ITERS {
let f = (1.0 - reg_lower_gamma(a, x)) - q; let fp = -gamma_pdf_scalar(a, x); if !fp.is_finite() || fp.abs() < 1e-300 {
break;
}
let mut dx = f / fp;
let max_step = x.max(1.0);
if dx.abs() > max_step {
dx = max_step * dx.signum();
}
let x_new = (x - dx).max(1e-300);
let tol = (REL_TOL * (1.0 + x_new.abs())).max(ABS_FLOOR);
if (x_new - x).abs() <= tol {
x = x_new;
break;
}
x = x_new;
}
let f = (1.0 - reg_lower_gamma(a, x)) - q;
let fp = -gamma_pdf_scalar(a, x);
if fp.is_finite() && fp.abs() > 0.0 && x.is_finite() && x > 0.0 {
let fpp = -fp * ((a - 1.0) / x - 1.0); let h = f / fp;
let denom = 1.0 - 0.5 * h * (fpp / fp);
if denom.is_finite() && denom != 0.0 {
let xh = (x - h / denom).max(1e-300);
if ((1.0 - reg_lower_gamma(a, xh)) - q).abs() <= f.abs() {
x = xh;
}
}
}
x
}
#[inline(always)]
pub fn reg_lower_gamma(a: f64, x: f64) -> f64 {
if !(a.is_finite() && x.is_finite()) {
return f64::NAN;
}
if x < 0.0 {
return f64::NAN;
}
if a < 0.0 {
return f64::NAN;
}
if a == 0.0 {
return 1.0;
} if x == 0.0 {
return 0.0;
}
if x <= 0.0 {
return 0.0;
}
if a <= 0.0 {
return 0.0;
}
if x < a + 1.0 {
let mut ap = a;
let mut sum = 1.0 / a;
let mut del = sum;
for _ in 0..100 {
ap += 1.0;
del *= x / ap;
sum += del;
if del.abs() < sum.abs() * 1e-14 {
break;
}
}
(-(x) + a * x.ln() - ln_gamma(a)).exp() * sum
} else {
let mut b = x + 1.0 - a;
let mut c = 1.0 / f64::MIN_POSITIVE;
let mut d = 1.0 / b;
let mut h = d;
for i in 1..100 {
let an = -i as f64 * (i as f64 - a);
b += 2.0;
d = an * d + b;
if d.abs() < 1e-30 {
d = 1e-30
}
c = b + an / c;
if c.abs() < 1e-30 {
c = 1e-30
}
d = 1.0 / d;
let delta = d * c;
h *= delta;
if (delta - 1.0).abs() < 1e-14 {
break;
}
}
1.0 - (-(x) + a * x.ln() - ln_gamma(a)).exp() * h
}
}
#[inline(always)]
pub fn gamma_pdf_scalar(a: f64, x: f64) -> f64 {
if x <= 0.0 {
return 0.0;
}
((a - 1.0) * x.ln() - x - ln_gamma(a)).exp()
}
#[inline(always)]
pub fn gamma_func(x: f64) -> f64 {
if x.is_nan() {
return f64::NAN;
}
if x == 0.0 {
return f64::INFINITY;
}
if x < 0.0 && (x.fract().abs() < 1e-14) {
return f64::NAN;
}
if x > 0.0 {
if x.fract().abs() < 1e-14 {
let n = x as usize;
if n >= 1 && n <= 171 {
return factorial_lookup((n - 1) as u64);
}
}
if (x - 0.5).fract().abs() < 1e-14 {
let n = (x - 0.5).round() as u64; if n <= 170 {
return half_integer_gamma(n);
}
}
return ln_gamma(x).exp();
}
let sin_pi_x = (std::f64::consts::PI * x).sin();
if sin_pi_x.abs() < f64::EPSILON {
return f64::NAN; }
let gamma_1_minus_x = gamma_func(1.0 - x);
std::f64::consts::PI / (sin_pi_x * gamma_1_minus_x)
}
#[inline(always)]
pub fn half_integer_gamma(n: u64) -> f64 {
let two_n_fact = factorial_lookup(2 * n);
let n_fact = factorial_lookup(n);
let pow4n = 2f64.powi(2 * n as i32); (PI.sqrt() * two_n_fact) / (pow4n * n_fact)
}
#[inline(always)]
pub fn incomplete_beta(a: f64, b: f64, x: f64) -> f64 {
if !(a.is_finite() && b.is_finite() && x.is_finite()) {
return f64::NAN;
}
if x <= 0.0 {
return 0.0;
}
if x >= 1.0 {
return 1.0;
}
if a == 0.0 {
return 1.0; }
if b == 0.0 {
return 0.0; }
if x > (a + 1.0) / (a + b + 2.0) {
return 1.0 - incomplete_beta(b, a, 1.0 - x);
};
let ln_beta = ln_gamma(a) + ln_gamma(b) - ln_gamma(a + b);
let front = (a * x.ln() + b * (1.0 - x).ln() - ln_beta).exp() / a;
const EPS: f64 = 1e-14;
const FPMIN: f64 = 1e-300;
const MAX_ITS: usize = 200;
let mut c = 1.0;
let mut d = 1.0 - (a + b) * x / (a + 1.0);
if d.abs() < FPMIN {
d = FPMIN;
}
d = 1.0 / d;
let mut h = d;
for m in 1..=MAX_ITS {
let m2 = 2 * m;
let aa = m as f64 * (b - m as f64) * x / ((a + m2 as f64 - 1.0) * (a + m2 as f64));
d = 1.0 + aa * d;
if d.abs() < FPMIN {
d = FPMIN;
}
c = 1.0 + aa / c;
if c.abs() < FPMIN {
c = FPMIN;
}
d = 1.0 / d;
h *= d * c;
let aa =
-(a + m as f64) * (a + b + m as f64) * x / ((a + m2 as f64) * (a + m2 as f64 + 1.0));
d = 1.0 + aa * d;
if d.abs() < FPMIN {
d = FPMIN;
}
c = 1.0 + aa / c;
if c.abs() < FPMIN {
c = FPMIN;
}
d = 1.0 / d;
let delta = d * c;
h *= delta;
if (delta - 1.0).abs() < EPS {
break;
}
}
front * h
}
#[inline(always)]
pub fn incomplete_beta_inv(a: f64, b: f64, p: f64) -> f64 {
if !(a.is_finite() && b.is_finite() && p.is_finite()) {
return f64::NAN;
}
if p < 0.0 || p > 1.0 {
return f64::NAN;
}
if p == 0.0 {
return 0.0;
}
if p == 1.0 {
return 1.0;
}
if a == 0.0 {
return 1.0;
} if b == 0.0 {
return 0.0;
}
let pp = if p < 0.5 { p } else { 1.0 - p };
let t = (-2.0 * pp.ln()).sqrt();
let mut x: f64;
if a > 1.0 && b > 1.0 {
let num = 2.30753 + 0.27061 * t;
let den = 1.0 + (0.99229 + 0.04481 * t) * t;
let mut xp = t - num / den;
if p < 0.5 {
xp = -xp;
}
let al = (xp * xp - 3.0) / 6.0;
let h = 2.0 / (1.0 / (2.0 * a - 1.0) + 1.0 / (2.0 * b - 1.0));
let w = xp * (al + h).sqrt() / h
- (1.0 / (2.0 * b - 1.0) - 1.0 / (2.0 * a - 1.0)) * (al + 5.0 / 6.0 - 2.0 / (3.0 * h));
x = a / (a + b * w.exp());
} else {
let r = a / (a + b);
let y = if p < r { p / r } else { (1.0 - p) / (1.0 - r) };
x = if p < r {
y.powf(1.0 / a)
} else {
1.0 - y.powf(1.0 / b)
};
}
const EPS_X: f64 = 1e-15;
if x <= 0.0 {
x = EPS_X;
}
if x >= 1.0 {
x = 1.0 - EPS_X;
}
const NEWTON_EPS: f64 = 1e-14;
const MAX_NEWTON: usize = 20;
let ln_beta = ln_gamma(a) + ln_gamma(b) - ln_gamma(a + b);
for _ in 0..MAX_NEWTON {
let f = incomplete_beta(a, b, x) - p;
let df = ((a - 1.0) * x.ln() + (b - 1.0) * (1.0 - x).ln() - ln_beta).exp();
if df == 0.0 || !df.is_finite() {
break;
}
let dx = f / df;
let mut x_new = x - dx;
if x_new <= 0.0 {
x_new = 0.5 * x;
}
if x_new >= 1.0 {
x_new = 0.5 * (1.0 + x);
}
if (dx).abs() < NEWTON_EPS * x_new.max(EPS_X) {
return x_new;
}
x = x_new;
}
let mut lo = 0.0;
let mut hi = 1.0;
if incomplete_beta(a, b, x) > p {
hi = x;
} else {
lo = x;
}
for _ in 0..64 {
let mid = 0.5 * (lo + hi);
let f_mid = incomplete_beta(a, b, mid);
if (f_mid - p).abs() < 1e-14 {
return mid;
}
if f_mid < p {
lo = mid;
} else {
hi = mid;
}
}
0.5 * (lo + hi)
}
#[inline(always)]
pub fn regularised_gamma_p(s: f64, x: f64) -> f64 {
if x <= 0.0 {
return 0.0;
}
if x < s + 1.0 {
let mut sum = 1.0 / s;
let mut value = sum;
let mut n = 1.0;
while value.abs() > 1e-15 * sum {
value *= x / (s + n);
sum += value;
n += 1.0;
if n > 200.0 {
break;
}
}
(-(x) + s * x.ln() - ln_gamma(s)).exp() * sum
} else {
let mut a = 1.0 - s;
let mut b = a + x + 1.0;
let mut c = 0.0;
let mut p1 = 1.0;
let mut q1 = x;
let mut p2 = x + 1.0;
let mut q2 = b * x;
let mut ans = p2 / q2;
let mut n = 1.0;
while (p2 - p1).abs() > 1e-15 * ans.abs() {
a += 1.0;
b += 2.0;
c += 1.0;
let an = a * c;
let p3 = b * p2 - an * p1;
let q3 = b * q2 - an * q1;
if q3.abs() < 1e-30 {
break;
}
p1 = p2;
q1 = q2;
p2 = p3;
q2 = q3;
if q2 != 0.0 {
ans = p2 / q2;
}
n += 1.0;
if n > 200.0 {
break;
}
}
1.0 - (-(x) + s * x.ln() - ln_gamma(s)).exp() * ans
}
}
#[inline(always)]
pub fn ln_choose_v(
n: core::simd::Simd<f64, W64>,
k: core::simd::Simd<f64, W64>,
) -> core::simd::Simd<f64, W64> {
ln_gamma_simd(n + core::simd::Simd::<f64, W64>::splat(1.0))
- ln_gamma_simd(k + core::simd::Simd::<f64, W64>::splat(1.0))
- ln_gamma_simd(n - k + core::simd::Simd::<f64, W64>::splat(1.0))
}
#[inline(always)]
pub fn ln_choose(n: u64, k: u64) -> f64 {
if k > n {
return f64::NEG_INFINITY;
} ln_gamma((n + 1) as f64) - ln_gamma((k + 1) as f64) - ln_gamma((n - k + 1) as f64)
}
#[inline(always)]
pub fn ln_choose_simd<const N: usize>(n: Simd<f64, N>, k: Simd<f64, N>) -> Simd<f64, N>
where
{
ln_gamma_simd(n + Simd::splat(1.0))
- ln_gamma_simd(k + Simd::splat(1.0))
- ln_gamma_simd(n - k + Simd::splat(1.0))
}
#[inline(always)]
pub fn binomial_quantile_cornish_fisher(pi: f64, n: u64, p_: f64) -> f64 {
if !pi.is_finite() || !p_.is_finite() || n > (i64::MAX as u64) {
return f64::NAN;
}
if pi <= 0.0 {
return -1.0;
}
if pi >= 1.0 {
return n as f64;
}
let mu = n as f64 * p_;
let sigma = (n as f64 * p_ * (1.0 - p_)).sqrt();
let z = normal_quantile_scalar(pi, 0.0, 1.0);
let skew = (1.0 - 2.0 * p_) / sigma;
let cf = z + skew * (z * z - 1.0) / 6.0;
let mut k = (mu + sigma * cf + 0.5).floor();
if k < 0.0 {
k = 0.0;
}
if k > n as f64 {
k = n as f64;
}
let mut k_int = k as i64;
let cdf = binomial_cdf_scalar(k_int, n, p_);
if cdf < pi {
while k_int < (n as i64) && binomial_cdf_scalar(k_int, n, p_) < pi {
k_int += 1;
}
} else {
while k_int > 0 && binomial_cdf_scalar(k_int - 1, n, p_) >= pi {
k_int -= 1;
}
}
k_int as f64
}
#[inline(always)]
pub fn binomial_cdf_scalar(k: i64, n: u64, p: f64) -> f64 {
if k < 0 {
return 0.0;
} if k as u64 >= n {
return 1.0;
} if p <= 0.0 {
return 1.0;
} if p >= 1.0 {
return 0.0;
}
let k = k as u64;
let mut cdf = 0.0_f64;
let mut prob = (1.0 - p).powf(n as f64); cdf += prob;
for i in 1..=k as usize {
prob *= ((n - i as u64 + 1) as f64) * p / ((i as f64) * (1.0 - p));
cdf += prob;
}
cdf
}
#[inline(always)]
pub fn inv_std_normal_core(p: f64) -> f64 {
debug_assert!(p > 0.0 && p <= 0.5);
if p > P_LOW {
let r = p - 0.5;
let s = r * r;
let num = (((((A[0] * s + A[1]) * s + A[2]) * s + A[3]) * s + A[4]) * s + A[5]) * r;
let den = ((((B[0] * s + B[1]) * s + B[2]) * s + B[3]) * s + B[4]) * s + 1.0;
num / den
} else {
let r = (-2.0 * p.ln()).sqrt();
let num = ((((C[0] * r + C[1]) * r + C[2]) * r + C[3]) * r + C[4]) * r + C[5];
let den = (((D[0] * r + D[1]) * r + D[2]) * r + D[3]) * r + 1.0;
num / den }
}
#[inline(always)]
pub fn inv_std_normal(p: f64) -> f64 {
if !(p > 0.0 && p < 1.0) {
return f64::NAN;
}
let (q, sign) = if p < 0.5 { (p, 1.0) } else { (1.0 - p, -1.0) };
let x = if q < 0.02425 {
let t = (-2.0 * q.ln()).sqrt();
(((((C[0] * t + C[1]) * t + C[2]) * t + C[3]) * t + C[4]) * t + C[5])
/ ((((D[0] * t + D[1]) * t + D[2]) * t + D[3]) * t + 1.0)
} else {
let t = q - 0.5;
let r = t * t;
(((((A[0] * r + A[1]) * r + A[2]) * r + A[3]) * r + A[4]) * r + A[5]) * t
/ (((((B[0] * r + B[1]) * r + B[2]) * r + B[3]) * r + B[4]) * r + 1.0)
};
sign * x
}
#[inline(always)]
pub fn chi2_newton_refine_extreme(mut x: f64, a: f64, p: f64) -> f64 {
let ln_norm = -a * std::f64::consts::LN_2 - ln_gamma(a);
let mut lo = 0.0_f64;
let mut hi = f64::INFINITY;
for _ in 0..16 {
if !x.is_finite() {
break;
}
let t = 0.5 * x;
let fx = reg_lower_gamma(a, t) - p;
if fx.abs() < 1e-15 {
break;
}
let log_pdf = ln_norm + (a - 1.0) * x.ln() - 0.5 * x;
let pdf = log_pdf.exp();
if pdf <= 0.0 || !pdf.is_finite() {
if fx > 0.0 {
hi = hi.min(x);
} else {
lo = lo.max(x);
}
x = if hi.is_finite() {
0.5 * (lo + hi)
} else {
(x + lo.max(0.0)) * 0.5
};
continue;
}
let mut step = fx / pdf;
let max_step = 0.5 * (x.max(1.0));
if step.abs() > max_step {
step = step.signum() * max_step;
}
let x_new = (x - step).max(0.0);
if fx > 0.0 {
hi = hi.min(x);
} else {
lo = lo.max(x);
}
if x_new < lo || x_new > hi {
x = if hi.is_finite() {
0.5 * (lo + hi)
} else {
(x + lo) * 0.5
};
} else {
x = x_new;
}
}
x
}
#[inline(always)]
pub fn chi2_newton_refine(mut x: f64, a: f64, p: f64) -> f64 {
let ln_norm = -a * std::f64::consts::LN_2 - ln_gamma(a); let mut lo = 0.0_f64;
let mut hi = f64::INFINITY;
for _ in 0..8 {
if !x.is_finite() {
break;
}
let t = 0.5 * x;
let fx = reg_lower_gamma(a, t) - p; if fx.abs() < 1e-14 {
break;
}
let log_pdf = ln_norm + (a - 1.0) * x.ln() - 0.5 * x;
let pdf = log_pdf.exp();
if pdf <= 0.0 || !pdf.is_finite() {
if fx > 0.0 {
hi = hi.min(x);
} else {
lo = lo.max(x);
}
x = if hi.is_finite() {
0.5 * (lo + hi)
} else {
(x + lo.max(0.0)) * 0.5
};
continue;
}
let mut step = fx / pdf;
let max_step = 0.5 * (x.max(1.0));
if step.abs() > max_step {
step = step.signum() * max_step;
}
let x_new = (x - step).max(0.0); if fx > 0.0 {
hi = hi.min(x);
} else {
lo = lo.max(x);
}
x = x_new;
}
x
}
#[inline(always)]
pub fn normal_cdf_scalar(z: f64) -> f64 {
if z < 0.0 {
0.5 * erfc(-z / SQRT_2)
} else {
1.0 - 0.5 * erfc(z / SQRT_2)
}
}
#[inline(always)]
pub fn normal_pdf_scalar(z: f64) -> f64 {
(-0.5 * z * z).exp() / (2.0 * PI).sqrt()
}
pub fn normal_quantile_scalar(q: f64, mean: f64, std: f64) -> f64 {
if !q.is_finite() || !mean.is_finite() || !std.is_finite() || std <= 0.0 {
return f64::NAN;
}
if q < 0.0 || q > 1.0 {
return f64::NAN;
}
if q == 0.0 {
return f64::NEG_INFINITY;
}
if q == 1.0 {
return f64::INFINITY;
}
if q == 0.5 {
return mean;
}
let (p_left, sign) = if q < 0.5 { (q, -1.0) } else { (1.0 - q, 1.0) };
const EPS_DBL: f64 = 1.110_223_024_625_156_5e-16;
if p_left < EPS_DBL {
let z_tail = -SQRT_2 * erfc_inv(2.0 * p_left);
return mean + std * sign * -z_tail; }
let mut z = inv_std_normal_core(p_left);
let pdf = normal_pdf_scalar(z);
let cdf = normal_cdf_scalar(z);
let f = cdf - p_left;
let u = f / pdf;
z -= u * (1.0 + 0.5 * z * u);
let z_final = sign * -z;
mean + std * z_final
}
#[inline(always)]
pub fn student_t_cdf_scalar(x: f64, df: f64) -> f64 {
if !x.is_finite() || !df.is_finite() || df < 1.0 {
return f64::NAN;
}
let xb = df / (df + x * x);
let ib = incomplete_beta(df / 2.0, 0.5, xb);
let sign = if x >= 0.0 { 1.0 } else { -1.0 };
0.5 + sign * 0.5 * (1.0 - ib)
}
#[inline(always)]
pub fn chi_square_cdf_scalar(x: f64, df: f64) -> f64 {
if !x.is_finite() || !df.is_finite() || df < 1.0 {
return f64::NAN;
}
if x <= 0.0 {
return 0.0;
}
reg_lower_gamma(df / 2.0, x / 2.0)
}
#[inline(always)]
pub fn f_cdf_scalar(x: f64, d1: f64, d2: f64) -> f64 {
if !x.is_finite() || !d1.is_finite() || !d2.is_finite() || d1 < 1.0 || d2 < 1.0 {
return f64::NAN;
}
if x <= 0.0 {
return 0.0;
}
let xt = d1 * x / (d1 * x + d2);
incomplete_beta(d1 / 2.0, d2 / 2.0, xt)
}
#[inline(always)]
pub fn student_t_quantile_scalar(q: f64, df: f64) -> f64 {
if !q.is_finite() || !df.is_finite() || df < 1.0 {
return f64::NAN;
}
if q < 0.0 || q > 1.0 {
return f64::NAN;
}
if q == 0.0 {
return f64::NEG_INFINITY;
}
if q == 1.0 {
return f64::INFINITY;
}
if df > 1e7 {
return normal_quantile_scalar(q, 0.0, 1.0);
}
let x = if q < 0.5 {
incomplete_beta_inv(0.5 * df, 0.5, 2.0 * q)
} else {
incomplete_beta_inv(0.5 * df, 0.5, 2.0 * (1.0 - q))
};
let t = ((df * (1.0 - x)) / x).sqrt();
if q < 0.5 { -t } else { t }
}
#[inline(always)]
pub fn normal_quantile_scalar_std(q: f64) -> f64 {
normal_quantile_scalar(q, 0.0, 1.0)
}
#[cfg(test)]
mod tests {
use super::*;
#[test]
fn test_ln_gamma() {
assert!((ln_gamma(1.0) - 0.0).abs() < 1e-14);
assert!((ln_gamma(5.0) - 3.1780538303479458).abs() < 1e-14);
assert!((ln_gamma(0.5) - 0.5723649429247).abs() < 1e-14);
assert!((ln_gamma(10.1) - 13.027526738633238).abs() < 1e-10);
assert!(ln_gamma(0.0).is_infinite() && ln_gamma(0.0).is_sign_positive());
assert!(ln_gamma(-1.0).is_infinite() && ln_gamma(-1.0).is_sign_positive());
assert!((ln_gamma(-0.5) - 1.2655121234846454).abs() < 1e-14);
assert!((ln_gamma(-10.1) - -13.020973271011497).abs() < 1e-12);
assert!(ln_gamma(f64::NAN).is_nan());
assert!((ln_gamma(1e-10) - 23.025850929882733).abs() < 1e-12);
assert!((ln_gamma(171.624) - 709.7807744366991).abs() < 1e-9);
}
#[test]
fn test_ln_gamma_plus1() {
assert!((ln_gamma_plus1(5.0) - 4.787491742782046).abs() < 1e-14);
}
#[test]
fn test_gamma_func() {
assert!((gamma_func(1.0) - 1.0).abs() < 1e-14);
assert!((gamma_func(5.0) - 24.0).abs() < 1e-14);
assert!((gamma_func(0.5) - 1.7724538509055159).abs() < 1e-14);
assert!((gamma_func(10.1) - 454760.7514415855).abs() < 1e-7);
assert!(gamma_func(0.0).is_infinite() && gamma_func(0.0).is_sign_positive());
assert!(gamma_func(-1.0).is_nan());
assert!((gamma_func(-0.5) + 3.5449077018110318).abs() < 1e-14);
assert!((gamma_func(-10.1) + 2.213416583085619e-6).abs() < 1e-14);
assert!((gamma_func(171.0) - 7.257415615308e+306).abs() < 1e292); assert!(gamma_func(f64::NAN).is_nan());
assert!(ln_gamma(1e308).is_infinite());
assert!(ln_gamma(f64::NEG_INFINITY).is_nan());
assert!(ln_gamma(f64::INFINITY).is_infinite());
}
#[test]
fn test_ln_choose() {
assert!((ln_choose(5, 2) - 2.302585092994046).abs() < 1e-14);
assert!(ln_choose(2, 3).is_infinite() && ln_choose(2, 3).is_sign_negative());
assert!((ln_choose(100, 3) - ln_choose(100, 97)).abs() < 1e-14);
assert!((ln_choose(1000, 10) - 53.927997037888275).abs() < 1e-10);
assert!((ln_choose(10, 0) - 0.0).abs() < 1e-14);
assert!((ln_choose(10, 10) - 0.0).abs() < 1e-14);
assert!((ln_choose(0, 0) - 0.0).abs() < 1e-14);
}
#[test]
fn test_incomplete_beta() {
assert!((incomplete_beta(2.0, 2.0, 0.5) - 0.5).abs() < 1e-14);
assert!((incomplete_beta(2.5, 1.5, 0.7) - 0.5843121477019746).abs() < 1e-12);
assert!((incomplete_beta(2.0, 2.0, 0.0) - 0.0).abs() < 1e-14);
assert!((incomplete_beta(2.0, 2.0, 1.0) - 1.0).abs() < 1e-14);
assert!((incomplete_beta(0.0, 2.0, 0.5) - 1.0).abs() < 1e-14);
assert!((incomplete_beta(2.0, 0.0, 0.5) - 0.0).abs() < 1e-14);
assert!(incomplete_beta(2.0, 2.0, f64::NAN).is_nan());
let a = 2.7;
let b = 5.3;
let x = 0.4;
let ix = incomplete_beta(a, b, x);
let ix2 = incomplete_beta(b, a, 1.0 - x);
assert!((ix + ix2 - 1.0).abs() < 1e-16);
assert!((incomplete_beta(50.0, 0.5, 0.01) - 7.998227417904836e-102).abs() < 1e-90);
assert!(incomplete_beta(1e8, 2.0, 0.5).is_finite());
assert!((incomplete_beta(2.0, 3.0, 1e-20) - 0.0).abs() < 1e-20);
assert!((incomplete_beta(2.0, 3.0, 1.0 - 1e-20) - 1.0).abs() < 1e-14);
assert!(incomplete_beta(2.0, 3.0, f64::MIN_POSITIVE).is_finite());
for &a in &[2.0, 5.0, 10.0] {
for &b in &[2.0, 5.0, 10.0] {
for &p in &[1e-15, 1e-6, 0.1, 0.5, 0.9, 1.0 - 1e-8, 1.0 - 1e-15] {
let x = incomplete_beta_inv(a, b, p);
let p2 = incomplete_beta(a, b, x);
assert!((p - p2).abs() < 1e-12);
}
}
}
}
#[test]
fn test_incomplete_beta_inv() {
assert!((incomplete_beta_inv(2.0, 2.0, 0.5) - 0.5).abs() < 1e-14);
assert!(
(incomplete_beta_inv(2.5, 1.5, 0.8433681004114891) - 0.8595961302031141).abs() < 1e-12
);
assert!((incomplete_beta_inv(2.0, 2.0, 0.0) - 0.0).abs() < 1e-14);
assert!((incomplete_beta_inv(2.0, 2.0, 1.0) - 1.0).abs() < 1e-14);
assert!(incomplete_beta_inv(2.0, 2.0, f64::NAN).is_nan());
assert!(incomplete_beta_inv(f64::NAN, 2.0, 0.5).is_nan());
let a = 7.0;
let b = 3.0;
let p = 1e-6;
let x = incomplete_beta_inv(a, b, p);
let p2 = incomplete_beta(a, b, x);
println!("{}", (p - p2).abs());
assert!((p - p2).abs() < 1e-16);
}
#[test]
fn test_reg_lower_gamma() {
assert!((reg_lower_gamma(2.0, 2.0) - 0.5939941502901616).abs() < 1e-12);
assert!((reg_lower_gamma(5.0, 1.0) - 0.003659846827343713).abs() < 1e-14);
assert!((reg_lower_gamma(2.0, 0.0) - 0.0).abs() < 1e-14);
assert!((reg_lower_gamma(0.0, 2.0) - 1.0).abs() < 1e-14);
assert!(reg_lower_gamma(2.0, -1.0).is_nan());
assert!(reg_lower_gamma(-1.0, 2.0).is_nan());
assert!(reg_lower_gamma(f64::NAN, 2.0).is_nan());
let a = 3.5;
let p = 0.123456;
let x = inv_reg_lower_gamma(a, p);
let p2 = reg_lower_gamma(a, x);
assert!((p - p2).abs() < 1e-10);
assert!(reg_lower_gamma(1e8, 1e8).is_finite());
assert!(reg_lower_gamma(1e-20, 1e20).is_finite());
}
#[test]
fn test_inv_reg_lower_gamma() {
assert!((inv_reg_lower_gamma(2.0, 0.5939941502901616) - 2.0).abs() < 1e-10);
assert!(
(inv_reg_lower_gamma(5.0, 0.003659846827343713) - 1.0000000000000002).abs() < 1e-10
);
assert!(inv_reg_lower_gamma(2.0, f64::NAN).is_nan());
assert!(inv_reg_lower_gamma(f64::NAN, 2.0).is_nan());
}
#[test]
fn test_binomial_cdf_scalar() {
assert!((binomial_cdf_scalar(2, 4, 0.5) - 0.6875).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 4, 0.5) - 0.0625).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 4, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(2, 4, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(4, 4, 1.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(2, 4, 1.0) - 0.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 4, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 0, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 4, 1.0) - 0.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(5, 4, 1.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 0, 0.5) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(1, 0, 0.5) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 0, 0.5) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(-1i64, 0, 0.5) - 0.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 0, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(0, 0, 1.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 0, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 0, 1.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 10, 0.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 10, 1.0) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(5, 10000, 1e-4) - 0.999406428148761).abs() < 1e-12);
for n in [0_i64, 1, 10, 100] {
for k in 0..=n {
assert!(
binomial_cdf_scalar(k, n.try_into().unwrap(), 0.3)
<= binomial_cdf_scalar(k + 1, n.try_into().unwrap(), 0.3)
);
}
}
assert!((binomial_cdf_scalar(0, 1000000, 0.5) - (0.5f64).powi(1000000)).abs() < 1e-14);
assert!((binomial_cdf_scalar(1000000, 1000000, 0.5) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 100, f64::MIN_POSITIVE) - 1.0).abs() < 1e-14);
assert!((binomial_cdf_scalar(10, 100, 1.0 - f64::EPSILON) - 0.0).abs() < 1e-14);
}
#[test]
fn test_binomial_quantile_cornish_fisher() {
assert!((binomial_quantile_cornish_fisher(0.5, 10, 0.5) - 5.0).abs() < 1e-14);
assert!((binomial_quantile_cornish_fisher(1e-10, 10, 0.5) - 0.0).abs() < 1e-8);
assert!((binomial_quantile_cornish_fisher(1.0 - 1e-10, 10, 0.5) - 10.0).abs() < 1e-8);
assert_eq!(binomial_quantile_cornish_fisher(0.0, 20, 0.3), -1.0);
assert_eq!(binomial_quantile_cornish_fisher(1.0, 20, 0.3), 20.0);
}
#[test]
fn test_ln_gamma_additional_edges() {
for i in 0..100 {
let x = -(i as f64);
assert!(ln_gamma(x).is_infinite());
}
for i in 0..10 {
let x = -(i as f64);
let just_above = x + 1e-12;
let just_below = x - 1e-12;
assert!(
ln_gamma(just_above).is_finite(),
"ln_gamma({}) not finite",
just_above
);
assert!(
ln_gamma(just_below).is_finite(),
"ln_gamma({}) not finite",
just_below
);
}
}
#[test]
fn test_student_t_cdf_scalar() {
assert!((student_t_cdf_scalar(0.0, 10.0) - 0.5).abs() < 1e-15);
assert!((student_t_cdf_scalar(1.0, 10.0) - 0.82955343384897).abs() < 1e-14);
assert!((student_t_cdf_scalar(-2.0, 5.0) - 0.05096973941492914).abs() < 1e-14);
assert!((student_t_cdf_scalar(3.0, 30.0) - 0.997305017967174).abs() < 1e-14);
assert!((student_t_cdf_scalar(0.5, 1.0) - 0.6475836176504333).abs() < 1e-14);
let cdf_pos = student_t_cdf_scalar(1.0, 10.0);
let cdf_neg = student_t_cdf_scalar(-1.0, 10.0);
assert!((cdf_pos + cdf_neg - 1.0).abs() < 1e-15);
assert!(student_t_cdf_scalar(f64::NAN, 10.0).is_nan());
assert!(student_t_cdf_scalar(1.0, f64::NAN).is_nan());
assert!(student_t_cdf_scalar(1.0, 0.5).is_nan()); }
#[test]
fn test_chi_square_cdf_scalar() {
assert!((chi_square_cdf_scalar(3.841459148, 1.0) - 0.9500000097600837).abs() < 1e-14);
assert!((chi_square_cdf_scalar(5.991464547, 2.0) - 0.9499999999973004).abs() < 1e-14);
assert!((chi_square_cdf_scalar(0.0, 5.0) - 0.0).abs() < 1e-15);
assert!((chi_square_cdf_scalar(10.0, 5.0) - 0.9247647538534878).abs() < 1e-14);
assert!((chi_square_cdf_scalar(1.0, 1.0) - 0.6826894921370859).abs() < 1e-14);
assert!((chi_square_cdf_scalar(20.0, 10.0) - 0.9707473119230389).abs() < 1e-14);
assert!((chi_square_cdf_scalar(-1.0, 5.0) - 0.0).abs() < 1e-15);
assert!(chi_square_cdf_scalar(f64::NAN, 5.0).is_nan());
assert!(chi_square_cdf_scalar(5.0, 0.5).is_nan()); }
#[test]
fn test_f_cdf_scalar() {
assert!((f_cdf_scalar(1.0, 5.0, 10.0) - 0.5348805734622).abs() < 1e-13);
assert!((f_cdf_scalar(3.0, 2.0, 20.0) - 0.9274618497135944).abs() < 1e-14);
assert!((f_cdf_scalar(0.0, 5.0, 10.0) - 0.0).abs() < 1e-15);
assert!((f_cdf_scalar(4.0, 3.0, 15.0) - 0.9718629601291221).abs() < 1e-14);
assert!((f_cdf_scalar(10.0, 1.0, 1.0) - 0.805017770957863).abs() < 1e-14);
assert!(f_cdf_scalar(f64::NAN, 5.0, 10.0).is_nan());
assert!(f_cdf_scalar(1.0, 0.5, 10.0).is_nan()); assert!(f_cdf_scalar(1.0, 5.0, 0.5).is_nan()); }
#[test]
fn test_student_t_quantile_scalar() {
assert!(
(student_t_quantile_scalar(0.975, 10.0) - 2.22813885198627481543e+00).abs() < 1e-14
);
assert!(
(student_t_quantile_scalar(0.025, 10.0) - -2.22813885198627481543e+00).abs() < 1e-14
);
assert!((student_t_quantile_scalar(0.99, 5.0) - 3.36492999890721877776e+00).abs() < 1e-14);
assert!((student_t_quantile_scalar(0.01, 5.0) - -3.36492999890721877776e+00).abs() < 1e-14);
assert!(student_t_quantile_scalar(0.5, 10.0).abs() < 1e-14);
assert!((student_t_quantile_scalar(0.9, 15.0) - 1.34060560785045557175e+00).abs() < 1e-14);
assert!((student_t_quantile_scalar(0.1, 15.0) - -1.34060560785045557175e+00).abs() < 1e-14);
assert!((student_t_quantile_scalar(0.75, 20.0) - 6.86954496448803464403e-01).abs() < 1e-14);
assert!((student_t_quantile_scalar(0.999, 3.0) - 1.02145318524073864808e+01).abs() < 1e-14);
assert!(
(student_t_quantile_scalar(0.001, 3.0) - -1.02145318524073864808e+01).abs() < 1e-14
);
assert!((student_t_quantile_scalar(0.95, 1.0) - 6.31375151467504291958e+00).abs() < 1e-14);
assert!(
(student_t_quantile_scalar(0.99, 100.0) - 2.36421736623848222081e+00).abs() < 1e-14
);
assert!((student_t_quantile_scalar(0.9, 2.0) - 1.88561808316412671260e+00).abs() < 1e-14);
assert!(
(student_t_quantile_scalar(0.999, 30.0) - 3.38518486682930497267e+00).abs() < 1e-14
);
for &df in &[1.0, 3.0, 10.0, 30.0, 100.0] {
for &q in &[0.01, 0.1, 0.25] {
let lo = student_t_quantile_scalar(q, df);
let hi = student_t_quantile_scalar(1.0 - q, df);
assert!(
(lo + hi).abs() < 1e-12,
"symmetry failed for q={}, df={}: lo={}, hi={}",
q,
df,
lo,
hi
);
}
}
assert!(
student_t_quantile_scalar(0.0, 10.0).is_infinite()
&& student_t_quantile_scalar(0.0, 10.0).is_sign_negative()
);
assert!(
student_t_quantile_scalar(1.0, 10.0).is_infinite()
&& student_t_quantile_scalar(1.0, 10.0).is_sign_positive()
);
assert!(student_t_quantile_scalar(f64::NAN, 10.0).is_nan());
assert!(student_t_quantile_scalar(0.5, 0.5).is_nan());
for &df in &[1.0, 2.0, 5.0, 10.0, 30.0, 50.0, 100.0, 500.0] {
for &q in &[
0.001, 0.01, 0.05, 0.1, 0.25, 0.5, 0.75, 0.9, 0.95, 0.975, 0.99, 0.999,
] {
let x = student_t_quantile_scalar(q, df);
let p = student_t_cdf_scalar(x, df);
assert!(
(p - q).abs() < 1e-14,
"roundtrip failed for q={}, df={}: got p={}, err={:.2e}",
q,
df,
p,
(p - q).abs()
);
}
}
}
#[test]
fn test_normal_quantile_scalar_std() {
assert!((normal_quantile_scalar_std(0.5) - 0.0).abs() < 1e-15);
assert!((normal_quantile_scalar_std(0.975) - 1.959963984540054).abs() < 1e-14);
assert!((normal_quantile_scalar_std(0.025) - -1.9599639845400545).abs() < 1e-14);
assert!((normal_quantile_scalar_std(0.01) - -2.3263478740408408).abs() < 1e-14);
assert!(normal_quantile_scalar_std(f64::NAN).is_nan());
}
}