STK++ 0.9.13
STK::Funct Namespace Reference

The namespace Funct enclose all usual and special functions. More...

Classes

class  Seriedl
 This Serie computes. More...
 
class  Serielog1p
 This Series computes. More...
 

Functions

Real betaRatio_sr (Real const &a, Real const &b, Real x, bool xm1, bool lower_tail)
 Compute the incomplete beta function ratio I_x(a,b) using its series representation.
 
Real betaRatio_up (Real const &a, Real const &b, Real const &x, bool xm1, bool lower_tail)
 Compute the incomplete beta function ratio I_x(a,b) using its recurrence formula and its asymptotic expansion.
 
Real betaRatio_cf (Real const &a, Real const &b, Real x, bool xm1, bool lower_tail=true)
 Compute the incomplete beta function ratio using the continued fraction method.
 
Real betaRatio_se (Real const &a, Real const &b, Real const &x, bool xm1, bool lower_tail)
 Compute the incomplete beta function ratio I_x(a,b) using the serie expansion method.
 
Real betaRatio_ae (Real const &a, Real const &b, Real const &x, bool xm1, bool lower_tail)
 Compute the incomplete beta function ratio I_x(a,b) using the asymptotic expansion method.
 
Real betaRatio (Real const &a, Real const &b, Real const &x, bool lower_tail=true)
 Compute the incomplete beta function ratio Compute the beta ratio function.
 
Real betaRatio_raw (Real const &a, Real const &b, Real const &x, bool lower_tail)
 
Real factorial (int n)
 This function computes $ n! $ for integer argument.
 
Real lfactorial (int n)
 This function computes $ \ln(n!) $ for integer argument.
 
Real factorial (Real const &z)
 This function computes $ z! $ when z is an integer in a Real format.
 
Real lfactorial (Real const &z)
 This function computes $ \ln(z!) $ when z is an integer in Real format.
 
Real gamma (Real const &z)
 This function computes the function $ \Gamma(z) $.
 
Real lgamma (Real const &z)
 This function computes the function $ \ln(\Gamma(z)) $.
 
Real factorial_raw (int)
 
Real lfactorial_raw (int)
 
Real factorial_raw (Real const &)
 
Real lfactorial_raw (Real const &)
 
Real gamma_raw (Real const &)
 
Real lgamma_raw (Real const &)
 
Real lanczosSerie (Real const &z)
 Compute the Lanzcos correction series for the gamma function with n = 21 terms.
 
Real gammaLanczos (Real const &z)
 Compute the gamma function using the Lanzcos expansion using n = 21 terms and r= 22.618910.
 
double stirlingSerie (Real const &z)
 Compute the Stirling's series for the lgamma function.
 
Real gammaStirling (Real const &z)
 This function computes the gamma function using the Stirling approximation.
 
Real lgammaStirling (Real const &z)
 This function computes the log gamma function using the Stirling approximation.
 
Real lgammaStirlingError (Real const &z)
 Compute the error when we compute $ \ln(\Gamma(z)) $ using the Stirling's formula.
 
Real lgammaStirlingError (int n)
 Compute the error when we compute $ \ln(\Gamma(z)) $ using the Stirling's formula and z is an integer.
 
Real gammaRatio_raw (Real const &a, Real const &x, bool lower_tail)
 
Real gammaRatioQ_raw (Real const &a, Real const &x)
 
Real gammaRatioP_raw (Real const &a, Real const &x)
 
Real gammaRatio (Real const &a, Real const &x, bool lower_tail)
 Compute the incomplete gamma functions ratio.
 
Real gammaRatioQ (Real const &a, Real const &x)
 Compute the incomplete gamma function ratio Q(a,x).
 
Real gammaRatioP (Real const &a, Real const &x)
 Compute the incomplete gamma function ratio P(a,x).
 
Real beta_pdf_raw (Real const &x, Real const &a, Real const &b)
 Compute the beta density function.
 
Real binomial_pdf_raw (Real const &x, Real const &n, Real const &p)
 Compute the generalized binomial probability mass function.
 
Real binomial_pdf_raw (int x, int n, Real const &p)
 Compute the binomial probability mass function.
 
Real binomial_lpdf_raw (Real const &x, Real const &n, Real const &p)
 Compute the generalized binomial log-probability mass function.
 
Real binomial_lpdf_raw (int x, int n, Real const &p)
 Compute the binomial log-probability mass function.
 
