mathru/special/
gamma.rs

1//! Provides [gamma](https://en.wikipedia.org/wiki/Beta_function) related functions
2
3use crate::algebra::abstr::cast::FromPrimitive;
4use crate::algebra::abstr::cast::ToPrimitive;
5use crate::algebra::abstr::Real;
6use crate::elementary::Power;
7use std::f64;
8
9/// Provides [gamma](https://en.wikipedia.org/wiki/Beta_function) related functions
10pub trait Gamma {
11    /// Gamma function
12    ///
13    /// ```math
14    /// \Gamma(z) = \int_0^\infty t^{z-1} {\mathrm e}^{-t} \mathrm dt
15    /// ```
16    ///
17    /// Fore more information:
18    /// <https://en.wikipedia.org/wiki/Gamma_function>
19    ///
20    /// # Arguments
21    ///
22    /// * `self` > 0.0
23    ///
24    /// # Panics
25    ///
26    /// *`self` <= 0.0
27    ///
28    /// # Example
29    ///
30    /// ```
31    /// use mathru::special::gamma::Gamma;
32    ///
33    /// let z: f64 = 0.3_f64;
34    /// let gamma: f64 = z.gamma();
35    /// ```
36    fn gamma(self) -> Self;
37
38    /// Log-gamma function
39    ///
40    /// ln&Gamma;(z)
41    //
42    /// Fore more information:
43    /// <https://en.wikipedia.org/wiki/Gamma_function#The_log-gamma_function>
44    ///
45    /// # Arguments
46    ///
47    /// * `self`
48    ///
49    /// ```
50    /// use mathru::special::gamma::Gamma;
51    ///
52    /// let x: f64 = 0.3_f64;
53    /// let ln_gamma: f64 = x.ln_gamma();
54    /// ```
55    fn ln_gamma(self) -> Self;
56
57    /// Digamma function
58    ///
59    /// ```math
60    /// \psi(x)=\frac{d}{dx}\ln\big(\Gamma(x)\big)=\frac{\Gamma'(x)}{\Gamma(x)}
61    /// ```
62    ///
63    /// Fore more information:
64    /// <https://en.wikipedia.org/wiki/Digamma_function>
65    ///
66    /// # Arguments
67    ///
68    /// * `self`
69    ///
70    /// # Example
71    ///
72    /// ```
73    /// use mathru::special::gamma::Gamma;
74    ///
75    /// let x: f64 = 0.3_f64;
76    /// let digamma: f64 = x.digamma();
77    /// ```
78    fn digamma(self) -> Self;
79
80    /// Upper incomplete gamma function
81    ///
82    /// ```math
83    /// \Gamma(a,x) = \int_x^{\infty} t^{a-1}\,\mathrm{e}^{-t}\,{\rm d}t
84    /// ```
85    ///
86    /// Fore more information:
87    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Upper_incomplete_Gamma_function>
88    ///
89    /// # Arguments
90    ///
91    /// * `self`
92    /// * `x`
93    ///
94    /// # Example
95    ///
96    /// ```
97    /// use mathru::special::gamma::Gamma;
98    ///
99    /// let a: f64 = 0.5_f64;
100    /// let x: f64 = 0.3_f64;
101    ///
102    /// let gamma_u: f64 = a.gamma_u(x);
103    /// ```
104    fn gamma_u(self, x: Self) -> Self;
105
106    /// Upper incomplete regularized gamma function
107    ///
108    /// ```math
109    /// Q(a,x)=\frac{\Gamma(a,x)}{\Gamma(a)}
110    /// ```
111    ///
112    /// Fore more information:
113    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
114    ///
115    /// # Arguments
116    ///
117    /// * `self`
118    /// * `x`
119    ///
120    /// # Example
121    ///
122    /// ```
123    /// use mathru::special::gamma::Gamma;
124    ///
125    /// let a: f64 = 0.5_f64;
126    /// let x: f64 = 0.3_f64;
127    /// let gamma_ur: f64 = a.gamma_ur(x);
128    /// ```
129    fn gamma_ur(self, x: Self) -> Self;
130
131    /// Lower incomplete gamma function
132    ///
133    /// ```math
134    /// \gamma(a,x) = \int_0^x t^{a-1}\,\mathrm{e}^{-t}\,{\rm d}t
135    /// ```
136    ///
137    /// Fore more information:
138    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
139    ///
140    /// # Arguments
141    ///
142    /// * `self`
143    /// * `x`
144    ///
145    /// # Example
146    ///
147    /// ```
148    /// use mathru::special::gamma::Gamma;
149    ///
150    /// let a: f64 = 0.5_f64;
151    /// let x: f64 = 0.3_f64;
152    /// let gamma_l: f64 = a.gamma_l(x);
153    /// ```
154    fn gamma_l(self, x: Self) -> Self;
155
156    /// Lower regularized incomplete gamma function
157    ///
158    /// ```math
159    /// P(a,x)=\frac{\gamma(a,x)}{\Gamma(a)}=1-Q(a,x)
160    /// ```
161    ///
162    /// Fore more information:
163    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
164    ///
165    /// # Arguments
166    ///
167    /// * `self`
168    /// * `x`
169    ///
170    /// # Example
171    ///
172    /// ```
173    /// use mathru::special::gamma::Gamma;
174    ///
175    /// let a: f64 = 0.5_f64;
176    /// let x: f64 = 0.3_f64;
177    /// let gamma_lr: f64 = a.gamma_lr(x);
178    /// ```
179    fn gamma_lr(self, x: Self) -> Self;
180
181    /// Inverse of the upper incomplete regularized gamma function
182    ///
183    /// ```math
184    /// Q^{-1}(q,x)
185    /// ```
186    ///
187    /// Fore more information:
188    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
189    ///
190    /// # Arguments
191    ///
192    /// * `q`
193    /// * `x`
194    ///
195    /// # Example
196    ///
197    /// ```
198    /// use mathru::special::gamma;
199    /// use mathru::assert_abs_diff_eq;
200    ///
201    /// let a: f64 = 0.5_f64;
202    /// let x: f64 = 0.3_f64;
203    /// let q = gamma::gamma_ur(a, x);
204    /// let x_s: f64 = gamma::gamma_ur_inv(a, q);
205    /// assert_abs_diff_eq!(x, x_s);
206    /// ```
207    fn gamma_ur_inv(self, p: Self) -> Self;
208
209    /// Inverse of the lower incomplete regularized gamma function
210    ///
211    /// ```math
212    /// P^{-1}(a,p)
213    /// ```
214    ///
215    /// Fore more information:
216    /// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
217    ///
218    /// # Arguments
219    ///
220    /// * `a`
221    /// * `x`
222    ///
223    /// # Example
224    ///
225    /// ```
226    /// use mathru::special::gamma;
227    /// use mathru::assert_abs_diff_eq;
228    ///
229    /// let a: f64 = 0.5_f64;
230    /// let x: f64 = 0.3_f64;
231    /// let p = gamma::gamma_lr(a, x);
232    /// let x_s: f64 = gamma::gamma_lr_inv(a, p);
233    /// assert_abs_diff_eq!(x, x_s);
234    /// ```
235    fn gamma_lr_inv(self, p: Self) -> Self;
236}
237
238macro_rules! impl_gamma {
239    ($t: ty, $PI: expr, $E: expr) => {
240        #[allow(clippy::excessive_precision)]
241        impl Gamma for $t {
242            fn gamma(self) -> Self {
243                if self == 0.0 {
244                    panic!("The gamma function is undefined for z == 0.0")
245                }
246
247                if self < 0.5 {
248                    return $PI / (($PI * self).sin() * (1.0 - self).gamma()); //reflection formula
249                }
250
251                let t: $t = self + 6.5;
252                let x: $t = 0.999_999_999_999_809_9 + 676.5203681218851 / self
253                    - 1259.1392167224028 / (self + 1.0)
254                    + 771.323_428_777_653_1 / (self + 2.0)
255                    - 176.615_029_162_140_6 / (self + 3.0)
256                    + 12.507343278686905 / (self + 4.0)
257                    - 0.13857109526572012 / (self + 5.0)
258                    + 9.984_369_578_019_572e-6 / (self + 6.0)
259                    + 1.5056327351493116e-7 / (self + 7.0);
260
261                return 2.0.sqrt() * $PI.sqrt() * t.pow((self - 0.5)) * (-t).exp() * x;
262            }
263
264            fn ln_gamma(self) -> Self {
265                // Auxiliary variable when evaluating the `gamma_ln` function
266                let gamma_r: Self = 10.900511;
267
268                // Polynomial coefficients for approximating the `gamma_ln` function
269                let gamma_dk: &[$t] = &[
270                    2.485_740_891_387_535_5e-5,
271                    1.051_423_785_817_219_7,
272                    -3.456_870_972_220_1625,
273                    4.512_277_094_668_948,
274                    -2.982_852_253_235_766_4,
275                    1.056_397_115_771_267,
276                    -1.954_287_731_916_458_7e-1,
277                    1.709_705_434_044_412e-2,
278                    -5.719_261_174_043_057e-4,
279                    4.633_994_733_599_057e-6,
280                    -2.719_949_084_886_077_2e-9,
281                ];
282
283                let x: Self = self;
284
285                if x < 0.5 {
286                    let s = gamma_dk
287                        .iter()
288                        .enumerate()
289                        .skip(1)
290                        .fold(gamma_dk[0], |s, t| s + *t.1 / ((t.0 as u64) as $t - x));
291
292                    $PI.ln()
293                        - ($PI * x).sin().ln()
294                        - s.ln()
295                        - (2.0 * ($E / $PI).pow(0.5)).ln()
296                        - (0.5 - x) * ((0.5 - x + gamma_r) / $E).ln()
297                } else {
298                    let s = gamma_dk
299                        .iter()
300                        .enumerate()
301                        .skip(1)
302                        .fold(gamma_dk[0], |s, t| {
303                            s + *t.1 / (x + (t.0 as u64) as $t - 1.0)
304                        });
305
306                    s.ln()
307                        + (2.0 * ($E / $PI).pow(0.5)).ln()
308                        + (x - 0.5) * ((x - 0.5 + gamma_r) / $E).ln()
309                }
310            }
311
312            //
313            fn digamma(self) -> Self {
314                let c: Self = 8.5;
315                let mut value: Self = 0.0;
316                let mut x2: Self = self;
317                //The comparison only compares the real part ot the number
318                while x2 < c {
319                    value -= 1.0 / x2;
320                    x2 += 1.0;
321                }
322                /*
323                  Use Stirling's (actually de Moivre's) expansion
324                */
325                let mut r: Self = 1.0 / x2;
326                value = value + x2.ln() - 0.5 * r;
327
328                r = r * r;
329
330                value -= r
331                    * (1.0 / 12.0
332                        - r * (1.0 / 120.0
333                            - r * (1.0 / 252.0
334                                - r * (1.0 / 240.0
335                                    - r * (1.0 / 132.0
336                                        - r * (691.0 / 32760.0 - r * (1.0 / 12.0)))))));
337
338                return value;
339            }
340
341            fn gamma_u(self, x: Self) -> Self {
342                return self.gamma_ur(x) * self.gamma();
343            }
344
345            fn gamma_ur(self, x: Self) -> Self {
346                return 1.0 - self.gamma_lr(x);
347            }
348
349            fn gamma_l(self, x: Self) -> Self {
350                return self.gamma_lr(x) * self.gamma();
351            }
352
353            fn gamma_lr(self, x: Self) -> Self {
354                if x <= 0.0 {
355                    panic!("Lower regularized gamma function is not defined for `x` <= 0.0");
356                }
357
358                let eps: Self = 0.000000000000001;
359                let big: Self = 4503599627370496.0;
360                let big_inv: Self = 2.220_446_049_250_313e-16;
361
362                if self == 0.0 {
363                    return 1.0;
364                }
365
366                if x == 0.0 {
367                    return 0.0;
368                }
369
370                let ax: Self = self * x.ln() - x - self.ln_gamma();
371
372                if ax < -709.782_712_893_384 {
373                    if self < x {
374                        return 1.0;
375                    }
376                    return 0.0;
377                }
378
379                if x <= 1.0 || x <= self {
380                    let mut r2: Self = self;
381                    let mut c2: Self = 1.0;
382                    let mut ans2: Self = 1.0;
383                    loop {
384                        r2 += 1.0;
385                        c2 *= x / r2;
386                        ans2 += c2;
387
388                        if c2 / ans2 <= eps {
389                            break;
390                        }
391                    }
392                    return ax.exp() * ans2 / self;
393                }
394
395                let mut y: Self = 1.0 - self;
396                let mut z: Self = x + y + 1.0;
397                let mut c: Self = 0.0;
398
399                let mut p3: Self = 1.0;
400                let mut q3: Self = x;
401                let mut p2: Self = x + 1.0;
402                let mut q2: Self = z * x;
403                let mut ans: Self = p2 / q2;
404
405                loop {
406                    y += 1.0;
407                    z += 2.0;
408                    c += 1.0;
409                    let yc: Self = y * c;
410
411                    let p = p2 * z - p3 * yc;
412                    let q = q2 * z - q3 * yc;
413
414                    p3 = p2;
415                    p2 = p;
416                    q3 = q2;
417                    q2 = q;
418
419                    if p.abs() > big {
420                        p3 *= big_inv;
421                        p2 *= big_inv;
422                        q3 *= big_inv;
423                        q2 *= big_inv;
424                    }
425
426                    if q != 0.0 {
427                        let nextans = p / q;
428                        let error = ((ans - nextans) / nextans).abs();
429                        ans = nextans;
430
431                        if error <= eps {
432                            break;
433                        }
434                    }
435                }
436
437                1.0 - ax.exp() * ans
438            }
439
440            /// Computation of the Incomplete Gamma Function Ratios and their Inverse
441            /// ARMIDO R. DIDONATO and ALFRED H. MORRIS, JR. U.S. Naval Surface Weapons Center
442            /// ACM Transactions on Mathematical Software, Vol. 12, No. 4, December 1986
443            fn gamma_ur_inv(self, q: Self) -> Self {
444                let a: Self = self;
445                let c: Self = 0.577_215_664_901_532_9;
446                let eps: Self = 1.0e-10;
447                let tau: Self = 1.0e-5;
448
449                if a == 1.0 {
450                    return -q.ln();
451                }
452
453                let b: Self = q * gamma(a);
454                let p: Self = 1.0 - q;
455                let x_0: Self;
456
457                if 0.6 < b || (0.45 <= b && a >= 0.3) {
458                    let u: Self = if b * q > 1.0e-8 {
459                        (p * gamma(a + 1.0)).pow(1.0 / a)
460                    } else {
461                        (-q / a - c).exp()
462                    };
463
464                    x_0 = u / (1.0 - u / (a + 1.0))
465                } else {
466                    if a < 0.3 && (0.35..=0.6).contains(&b) {
467                        let t: Self = (-c - b).exp();
468                        let u: Self = t * t.exp();
469                        x_0 = t * u.exp()
470                    } else {
471                        if (0.15..0.35).contains(&b) || (0.35 <= a && (0.15..0.45).contains(&b)) {
472                            let y: Self = -b.ln();
473                            let v: Self = y - (1.0 - a) * y.ln();
474                            x_0 = y - (1.0 - a) * v.ln() - (1.0 + (1.0 - a) / (1.0 + v)).ln();
475                        } else {
476                            if 0.01 < b && b < 0.15 {
477                                let y: Self = -b.ln();
478                                let v: Self = y - (1.0 - a) * y.ln();
479                                x_0 = y
480                                    - (1.0 - a) * v.ln()
481                                    - ((v * v + 2.0 * (3.0 - a) * v + (2.0 - a) * (3.0 - a))
482                                        / (v * v + (5.0 - a) * v + 2.0))
483                                        .ln()
484                            } else {
485                                let y: Self = -b.ln();
486                                let c_1: Self = (a - 1.0) * y.ln();
487
488                                let c_2: Self = (a - 1.0) * (1.0 + c_1);
489
490                                let c_1_2: Self = c_1 * c_1;
491                                let c_3: Self = (a - 1.0)
492                                    * (-c_1_2 / 2.0 + (a - 2.0) * c_1 + (3.0 * a - 5.0) / 2.0);
493
494                                let a_2: Self = a * a;
495
496                                let c_1_3: Self = c_1_2 * c_1;
497                                let c_4: Self = (a - 1.0)
498                                    * (c_1_3 / 3.0 - (3.0 * a - 5.0) / 2.0 * c_1_2
499                                        + (a_2 + -6.0 * a + 7.0) * c_1
500                                        + (11.0 * a_2 - 46.0 * a + 47.0) / 6.0);
501
502                                let c_1_4: Self = c_1_3 * c_1;
503
504                                let a_3: Self = a_2 * a;
505
506                                let c_5: Self = (a - 1.0)
507                                    * (-c_1_4 / 4.0
508                                        + (11.0 * a - 17.0) / 6.0 * c_1_3
509                                        + (-3.0 * a_2 + 13.0 * a - 13.0) * c_1_2
510                                        + (2.0 * a_3 - 25.0 * a_2 + 72.0 * a - 61.0) / 2.0 * c_1
511                                        + (25.0 * a_3 - 195.0 * a_2 + 477.0 * a - 379.0) / 12.0);
512
513                                let y_2: Self = y * y;
514                                let y_3: Self = y_2 * y;
515                                let y_4: Self = y_3 * y;
516
517                                x_0 = y + c_1 + c_2 / y + c_3 / y_2 + c_4 / y_3 + c_5 / y_4;
518                            }
519                        }
520                    }
521                }
522
523                let mut x_n: Self = x_0;
524
525                for _i in 0..20 {
526                    let r_x: Self = Self::from_f64(r(a.to_f64(), x_n.to_f64()));
527                    let w_n: Self = (a - 1.0 - x_n) / 2.0;
528
529                    let t_n: Self = if p <= 0.5 {
530                        (gamma_lr(a, x_n) - p) / r_x
531                    } else {
532                        (q - gamma_ur(a, x_n)) / r_x
533                    };
534
535                    let x_n_1: Self = if t_n.abs() <= 0.1 && (w_n * t_n).abs() <= 0.1 {
536                        // Schröder method
537                        let k: Self = w_n * t_n * t_n;
538                        let h_n: Self = t_n + k;
539                        let x_n_1: Self = x_n * (1.0 - h_n);
540
541                        if w_n.abs() >= 1.0 && k.abs() <= eps {
542                            return x_n_1;
543                        }
544
545                        x_n_1
546                    } else {
547                        // Newton-Raphson
548                        let h_n: Self = t_n;
549                        let x_n_1: Self = x_n * (1.0 - h_n);
550
551                        if h_n.abs() < eps {
552                            return x_n_1;
553                        }
554
555                        if h_n.abs() <= tau {
556                            if p <= 0.5 {
557                                if (gamma_lr(a, x_n) - p).abs() < tau * p {
558                                    return x_n_1;
559                                }
560                            } else {
561                                if (q - gamma_ur(a, x_n)).abs() <= tau * q {
562                                    return x_n_1;
563                                }
564                            }
565                        }
566
567                        x_n_1
568                    };
569
570                    x_n = x_n_1;
571                }
572
573                return x_n;
574
575                fn r(a: f64, x: f64) -> f64 {
576                    if a <= 0.0 {
577                        panic!();
578                    }
579                    if x < 0.0 {
580                        panic!();
581                    }
582
583                    x.powf(a) * (-x).exp() / gamma(a)
584                }
585            }
586
587            fn gamma_lr_inv(self, q: Self) -> Self {
588                self.gamma_ur_inv(1.0 - q)
589            }
590        }
591    };
592}
593
594impl_gamma!(f64, std::f64::consts::PI, std::f64::consts::E);
595impl_gamma!(f32, std::f32::consts::PI, std::f32::consts::E);
596
597/// Gamma function
598///
599/// ```math
600/// \Gamma(z) = \int_0^\infty t^{z-1} {\mathrm e}^{-t} \mathrm dt
601/// ```
602///
603/// Fore more information:
604/// <https://en.wikipedia.org/wiki/Gamma_function>
605///
606/// # Arguments
607///
608/// * `z`
609///
610/// # Panics
611///
612/// *`z` == 0.0
613///
614/// # Example
615///
616/// ```
617/// use mathru::special::gamma;
618///
619/// let z: f64 = 0.3_f64;
620/// let gamma: f64 = gamma::gamma(z);
621/// ```
622/// The following approximation is implemented
623/// <https://en.wikipedia.org/wiki/Lanczos_approximation>
624pub fn gamma<T>(z: T) -> T
625where
626    T: Real + Gamma,
627{
628    z.gamma()
629}
630
631/// Log-gamma function
632///
633/// ln&Gamma;(z)
634///
635/// Fore more information:
636/// <https://en.wikipedia.org/wiki/Gamma_function#The_log-gamma_function>
637///
638/// # Arguments
639///
640/// * `z`
641///
642/// # Example
643///
644/// ```
645/// use mathru::special::gamma;
646///
647/// let x: f64 = 0.3_f64;
648/// let ln_gamma: f64 = gamma::ln_gamma(x);
649/// ```
650pub fn ln_gamma<T>(x: T) -> T
651where
652    T: Gamma,
653{
654    x.gamma()
655}
656
657/// Digamma function
658///
659/// ```math
660/// \psi(x)=\frac{d}{dx}\ln\big(\Gamma(x)\big)=\frac{\Gamma'(x)}{\Gamma(x)}
661/// ```
662///
663/// Fore more information:
664/// <https://en.wikipedia.org/wiki/Digamma_function>
665///
666/// # Arguments
667///
668/// * `x`
669///
670/// # Example
671///
672/// ```
673/// use mathru::special::gamma;
674///
675/// let x: f64 = 0.3_f64;
676/// let digamma: f64 = gamma::digamma(x);
677/// ```
678pub fn digamma<T>(x: T) -> T
679where
680    T: Gamma,
681{
682    x.digamma()
683}
684
685/// Upper incomplete gamma function
686///
687/// ```math
688/// \Gamma(a,x) = \int_x^{\infty} t^{a-1}\,\mathrm{e}^{-t}\,{\rm d}t
689/// ```
690///
691/// Fore more information:
692/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Upper_incomplete_Gamma_function>
693///
694/// # Arguments
695///
696/// * `a`
697/// * `x`
698///
699/// # Example
700///
701/// ```
702/// use mathru::special::gamma;
703///
704/// let a: f64 = 0.5_f64;
705/// let x: f64 = 0.3_f64;
706///
707/// let gamma_u: f64 = gamma::gamma_u(a, x);
708/// ```
709pub fn gamma_u<T>(a: T, x: T) -> T
710where
711    T: Real + Gamma,
712{
713    a.gamma_u(x)
714}
715
716/// Upper incomplete regularized gamma function
717///
718/// ```math
719/// Q(a,x)=\frac{\Gamma(a,x)}{\Gamma(a)}
720/// ```
721///
722/// Fore more information:
723/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
724///
725/// # Arguments
726///
727/// * `a`
728/// * `x`
729///
730/// # Example
731///
732/// ```
733/// use mathru::special::gamma;
734///
735/// let a: f64 = 0.5_f64;
736/// let x: f64 = 0.3_f64;
737/// let gamma_u: f64 = gamma::gamma_ur(a, x);
738/// ```
739pub fn gamma_ur<T>(a: T, x: T) -> T
740where
741    T: Real + Gamma,
742{
743    a.gamma_ur(x)
744}
745
746/// Inverse of the upper incomplete regularized gamma function
747///
748/// ```math
749/// Q^{-1}(a,q)
750/// ```
751///
752/// Fore more information:
753/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
754///
755/// # Arguments
756///
757/// * `a`
758/// * `q`
759///
760/// # Example
761///
762/// ```
763/// use mathru::special::gamma;
764///
765/// let a: f64 = 0.5_f64;
766/// let x: f64 = 0.3_f64;
767/// let q = gamma::gamma_ur(a, x);
768/// let x_s: f64 = gamma::gamma_ur_inv(a, q);
769/// ```
770pub fn gamma_ur_inv<T>(a: T, q: T) -> T
771where
772    T: Real + Gamma,
773{
774    a.gamma_ur_inv(q)
775}
776
777/// Lower incomplete gamma function
778///
779/// ```math
780/// \gamma(a,x) = \int_0^x t^{a-1}\,\mathrm{e}^{-t}\,{\rm d}t
781/// ```
782///
783/// Fore more information:
784/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
785///
786/// # Arguments
787///
788/// * `a`
789/// * `x`
790///
791/// # Example
792///
793/// ```
794/// use mathru::special::gamma;
795///
796/// let a: f64 = 0.5_f64;
797/// let x: f64 = 0.3_f64;
798/// let gamma_l: f64 = gamma::gamma_l(a, x);
799/// ```
800pub fn gamma_l<T>(a: T, x: T) -> T
801where
802    T: Real + Gamma,
803{
804    a.gamma_l(x)
805}
806
807/// Lower regularized incomplete gamma function
808///
809/// ```math
810/// P(a,x)=\frac{\gamma(a,x)}{\Gamma(a)}=1-Q(a,x)
811/// ```
812///
813/// Fore more information:
814/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
815///
816/// <https://people.sc.fsu.edu/~jburkardt/c_src/asa239/asa239.c>
817/// # Arguments
818///
819/// * `a`
820/// * `x`
821///
822/// # Example
823///
824/// ```
825/// use mathru::special::gamma;
826///
827/// let a: f64 = 0.5_f64;
828/// let x: f64 = 0.3_f64;
829/// let gamma_lr: f64 = gamma::gamma_lr(a, x);
830/// ```
831pub fn gamma_lr<T>(a: T, x: T) -> T
832where
833    T: Real + Gamma,
834{
835    a.gamma_lr(x)
836}
837
838/// Inverse of the lower incomplete regularized gamma function
839///
840/// ```math
841/// P^{-1}(a,p)
842/// ```
843///
844/// Fore more information:
845/// <https://en.wikipedia.org/wiki/Incomplete_gamma_function#Regularized_Gamma_functions_and_Poisson_random_variables>
846///
847/// # Arguments
848///
849/// * `a`
850/// * `p`
851///
852/// # Example
853///
854/// ```
855/// use mathru::special::gamma;
856///
857/// let a: f64 = 0.5_f64;
858/// let x: f64 = 0.3_f64;
859/// let p = gamma::gamma_lr(a, x);
860/// let x_s: f64 = gamma::gamma_lr_inv(a, p);
861/// ```
862pub fn gamma_lr_inv<T>(a: T, p: T) -> T
863where
864    T: Real + Gamma,
865{
866    a.gamma_lr_inv(p)
867}