Skip to main content

special_fun/
lib.rs

1#![no_std]
2
3use core::ops::{Add, Sub};
4
5/// Low-level bindings to double-precision Cephes functions.
6///
7/// This provides access to the C function implemented in Cephes for maximum
8/// flexibility.
9#[allow(dead_code)]
10pub mod unsafe_cephes_double {
11    extern "C" {
12        // Floating point numeric utilities
13        /// Round to nearest or event integer valued f64.
14        pub fn round(x: f64) -> f64;
15        /// Largest integer less than or equal to x.
16        pub fn floor(x: f64) -> f64;
17        /// Smallest integer greater than or equal to x.
18        pub fn ceil(x: f64) -> f64;
19        /// Return the significand between 0.5 and 1. Write exponent to expnt.
20        /// x = y * 2**expn
21        pub fn frexp(x: f64, expnt: &mut i32) -> f64;
22        /// Multiply x by 2**n.
23        pub fn ldexp(x: f64, n: i32) -> f64;
24        /// Absolute value.
25        pub fn fabs(x: f64) -> f64;
26        /// Return 1 if the sign bit of x is 1, else 0.
27        pub fn signbit(x: f64) -> i32;
28        /// Return 1 if x is NaN, else 0.
29        pub fn isnan(x: f64) -> i32;
30        /// Return 1 if x is finite, else 0.
31        pub fn isfinite(x: f64) -> i32;
32
33        // Roots
34        /// Cube root.
35        pub fn cbrt(x: f64) -> f64;
36        /// Square root.
37        pub fn sqrt(x: f64) -> f64;
38        /// Integer square root.
39        pub fn lsqrt(x: i64) -> i64;
40
41        // Exponential functions
42        /// Exponential function.
43        pub fn exp(x: f64) -> f64;
44        /// Base 10 exponential function.
45        pub fn exp10(x: f64) -> f64;
46        /// Base 2 exponential function.
47        pub fn exp2(x: f64) -> f64;
48        /// Compute accurately exponential of squared argument.
49        pub fn expm1(x: f64) -> f64;
50        /// Compute accurately exp(x) - 1 for x close to 0.
51        pub fn expx2(x: f64, sign: i32) -> f64;
52        /// Exponential integral.
53        pub fn ei(x: f64) -> f64;
54        /// Error function.
55        pub fn erf(x: f64) -> f64;
56        /// Complementary error function.
57        pub fn erfc(x: f64) -> f64;
58        /// Power function.
59        pub fn pow(x: f64, y: f64) -> f64;
60        /// Integer power function.
61        pub fn powi(x: f64, n: i32) -> f64;
62
63        // Logarithmic functions
64        /// Natural logarithm.
65        pub fn log(x: f64) -> f64;
66        /// Common logarithm.
67        pub fn log10(x: f64) -> f64;
68        /// Base 2 logarithm.
69        pub fn log2(x: f64) -> f64;
70        /// Compute accurately log(1 + x) for x close to 0.
71        pub fn log1p(x: f64) -> f64;
72        /// Dilogarithm (Spence's function).
73        pub fn spence(x: f64) -> f64;
74
75        // Trigonometric functions
76        /// Circular sine.
77        pub fn sin(x: f64) -> f64;
78        /// Circular cosine.
79        pub fn cos(x: f64) -> f64;
80        /// Circular tangent.
81        pub fn tan(x: f64) -> f64;
82        /// Inverse circular sine.
83        pub fn asin(x: f64) -> f64;
84        /// Inverse circular cosine.
85        pub fn acos(x: f64) -> f64;
86        /// Inverse circular tangent.
87        pub fn atan(x: f64) -> f64;
88        /// Quadrant-correct inverse circular tangent.
89        pub fn atan2(y: f64, x: f64) -> f64;
90        /// Compute accurately cos(x) - 1 for x close to 0.
91        pub fn cosm1(x: f64) -> f64;
92        /// Sine and cosine integrals.
93        pub fn sici(x: f64, si: &mut f64, ci: &mut f64) -> i32;
94
95        // Hyperbolic functions
96        /// Hyperbolic sine.
97        pub fn sinh(x: f64) -> f64;
98        /// Hyperbolic cosine.
99        pub fn cosh(x: f64) -> f64;
100        /// Hyperbolic tangent.
101        pub fn tanh(x: f64) -> f64;
102        /// Inverse hyperbolic sine.
103        pub fn asinh(x: f64) -> f64;
104        /// Inverse hyperbolic cosine.
105        pub fn acosh(x: f64) -> f64;
106        /// Inverse hyperbolic tangent.
107        pub fn atanh(x: f64) -> f64;
108        /// Hyperbolic sine and cosine integrals.
109        pub fn shichi(x: f64, chi: &mut f64, shi: &mut f64) -> i32;
110
111        // Beta functions
112        /// Beta function.
113        pub fn beta(a: f64, b: f64) -> f64;
114        /// Regularized incomplete beta function.
115        pub fn incbet(a: f64, b: f64, x: f64) -> f64;
116        /// Inverse of incomplete beta integral.
117        pub fn incbi(a: f64, b: f64, y: f64) -> f64;
118
119        // Gamma functions
120        /// Gamma function.
121        pub fn gamma(x: f64) -> f64;
122        /// Reciprocal gamma function.
123        pub fn rgamma(x: f64) -> f64;
124        /// Natural logarithm of gamma function.
125        pub fn lgam(x: f64) -> f64;
126        /// Regularized incomplete gamma integral.
127        pub fn igam(a: f64, x: f64) -> f64;
128        /// Complemented incomplete gamma integral.
129        pub fn igamc(a: f64, x: f64) -> f64;
130        /// Inverse of complemented incomplete gamma integral.
131        pub fn igami(a: f64, p: f64) -> f64;
132        /// Psi (digamma) function.
133        pub fn psi(x: f64) -> f64;
134        /// Factorial function.
135        pub fn fac(i: i32) -> f64;
136
137        // Bessel functions
138        /// Bessel function of order zero.
139        pub fn j0(x: f64) -> f64;
140        /// Bessel function of order one.
141        pub fn j1(x: f64) -> f64;
142        /// Bessel function of integer order.
143        pub fn jn(n: i32, x: f64) -> f64;
144        /// Bessel function of real order.
145        pub fn jv(n: f64, x: f64) -> f64;
146
147        /// Bessel function of the second kind, order zero.
148        pub fn y0(x: f64) -> f64;
149        /// Bessel function of the second kind, order one.
150        pub fn y1(x: f64) -> f64;
151        /// Bessel function of the second kind, integer order.
152        pub fn yn(n: i32, x: f64) -> f64;
153        /// Bessel function of the second kind, real order.
154        pub fn yv(v: f64, x: f64) -> f64;
155
156        /// Modified Bessel function of order zero.
157        pub fn i0(x: f64) -> f64;
158        /// Modified Bessel function of order zero, exponentially scaled.
159        pub fn i0e(x: f64) -> f64;
160        /// Modified Bessel function of order one.
161        pub fn i1(x: f64) -> f64;
162        /// Modified Bessel function of order one, exponentially scaled.
163        pub fn i1e(x: f64) -> f64;
164        /// Modified Bessel function of real order.
165        pub fn iv(v: f64, x: f64) -> f64;
166
167        /// Modified Bessel function of the third kind, order zero.
168        pub fn k0(x: f64) -> f64;
169        /// Modified Bessel function of the third kind, order zero,
170        /// exponentially scaled.
171        pub fn k0e(x: f64) -> f64;
172        /// Modified Bessel function of the third kind, order one.
173        pub fn k1(x: f64) -> f64;
174        /// Modified Bessel function of the third kind, order one,
175        /// exponentially scaled.
176        pub fn k1e(x: f64) -> f64;
177        /// Modified Bessel function of the third kind, integer order.
178        pub fn kn(n: i32, x: f64) -> f64;
179
180        // Elliptic functions
181        /// Incomplete elliptic integral of the first kind.
182        pub fn ellik(phi: f64, m: f64) -> f64;
183        /// Incomplete elliptic integral of the second kind.
184        pub fn ellie(phi: f64, m: f64) -> f64;
185        /// Complete elliptic integral of the first kind.
186        pub fn ellpk(m1: f64) -> f64;
187        /// Complete elliptic integral of the second kind.
188        pub fn ellpe(m1: f64) -> f64;
189        /// Jacobian elliptic function.
190        pub fn ellpj(u: f64, m: f64, sn: &mut f64, cn: &mut f64, dn: &mut f64, phi: &mut f64) -> i32;
191
192        // Hypergeometric functions
193        /// Confluent hypergeometric function 1F1.
194        pub fn hyperg(a: f64, b: f64, x: f64) -> f64;
195        /// Hypergeometric function 1F2.
196        pub fn onef2(a: f64, b: f64, c: f64, x: f64, err: &mut f64) -> f64;
197        /// Gauss hypergeometric function 2F1.
198        pub fn hyp2f1(a: f64, b: f64, c: f64, x: f64) -> f64;
199        /// Hypergeometric function 3F0.
200        pub fn threef0(a: f64, b: f64, c: f64, x: f64, err: &mut f64) -> f64;
201
202        // Distributions
203        /// Binomial distribution.
204        pub fn bdtr(k: i32, n: i32, p: f64) -> f64;
205        /// Complemented binomial distribution.
206        pub fn bdtrc(k: i32, n: i32, p: f64) -> f64;
207        /// Inverse of binomial distribution.
208        pub fn bdtri(k: i32, n: i32, y: f64) -> f64;
209
210        /// Negative binomial distribution.
211        pub fn nbdtr(k: i32, n: i32, p: f64) -> f64;
212        /// Complemented negative binomial distribution.
213        pub fn nbdtrc(k: i32, n: i32, p: f64) -> f64;
214        /// Inverse of negative binomial distribution.
215        pub fn nbdtri(k: i32, n: i32, p: f64) -> f64;
216
217        /// Beta distribution.
218        pub fn btdtr(a: f64, b: f64, x: f64) -> f64;
219
220        /// Chi-square distribution.
221        pub fn chdtr(df: f64, x: f64) -> f64;
222        /// Complemented chi-square distribution.
223        pub fn chdtrc(v: f64, x: f64) -> f64;
224        /// Inverse of complemented chi-square distribution.
225        pub fn chdtri(df: f64, y: f64) -> f64;
226
227        /// F distribution.
228        pub fn fdtr(df1: i32, df2: i32, x: f64) -> f64;
229        /// Complemented F distribution.
230        pub fn fdtrc(df1: i32, df2: i32, x: f64) -> f64;
231        /// Inverse of complemented F distribution.
232        pub fn fdtri(df1: i32, df2: i32, p: f64) -> f64;
233
234        /// Gamma distribution.
235        pub fn gdtr(a: f64, b: f64, x: f64) -> f64;
236        /// Complemented gamma distribution.
237        pub fn gdtrc(a: f64, b: f64, x: f64) -> f64;
238
239        /// Normal distribution.
240        pub fn ndtr(x: f64) -> f64;
241        /// Inverse of normal distribution.
242        pub fn ndtri(y: f64) -> f64;
243
244        /// Poisson distribution.
245        pub fn pdtr(k: i32, m: f64) -> f64;
246        /// Complemented Poisson distribution.
247        pub fn pdtrc(k: i32, m: f64) -> f64;
248        /// Inverse of Poisson distribution.
249        pub fn pdtri(k: i32, y: f64) -> f64;
250
251        /// Student's t distribution.
252        pub fn stdtr(k: i16, t: f64) -> f64;
253        /// Inverse of Student's t distribution.
254        pub fn stdtri(k: i32, p: f64) -> f64;
255
256        // Misc special functions
257        /// Airy function.
258        pub fn airy(x: f64, ai: &mut f64, aip: &mut f64, bi: &mut f64, bip: &mut f64) -> i32;
259        /// Dawson's integral.
260        pub fn dawsn(x: f64) -> f64;
261        /// Fresnel integral.
262        pub fn fresnl(x: f64, s: &mut f64, c: &mut f64);
263        /// Integral of Planck's black body radiation formula.
264        pub fn plancki(lambda: f64, temperature: f64) -> f64;
265        /// Struve function.
266        pub fn struve(v: f64, x: f64) -> f64;
267        /// Riemann zeta function.
268        pub fn zetac(x: f64) -> f64;
269        /// Riemann zeta function of two arguments.
270        pub fn zeta(x: f64, q: f64) -> f64;
271    }
272}
273
274/// High-level bindings to double-precision Cephes functions.
275///
276/// This provides safe access to a subset of the Cephes functions.
277pub mod cephes_double {
278    use unsafe_cephes_double;
279
280    /// Function to compute the base-10 exponential of x
281    ///
282    /// # Original Description from Stephen L. Moshier
283    /// Returns 10 raised to the x power.
284    ///
285    /// Range reduction is accomplished by expressing the argument
286    /// as 10^x = 2^n 10^f, with |f| < 0.5 log10(2).
287    /// The Pade' form
288    ///
289    /// 10^x = 1 + 2x P( x^2)/( Q( x^2 ) - P( x^2 ) )
290    ///
291    /// is used to approximate 10^f.
292    ///
293    /// ACCURACY (Relative error):
294    ///
295    /// | arithmetic | domain | # trials | peak | rms |
296    /// | :--------: | :----: | :-------: | :---: | :--: |
297    /// | IEEE       | -307,+307 | 30000 | 2.2e-16 | 5.5e-17 |
298    ///
299    /// ERROR MESSAGES:
300    ///
301    /// | message | condition | value returned |
302    /// | :-----: | :-------: | :------------: |
303    /// | exp10 underflow | x < -MAXL10 | 0.0  |
304    /// | exp10 overflow | x > MAXL10 | MAXNUM |
305    ///
306    /// IEEE arithmetic: MAXL10 = 308.2547155599167.
307    ///
308    ///
309    /// # Arguments
310    /// * `x`: A double precision number
311    ///
312    /// # Example
313    /// ```
314    /// use special_fun::cephes_double::exp10;
315    /// let x = 1.0f64;
316    /// exp10(x);
317    /// ```
318    pub fn exp10(x: f64) -> f64 {
319        unsafe { unsafe_cephes_double::exp10(x) }
320    }
321
322    /// Function to accurately compute `exp(x) - 1` for small x
323    ///
324    /// # Arguments
325    /// * `x`: A double precision number
326    ///
327    /// # Example
328    /// ```
329    /// use special_fun::cephes_double::expm1;
330    /// let x = 1.0f64;
331    /// expm1(x);
332    /// ```
333    pub fn expm1(x: f64) -> f64 {
334        unsafe { unsafe_cephes_double::expm1(x) }
335    }
336
337    /// Function to accurately compute the exponential of a squared argument
338    ///
339    /// # Original Description from Stephen L. Moshier
340    ///
341    /// Computes y = exp(x*x) while suppressing error amplification
342    /// that would ordinarily arise from the inexactness of the
343    /// exponential argument x*x.
344    ///
345    /// If sign < 0, the result is inverted; i.e., y = exp(-x*x) .
346    ///
347    /// ACCURACY (Relative error):
348    ///
349    /// | arithmetic | domain | # trials | peak | rms |
350    /// | :--------: | :----: | :------: | :--: | :--:|
351    /// | IEEE | -26.6, 26.6 | 10^7 | 3.9e-16 | 8.9e-17 |
352    ///
353    /// # Arguments
354    /// * `x`: A double precision number
355    /// * 'sign': An integer representing the sign of the exponent
356    ///
357    /// # Example
358    /// ```
359    /// use special_fun::cephes_double::expx2;
360    /// let x = 1.0f64;
361    /// expx2(x, -1);
362    /// ```
363    pub fn expx2(x: f64, sign: i32) -> f64 {
364        unsafe { unsafe_cephes_double::expx2(x, sign) }
365    }
366
367    /// Function to accurately compute the exponential integral of x
368    ///
369    /// The exponential integral is given by:
370    ///
371    /// $Ei(x) = -\int^{\infty}_{-x} \frac{e^{-t}}{t} dt$
372    ///
373    /// # Original Description from Stephen L. Moshier
374    ///
375    /// ACCURACY (Relative error):
376    ///
377    /// | arithmetic | domain | # trials | peak | rms |
378    /// | :--------: | :----: | :------: | :--: | :--:|
379    /// | IEEE | 0, 100 | 50000 | 8.6e-16 | 1.3e-16 |
380    ///
381    /// # Arguments
382    /// * `x`: A double precision number
383    ///
384    /// # Example
385    /// ```
386    /// use special_fun::cephes_double::ei;
387    /// let x = 1.0f64;
388    /// ei(x);
389    /// ```
390    pub fn ei(x: f64) -> f64 {
391        unsafe { unsafe_cephes_double::ei(x) }
392    }
393
394    /// Function to accurately compute the error function of x
395    ///
396    /// The error function is given by:
397    ///
398    /// $Erf(x) = \frac{2}{\sqrt{pi}} -\int^{x}_{0} e^{-t^2} dt$
399    ///
400    /// # Arguments
401    /// * `x`: A double precision number
402    ///
403    /// # Example
404    /// ```
405    /// use special_fun::cephes_double::erf;
406    /// let x = 1.0f64;
407    /// erf(x);
408    /// ```
409    pub fn erf(x: f64) -> f64 {
410        unsafe { unsafe_cephes_double::erf(x) }
411    }
412
413    /// Function to accurately compute the complementary error function of x
414    ///
415    /// The error function is given by:
416    ///
417    /// $Erf(x) = 1 + \frac{2}{\sqrt{pi}} -\int^{x}_{0} e^{-t^2} dt$
418    ///
419    /// # Arguments
420    /// * `x`: A double precision number
421    ///
422    /// # Example
423    /// ```
424    /// use special_fun::cephes_double::erfc;
425    /// let x = 1.0f64;
426    /// erfc(x);
427    /// ```
428    pub fn erfc(x: f64) -> f64 {
429        unsafe { unsafe_cephes_double::erfc(x) }
430    }
431
432    /// Function to accurately compute `log(1 + x)`
433    ///
434    /// # Arguments
435    /// * `x`: A double precision number
436    ///
437    /// # Example
438    /// ```
439    /// use special_fun::cephes_double::log1p;
440    /// let x = 0.1f64;
441    /// log1p(x);
442    /// ```
443    pub fn log1p(x: f64) -> f64 {
444        unsafe { unsafe_cephes_double::log1p(x) }
445    }
446
447    /// Function to accurately compute the dilogarithm function
448    ///
449    /// Spence's function is a special case of the polylogarithm function, and
450    /// is defined as:
451    ///
452    /// $ Li_2(x) = \int^x_1 \frac{ln(t)}{t - 1} du
453    ///
454    /// Note, this implies that the domain has been shifted by 1 from the
455    /// standard definition (as per wikipedia) of:
456    ///
457    /// $ Li_2(x) = - \int^x_0 \frac{ ln(1 - t) }{ t } dt
458    ///
459    /// # Original Description from Stephen L. Moshier
460    /// Computes the integral for x >= 0.  A rational approximation gives the
461    /// integral in the interval (0.5, 1.5).  Transformation formulas for 1/x
462    /// and 1-x are employed outside the basic expansion range.
463    ///
464    /// ACCURACY(Relative error):
465    ///
466    /// | arithmetic | domain | # trials | peak | rms |
467    /// | :---: | :--: | :--: | :--: | :--: | :--: |
468    /// | IEEE  | 0,4 | 30000 | 3.9e-15 | 5.4e-16 |
469    ///
470    /// # Arguments
471    /// * `x`: A double precision number
472    ///
473    /// # Example
474    /// ```
475    /// use special_fun::cephes_double::spence;
476    /// let x = 0.1f64;
477    /// spence(x);
478    /// ```
479    pub fn spence(x: f64) -> f64 {
480        unsafe { unsafe_cephes_double::spence(x) }
481    }
482
483    /// Function to accurately compute `cos(x) - 1` for small x
484    ///
485    /// # Arguments
486    /// * `x`: A double precision number
487    ///
488    /// # Example
489    /// ```
490    /// use special_fun::cephes_double::cosm1;
491    /// let x = 0.1f64;
492    /// cosm1(x);
493    /// ```
494    pub fn cosm1(x: f64) -> f64 {
495        unsafe { unsafe_cephes_double::cosm1(x) }
496    }
497
498    /// Function to accurately compute sine and cosine integrals
499    ///
500    /// The sine integral is defined as:
501    ///
502    /// $ Si(x) = \int_0^{\infty} \frac{ sin(t) }{t} dt $
503    ///
504    /// and the cosine integral is defined as:
505    ///
506    /// $ Ci(x) = -\int_x^{\infty} \frac{ cos(t) }{t} dt $
507    ///
508    /// which reduces to:
509    ///
510    /// $ Ci(x) = \gamma + ln(x) + \int_0^x \frac{cos(t) - 1}{t} dt $
511    ///
512    /// where $\gamma$ is the euler constant.
513    ///
514    /// # Original Description from Stephen L. Moshier
515    /// The integrals are approximated by rational functions.
516    /// For x > 8 auxiliary functions f(x) and g(x) are employed
517    /// such that
518    ///
519    /// Ci(x) = f(x) sin(x) - g(x) cos(x)
520    /// Si(x) = pi/2 - f(x) cos(x) - g(x) sin(x)
521    ///
522    /// ACCURACY(Absolute error, except relative when > 1):
523    ///
524    /// | arithmetic | Function | domain | # trials | peak | rms |
525    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
526    /// | IEEE | Si | 0,50 | 30000 | 4.4e-16 | 7.3e-17 |
527    /// | IEEE | Ci | 0,50 | 30000 | 6.9e-16 | 5.1e-17 |
528    ///
529    /// # Arguments
530    /// * `x`: A double precision number
531    /// * `Si`: A mutable double precision number that will return the sine
532    ///         integral value
533    /// * `Ci`: A mutable double precision number that will return the cosine
534    ///         integral value
535    ///
536    /// # Example
537    /// ```
538    /// use special_fun::cephes_double::sici;
539    /// let x = 0.1f64;
540    /// let mut si = 0.0_f64;
541    /// let mut ci = 0.0_f64;
542    /// let mut code = 0.0_f64;
543    /// sici(0.5_f64, &mut si, &mut ci);
544    /// ```
545    pub fn sici(x: f64, si: &mut f64, ci: &mut f64) {
546        let code = unsafe { unsafe_cephes_double::sici(x, si, ci) };
547        assert_eq!(code, 0);
548    }
549
550    /// Function to accurately compute hyperbolic sine and cosine integrals
551    ///
552    /// The hyperbolic sine integral is defined as:
553    ///
554    /// $ Shi(x) = \int_0^{\infty} \frac{ sinh(t) }{t} dt $
555    ///
556    /// and the hyperbolic cosine integral is defined as:
557    ///
558    /// $ Chi(x) = -\int_x^{\infty} \frac{ cosh(t) }{t} dt $
559    ///
560    /// which reduces to:
561    ///
562    /// $ Chi(x) = \gamma + ln(x) + \int_0^x \frac{cosh(t) - 1}{t} dt $
563    ///
564    /// where $\gamma$ is the euler constant.
565    ///
566    /// # Original Description from Stephen L. Moshier
567    /// The integrals are approximated by rational functions.
568    /// For x > 8 auxiliary functions f(x) and g(x) are employed
569    /// such that
570    ///
571    /// Ci(x) = f(x) sin(x) - g(x) cos(x)
572    /// Si(x) = pi/2 - f(x) cos(x) - g(x) sin(x)
573    ///
574    /// ACCURACY(Relative error)):
575    ///
576    /// | arithmetic | Function | domain | # trials | peak | rms |
577    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
578    /// | IEEE | Shi | 0,88 | 30000 | 6.9-16 | 1.6e-16 |
579    ///
580    /// ACCURACY(Absolute error, except relative when |Chi| > 1):
581    ///
582    /// | arithmetic | Function | domain | # trials | peak | rms |
583    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
584    /// | IEEE | Chi | 0,88 | 30000 | 8.4e-16 | 1.4e-16 |
585    ///
586    /// # Arguments
587    /// * `x`: A double precision number
588    /// * `Shi`: A mutable double precision number that will return the
589    ///          hyperbolic sine integral value
590    /// * `Chi`: A mutable double precision number that will return the
591    ///          hyperbolic cosine integral value
592    ///
593    /// # Example
594    /// ```
595    /// use special_fun::cephes_double::shichi;
596    /// let x = 0.1f64;
597    /// let mut shi = 0.0_f64;
598    /// let mut chi = 0.0_f64;
599    /// shichi(0.5_f64, &mut shi, &mut chi);
600    /// ```
601    pub fn shichi(x: f64, chi: &mut f64, shi: &mut f64){
602        let code = unsafe { unsafe_cephes_double::shichi(x, chi, shi) };
603        assert_eq!(code, 0);
604    }
605
606    /// Function to compute the beta function.
607    ///
608    /// The beta function is given by:
609    ///
610    /// $B(a, b) = \int_0^1 t^{a - 1} ( 1 - t )^{b - 1} dt $
611    ///
612    /// which simplifies to
613    ///
614    /// $ B(a, b) = \frac{ \Gamma(a) \Gamma(b)}{ \Gamma(a + b) } $
615    ///
616    /// # Original Description from Stephen L. Moshier
617    /// For large arguments the logarithm of the function is
618    /// evaluated using lgam(), then exponentiated.
619    ///
620    /// ACCURACY (Relative error):
621    ///
622    /// | arithmetic | domain | # trials | peak | rms |
623    /// | :--------: | :----: | :------: | :--: | :-: |
624    /// | IEEE       | 0,30   | 30000    | 8.1e-14 | 1.1e-14 |
625    ///
626    /// ERROR MESSAGES:
627    ///
628    /// | message   | condition | value returned |
629    /// | :-------: | :--------: | :----------: |
630    /// | beta overflow | log(beta) > MAXLOG | 0.0 |
631    /// | beta overflow | a or b <0 integer      | 0.0 |
632    ///
633    /// # Arguments
634    /// * `a`: A double precision parameter, assumed to be greater than 0
635    /// * `b`: A double precision parameter, assumed to be greater than 0
636    ///
637    /// # Example
638    /// ```
639    /// use special_fun::cephes_double::beta;
640    /// let a = 0.1f64;
641    /// let b = 0.2f64;
642    /// beta(a, b);
643    /// ```
644    pub fn beta(a: f64, b: f64) -> f64 {
645        unsafe { unsafe_cephes_double::beta(a, b) }
646    }
647
648    /// Function to compute the regularized incomplete beta function.
649    ///
650    /// The regularized incomplete beta function is given by:
651    ///
652    /// $I_x(a, b) = \frac{1}{B(a, b)} \int_0^x t^{a - 1} ( 1 - t )^{b - 1} dt $
653    ///
654    /// and for $x = 1$, the regularized incomplete beta function evaluates to
655    /// exactly 1. Note, $ 1 - I_x(a, b) = I_{1 - x}(a, b)$.
656    ///
657    /// # Original Description from Stephen L. Moshier
658    /// The integral is evaluated by a continued fraction expansion
659    /// or, when b*x is small, by a power series.
660    ///
661    /// ACCURACY (Relative error):
662    ///
663    /// | arithmetic | domain | # trials | peak | rms |
664    /// | :--------: | :----: | :------: | :--: | :-: |
665    /// | IEEE       | 0,5    | 10000    | 6.9e-15 | 4.5e-16 |
666    /// | IEEE | 0,85 | 250000 | 2.2e-13 | 1.7e-14 |
667    /// | IEEE | 0,1000   |  30000  |  5.3e-12 | 6.3e-13 |
668    /// | IEEE | 0,10000  | 250000  |  9.3e-11 | 7.1e-12 |
669    /// | IEEE | 0,100000 |  10000  |  8.7e-10 | 4.8e-11 |
670    ///
671    /// ERROR MESSAGES:
672    ///
673    /// | message   | condition | value returned |
674    /// | :-------: | :--------: | :----------: |
675    /// | incbet domain | x<0, x>1 | 0.0 |
676    /// | incbet underflow | - | 0.0 |
677    ///
678    /// # Arguments
679    /// * `a`: A double precision parameter, assumed to be greater than 0
680    /// * `b`: A double precision parameter, assumed to be greater than 0
681    /// * `x`: A double precision parameter, assumed to be between 0 and 1
682    ///
683    /// # Example
684    /// ```
685    /// use special_fun::cephes_double::incbet;
686    /// let a = 0.1f64;
687    /// let b = 0.2f64;
688    /// incbet(a, b, 0.5f64);
689    /// ```
690    pub fn incbet(a: f64, b: f64, x: f64) -> f64 {
691        unsafe { unsafe_cephes_double::incbet(a, b, x) }
692    }
693
694    /// Inverse of incomplete beta integral.
695    pub fn incbi(a: f64, b: f64, y: f64) -> f64 {
696        unsafe { unsafe_cephes_double::incbi(a, b, y) }
697    }
698
699    /// Gamma function.
700    pub fn gamma(x: f64) -> f64 {
701        unsafe { unsafe_cephes_double::gamma(x) }
702    }
703
704    /// Reciprocal gamma function.
705    pub fn rgamma(x: f64) -> f64 {
706        unsafe { unsafe_cephes_double::rgamma(x) }
707    }
708
709    /// Natural logarithm of gamma function.
710    pub fn lgam(x: f64) -> f64 {
711        unsafe { unsafe_cephes_double::lgam(x) }
712    }
713
714    /// Regularized incomplete gamma integral.
715    pub fn igam(a: f64, x: f64) -> f64 {
716        unsafe { unsafe_cephes_double::igam(a, x) }
717    }
718
719    /// Complemented incomplete gamma integral.
720    pub fn igamc(a: f64, x: f64) -> f64 {
721        unsafe { unsafe_cephes_double::igamc(a, x) }
722    }
723
724    /// Inverse of complemented incomplete gamma integral.
725    pub fn igami(a: f64, p: f64) -> f64 {
726        unsafe { unsafe_cephes_double::igami(a, p) }
727    }
728
729    /// Psi (digamma) function.
730    pub fn psi(x: f64) -> f64 {
731        unsafe { unsafe_cephes_double::psi(x) }
732    }
733
734    /// Factorial function.
735    pub fn fac(i: i32) -> f64 {
736        unsafe { unsafe_cephes_double::fac(i) }
737    }
738}
739
740/// Low-level bindings to single-precision Cephes functions.
741///
742/// This provides access to the C function implemented in Cephes for maximum
743/// flexibility. Note that Cephes implements less functions for single than for
744/// double precision.
745#[allow(dead_code)]
746pub mod unsafe_cephes_single {
747    extern "C" {
748        // Floating point numeric utilities
749        /// Round to nearest or event integer valued f32.
750        pub fn roundf(x: f32) -> f32;
751        /// Largest integer less than or equal to x.
752        pub fn floorf(x: f32) -> f32;
753        /// Smallest integer greater than or equal to x.
754        pub fn ceilf(x: f32) -> f32;
755        /// Return the significand between 0.5 and 1. Write exponent to expnt.
756        /// x = y * 2**expn
757        pub fn frexpf(x: f32, expnt: &mut i32) -> f32;
758        /// Multiply x by 2**n.
759        pub fn ldexpf(x: f32, n: i32) -> f32;
760        /// Absolute value.
761        pub fn fabsf(x: f32) -> f32;
762        /// Return 1 if the sign bit of x is 1, else 0.
763        pub fn signbitf(x: f32) -> i32;
764        /// Return 1 if x is NaN, else 0.
765        pub fn isnanf(x: f32) -> i32;
766        /// Return 1 if x is finite, else 0.
767        pub fn isfinitef(x: f32) -> i32;
768
769        // Roots
770        /// Cube root.
771        pub fn cbrtf(x: f32) -> f32;
772        /// Square root.
773        pub fn sqrtf(x: f32) -> f32;
774        /// Integer square root.
775        pub fn lsqrtf(x: i64) -> i64;
776
777        // Exponential functions
778        /// Exponential function.
779        pub fn expf(x: f32) -> f32;
780        /// Base 10 exponential function.
781        pub fn exp10f(x: f32) -> f32;
782        /// Base 2 exponential function.
783        pub fn exp2f(x: f32) -> f32;
784        /// Compute accurately exponential of squared argument.
785        pub fn expm1f(x: f32) -> f32;
786        /// Compute accurately exp(x) - 1 for x close to 0.
787        pub fn expx2f(x: f32, sign: i32) -> f32;
788        /// Error function.
789        pub fn erff(x: f32) -> f32;
790        /// Complementary error function.
791        pub fn erfcf(x: f32) -> f32;
792        /// Power function.
793        pub fn powf(x: f32, y: f32) -> f32;
794        /// Integer power function.
795        pub fn powif(x: f32, n: i32) -> f32;
796
797        // Logarithmic functions
798        /// Natural logarithm.
799        pub fn logf(x: f32) -> f32;
800        /// Common logarithm.
801        pub fn log10f(x: f32) -> f32;
802        /// Base 2 logarithm.
803        pub fn log2f(x: f32) -> f32;
804        /// Compute accurately log(1 + x) for x close to 0.
805        pub fn log1pf(x: f32) -> f32;
806        /// Dilogarithm (Spence's function).
807        pub fn spencef(x: f32) -> f32;
808
809        // Trigonometric functions
810        /// Circular sine.
811        pub fn sinf(x: f32) -> f32;
812        /// Circular cosine.
813        pub fn cosf(x: f32) -> f32;
814        /// Circular tangent.
815        pub fn tanf(x: f32) -> f32;
816        /// Inverse circular sine.
817        pub fn asinf(x: f32) -> f32;
818        /// Inverse circular cosine.
819        pub fn acosf(x: f32) -> f32;
820        /// Inverse circular tangent.
821        pub fn atanf(x: f32) -> f32;
822        /// Quadrant-correct inverse circular tangent.
823        pub fn atan2f(y: f32, x: f32) -> f32;
824        /// Sine and cosine integrals.
825        pub fn sicif(x: f32, si: &mut f32, ci: &mut f32) -> i32;
826
827        // Hyperbolic functions
828        /// Hyperbolic sine.
829        pub fn sinhf(x: f32) -> f32;
830        /// Hyperbolic cosine.
831        pub fn coshf(x: f32) -> f32;
832        /// Hyperbolic tangent.
833        pub fn tanhf(x: f32) -> f32;
834        /// Inverse hyperbolic sine.
835        pub fn asinhf(x: f32) -> f32;
836        /// Inverse hyperbolic cosine.
837        pub fn acoshf(x: f32) -> f32;
838        /// Inverse hyperbolic tangent.
839        pub fn atanhf(x: f32) -> f32;
840        /// Hyperbolic sine and cosine integrals.
841        pub fn shichif(x: f32, chi: &mut f32, shi: &mut f32) -> i32;
842
843        // Beta functions
844        /// Beta function.
845        pub fn betaf(a: f32, b: f32) -> f32;
846        /// Regularized incomplete beta function.
847        pub fn incbetf(a: f32, b: f32, x: f32) -> f32;
848        /// Inverse of incomplete beta integral.
849        pub fn incbif(a: f32, b: f32, y: f32) -> f32;
850
851        // Gamma functions
852        /// Gamma function.
853        pub fn gammaf(x: f32) -> f32;
854        /// Reciprocal gamma function.
855        pub fn rgammaf(x: f32) -> f32;
856        /// Natural logarithm of gamma function.
857        pub fn lgamf(x: f32) -> f32;
858        /// Regularized incomplete gamma integral.
859        pub fn igamf(a: f32, x: f32) -> f32;
860        /// Complemented incomplete gamma integral.
861        pub fn igamcf(a: f32, x: f32) -> f32;
862        /// Inverse of complemented incomplete gamma integral.
863        pub fn igamif(a: f32, p: f32) -> f32;
864        /// Psi (digamma) function.
865        pub fn psif(x: f32) -> f32;
866        /// Factorial function.
867        pub fn facf(i: i32) -> f32;
868
869        // Bessel functions
870        /// Bessel function of order zero.
871        pub fn j0f(x: f32) -> f32;
872        /// Bessel function of order one.
873        pub fn j1f(x: f32) -> f32;
874        /// Bessel function of integer order.
875        pub fn jnf(n: i32, x: f32) -> f32;
876        /// Bessel function of real order.
877        pub fn jvf(n: f32, x: f32) -> f32;
878
879        /// Bessel function of the second kind, order zero.
880        pub fn y0f(x: f32) -> f32;
881        /// Bessel function of the second kind, order one.
882        pub fn y1f(x: f32) -> f32;
883        /// Bessel function of the second kind, integer order.
884        pub fn ynf(n: i32, x: f32) -> f32;
885        /// Bessel function of the second kind, real order.
886        pub fn yvf(v: f32, x: f32) -> f32;
887
888        /// Modified Bessel function of order zero.
889        pub fn i0f(x: f32) -> f32;
890        /// Modified Bessel function of order zero, exponentially scaled.
891        pub fn i0ef(x: f32) -> f32;
892        /// Modified Bessel function of order one.
893        pub fn i1f(x: f32) -> f32;
894        /// Modified Bessel function of order one, exponentially scaled.
895        pub fn i1ef(x: f32) -> f32;
896        /// Modified Bessel function of real order.
897        pub fn ivf(v: f32, x: f32) -> f32;
898
899        /// Modified Bessel function of the third kind, order zero.
900        pub fn k0f(x: f32) -> f32;
901        /// Modified Bessel function of the third kind, order zero,
902        /// exponentially scaled.
903        pub fn k0ef(x: f32) -> f32;
904        /// Modified Bessel function of the third kind, order one.
905        pub fn k1f(x: f32) -> f32;
906        /// Modified Bessel function of the third kind, order one,
907        /// exponentially scaled.
908        pub fn k1ef(x: f32) -> f32;
909        /// Modified Bessel function of the third kind, integer order.
910        pub fn knf(n: i32, x: f32) -> f32;
911
912        // Elliptic functions
913        /// Incomplete elliptic integral of the first kind.
914        pub fn ellikf(phi: f32, m: f32) -> f32;
915        /// Incomplete elliptic integral of the second kind.
916        pub fn ellief(phi: f32, m: f32) -> f32;
917        /// Complete elliptic integral of the first kind.
918        pub fn ellpkf(m1: f32) -> f32;
919        /// Complete elliptic integral of the second kind.
920        pub fn ellpef(m1: f32) -> f32;
921        /// Jacobian elliptic function.
922        pub fn ellpjf(u: f32, m: f32, sn: &mut f32, cn: &mut f32, dn: &mut f32, phi: &mut f32) -> i32;
923
924        // Hypergeometric functions
925        /// Confluent hypergeometric function 1F1.
926        pub fn hypergf(a: f32, b: f32, x: f32) -> f32;
927        /// Hypergeometric function 1F2.
928        pub fn onef2f(a: f32, b: f32, c: f32, x: f32, err: &mut f32) -> f32;
929        /// Gauss hypergeometric function 2F1.
930        pub fn hyp2f1f(a: f32, b: f32, c: f32, x: f32) -> f32;
931        /// Hypergeometric function 3F0.
932        pub fn threef0f(a: f32, b: f32, c: f32, x: f32, err: &mut f32) -> f32;
933
934        // Distributions
935        /// Binomial distribution.
936        pub fn bdtrf(k: i32, n: i32, p: f32) -> f32;
937        /// Complemented binomial distribution.
938        pub fn bdtrcf(k: i32, n: i32, p: f32) -> f32;
939        /// Inverse of binomial distribution.
940        pub fn bdtrif(k: i32, n: i32, y: f32) -> f32;
941
942        /// Negative binomial distribution.
943        pub fn nbdtrf(k: i32, n: i32, p: f32) -> f32;
944        /// Complemented negative binomial distribution.
945        pub fn nbdtrcf(k: i32, n: i32, p: f32) -> f32;
946        /// Inverse of negative binomial distribution.
947        pub fn nbdtrif(k: i32, n: i32, p: f32) -> f32;
948
949        /// Beta distribution.
950        pub fn btdtrf(a: f32, b: f32, x: f32) -> f32;
951
952        /// Chi-square distribution.
953        pub fn chdtrf(df: f32, x: f32) -> f32;
954        /// Complemented chi-square distribution.
955        pub fn chdtrcf(v: f32, x: f32) -> f32;
956        /// Inverse of complemented chi-square distribution.
957        pub fn chdtrif(df: f32, y: f32) -> f32;
958
959        /// F distribution.
960        pub fn fdtrf(df1: i32, df2: i32, x: f32) -> f32;
961        /// Complemented F distribution.
962        pub fn fdtrcf(df1: i32, df2: i32, x: f32) -> f32;
963        /// Inverse of complemented F distribution.
964        pub fn fdtrif(df1: i32, df2: i32, p: f32) -> f32;
965
966        /// Gamma distribution.
967        pub fn gdtrf(a: f32, b: f32, x: f32) -> f32;
968        /// Complemented gamma distribution.
969        pub fn gdtrcf(a: f32, b: f32, x: f32) -> f32;
970
971        /// Normal distribution.
972        pub fn ndtrf(x: f32) -> f32;
973        /// Inverse of normal distribution.
974        pub fn ndtrif(y: f32) -> f32;
975
976        /// Poisson distribution.
977        pub fn pdtrf(k: i32, m: f32) -> f32;
978        /// Complemented Poisson distribution.
979        pub fn pdtrcf(k: i32, m: f32) -> f32;
980        /// Inverse of Poisson distribution.
981        pub fn pdtrif(k: i32, y: f32) -> f32;
982
983        /// Student's t distribution.
984        pub fn stdtrf(k: i16, t: f32) -> f32;
985        /// Inverse of Student's t distribution.
986        pub fn stdtrif(k: i32, p: f32) -> f32;
987
988        // Misc special functions
989        /// Airy function.
990        pub fn airyf(x: f32, ai: &mut f32, aip: &mut f32, bi: &mut f32, bip: &mut f32) -> i32;
991        /// Dawson's integral.
992        pub fn dawsnf(x: f32) -> f32;
993        /// Fresnel integral.
994        pub fn fresnlf(x: f32, s: &mut f32, c: &mut f32);
995        /// Integral of Planck's black body radiation formula.
996        pub fn planckif(lambda: f32, temperature: f32) -> f32;
997        /// Struve function.
998        pub fn struvef(v: f32, x: f32) -> f32;
999        /// Riemann zeta function.
1000        pub fn zetacf(x: f32) -> f32;
1001        /// Riemann zeta function of two arguments.
1002        pub fn zetaf(x: f32, q: f32) -> f32;
1003    }
1004}
1005
1006/// High-level bindings to single-precision Cephes functions.
1007///
1008/// This provides safe access to a subset of the Cephes functions. Note that
1009/// Cephes implements less functions for single than for double precision.
1010pub mod cephes_single {
1011    use unsafe_cephes_single;
1012
1013    /// Function to compute the base-10 exponential of x
1014    ///
1015    /// # Original Description from Stephen L. Moshier
1016    /// Returns 10 raised to the x power.
1017    ///
1018    /// Range reduction is accomplished by expressing the argument
1019    /// as 10^x = 2^n 10^f, with |f| < 0.5 log10(2).
1020    /// The Pade' form
1021    ///
1022    /// 10^x = 1 + 2x P( x^2)/( Q( x^2 ) - P( x^2 ) )
1023    ///
1024    /// is used to approximate 10^f.
1025    ///
1026    /// ACCURACY (Relative error):
1027    ///
1028    /// | arithmetic | domain | # trials | peak | rms |
1029    /// | :--------: | :----: | :-------: | :---: | :--: |
1030    /// | IEEE       | -38,+38 | 100000 | 9.8e-8 | 2.8e-8 |
1031    ///
1032    /// ERROR MESSAGES:
1033    ///
1034    /// | message | condition | value returned |
1035    /// | :-----: | :-------: | :------------: |
1036    /// | exp10 underflow | x < -MAXL10 | 0.0  |
1037    /// | exp10 overflow | x > MAXL10 | MAXNUM |
1038    ///
1039    /// IEEE single arithmetic: MAXL10 = 38.230809449325611792.
1040    ///
1041    ///
1042    /// # Arguments
1043    /// * `x`: A single precision number
1044    ///
1045    /// # Example
1046    /// ```
1047    /// use special_fun::cephes_single::exp10;
1048    /// let x = 1.0f32;
1049    /// exp10(x);
1050    /// ```
1051    pub fn exp10(x: f32) -> f32 {
1052        unsafe { unsafe_cephes_single::exp10f(x) }
1053    }
1054
1055    /// Function to accurately compute `exp(x) - 1` for small x
1056    ///
1057    /// # Arguments
1058    /// * `x`: A single precision number
1059    ///
1060    /// # Example
1061    /// ```
1062    /// use special_fun::cephes_single::expm1;
1063    /// let x = 1.0f32;
1064    /// expm1(x);
1065    /// ```
1066    pub fn expm1(x: f32) -> f32 {
1067        unsafe { unsafe_cephes_single::expm1f(x) }
1068    }
1069
1070    /// Function to accurately compute the exponential of a squared argument
1071    ///
1072    /// # Original Description from Stephen L. Moshier
1073    ///
1074    /// Computes y = exp(x*x) while suppressing error amplification
1075    /// that would ordinarily arise from the inexactness of the
1076    /// exponential argument x*x.
1077    ///
1078    /// If sign < 0, the result is inverted; i.e., y = exp(-x*x) .
1079    ///
1080    /// ACCURACY (Relative error):
1081    ///
1082    /// | arithmetic | domain | # trials | peak | rms |
1083    /// | :--------: | :----: | :------: | :--: | :--:|
1084    /// | IEEE | -9.4, 9.4 | 10^7 | 1.7e-7 | 4.7e-8 |
1085    ///
1086    /// # Arguments
1087    /// * `x`: A single precision number
1088    /// * `sign`: An integer representing the sign of the exponent
1089    ///
1090    /// # Example
1091    /// ```
1092    /// use special_fun::cephes_single::expx2;
1093    /// let x = 1.0f32;
1094    /// expx2(x, -1);
1095    /// ```
1096    pub fn expx2(x: f32, sign: i32) -> f32 {
1097        unsafe { unsafe_cephes_single::expx2f(x, sign) }
1098    }
1099
1100    /// Function to accurately compute the error function of x
1101    ///
1102    /// The error function is given by:
1103    ///
1104    /// $Erf(x) = \frac{2}{\sqrt{pi}} -\int^{x}_{0} e^{-t^2} dt$
1105    ///
1106    /// # Arguments
1107    /// * `x`: A single precision number
1108    ///
1109    /// # Example
1110    /// ```
1111    /// use special_fun::cephes_single::erf;
1112    /// let x = 1.0f32;
1113    /// erf(x);
1114    /// ```
1115    pub fn erf(x: f32) -> f32 {
1116        unsafe { unsafe_cephes_single::erff(x) }
1117    }
1118
1119    /// Function to accurately compute the complementary error function of x
1120    ///
1121    /// The error function is given by:
1122    ///
1123    /// $Erf(x) = 1 + \frac{2}{\sqrt{pi}} -\int^{x}_{0} e^{-t^2} dt$
1124    ///
1125    /// # Arguments
1126    /// * `x`: A single precision number
1127    ///
1128    /// # Example
1129    /// ```
1130    /// use special_fun::cephes_single::erfc;
1131    /// let x = 1.0f32;
1132    /// erfc(x);
1133    /// ```
1134    pub fn erfc(x: f32) -> f32 {
1135        unsafe { unsafe_cephes_single::erfcf(x) }
1136    }
1137
1138    /// Function to accurately compute `log(1 + x)`
1139    ///
1140    /// # Arguments
1141    /// * `x`: A single precision number
1142    ///
1143    /// # Example
1144    /// ```
1145    /// use special_fun::cephes_single::log1p;
1146    /// let x = 0.1f32;
1147    /// log1p(x);
1148    /// ```
1149    pub fn log1p(x: f32) -> f32 {
1150        unsafe { unsafe_cephes_single::log1pf(x) }
1151    }
1152
1153    /// Function to accurately compute the dilogarithm function
1154    ///
1155    /// Spence's function is a special case of the polylogarithm function, and
1156    /// is defined as:
1157    ///
1158    /// $ Li_2(x) = \int^x_1 \frac{ln(t)}{t - 1} du
1159    ///
1160    /// Note, this implies that the domain has been shifted by 1 from the
1161    /// standard definition (as per wikipedia) of:
1162    ///
1163    /// $ Li_2(x) = - \int^x_0 \frac{ ln(1 - t) }{ t } dt
1164    ///
1165    /// # Original Description from Stephen L. Moshier
1166    /// Computes the integral for x >= 0.  A rational approximation gives the
1167    /// integral in the interval (0.5, 1.5).  Transformation formulas for 1/x
1168    /// and 1-x are employed outside the basic expansion range.
1169    ///
1170    /// ACCURACY(Relative error):
1171    ///
1172    /// | arithmetic | domain | # trials | peak | rms |
1173    /// | :---: | :--: | :--: | :--: | :--: | :--: |
1174    /// | IEEE  | 0,4 | 30000 | 4.4e-7 | 6.3e-8 |
1175    ///
1176    /// # Arguments
1177    /// * `x`: A single precision number
1178    ///
1179    /// # Example
1180    /// ```
1181    /// use special_fun::cephes_single::spence;
1182    /// let x = 0.1f32;
1183    /// spence(x);
1184    /// ```
1185    pub fn spence(x: f32) -> f32 {
1186        unsafe { unsafe_cephes_single::spencef(x) }
1187    }
1188
1189    /// Function to accurately compute sine and cosine integrals
1190    ///
1191    /// The sine integral is defined as:
1192    ///
1193    /// $ Si(x) = \int_0^{\infty} \frac{ sin(t) }{t} dt $
1194    ///
1195    /// and the cosine integral is defined as:
1196    ///
1197    /// $ Ci(x) = -\int_x^{\infty} \frac{ cos(t) }{t} dt $
1198    ///
1199    /// which reduces to:
1200    ///
1201    /// $ Ci(x) = \gamma + ln(x) + \int_0^x \frac{cos(t) - 1}{t} dt $
1202    ///
1203    /// where $\gamma$ is the euler constant.
1204    ///
1205    /// # Original Description from Stephen L. Moshier
1206    /// The integrals are approximated by rational functions.
1207    /// For x > 8 auxiliary functions f(x) and g(x) are employed
1208    /// such that
1209    ///
1210    /// Ci(x) = f(x) sin(x) - g(x) cos(x)
1211    /// Si(x) = pi/2 - f(x) cos(x) - g(x) sin(x)
1212    ///
1213    /// ACCURACY(Absolute error, except relative when > 1):
1214    ///
1215    /// | arithmetic | Function | domain | # trials | peak | rms |
1216    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
1217    /// | IEEE | Si | 0,50 | 30000 | 2.1e-7 | 4.3e-8 |
1218    /// | IEEE | Ci | 0,50 | 30000 | 3.9e-7 | 2.2e-8 |
1219    ///
1220    /// # Arguments
1221    /// * `x`: A single precision number
1222    /// * `Si`: A mutable single precision number that will return the sine
1223    ///         integral value
1224    /// * `Ci`: A mutable single precision number that will return the cosine
1225    ///         integral value
1226    ///
1227    /// # Example
1228    /// ```
1229    /// use special_fun::cephes_single::sici;
1230    /// let x = 0.1f32;
1231    /// let mut si = 0.0_f32;
1232    /// let mut ci = 0.0_f32;
1233    /// sici(0.5_f32, &mut si, &mut ci);
1234    /// ```
1235    pub fn sici(x: f32, si: &mut f32, ci: &mut f32) {
1236        let code = unsafe { unsafe_cephes_single::sicif(x, si, ci) };
1237        assert_eq!(code, 0);
1238    }
1239
1240    /// Function to accurately compute hyperbolic sine and cosine integrals
1241    ///
1242    /// The hyperbolic sine integral is defined as:
1243    ///
1244    /// $ Shi(x) = \int_0^{\infty} \frac{ sinh(t) }{t} dt $
1245    ///
1246    /// and the hyperbolic cosine integral is defined as:
1247    ///
1248    /// $ Chi(x) = -\int_x^{\infty} \frac{ cosh(t) }{t} dt $
1249    ///
1250    /// which reduces to:
1251    ///
1252    /// $ Chi(x) = \gamma + ln(x) + \int_0^x \frac{cosh(t) - 1}{t} dt $
1253    ///
1254    /// where $\gamma$ is the euler constant.
1255    ///
1256    /// # Original Description from Stephen L. Moshier
1257    /// The integrals are approximated by rational functions.
1258    /// For x > 8 auxiliary functions f(x) and g(x) are employed
1259    /// such that
1260    ///
1261    /// Ci(x) = f(x) sin(x) - g(x) cos(x)
1262    /// Si(x) = pi/2 - f(x) cos(x) - g(x) sin(x)
1263    ///
1264    /// ACCURACY(Relative error)):
1265    ///
1266    /// | arithmetic | Function | domain | # trials | peak | rms |
1267    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
1268    /// | IEEE | Shi | 0,88 | 30000 | 3.5e-7 | 7.0e-8 |
1269    ///
1270    /// ACCURACY(Absolute error, except relative when |Chi| > 1):
1271    ///
1272    /// | arithmetic | Function | domain | # trials | peak | rms |
1273    /// | :---: | :--: | :--: | :--: | :--: | :--: | :--: |
1274    /// | IEEE | Chi | 0,88 | 30000 | 3.8e-7 | 7.6e-8 |
1275    ///
1276    /// # Arguments
1277    /// * `x`: A single precision number
1278    /// * `Shi`: A mutable single precision number that will return the
1279    ///          hyperbolic sine integral value
1280    /// * `Chi`: A mutable single precision number that will return the
1281    ///          hyperbolic cosine integral value
1282    ///
1283    /// # Example
1284    /// ```
1285    /// use special_fun::cephes_single::shichi;
1286    /// let x = 0.1f32;
1287    /// let mut shi = 0.0_f32;
1288    /// let mut chi = 0.0_f32;
1289    /// shichi(0.5_f32, &mut shi, &mut chi);
1290    /// ```
1291    pub fn shichi(x: f32, chi: &mut f32, shi: &mut f32){
1292        let code = unsafe { unsafe_cephes_single::shichif(x, chi, shi) };
1293        assert_eq!(code, 0);
1294    }
1295
1296    /// Function to compute the beta function.
1297    ///
1298    /// The beta function is given by:
1299    ///
1300    /// $B(a, b) = \int_0^1 t^{a - 1} ( 1 - t )^{b - 1} dt $
1301    ///
1302    /// which simplifies to
1303    ///
1304    /// $ B(a, b) = \frac{ \Gamma(a) \Gamma(b)}{ \Gamma(a + b) } $
1305    ///
1306    /// # Original Description from Stephen L. Moshier
1307    /// For large arguments the logarithm of the function is
1308    /// evaluated using lgam(), then exponentiated.
1309    ///
1310    /// ACCURACY (Relative error):
1311    ///
1312    /// | arithmetic | domain | # trials | peak | rms |
1313    /// | :--------: | :----: | :------: | :--: | :-: |
1314    /// | IEEE       | 0,30   | 10000    | 4.0e-5 | 6.0e-6 |
1315    /// | IEEE       | -20,0  | 10000    | 4.9e-3 | 5.4e-5 |
1316    ///
1317    /// ERROR MESSAGES:
1318    ///
1319    /// | message   | condition | value returned |
1320    /// | :-------: | :--------: | :----------: |
1321    /// | beta overflow | log(beta) > MAXLOG | 0.0 |
1322    /// | beta overflow | a or b <0 integer      | 0.0 |
1323    ///
1324    /// # Arguments
1325    /// * `a`: A single precision parameter, assumed to be greater than 0
1326    /// * `b`: A single precision parameter, assumed to be greater than 0
1327    ///
1328    /// # Example
1329    /// ```
1330    /// use special_fun::cephes_single::beta;
1331    /// let a = 0.1f32;
1332    /// let b = 0.2f32;
1333    /// beta(a, b);
1334    /// ```
1335    pub fn beta(a: f32, b: f32) -> f32 {
1336        unsafe { unsafe_cephes_single::betaf(a, b) }
1337    }
1338
1339    /// Function to compute the regularized incomplete beta function.
1340    ///
1341    /// The regularized incomplete beta function is given by:
1342    ///
1343    /// $I_x(a, b) = \frac{1}{B(a, b)} \int_0^x t^{a - 1} ( 1 - t )^{b - 1} dt $
1344    ///
1345    /// and for $x = 1$, the regularized incomplete beta function evaluates to
1346    /// exactly 1. Note, $ 1 - I_x(a, b) = I_{1 - x}(a, b)$.
1347    ///
1348    /// # Original Description from Stephen L. Moshier
1349    /// The integral is evaluated by a continued fraction expansion
1350    /// or, when b*x is small, by a power series.
1351    ///
1352    /// ACCURACY (Relative error):
1353    ///
1354    /// | arithmetic | domain | # trials | peak | rms |
1355    /// | :--------: | :----: | :------: | :--: | :-: |
1356    /// | IEEE       | 0,30   | 10000    | 3.7e-5 | 5.1e-6 |
1357    /// | IEEE       | 0,100  | 10000    | 1.7e-4 | 2.5e-5 |
1358    ///
1359    /// The useful domain for relative error is limited by underflow of the
1360    /// single precision exponential function.
1361    ///
1362    /// ACCURACY (Absolute error):
1363    ///
1364    /// | arithmetic | domain | # trials | peak | rms |
1365    /// | :--------: | :----: | :------: | :--: | :-: |
1366    /// | IEEE       | 0,30   | 100000   | 2.2e-5 | 5.1e-6 |
1367    /// | IEEE       | 0,100  | 10000    | 6.5e-5 | 3.7e-6 |
1368    ///
1369    /// Larger errors may occur for extreme ratios of `a` and `b`.
1370    ///
1371    /// ERROR MESSAGES:
1372    ///
1373    /// | message   | condition | value returned |
1374    /// | :-------: | :--------: | :----------: |
1375    /// | incbet domain | x<0, x>1 | 0.0 |
1376    ///
1377    /// # Arguments
1378    /// * `a`: A single precision parameter, assumed to be greater than 0
1379    /// * `b`: A single precision parameter, assumed to be greater than 0
1380    /// * `x`: A single precision parameter, assumed to be between 0 and 1
1381    ///
1382    /// # Example
1383    /// ```
1384    /// use special_fun::cephes_single::incbet;
1385    /// let a = 0.1f32;
1386    /// let b = 0.2f32;
1387    /// incbet(a, b, 0.5f32);
1388    /// ```
1389    pub fn incbet(a: f32, b: f32, x: f32) -> f32 {
1390        unsafe { unsafe_cephes_single::incbetf(a, b, x) }
1391    }
1392
1393    /// Inverse of incomplete beta integral.
1394    pub fn incbi(a: f32, b: f32, y: f32) -> f32 {
1395        unsafe { unsafe_cephes_single::incbif(a, b, y) }
1396    }
1397
1398    /// Gamma function.
1399    pub fn gamma(x: f32) -> f32 {
1400        unsafe { unsafe_cephes_single::gammaf(x) }
1401    }
1402
1403    /// Reciprocal gamma function.
1404    pub fn rgamma(x: f32) -> f32 {
1405        unsafe { unsafe_cephes_single::rgammaf(x) }
1406    }
1407
1408    /// Natural logarithm of gamma function.
1409    pub fn lgam(x: f32) -> f32 {
1410        unsafe { unsafe_cephes_single::lgamf(x) }
1411    }
1412
1413    /// Regularized incomplete gamma integral.
1414    pub fn igam(a: f32, x: f32) -> f32 {
1415        unsafe { unsafe_cephes_single::igamf(a, x) }
1416    }
1417
1418    /// Complemented incomplete gamma integral.
1419    pub fn igamc(a: f32, x: f32) -> f32 {
1420        unsafe { unsafe_cephes_single::igamcf(a, x) }
1421    }
1422
1423    /// Inverse of complemented incomplete gamma integral.
1424    pub fn igami(a: f32, p: f32) -> f32 {
1425        unsafe { unsafe_cephes_single::igamif(a, p) }
1426    }
1427
1428    /// Psi (digamma) function.
1429    pub fn psi(x: f32) -> f32 {
1430        unsafe { unsafe_cephes_single::psif(x) }
1431    }
1432
1433    /// Factorial function.
1434    pub fn fac(i: i32) -> f32 {
1435        unsafe { unsafe_cephes_single::facf(i) }
1436    }
1437}
1438
1439/// Special functions on primitive floating point numbers.
1440///
1441/// This provides safe access to a subset of the Cephes functions implemented
1442/// for double and single precision.
1443pub trait FloatSpecial: Copy + Add<Output=Self> + Sub<Output=Self> {
1444    /// Beta function.
1445    fn beta(self, b: Self) -> Self;
1446    /// Logarithm of beta function.
1447    fn logbeta(self, b: Self) -> Self {
1448        self.loggamma() + b.loggamma() - (self + b).loggamma()
1449    }
1450    /// Regularized incomplete beta function.
1451    fn betainc(self, a: Self, b: Self) -> Self;
1452    /// Inverse of incomplete beta integral.
1453    fn betainc_inv(self, a: Self, b: Self) -> Self;
1454
1455    /// Factorial.
1456    fn factorial(self) -> Self;
1457    /// Gamma function.
1458    fn gamma(self) -> Self;
1459    /// Reciprocal gamma function.
1460    fn rgamma(self) -> Self;
1461    /// Logarithm of gamma function.
1462    fn loggamma(self) -> Self;
1463
1464    /// Regularized incomplete gamma integral.
1465    fn gammainc(self, a: Self) -> Self;
1466    /// Complemented incomplete gamma integral.
1467    fn gammac(self, a: Self) -> Self;
1468    /// Inverse of complemented incomplete gamma integral.
1469    fn gammac_inv(self, a: Self) -> Self;
1470
1471    /// Digamma function.
1472    fn digamma(self) -> Self;
1473
1474    /// Error function.
1475    fn erf(self) -> Self;
1476    /// Complementary error function.
1477    fn erfc(self) -> Self;
1478
1479    /// Confluent hypergeometric function 1F1.
1480    fn hyp1f1(self, a: Self, b: Self) -> Self;
1481    /// Hypergeometric function 1F2.
1482    fn hyp1f2(self, a: Self, b: Self, c: Self) -> Self;
1483    /// Gauss hypergeometric function 2F1.
1484    fn hyp2f1(self, a: Self, b: Self, c: Self) -> Self;
1485    /// Hypergeometric function 3F0.
1486    fn hyp3f0(self, a: Self, b: Self, c: Self) -> Self;
1487
1488    /// Normal distribution function.
1489    fn norm(self) -> Self;
1490    /// Inverse of Normal distribution function.
1491    fn norm_inv(self) -> Self;
1492
1493    /// Bessel function of real order of the first kind.
1494    fn besselj(self, v: Self) -> Self;
1495    /// Bessel function of real order of the second kind.
1496    fn bessely(self, v: Self) -> Self;
1497    /// Modified bessel function of real order of the first kind.
1498    fn besseli(self, v: Self) -> Self;
1499    /// Modified bessel function of integer order of the second kind.
1500    fn besselk(self, v: i32) -> Self;
1501
1502    /// Riemann zeta function.
1503    fn riemann_zeta(self) -> Self;
1504    /// Hurwitz zeta function.
1505    fn hurwitz_zeta(self, q: Self) -> Self;
1506}
1507
1508impl FloatSpecial for f64 {
1509    fn beta(self, b: f64) -> f64 {
1510        unsafe { unsafe_cephes_double::beta(self, b) }
1511    }
1512    fn betainc(self, a: f64, b: f64) -> f64 {
1513        unsafe { unsafe_cephes_double::incbet(a, b, self) }
1514    }
1515    fn betainc_inv(self, a: f64, b: f64) -> f64 {
1516        unsafe { unsafe_cephes_double::incbi(a, b, self) }
1517    }
1518
1519    fn factorial(self) -> f64 {
1520        unsafe { unsafe_cephes_double::gamma(self + 1.0) }
1521    }
1522    fn gamma(self) -> f64 {
1523        unsafe { unsafe_cephes_double::gamma(self) }
1524    }
1525    fn rgamma(self) -> f64 {
1526        unsafe { unsafe_cephes_double::rgamma(self) }
1527    }
1528    fn loggamma(self) -> f64 {
1529        unsafe { unsafe_cephes_double::lgam(self) }
1530    }
1531
1532    fn gammainc(self, a: f64) -> f64 {
1533        unsafe { unsafe_cephes_double::igam(a, self) }
1534    }
1535    fn gammac(self, a: f64) -> f64 {
1536        unsafe { unsafe_cephes_double::igamc(a, self) }
1537    }
1538    fn gammac_inv(self, a: f64) -> f64 {
1539        unsafe { unsafe_cephes_double::igami(a, self) }
1540    }
1541
1542    fn digamma(self) -> f64 {
1543        unsafe { unsafe_cephes_double::psi(self) }
1544    }
1545
1546    fn erf(self) -> f64 {
1547        unsafe { unsafe_cephes_double::erf(self) }
1548    }
1549    fn erfc(self) -> f64 {
1550        unsafe { unsafe_cephes_double::erfc(self) }
1551    }
1552
1553    fn hyp1f1(self, a: f64, b: f64) -> f64 {
1554        unsafe { unsafe_cephes_double::hyperg(a, b, self) }
1555    }
1556    fn hyp1f2(self, a: f64, b: f64, c: f64) -> f64 {
1557        let mut err = 0.0;
1558        unsafe { unsafe_cephes_double::onef2(a, b, c, self, &mut err) }
1559    }
1560    fn hyp2f1(self, a: f64, b: f64, c: f64) -> f64 {
1561        unsafe { unsafe_cephes_double::hyp2f1(a, b, c, self) }
1562    }
1563    fn hyp3f0(self, a: f64, b: f64, c: f64) -> f64 {
1564        let mut err = 0.0;
1565        unsafe { unsafe_cephes_double::threef0(a, b, c, self, &mut err) }
1566    }
1567
1568    fn norm(self) -> f64 {
1569        unsafe { unsafe_cephes_double::ndtr(self) }
1570    }
1571    fn norm_inv(self) -> f64 {
1572        unsafe { unsafe_cephes_double::ndtri(self) }
1573    }
1574
1575    fn besselj(self, v: f64) -> f64 {
1576        unsafe { unsafe_cephes_double::jv(v, self) }
1577    }
1578    fn bessely(self, v: f64) -> f64 {
1579        unsafe { unsafe_cephes_double::yv(v, self) }
1580    }
1581    fn besseli(self, v: f64) -> f64 {
1582        unsafe { unsafe_cephes_double::iv(v, self) }
1583    }
1584    fn besselk(self, n: i32) -> f64 {
1585        unsafe { unsafe_cephes_double::kn(n, self) }
1586    }
1587
1588    fn riemann_zeta(self) -> f64 {
1589        unsafe { 1.0 + unsafe_cephes_double::zetac(self) }
1590    }
1591    fn hurwitz_zeta(self, q: f64) -> f64 {
1592        unsafe { unsafe_cephes_double::zeta(self, q) }
1593    }
1594}
1595
1596impl FloatSpecial for f32 {
1597    fn beta(self, b: f32) -> f32 {
1598        unsafe { unsafe_cephes_single::betaf(self, b) }
1599    }
1600    fn betainc(self, a: f32, b: f32) -> f32 {
1601        unsafe { unsafe_cephes_single::incbetf(a, b, self) }
1602    }
1603    fn betainc_inv(self, a: f32, b: f32) -> f32 {
1604        unsafe { unsafe_cephes_single::incbif(a, b, self) }
1605    }
1606
1607    fn factorial(self) -> f32 {
1608        unsafe { unsafe_cephes_single::gammaf(self + 1.0) }
1609    }
1610    fn gamma(self) -> f32 {
1611        unsafe { unsafe_cephes_single::gammaf(self) }
1612    }
1613    fn rgamma(self) -> f32 {
1614        unsafe { unsafe_cephes_single::rgammaf(self) }
1615    }
1616    fn loggamma(self) -> f32 {
1617        unsafe { unsafe_cephes_single::lgamf(self) }
1618    }
1619
1620    fn gammainc(self, a: f32) -> f32 {
1621        unsafe { unsafe_cephes_single::igamf(a, self) }
1622    }
1623    fn gammac(self, a: f32) -> f32 {
1624        unsafe { unsafe_cephes_single::igamcf(a, self) }
1625    }
1626    fn gammac_inv(self, a: f32) -> f32 {
1627        unsafe { unsafe_cephes_single::igamif(a, self) }
1628    }
1629
1630    fn digamma(self) -> f32 {
1631        unsafe { unsafe_cephes_single::psif(self) }
1632    }
1633
1634    fn erf(self) -> f32 {
1635        unsafe { unsafe_cephes_single::erff(self) }
1636    }
1637    fn erfc(self) -> f32 {
1638        unsafe { unsafe_cephes_single::erfcf(self) }
1639    }
1640
1641    fn hyp1f1(self, a: f32, b: f32) -> f32 {
1642        unsafe { unsafe_cephes_single::hypergf(a, b, self) }
1643    }
1644    fn hyp1f2(self, a: f32, b: f32, c: f32) -> f32 {
1645        let mut err = 0.0;
1646        unsafe { unsafe_cephes_single::onef2f(a, b, c, self, &mut err) }
1647    }
1648    fn hyp2f1(self, a: f32, b: f32, c: f32) -> f32 {
1649        unsafe { unsafe_cephes_single::hyp2f1f(a, b, c, self) }
1650    }
1651    fn hyp3f0(self, a: f32, b: f32, c: f32) -> f32 {
1652        let mut err = 0.0;
1653        unsafe { unsafe_cephes_single::threef0f(a, b, c, self, &mut err) }
1654    }
1655
1656    fn norm(self) -> f32 {
1657        unsafe { unsafe_cephes_single::ndtrf(self) }
1658    }
1659    fn norm_inv(self) -> f32 {
1660        unsafe { unsafe_cephes_single::ndtrif(self) }
1661    }
1662
1663    fn besselj(self, v: f32) -> f32 {
1664        unsafe { unsafe_cephes_single::jvf(v, self) }
1665    }
1666    fn bessely(self, v: f32) -> f32 {
1667        unsafe { unsafe_cephes_single::yvf(v, self) }
1668    }
1669    fn besseli(self, v: f32) -> f32 {
1670        unsafe { unsafe_cephes_single::ivf(v, self) }
1671    }
1672    fn besselk(self, n: i32) -> f32 {
1673        unsafe { unsafe_cephes_single::knf(n, self) }
1674    }
1675
1676    fn riemann_zeta(self) -> f32 {
1677        unsafe { 1.0 + unsafe_cephes_single::zetacf(self) }
1678    }
1679    fn hurwitz_zeta(self, q: f32) -> f32 {
1680        unsafe { unsafe_cephes_single::zetaf(self, q) }
1681    }
1682}
1683