Real poisson_pdf_raw (Real const &x, Real const &lambda)
 Compute the Poisson density.
 
Real poisson_pdf_raw (int x, Real const &lambda)
 Compute the poisson density with integer value.
 
Real poisson_lpdf_raw (Real const &x, Real const &lambda)
 Compute the log-poisson density.
 
Real poisson_lpdf_raw (int x, Real const &lambda)
 Compute the log-poisson density with integer value.
 
Real erf_raw (Real const &a)
 Compute the error function erf(a) Compute the function.
 
Real erfc_raw (Real const &a)
 Compute the complementary error function erfc(a) Compute the function.
 
Real normal_cdf_raw (Real const &x)
 Compute the cumulative distribution function of the normal density.
 
Real normal_pdf_raw (Real const &x)
 compute the probability distribution function of the normal density.
 
Real psi_raw (Real x)
 Compute the psi function.
 
Real dev0 (Real const &a, Real const &b)
 compute the partial deviance $ d_0(a,b) = a\log(a/b)+ b - a $.
 
Real b1 (Real const &a, Real const &b, Real const &x, bool xm1)
 Compute the function.
 
Real g0 (Real const &x)
 compute the partial deviance $g_0(x) = x\log(x)+ 1 - x$.
 
Real log1p (Real const &x)
 compute the function $ \log(1+x) $.
 
Real expm1 (Real const &x)
 compute the function $ \exp(x)-1 $.
 
template<int N>
Real evalPolynomial (Real x, const Real *P)
 Polynomial evaluator.
 
template<>
Real evalPolynomial< 0 > (Real x, const Real *P)
 
template<int N>
Real evalPolynomial1 (Real x, const Real *P)
 Polynomial evaluator.
 
template<>
Real evalPolynomial1< 0 > (Real x, const Real *P)
 

Detailed Description

The namespace Funct enclose all usual and special functions.

The namespace Funct is the domain space of the special function like gamma function, beta function, incomplete gamma function, incomplete beta function... It include also some useful raw functions like log1p...

Function Documentation

◆ betaRatio_raw()

Real STK::Funct::betaRatio_raw ( Real const a,
Real const b,
Real const x,
bool  lower_tail 
)

Definition at line 556 of file STK_Funct_betaRatio.cpp.

557{
558 // trivial case
559 if (x<=0) return lower_tail ? 0. : 1.;
560 if (x>=1) return lower_tail ? 1. : 0.;
561 // parameters
562 Real p=a/(a+b);
563 // small a,b max(a,b) <= 15
564 if (std::max(a,b)<=15)
565 {
566 if ((a<=1 && x <=0.5)||(b<=1 && x>=0.5))
567 {
568 return (x <=0.5) ? betaRatio_sr(a, b, x, false, lower_tail)
569 : betaRatio_sr(b, a, x, true, !lower_tail);
570 }
571 // (a<=1) and (x>0.5), use serie_up for b
572 if (a<=1) return betaRatio_up(b, a, x, true, !lower_tail);
573 // (b<=1) and (x<0.5), use serie_up for a
574 if (b<=1) return betaRatio_up(a, b, x, false, lower_tail);
575 // min(a,b)>1, use serie_up for a or b
576 return (x <=0.5) ? betaRatio_up(a, b, x, false, lower_tail)
577 : betaRatio_up(b, a, x, true, !lower_tail);
578 }
579 // min(a,b) < 1 and max(a,b) > 15
580 if (std::min(a,b)<=1)
581 {
582// std::cout << "case 1" << std::endl;
583 // case b<=1 and a>=15
584 if (a>=15) return betaRatio_ae(a, b, x, false, lower_tail);
585 // case a<=1 and b>=15
586 if (b>=15) return betaRatio_ae(b, a, x, true, !lower_tail);
587 }
588
589 // ((1<a<=60)) or ((1<b<=60))
590 if ((std::min(a,b)<=60))
591 {
592// std::cout << "case 3.2" << std::endl;
593 return (a<b) ? betaRatio_up(a, b, x, false, lower_tail)
594 : betaRatio_up(b, a, x, true, !lower_tail);
595 }
596 // general case
597 if ((x<0.97*p)||(x>1.03*p))
598 {
599// std::cout << "case 2.1" << std::endl;
600 return (x < p) ? betaRatio_cf(a, b, x, false, lower_tail)
601 : betaRatio_cf(b, a, x, true, !lower_tail);
602 }
603// std::cout << "case 2.2" << std::endl;
604 return (a<b) ? betaRatio_se(a-1, b-1, x, false, lower_tail)
605 : betaRatio_se(b-1, a-1, x, true, !lower_tail);
606}
The MultidimRegression class allows to regress a multidimensional output variable among a multivariat...
Real betaRatio_se(Real const &a, Real const &b, Real const &x, bool xm1, bool lower_tail)
Compute the incomplete beta function ratio I_x(a,b) using the serie expansion method.
Real betaRatio_cf(Real const &a, Real const &b, Real x, bool xm1, bool lower_tail=true)
Compute the incomplete beta function ratio using the continued fraction method.
Real betaRatio_up(Real const &a, Real const &b, Real const &x, bool xm1, bool lower_tail)
Compute the incomplete beta function ratio I_x(a,b) using its recurrence formula and its asymptotic e...
Real betaRatio_sr(Real const &a, Real const &b, Real x, bool xm1, bool lower_tail)
Compute the incomplete beta function ratio I_x(a,b) using its series representation.
double Real
STK fundamental type of Real values.

