36#ifndef IS_RTKPP_LIBRARY
39#include "../include/STK_Funct_gammaRatio.h"
40#include "../include/STK_Funct_raw.h"
65 if (b==0.)
return( 0. );
146static Real gammaRatio_cf(
Real const& a,
Real const& x,
bool lower_tail)
151 Real Nold = 0, Ncur = 1.;
153 Real Dold = 1, Dcur = bn;
156 if (std::abs(Dcur) > 4294967296.0)
171 Real an = n*(--aminusn);
178 Real anew = bn * Ncur + an * Nold;
181 anew = bn * Dcur + an * Dold;
185 if (std::abs(Dcur) > 4294967296.0)
193 if (std::abs(Ncur) > 4294967296.0)
195 Ncur /= 4294967296.0;
196 Nold /= 4294967296.0;
197 Dold /= 4294967296.0;
198 Dcur /= 4294967296.0;
201 if (std::abs(Dcur) != 0)
206 if (std::abs(cf - cfold) < std::abs(cf)*Arithmetic<Real>::epsilon())
211 return (lower_tail) ? (1 - apois(a, x)*cf)
229static Real gammaRatio_sr(
Real const& a,
Real const& x,
bool lower_tail)
231 Real y = a+1., term = x / (a+1.),
sum = term;
236 sum += (term *= x / y);
238 while (term >
sum * Arithmetic<Real>::epsilon());
267static Real gammaRatio_ae(
Real const& a,
Real const& x,
bool lower_tail)
271 while ((b > 1.) && (term >
sum * Arithmetic<Real>::epsilon()))
273 sum += (term *= b / x);
303static inline Real poisson_ae(
Real const&
a1,
Real const& apd,
bool lower_tail =
true)
306 static const Real coefs_i[8] =
318 static const Real coefs_s[8] =
326 5246819/75246796800.,
327 -534703531/902961561600.
332 Real sqrt2D = sqrt (2. * D);
337 if (apd -
a1 < 0) sqrt2D = -sqrt2D;
341 Real sum1, num1_term, sum2, num2_term;
342 num1_term = sum1 = sqrt(
a1);
343 num2_term = sum2 = sqrt2D;
349 for (
int i = 1; i < 8; i++)
352 num += num1_term * coefs_i[i];
353 num += num2_term * coefs_s[i];
356 num1_term = num1_term/
a1 + sum1;
358 sum2 *= (D_x / (i + 0.5));
359 num2_term = num2_term/
a1 + sum2;
361 den += den_term * coefs_s[i];
395 if (x <= 1.)
return gammaRatio_dl(a, x,
lower_tail);
397 if (0.75*a> x)
return gammaRatio_sr(a, x,
lower_tail);
401 if (a>1)
return gammaRatio_ae(a, x,
lower_tail);
408 if (a>x)
return gammaRatio_sr(a, x,
lower_tail);
411 if (a>1)
return gammaRatio_ae(a, x,
lower_tail);
#define STKDOMAIN_ERROR_2ARG(Where, Arg1, Arg2, Error)
Real firstImpl() const
compute the first term of the serie
Real nextImpl() const
compute the next terms
Seriedl(Real const &a, Real const &x)
constructor
Interface base class for Series.
The MultidimRegression class allows to regress a multidimensional output variable among a multivariat...
Real dev0(Real const &a, Real const &b)
compute the partial deviance .
Real expm1(Real const &x)
compute the function .
Real gammaRatio(Real const &a, Real const &x, bool lower_tail)
Compute the incomplete gamma functions ratio.
Real normal_pdf_raw(Real const &x)
compute the probability distribution function of the normal density.
Real sumAlternateSerie(const ISerie< Serie > &f, int const &n=0)
Sum an alternating series using the Chebichev polynomials.
Real poisson_pdf_raw(Real const &x, Real const &lambda)
Compute the Poisson density.
Real lgammaStirlingError(Real const &z)
Compute the error when we compute using the Stirling's formula.
Real normal_cdf_raw(Real const &x)
Compute the cumulative distribution function of the normal density.
Real lgamma(Real const &)
This function computes the function .
Real gammaRatioP(Real const &a, Real const &x)
Compute the incomplete gamma function ratio P(a,x).
bool isNA(Type const &x)
utility method allowing to know if a value is a NA (Not Available) value
Arrays::SumOp< Lhs, Rhs >::result_type sum(Lhs const &lhs, Rhs const &rhs)
convenience function for summing two arrays
double Real
STK fundamental type of Real values.
Real gammaRatio_raw(Real const &a, Real const &x, bool lower_tail)
The namespace STK is the main domain space of the Statistical ToolKit project.
Arithmetic properties of STK fundamental types.