openGPMP
Open Source Mathematics Package
Static Public Member Functions | List of all members
gpmp::stats::PDF Class Reference

Class representing Probability Distribution Functions (PDFs) More...

#include <pdfs.hpp>

Static Public Member Functions

static double bernoulli (double x, double p)
 Calculates the probability of success in a Bernoulli trial. More...
 
static double beta (double x, double alpha, double beta)
 Calculates the probability density function (PDF) of the Beta distribution. More...
 
static double binomial (int k, int n, double p)
 Calculates the probability of observing k successes in n independent Bernoulli trials. More...
 
static double cauchy (double x, double x0, double gamma)
 Calculates the probability density function (PDF) of the Cauchy distribution. More...
 
static double chi_squared (double x, int k)
 Calculates the probability density function (PDF) of the chi-squared distribution. More...
 
static double exponential (double x, double lambda)
 Calculates the probability density function (PDF) of the exponential distribution. More...
 
static double f_dist (double x, int df1, int df2)
 Calculates the probability density function (PDF) of the F distribution. More...
 
static double gamma (double x, double alpha, double beta)
 Calculates the probability density function (PDF) of the gamma distribution. More...
 
static double inverse_gamma (double x, double alpha, double beta)
 Calculates the probability density function (PDF) of the inverse gamma distribution. More...
 
static double inverse_gaussian (double x, double mu, double lambda)
 Calculates the probability density function (PDF) of the inverse Gaussian distribution. More...
 
static double laplace (double x, double mu, double b)
 Calculates the probability density function (PDF) of the Laplace distribution. More...
 
static double logistic (double x, double mu, double s)
 Calculates the probability density function (PDF) of the logistic distribution. More...
 
static double log_normal (double x, double mu, double sigma)
 Calculates the probability density function (PDF) of the log-normal distribution. More...
 
static double gaussian (double x, double mu, double sigma)
 Calculates the probability density function (PDF) of the Gaussian (normal) distribution. More...
 
static double poisson (int k, double lambda)
 Calculates the probability density function (PDF) of the Poisson distribution. More...
 
static double rademacher (int k)
 Calculates the probability density function (PDF) of the Rademacher distribution. More...
 
static double student_t (double x, int df)
 Calculates the probability density function (PDF) of Student's t distribution. More...
 
static double uniform (double x, double a, double b)
 Calculates the probability density function (PDF) of the uniform distribution. More...
 
static double weibull (double x, double k, double lambda)
 Calculates the probability density function (PDF) of the Weibull distribution. More...
 
static double binomial_coefficient (int n, int k)
 Calculates the binomial coefficient "n choose k". More...
 
static double beta_function (double alpha, double beta)
 Calculates the beta function. More...
 
static double factorial (int n)
 Calculates the factorial of an integer. More...
 

Detailed Description

Class representing Probability Distribution Functions (PDFs)

Definition at line 41 of file pdfs.hpp.

Member Function Documentation

◆ bernoulli()

double gpmp::stats::PDF::bernoulli ( double  x,
double  p 
)
static

Calculates the probability of success in a Bernoulli trial.

Parameters
xThe outcome (0 or 1)
pThe probability of success
Returns
The probability of observing the outcome

Definition at line 44 of file pdfs.cpp.

44  {
45  if (std::abs(x - 0) < std::numeric_limits<double>::epsilon())
46  return 1 - p;
47  else if (std::abs(x - 1) < std::numeric_limits<double>::epsilon())
48  return p;
49  else
50  return 0;
51 }

Referenced by main().

◆ beta()

double gpmp::stats::PDF::beta ( double  x,
double  alpha,
double  beta 
)
static

Calculates the probability density function (PDF) of the Beta distribution.

Parameters
xThe value at which to evaluate the PDF
alphaThe shape parameter alpha (> 0)
betaThe shape parameter beta (> 0)
Returns
The value of the PDF at x

Definition at line 54 of file pdfs.cpp.