References betaRatio_ae(), betaRatio_cf(), betaRatio_se(), betaRatio_sr(), and betaRatio_up().

Referenced by betaRatio().

◆ evalPolynomial1< 0 >()

template<>
Real STK::Funct::evalPolynomial1< 0 > ( Real  x,
const Real P 
)
inline

Definition at line 93 of file STK_Funct_Util.h.

94{ return 1.;}

◆ evalPolynomial< 0 >()

template<>
Real STK::Funct::evalPolynomial< 0 > ( Real  x,
const Real P 
)
inline

Definition at line 76 of file STK_Funct_Util.h.

77{ return P[0];}

◆ factorial_raw() [1/2]

Real STK::Funct::factorial_raw ( int  n)
inline

Definition at line 193 of file STK_Funct_gamma.h.

194{ return (n < 51) ? Const::factorialArray[n] : Funct::gamma_raw(n + 1.);}

References STK::Const::factorialArray, and gamma_raw().

Referenced by factorial(), and factorial().

◆ factorial_raw() [2/2]

Real STK::Funct::factorial_raw ( Real const z)
inline

Definition at line 212 of file STK_Funct_gamma.h.

213{
214 const int n = std::floor(z);
215 return (n < 51) ? Const::factorialArray[n] : Funct::gamma_raw(z + 1.);
216}

References STK::Const::factorialArray, and gamma_raw().

◆ gamma_raw()

Real STK::Funct::gamma_raw ( Real const z)
inline

Definition at line 282 of file STK_Funct_gamma.h.

283{
284 int n = (int)std::floor(z);
285 // If z is an integer (z integer and z<0 is not possible)
286 if (z == n)
287 { return (z < 51) ? Const::factorialArray[(int)z-1] : gammaStirling(z);}
288 // compute the absolute value of x -> y and compute the sign of the gamma function
289 Real y = std::abs(z);
290 int ny = std::floor(y);
291 Real signgam = (z<0 && isEven(ny)) ? -1 : 1, value;
292 // if y is an integer and a half, use reflection formula
293 if (y == ny+0.5)
294 { value = (ny<50) ? Const::factorialHalvesArray[ny] : gammaStirling(y);}
295 else // normal case
296 {
297 if (y <= 8)
298 {
299 Real r = y-ny; // r in [0,1]
300 value = gammaLanczos(r); // use Lanzcos approximation
301 // scale result
302 for (int i=0; i<ny; i++) value *= (r+i);
303 }
304 else // shift value
305 {
306 if (y <=16)
307 {
308 value = gammaStirling(y+6.);
309 for (int i=5; i>=0; --i) value /= (y+i);
310 }
311 else
312 value = gammaStirling(y);
313 }
314 }
315 // z >0 terminate
316 if (z>0) return value;
317 // z <0 -> use reflection formula and check overflow
318 Real den = y*sin(Const::_PI_ *y)*value;
319 return (den == 0.0) ? signgam * Arithmetic<Real>::infinity() : -Const::_PI_/den;
320}
Real gammaStirling(Real const &z)
This function computes the gamma function using the Stirling approximation.
bool isEven(int const &x)
is x an even number ?
Definition STK_Misc.h:83

References STK::Const::factorialArray, gammaLanczos(), gammaStirling(), and STK::isEven().

Referenced by factorial_raw(), factorial_raw(), gamma(), and lgamma_raw().

◆ gammaRatio_raw()

