37#ifndef ROOT_Math_PdfFuncMathCore
38#define ROOT_Math_PdfFuncMathCore
85 if (x < 0 || x > 1.0)
return 0;
88 if (
a < 1)
return std::numeric_limits<double>::infinity();
89 else if (
a > 1)
return 0;
90 else if (
a == 1)
return b;
94 if (
b < 1)
return std::numeric_limits<double>::infinity();
95 else if (
b > 1)
return 0;
96 else if (
b == 1)
return a;
151 if (
n < 0)
return 0.0;
152 if (p < 0 || p > 1.0)
return 0.0;
177 double gammahalf =
gamma/2.0;
178 return gammahalf/(
M_PI * ((
x-x0)*(
x-x0) + gammahalf*gammahalf));
233 if (
x == x0 &&
a == 0)
return 0.5;
258 if (
sigma < 0.)
return 0.;
259 double z = (
x - mean)/
sigma;
260 if (alpha < 0) z = -z;
261 double abs_alpha = std::abs(alpha);
268 double nDivAlpha =
n/abs_alpha;
269 double AA =
std::exp(-0.5*abs_alpha*abs_alpha);
270 double B = nDivAlpha -abs_alpha;
271 double arg = nDivAlpha/(
B-z);
282 if (
sigma < 0.)
return 0.;
283 if (
n <= 1)
return std::numeric_limits<double>::quiet_NaN();
284 double abs_alpha = std::abs(alpha);
285 double C =
n/abs_alpha * 1./(
n-1.) *
std::exp(-alpha*alpha/2.);
311 return lambda *
std::exp (-lambda * (
x-x0));
337 return std::numeric_limits<double>::quiet_NaN();
363 inline double gamma_pdf(
double x,
double alpha,
double theta,
double x0 = 0) {
368 }
else if ((
x-x0) == 0) {
376 }
else if (alpha == 1) {
404 double tmp = (
x-x0)/
sigma;
425 inline double bigaussian_pdf(
double x,
double y,
double sigmax = 1,
double sigmay = 1,
double rho = 0,
double x0 = 0,
double y0 = 0) {
426 double u = (
x-x0)/sigmax;
427 double v = (
y-y0)/sigmay;
428 double c = 1. - rho*rho;
429 double z = u*u - 2.*rho*u*
v +
v*
v;
455 double landau_pdf(
double x,
double xi = 1,
double x0 = 0);
504 double tmp = (
x-x0)/
sigma;
585 if ((
x-x0) <
b && (
x-x0) >=
a)
double pow(double, double)
double uniform_pdf(double x, double a, double b, double x0=0)
Probability density function of the uniform (flat) distribution.
double bigaussian_pdf(double x, double y, double sigmax=1, double sigmay=1, double rho=0, double x0=0, double y0=0)
Probability density function of the bi-dimensional (Gaussian) distribution.
double normal_pdf(double x, double sigma=1, double x0=0)
Probability density function of the normal (Gaussian) distribution.
double lognormal_pdf(double x, double m, double s, double x0=0)
Probability density function of the lognormal distribution.
double binomial_pdf(unsigned int k, double p, unsigned int n)
Probability density function of the binomial distribution.
double fdistribution_pdf(double x, double n, double m, double x0=0)
Probability density function of the F-distribution.
double negative_binomial_pdf(unsigned int k, double p, double n)
Probability density function of the negative binomial distribution.
double gamma_pdf(double x, double alpha, double theta, double x0=0)
Probability density function of the gamma distribution.
double landau_pdf(double x, double xi=1, double x0=0)
Probability density function of the Landau distribution:
double exponential_pdf(double x, double lambda, double x0=0)
Probability density function of the exponential distribution.
double gaussian_pdf(double x, double sigma=1, double x0=0)
Probability density function of the normal (Gaussian) distribution.
double chisquared_pdf(double x, double r, double x0=0)
Probability density function of the distribution with degrees of freedom.
double breitwigner_pdf(double x, double gamma, double x0=0)
Probability density function of Breit-Wigner distribution, which is similar, just a different definit...
double crystalball_function(double x, double alpha, double n, double sigma, double mean=0)
Crystal ball function.
double crystalball_pdf(double x, double alpha, double n, double sigma, double mean=0)
pdf definition of the crystal_ball which is defined only for n > 1 otherwise integral is diverging
double beta_pdf(double x, double a, double b)
Probability density function of the beta distribution.
double poisson_pdf(unsigned int n, double mu)
Probability density function of the Poisson distribution.
double cauchy_pdf(double x, double b=1, double x0=0)
Probability density function of the Cauchy distribution which is also called Lorentzian distribution.
double tdistribution_pdf(double x, double r, double x0=0)
Probability density function of Student's t-distribution.
double lgamma(double x)
Calculates the logarithm of the gamma function.
double erf(double x)
Error function encountered in integrating the normal distribution.
Namespace for new Math classes and functions.
double log1p(double x)
declarations for functions which are not implemented by some compilers
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
static constexpr double s