54  {
55  if (x < 0 || x > 1)
56  return 0;
57  else
58  return (std::pow(x, alpha - 1) * std::pow(1 - x, beta - 1)) /
59  beta_function(alpha, beta);
60 }
static double beta(double x, double alpha, double beta)
Calculates the probability density function (PDF) of the Beta distribution.
Definition: pdfs.cpp:54
static double beta_function(double alpha, double beta)
Calculates the beta function.
Definition: pdfs.cpp:212

Referenced by main().

◆ beta_function()

double gpmp::stats::PDF::beta_function ( double  alpha,
double  beta 
)
static

Calculates the beta function.

Parameters
alphaThe first shape parameter (> 0)
betaThe second shape parameter (> 0)
Returns
The value of the beta function

Definition at line 212 of file pdfs.cpp.

212  {
213  return std::tgamma(alpha) * std::tgamma(beta) / std::tgamma(alpha + beta);
214 }

◆ binomial()

double gpmp::stats::PDF::binomial ( int  k,
int  n,
double  p 
)
static

Calculates the probability of observing k successes in n independent Bernoulli trials.

Parameters
kThe number of successes
nThe number of trials
pThe probability of success in each trial
Returns
The probability of observing k successes

Definition at line 63 of file pdfs.cpp.

63  {
64  if (k < 0 || k > n)
65  return 0;
66  else
67  return binomial_coefficient(n, k) * std::pow(p, k) *
68  std::pow(1 - p, n - k);
69 }
static double binomial_coefficient(int n, int k)
Calculates the binomial coefficient "n choose k".
Definition: pdfs.cpp:203

Referenced by main().

◆ binomial_coefficient()

double gpmp::stats::PDF::binomial_coefficient ( int  n,
int  k 
)
static

Calculates the binomial coefficient "n choose k".

Parameters
nThe total number of items
kThe number of items to choose
Returns
The binomial coefficient

Definition at line 203 of file pdfs.cpp.

203  {
204  double result = 1.0;
205  for (int i = 1; i <= k; ++i) {
206  result *= static_cast<double>(n - (k - i)) / i;
207  }
208  return result;
209 }

◆ cauchy()

double gpmp::stats::PDF::cauchy ( double  x,
double  x0,
double  gamma 
)
static

Calculates the probability density function (PDF) of the Cauchy distribution.

Parameters
xThe value at which to evaluate the PDF
x0The location parameter
gammaThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 72 of file pdfs.cpp.

72  {
73  return (1 / (M_PI * gamma)) *
74  (gamma * gamma / ((x - x0) * (x - x0) + gamma * gamma));
75 }
static double gamma(double x, double alpha, double beta)
Calculates the probability density function (PDF) of the gamma distribution.
Definition: pdfs.cpp:106

Referenced by main().

◆ chi_squared()

double gpmp::stats::PDF::chi_squared ( double  x,
int  k 
)
static

Calculates the probability density function (PDF) of the chi-squared distribution.

Parameters
xThe value at which to evaluate the PDF
kThe degrees of freedom (> 0)
Returns
The value of the PDF at x

Definition at line 78 of file pdfs.cpp.

78  {
79  if (x < 0)
80  return 0;
81  else
82  return (std::pow(x, k / 2.0 - 1) * std::exp(-x / 2.0)) /
83  (std::pow(2, k / 2.0) * std::tgamma(k / 2.0));
84 }

Referenced by main().

◆ exponential()

double gpmp::stats::PDF::exponential ( double  x,
double  lambda 
)
static

Calculates the probability density function (PDF) of the exponential distribution.

Parameters
xThe value at which to evaluate the PDF
lambdaThe rate parameter (> 0)
Returns
The value of the PDF at x

Definition at line 87 of file pdfs.cpp.

87  {
88  if (x < 0)
89  return 0;
90  else
91  return lambda * std::exp(-lambda * x);
92 }

Referenced by main().

◆ f_dist()

double gpmp::stats::PDF::f_dist ( double  x,
int  df1,
int  df2 
)
static