Real STK::Funct::gammaRatio_raw ( Real const a,
Real const x,
bool  lower_tail 
)

Definition at line 389 of file STK_Funct_gammaRatio.cpp.

390{
391 // trivial case
392 if (x<=0) return lower_tail ? 0 : 1.;
393 if (Arithmetic<Real>::isInfinite(x)) return lower_tail ? 1. : 0.;
394 // small values of x
395 if (x <= 1.) return gammaRatio_dl(a, x, lower_tail);
396 // large a compared to x
397 if (0.75*a> x) return gammaRatio_sr(a, x, lower_tail);
398 // large x compared to a
399 if (0.75*x > a)
400 { // if a>1 we can use the a.e. before computing
401 if (a>1) return gammaRatio_ae(a, x, lower_tail);
402 else return gammaRatio_cf(a, x, lower_tail);
403 }
404 // (a<100) and (x ~ a)
405 if (a<100.)
406 {
407 // a greater than x
408 if (a>x) return gammaRatio_sr(a, x, lower_tail);
409 else
410 { // if a>1 we can use the a.e. before computing
411 if (a>1) return gammaRatio_ae(a, x, lower_tail);
412 else return gammaRatio_cf(a, x, lower_tail);}
413 }
414 // (x ~ a) and (a>=100)
415 return poisson_ae(a-1, x, lower_tail);
416}
Arithmetic properties of STK fundamental types.

Referenced by gammaRatio(), gammaRatioP_raw(), and gammaRatioQ_raw().

◆ gammaRatioP_raw()

Real STK::Funct::gammaRatioP_raw ( Real const a,
Real const x 
)
inline

Definition at line 59 of file STK_Funct_gammaRatio.h.

60{ return gammaRatio_raw(a, x, true);}
Real gammaRatio_raw(Real const &a, Real const &x, bool lower_tail)

References gammaRatio_raw().

◆ gammaRatioQ_raw()

Real STK::Funct::gammaRatioQ_raw ( Real const a,
Real const x 
)
inline

Definition at line 56 of file STK_Funct_gammaRatio.h.

57{ return gammaRatio_raw(a, x, false);}

References gammaRatio_raw().

◆ lfactorial_raw() [1/2]

Real STK::Funct::lfactorial_raw ( int  n)
inline

Definition at line 234 of file STK_Funct_gamma.h.

235{ return (n < 51) ? Const::factorialLnArray[n] : lgamma( n + 1.);}
Real lgamma(Real const &)
This function computes the function .

References lgamma().

Referenced by lfactorial(), and lfactorial().

◆ lfactorial_raw() [2/2]

Real STK::Funct::lfactorial_raw ( Real const z)
inline

Definition at line 255 of file STK_Funct_gamma.h.

256{
257 const int n = (int)std::floor(z);
258 return (n < 51) ? Const::factorialLnArray[n] : lgamma(z + 1.);
259}

References lgamma().

◆ lgamma_raw()

Real STK::Funct::lgamma_raw ( Real const z)
inline

Definition at line 339 of file STK_Funct_gamma.h.

340{
341 // compute the absolute value of x -> y
342 Real y = std::abs(z), value;
343 int ny = std::floor(y);
344 // If x is an integer
345 if (y == ny)
346 { value = (ny<=50) ? Const::factorialLnArray[ny-1] : lgammaStirling(y);}
347 else
348 {
349 // If x is an integer and a half
350 if (y == ny+0.5)
351 { value = (ny<50) ? Const::factorialLnHalvesArray[ny] : lgammaStirling(y);}
352 else
353 {
354 // small values -> use gamma function
355 if (y <= 16) return log(std::abs(Funct::gamma_raw(z)));
356 else value = lgammaStirling(y);
357 }
358 }
359 // z >0 terminate
360 if (z>0) return value;
361 // z <0 -> use reflectiong formula
362 Real sinpiy = std::abs(sin(Const::_PI_ *y));
363 // overflow
364 return (sinpiy == 0.0) ? -Arithmetic<Real>::infinity()
365 : Const::_LNSQRTPI_2_+(z-0.5)*log(y)-z-log(sinpiy)+stirlingSerie(z);
366}
Real lgammaStirling(Real const &z)
This function computes the log gamma function using the Stirling approximation.
double stirlingSerie(Real const &z)
Compute the Stirling's series for the lgamma function.

References gamma_raw(), lgammaStirling(), and stirlingSerie().

Referenced by lgamma().