Calculates the probability density function (PDF) of the F distribution.

Parameters
xThe value at which to evaluate the PDF
df1The degrees of freedom of the numerator (> 0)
df2The degrees of freedom of the denominator (> 0)
Returns
The value of the PDF at x

Definition at line 95 of file pdfs.cpp.

95  {
96  if (x < 0)
97  return 0;
98  else
99  return std::pow(df1, df1 / 2.0) * std::pow(df2, df2 / 2.0) *
100  std::pow(x, df1 / 2.0 - 1) /
101  (std::pow(df2 + df1 * x, (df1 + df2) / 2.0) *
102  std::tgamma(df1 / 2.0) * std::tgamma(df2 / 2.0));
103 }

Referenced by main().

◆ factorial()

static double gpmp::stats::PDF::factorial ( int  n)
static

Calculates the factorial of an integer.

Parameters
nThe integer
Returns
The factorial of n

◆ gamma()

double gpmp::stats::PDF::gamma ( double  x,
double  alpha,
double  beta 
)
static

Calculates the probability density function (PDF) of the gamma distribution.

Parameters
xThe value at which to evaluate the PDF
alphaThe shape parameter (> 0)
betaThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 106 of file pdfs.cpp.

106  {
107  if (x < 0)
108  return 0;
109  else
110  return (std::pow(x, alpha - 1) * std::exp(-x / beta)) /
111  (std::pow(beta, alpha) * std::tgamma(alpha));
112 }

Referenced by main().

◆ gaussian()

double gpmp::stats::PDF::gaussian ( double  x,
double  mu,
double  sigma 
)
static

Calculates the probability density function (PDF) of the Gaussian (normal) distribution.

Parameters
xThe value at which to evaluate the PDF
muThe mean of the distribution
sigmaThe standard deviation of the distribution (> 0)
Returns
The value of the PDF at x

Definition at line 154 of file pdfs.cpp.

154  {
155  return (1 / (sigma * std::sqrt(2 * M_PI))) *
156  std::exp(-0.5 * std::pow((x - mu) / sigma, 2));
157 }

Referenced by main().

◆ inverse_gamma()

double gpmp::stats::PDF::inverse_gamma ( double  x,
double  alpha,
double  beta 
)
static

Calculates the probability density function (PDF) of the inverse gamma distribution.

Parameters
xThe value at which to evaluate the PDF
alphaThe shape parameter (> 0)
betaThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 115 of file pdfs.cpp.

115  {
116  if (x <= 0)
117  return 0;
118  else
119  return (std::pow(beta, alpha) * std::pow(x, -alpha - 1) *
120  std::exp(-beta / x)) /
121  std::tgamma(alpha);
122 }

Referenced by main().

◆ inverse_gaussian()

double gpmp::stats::PDF::inverse_gaussian ( double  x,
double  mu,
double  lambda 
)
static

Calculates the probability density function (PDF) of the inverse Gaussian distribution.

Parameters
xThe value at which to evaluate the PDF
muThe mean parameter (> 0)
lambdaThe shape parameter (> 0)
Returns
The value of the PDF at x

Definition at line 125 of file pdfs.cpp.

125  {
126  if (x <= 0)
127  return 0;
128  else
129  return std::sqrt(lambda / (2 * M_PI * x * x * x)) *
130  std::exp(-lambda * (x - mu) * (x - mu) / (2 * mu * mu * x));
131 }

Referenced by main().

◆ laplace()

double gpmp::stats::PDF::laplace ( double  x,
double  mu,
double  b 
)
static

Calculates the probability density function (PDF) of the Laplace distribution.

Parameters
xThe value at which to evaluate the PDF
muThe location parameter
bThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 134 of file pdfs.cpp.

134  {
135  return 0.5 * std::exp(-std::abs(x - mu) / b) / b;
136 }

Referenced by main().

◆ log_normal()

double gpmp::stats::PDF::log_normal ( double  x,
double  mu,
double  sigma 
)
static

Calculates the probability density function (PDF) of the log-normal distribution.

Parameters
xThe value at which to evaluate the PDF
muThe mean of the underlying normal distribution
sigmaThe standard deviation of the underlying normal distribution (> 0)
Returns
The value of the PDF at x

Definition at line 145 of file pdfs.cpp.

145  {
146  if (x <= 0)
147  return 0;
148  else
149  return (1 / (x * sigma * std::sqrt(2 * M_PI))) *
150  std::exp(-0.5 * std::pow((std::log(x) - mu) / sigma, 2));
151 }

Referenced by main().

◆ logistic()

double gpmp::stats::PDF::logistic ( double  x,
double  mu,
double  s 
)
static

Calculates the probability density function (PDF) of the logistic distribution.

Parameters
xThe value at which to evaluate the PDF
muThe location parameter
sThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 139 of file pdfs.cpp.

139  {
140  double z = (x - mu) / s;
141  return std::exp(-z) / (s * std::pow(1 + std::exp(-z), 2));
142 }

Referenced by main().

◆ poisson()

double gpmp::stats::PDF::poisson ( int  k,
double  lambda 
)
static

Calculates the probability density function (PDF) of the Poisson distribution.

Parameters
kThe number of events
lambdaThe rate parameter (> 0)
Returns
The value of the PDF at k

Definition at line 160 of file pdfs.cpp.

160  {
161  if (k < 0)
162  return 0;
163  else
164  return std::exp(-lambda) * std::pow(lambda, k) /
166 }
static T factorial(T n)
Calculates the factorial of the given integer or double.
Definition: utils.hpp:211

References gpmp::core::Misc::factorial().

Referenced by main().

◆ rademacher()

double gpmp::stats::PDF::rademacher ( int  k)
static

Calculates the probability density function (PDF) of the Rademacher distribution.

Parameters
kThe value (1 or -1)
Returns
The value of the PDF at k

Definition at line 169 of file pdfs.cpp.

169  {
170  if (k == -1)
171  return 0.5;
172  else if (k == 1)
173  return 0.5;
174  else
175  return 0;
176 }

Referenced by main().

◆ student_t()

double gpmp::stats::PDF::student_t ( double  x,
int  df 
)
static

Calculates the probability density function (PDF) of Student's t distribution.

Parameters
xThe value at which to evaluate the PDF
dfThe degrees of freedom (> 0)
Returns
The value of the PDF at x

Definition at line 179 of file pdfs.cpp.

179  {
180  double numerator = std::tgamma((df + 1) / 2.0);
181  double denominator = std::sqrt(df * M_PI) * std::tgamma(df / 2.0);
182  return std::pow(1 + x * x / df, -(df + 1) / 2.0) * numerator / denominator;
183 }

Referenced by main().

◆ uniform()

double gpmp::stats::PDF::uniform ( double  x,
double  a,
double  b 
)
static

Calculates the probability density function (PDF) of the uniform distribution.

Parameters
xThe value at which to evaluate the PDF
aThe minimum value of the interval
bThe maximum value of the interval
Returns
The value of the PDF at x

Definition at line 186 of file pdfs.cpp.

186  {
187  if (x >= a && x <= b)
188  return 1 / (b - a);
189  else
190  return 0;
191 }

Referenced by main().

◆ weibull()

double gpmp::stats::PDF::weibull ( double  x,
double  k,
double  lambda 
)
static

Calculates the probability density function (PDF) of the Weibull distribution.

Parameters
xThe value at which to evaluate the PDF
kThe shape parameter (> 0)
lambdaThe scale parameter (> 0)
Returns
The value of the PDF at x

Definition at line 194 of file pdfs.cpp.

194  {
195  if (x < 0)
196  return 0;
197  else
198  return (k / lambda) * std::pow(x / lambda, k - 1) *
199  std::exp(-std::pow(x / lambda, k));
200 }

Referenced by main().


The documentation for this class was generated from the following files: