Univariate Distributions
Univariate distributions are the distributions whose variate forms are Univariate
(i.e each sample is a scalar). Abstract types for univariate distributions:
const UnivariateDistribution{S<:ValueSupport} = Distribution{Univariate,S}
const DiscreteUnivariateDistribution = Distribution{Univariate, Discrete}
const ContinuousUnivariateDistribution = Distribution{Univariate, Continuous}
Common Interface
A series of methods are implemented for each univariate distribution, which provide useful functionalities such as moment computation, pdf evaluation, and sampling (i.e. random number generation).
Parameter Retrieval
Note: params
are defined for all univariate distributions, while other parameter retrieval methods are only defined for those distributions for which these parameters make sense. See below for details.
StatsBase.params
— Methodparams(d::UnivariateDistribution)
Return a tuple of parameters. Let d
be a distribution of type D
, then D(params(d)...)
will construct exactly the same distribution as $d$.
Distributions.scale
— Methodscale(d::UnivariateDistribution)
Get the scale parameter.
Distributions.location
— Methodlocation(d::UnivariateDistribution)
Get the location parameter.
Distributions.shape
— Methodshape(d::UnivariateDistribution)
Get the shape parameter.
Distributions.rate
— Methodrate(d::UnivariateDistribution)
Get the rate parameter.
Distributions.ncategories
— Methodncategories(d::UnivariateDistribution)
Get the number of categories.
Distributions.ntrials
— Methodntrials(d::UnivariateDistribution)
Get the number of trials.
StatsBase.dof
— Methoddof(d::UnivariateDistribution)
Get the degrees of freedom.
For distributions for which success and failure have a meaning, the following methods are defined:
Distributions.succprob
— Methodsuccprob(d::DiscreteUnivariateDistribution)
Get the probability of success.
Distributions.failprob
— Methodfailprob(d::DiscreteUnivariateDistribution)
Get the probability of failure.
Computation of statistics
Base.maximum
— Methodmaximum(d::UnivariateDistribution)
Return the maximum of the support of d
.
Base.minimum
— Methodminimum(d::UnivariateDistribution)
Return the minimum of the support of d
.
Base.extrema
— Methodextrema(d::UnivariateDistribution)
Return the minimum and maximum of the support of d
as a 2-tuple.
Statistics.mean
— Methodmean(d::UnivariateDistribution)
Compute the expectation.
Statistics.var
— Methodvar(d::UnivariateDistribution)
Compute the variance. (A generic std is provided as std(d) = sqrt(var(d))
)
Statistics.std
— Methodstd(d::UnivariateDistribution)
Return the standard deviation of distribution d
, i.e. sqrt(var(d))
.
Statistics.median
— Methodmedian(d::UnivariateDistribution)
Return the median value of distribution d
.
StatsBase.modes
— Methodmodes(d::UnivariateDistribution)
Get all modes (if this makes sense).
StatsBase.mode
— Methodmode(d::UnivariateDistribution)
Returns the first mode.
StatsBase.skewness
— Methodskewness(d::UnivariateDistribution)
Compute the skewness.
StatsBase.kurtosis
— Methodkurtosis(d::UnivariateDistribution)
Compute the excessive kurtosis.
StatsBase.kurtosis
— Methodkurtosis(d::Distribution, correction::Bool)
Computes excess kurtosis by default. Proper kurtosis can be returned with correction=false
Distributions.isplatykurtic
— Methodisplatykurtic(d)
Return whether d
is platykurtic (i.e kurtosis(d) < 0
).
Distributions.isleptokurtic
— Methodisleptokurtic(d)
Return whether d
is leptokurtic (i.e kurtosis(d) > 0
).
Distributions.ismesokurtic
— Methodismesokurtic(d)
Return whether d
is mesokurtic (i.e kurtosis(d) == 0
).
StatsBase.entropy
— Methodentropy(d::UnivariateDistribution)
Compute the entropy value of distribution d
.
StatsBase.entropy
— Methodentropy(d::UnivariateDistribution, b::Real)
Compute the entropy value of distribution d
, w.r.t. a given base.
StatsBase.entropy
— Methodentropy(d::UnivariateDistribution, b::Real)
Compute the entropy value of distribution d
, w.r.t. a given base.
Distributions.mgf
— Methodmgf(d::UnivariateDistribution, t)
Evaluate the moment generating function of distribution d
.
Distributions.cf
— Methodcf(d::UnivariateDistribution, t)
Evaluate the characteristic function of distribution d
.
Probability Evaluation
Distributions.insupport
— Methodinsupport(d::UnivariateDistribution, x::Any)
When x
is a scalar, it returns whether x is within the support of d
(e.g., insupport(d, x) = minimum(d) <= x <= maximum(d)
). When x
is an array, it returns whether every element in x is within the support of d
.
Generic fallback methods are provided, but it is often the case that insupport
can be done more efficiently, and a specialized insupport
is thus desirable. You should also override this function if the support is composed of multiple disjoint intervals.
Distributions.pdf
— Methodpdf(d::UnivariateDistribution, x::Real)
Evaluate the probability density (mass) at x
.
See also: logpdf
.
Distributions.logpdf
— Methodlogpdf(d::UnivariateDistribution, x::Real)
Evaluate the logarithm of probability density (mass) at x
.
See also: pdf
.
StatsBase.loglikelihood
— Methodloglikelihood(d::UnivariateDistribution, x::Union{Real,AbstractArray})
The log-likelihood of distribution d
with respect to all samples contained in x
.
Here x
can be a single scalar sample or an array of samples.
Distributions.cdf
— Methodcdf(d::UnivariateDistribution, x::Real)
Evaluate the cumulative probability at x
.
Distributions.logcdf
— Methodlogcdf(d::UnivariateDistribution, x::Real)
The logarithm of the cumulative function value(s) evaluated at x
, i.e. log(cdf(x))
.
Distributions.logdiffcdf
— Methodlogdiffcdf(d::UnivariateDistribution, x::Real, y::Real)
The natural logarithm of the difference between the cumulative density function at x
and y
, i.e. log(cdf(x) - cdf(y))
.
Distributions.ccdf
— Methodccdf(d::UnivariateDistribution, x::Real)
The complementary cumulative function evaluated at x
, i.e. 1 - cdf(d, x)
.
Distributions.logccdf
— Methodlogccdf(d::UnivariateDistribution, x::Real)
The logarithm of the complementary cumulative function values evaluated at x, i.e. log(ccdf(x))
.
Statistics.quantile
— Methodquantile(d::UnivariateDistribution, q::Real)
Evaluate the inverse cumulative distribution function at q
.
See also: cquantile
, invlogcdf
, and invlogccdf
.
Distributions.cquantile
— Methodcquantile(d::UnivariateDistribution, q::Real)
The complementary quantile value, i.e. quantile(d, 1-q)
.
Distributions.invlogcdf
— Methodinvlogcdf(d::UnivariateDistribution, lp::Real)
The inverse function of logcdf.
Distributions.invlogccdf
— Methodinvlogccdf(d::UnivariateDistribution, lp::Real)
The inverse function of logccdf.
Sampling (Random number generation)
Base.rand
— Methodrand(rng::AbstractRNG, d::UnivariateDistribution)
Generate a scalar sample from d
. The general fallback is quantile(d, rand())
.
Random.rand!
— Methodrand!(rng::AbstractRNG, ::UnivariateDistribution, ::AbstractArray)
Sample a univariate distribution and store the results in the provided array.
Continuous Distributions
Distributions.Arcsine
— TypeArcsine(a,b)
The Arcsine distribution has probability density function
\[f(x) = \frac{1}{\pi \sqrt{(x - a) (b - x)}}, \quad x \in [a, b]\]
Arcsine() # Arcsine distribution with support [0, 1]
Arcsine(b) # Arcsine distribution with support [0, b]
Arcsine(a, b) # Arcsine distribution with support [a, b]
params(d) # Get the parameters, i.e. (a, b)
minimum(d) # Get the lower bound, i.e. a
maximum(d) # Get the upper bound, i.e. b
location(d) # Get the left bound, i.e. a
scale(d) # Get the span of the support, i.e. b - a
External links
Use Arcsine(a, b, check_args=false)
to bypass argument checks.
Distributions.Beta
— TypeBeta(α,β)
The Beta distribution has probability density function
\[f(x; \alpha, \beta) = \frac{1}{B(\alpha, \beta)} x^{\alpha - 1} (1 - x)^{\beta - 1}, \quad x \in [0, 1]\]
The Beta distribution is related to the Gamma
distribution via the property that if $X \sim \operatorname{Gamma}(\alpha)$ and $Y \sim \operatorname{Gamma}(\beta)$ independently, then $X / (X + Y) \sim Beta(\alpha, \beta)$.
Beta() # equivalent to Beta(1, 1)
Beta(a) # equivalent to Beta(a, a)
Beta(a, b) # Beta distribution with shape parameters a and b
params(d) # Get the parameters, i.e. (a, b)
External links
Distributions.BetaPrime
— TypeBetaPrime(α,β)
The Beta prime distribution has probability density function
\[f(x; \alpha, \beta) = \frac{1}{B(\alpha, \beta)} x^{\alpha - 1} (1 + x)^{- (\alpha + \beta)}, \quad x > 0\]
The Beta prime distribution is related to the Beta
distribution via the relation ship that if $X \sim \operatorname{Beta}(\alpha, \beta)$ then $\frac{X}{1 - X} \sim \operatorname{BetaPrime}(\alpha, \beta)$
BetaPrime() # equivalent to BetaPrime(1, 1)
BetaPrime(a) # equivalent to BetaPrime(a, a)
BetaPrime(a, b) # Beta prime distribution with shape parameters a and b
params(d) # Get the parameters, i.e. (a, b)
External links
Distributions.Biweight
— TypeBiweight(μ, σ)
Distributions.Cauchy
— TypeCauchy(μ, σ)
The Cauchy distribution with location μ
and scale σ
has probability density function
\[f(x; \mu, \sigma) = \frac{1}{\pi \sigma \left(1 + \left(\frac{x - \mu}{\sigma} \right)^2 \right)}\]
Cauchy() # Standard Cauchy distribution, i.e. Cauchy(0, 1)
Cauchy(u) # Cauchy distribution with location u and unit scale, i.e. Cauchy(u, 1)
Cauchy(u, b) # Cauchy distribution with location u and scale b
params(d) # Get the parameters, i.e. (u, b)
location(d) # Get the location parameter, i.e. u
scale(d) # Get the scale parameter, i.e. b
External links
Distributions.Chernoff
— TypeChernoff()
The Chernoff distribution is the distribution of the random variable
\[\underset{t \in (-\infty,\infty)}{\arg\max} ( G(t) - t^2 ),\]
where $G$ is standard two–sided Brownian motion.
The distribution arises as the limit distribution of various cube-root-n consistent estimators, including the isotonic regression estimator of Brunk, the isotonic density estimator of Grenander, the least median of squares estimator of Rousseeuw, and the maximum score estimator of Manski.
For theoretical results, see e.g. Kim and Pollard, Annals of Statistics, 1990. The code for the computation of pdf and cdf is based on the algorithm described in Groeneboom and Wellner, Journal of Computational and Graphical Statistics, 2001.
Chernoff()
pdf(Chernoff(),x::Real)
cdf(Chernoff(),x::Real)
logpdf(Chernoff(),x::Real)
survivor(Chernoff(),x::Real)
mean(Chernoff())
var(Chernoff())
skewness(Chernoff())
kurtosis(Chernoff())
kurtosis(Chernoff(), excess::Bool)
mode(Chernoff())
entropy(Chernoff())
rand(Chernoff())
rand(rng, Chernoff()
cdf(Chernoff(),-x) #For tail probabilities, use this instead of 1-cdf(Chernoff(),x)
Distributions.Chi
— TypeChi(ν)
The Chi distribution ν
degrees of freedom has probability density function
\[f(x; k) = \frac{1}{\Gamma(k/2)} 2^{1 - k/2} x^{k-1} e^{-x^2/2}, \quad x > 0\]
It is the distribution of the square-root of a Chisq
variate.
Chi(k) # Chi distribution with k degrees of freedom
params(d) # Get the parameters, i.e. (k,)
dof(d) # Get the degrees of freedom, i.e. k
External links
Distributions.Chisq
— TypeChisq(ν)
The Chi squared distribution (typically written χ²) with ν
degrees of freedom has the probability density function
\[f(x; k) = \frac{x^{k/2 - 1} e^{-x/2}}{2^{k/2} \Gamma(k/2)}, \quad x > 0.\]
If ν
is an integer, then it is the distribution of the sum of squares of ν
independent standard Normal
variates.
Chisq(k) # Chi-squared distribution with k degrees of freedom
params(d) # Get the parameters, i.e. (k,)
dof(d) # Get the degrees of freedom, i.e. k
External links
Distributions.Cosine
— TypeDistributions.Epanechnikov
— TypeEpanechnikov(μ, σ)
Distributions.Erlang
— TypeErlang(α,θ)
The Erlang distribution is a special case of a Gamma
distribution with integer shape parameter.
Erlang() # Erlang distribution with unit shape and unit scale, i.e. Erlang(1, 1)
Erlang(a) # Erlang distribution with shape parameter a and unit scale, i.e. Erlang(a, 1)
Erlang(a, s) # Erlang distribution with shape parameter a and scale s
External links
Distributions.Exponential
— TypeExponential(θ)
The Exponential distribution with scale parameter θ
has probability density function
\[f(x; \theta) = \frac{1}{\theta} e^{-\frac{x}{\theta}}, \quad x > 0\]
Exponential() # Exponential distribution with unit scale, i.e. Exponential(1)
Exponential(b) # Exponential distribution with scale b
params(d) # Get the parameters, i.e. (b,)
scale(d) # Get the scale parameter, i.e. b
rate(d) # Get the rate parameter, i.e. 1 / b
External links
Distributions.FDist
— TypeFDist(ν1, ν2)
The F distribution has probability density function
\[f(x; \nu_1, \nu_2) = \frac{1}{x B(\nu_1/2, \nu_2/2)} \sqrt{\frac{(\nu_1 x)^{\nu_1} \cdot \nu_2^{\nu_2}}{(\nu_1 x + \nu_2)^{\nu_1 + \nu_2}}}, \quad x>0\]
It is related to the Chisq
distribution via the property that if $X_1 \sim \operatorname{Chisq}(\nu_1)$ and $X_2 \sim \operatorname{Chisq}(\nu_2)$, then $(X_1/\nu_1) / (X_2 / \nu_2) \sim \operatorname{FDist}(\nu_1, \nu_2)$.
FDist(ν1, ν2) # F-Distribution with parameters ν1 and ν2
params(d) # Get the parameters, i.e. (ν1, ν2)
External links
Distributions.Frechet
— TypeFrechet(α,θ)
The Fréchet distribution with shape α
and scale θ
has probability density function
\[f(x; \alpha, \theta) = \frac{\alpha}{\theta} \left( \frac{x}{\theta} \right)^{-\alpha-1} e^{-(x/\theta)^{-\alpha}}, \quad x > 0\]
Frechet() # Fréchet distribution with unit shape and unit scale, i.e. Frechet(1, 1)
Frechet(a) # Fréchet distribution with shape a and unit scale, i.e. Frechet(a, 1)
Frechet(a, b) # Fréchet distribution with shape a and scale b
params(d) # Get the parameters, i.e. (a, b)
shape(d) # Get the shape parameter, i.e. a
scale(d) # Get the scale parameter, i.e. b
External links
Distributions.Gamma
— TypeGamma(α,θ)
The Gamma distribution with shape parameter α
and scale θ
has probability density function
\[f(x; \alpha, \theta) = \frac{x^{\alpha-1} e^{-x/\theta}}{\Gamma(\alpha) \theta^\alpha}, \quad x > 0\]
Gamma() # Gamma distribution with unit shape and unit scale, i.e. Gamma(1, 1)
Gamma(α) # Gamma distribution with shape α and unit scale, i.e. Gamma(α, 1)
Gamma(α, θ) # Gamma distribution with shape α and scale θ
params(d) # Get the parameters, i.e. (α, θ)
shape(d) # Get the shape parameter, i.e. α
scale(d) # Get the scale parameter, i.e. θ
External links
Distributions.GeneralizedExtremeValue
— TypeGeneralizedExtremeValue(μ, σ, ξ)
The Generalized extreme value distribution with shape parameter ξ
, scale σ
and location μ
has probability density function
\[f(x; \xi, \sigma, \mu) = \begin{cases} \frac{1}{\sigma} \left[ 1+\left(\frac{x-\mu}{\sigma}\right)\xi\right]^{-1/\xi-1} \exp\left\{-\left[ 1+ \left(\frac{x-\mu}{\sigma}\right)\xi\right]^{-1/\xi} \right\} & \text{for } \xi \neq 0 \\\ \frac{1}{\sigma} \exp\left\{-\frac{x-\mu}{\sigma}\right\} \exp\left\{-\exp\left[-\frac{x-\mu}{\sigma}\right]\right\} & \text{for } \xi = 0 \\ \end{cases}\]
for
\[x \in \begin{cases} \left[ \mu - \frac{\sigma}{\xi}, + \infty \right) & \text{for } \xi > 0 \\ \left( - \infty, + \infty \right) & \text{for } \xi = 0 \\ \left( - \infty, \mu - \frac{\sigma}{\xi} \right] & \text{for } \xi < 0 \end{cases}\]
GeneralizedExtremeValue(m, s, k) # Generalized Pareto distribution with shape k, scale s and location m.
params(d) # Get the parameters, i.e. (m, s, k)
location(d) # Get the location parameter, i.e. m
scale(d) # Get the scale parameter, i.e. s
shape(d) # Get the shape parameter, i.e. k (sometimes called c)
External links
Distributions.GeneralizedPareto
— TypeGeneralizedPareto(μ, σ, ξ)
The Generalized Pareto distribution with shape parameter ξ
, scale σ
and location μ
has probability density function
\[f(x; \mu, \sigma, \xi) = \begin{cases} \frac{1}{\sigma}(1 + \xi \frac{x - \mu}{\sigma} )^{-\frac{1}{\xi} - 1} & \text{for } \xi \neq 0 \\ \frac{1}{\sigma} e^{-\frac{\left( x - \mu \right) }{\sigma}} & \text{for } \xi = 0 \end{cases}~, \quad x \in \begin{cases} \left[ \mu, \infty \right] & \text{for } \xi \geq 0 \\ \left[ \mu, \mu - \sigma / \xi \right] & \text{for } \xi < 0 \end{cases}\]
GeneralizedPareto() # Generalized Pareto distribution with unit shape and unit scale, i.e. GeneralizedPareto(0, 1, 1)
GeneralizedPareto(k, s) # Generalized Pareto distribution with shape k and scale s, i.e. GeneralizedPareto(0, k, s)
GeneralizedPareto(m, k, s) # Generalized Pareto distribution with shape k, scale s and location m.
params(d) # Get the parameters, i.e. (m, s, k)
location(d) # Get the location parameter, i.e. m
scale(d) # Get the scale parameter, i.e. s
shape(d) # Get the shape parameter, i.e. k
External links
Distributions.Gumbel
— TypeGumbel(μ, θ)
The Gumbel (maxima) distribution with location μ
and scale θ
has probability density function
\[f(x; \mu, \theta) = \frac{1}{\theta} e^{-(z + e^{-z})}, \quad \text{ with } z = \frac{x - \mu}{\theta}\]
Gumbel() # Gumbel distribution with zero location and unit scale, i.e. Gumbel(0, 1)
Gumbel(u) # Gumbel distribution with location u and unit scale, i.e. Gumbel(u, 1)
Gumbel(u, b) # Gumbel distribution with location u and scale b
params(d) # Get the parameters, i.e. (u, b)
location(d) # Get the location parameter, i.e. u
scale(d) # Get the scale parameter, i.e. b
External links
Distributions.InverseGamma
— TypeInverseGamma(α, θ)
The inverse Gamma distribution with shape parameter α
and scale θ
has probability density function
\[f(x; \alpha, \theta) = \frac{\theta^\alpha x^{-(\alpha + 1)}}{\Gamma(\alpha)} e^{-\frac{\theta}{x}}, \quad x > 0\]
It is related to the Gamma
distribution: if $X \sim \operatorname{Gamma}(\alpha, \beta)$, then $1 / X \sim \operatorname{InverseGamma}(\alpha, \beta^{-1})$.
InverseGamma() # Inverse Gamma distribution with unit shape and unit scale, i.e. InverseGamma(1, 1)
InverseGamma(α) # Inverse Gamma distribution with shape α and unit scale, i.e. InverseGamma(α, 1)
InverseGamma(α, θ) # Inverse Gamma distribution with shape α and scale θ
params(d) # Get the parameters, i.e. (α, θ)
shape(d) # Get the shape parameter, i.e. α
scale(d) # Get the scale parameter, i.e. θ
External links
Distributions.InverseGaussian
— TypeInverseGaussian(μ,λ)
The inverse Gaussian distribution with mean μ
and shape λ
has probability density function
\[f(x; \mu, \lambda) = \sqrt{\frac{\lambda}{2\pi x^3}} \exp\!\left(\frac{-\lambda(x-\mu)^2}{2\mu^2x}\right), \quad x > 0\]
InverseGaussian() # Inverse Gaussian distribution with unit mean and unit shape, i.e. InverseGaussian(1, 1)
InverseGaussian(mu), # Inverse Gaussian distribution with mean mu and unit shape, i.e. InverseGaussian(u, 1)
InverseGaussian(mu, lambda) # Inverse Gaussian distribution with mean mu and shape lambda
params(d) # Get the parameters, i.e. (mu, lambda)
mean(d) # Get the mean parameter, i.e. mu
shape(d) # Get the shape parameter, i.e. lambda
External links
Distributions.Kolmogorov
— TypeKolmogorov()
Kolmogorov distribution defined as
\[\sup_{t \in [0,1]} |B(t)|\]
where $B(t)$ is a Brownian bridge used in the Kolmogorov–Smirnov test for large n.
Distributions.KSDist
— TypeKSDist(n)
Distribution of the (two-sided) Kolmogorov-Smirnoff statistic
\[D_n = \sup_x | \hat{F}_n(x) -F(x)| \sqrt(n)\]
$D_n$ converges a.s. to the Kolmogorov distribution.
Distributions.KSOneSided
— TypeKSOneSided(n)
Distribution of the one-sided Kolmogorov-Smirnov test statistic:
\[D^+_n = \sup_x (\hat{F}_n(x) -F(x))\]
Distributions.Laplace
— TypeLaplace(μ,β)
The Laplace distribution with location μ
and scale β
has probability density function
\[f(x; \mu, \beta) = \frac{1}{2 \beta} \exp \left(- \frac{|x - \mu|}{\beta} \right)\]
Laplace() # Laplace distribution with zero location and unit scale, i.e. Laplace(0, 1)
Laplace(μ) # Laplace distribution with location μ and unit scale, i.e. Laplace(μ, 1)
Laplace(μ, β) # Laplace distribution with location μ and scale β
params(d) # Get the parameters, i.e., (μ, β)
location(d) # Get the location parameter, i.e. μ
scale(d) # Get the scale parameter, i.e. β
External links
Distributions.Levy
— TypeLevy(μ, σ)
The Lévy distribution with location μ
and scale σ
has probability density function
\[f(x; \mu, \sigma) = \sqrt{\frac{\sigma}{2 \pi (x - \mu)^3}} \exp \left( - \frac{\sigma}{2 (x - \mu)} \right), \quad x > \mu\]
Levy() # Levy distribution with zero location and unit scale, i.e. Levy(0, 1)
Levy(u) # Levy distribution with location u and unit scale, i.e. Levy(u, 1)
Levy(u, c) # Levy distribution with location u ans scale c
params(d) # Get the parameters, i.e. (u, c)
location(d) # Get the location parameter, i.e. u
External links
Distributions.LocationScale
— TypeLocationScale(μ,σ,ρ)
A location-scale transformed distribution with location parameter μ
, scale parameter σ
, and given distribution ρ
.
\[f(x) = \frac{1}{σ} ρ \! \left( \frac{x-μ}{σ} \right)\]
LocationScale(μ,σ,ρ) # location-scale transformed distribution
params(d) # Get the parameters, i.e. (μ, σ, and the base distribution)
location(d) # Get the location parameter
scale(d) # Get the scale parameter
External links Location-Scale family on Wikipedia
Distributions.Logistic
— TypeLogistic(μ,θ)
The Logistic distribution with location μ
and scale θ
has probability density function
\[f(x; \mu, \theta) = \frac{1}{4 \theta} \mathrm{sech}^2 \left( \frac{x - \mu}{2 \theta} \right)\]
Logistic() # Logistic distribution with zero location and unit scale, i.e. Logistic(0, 1)
Logistic(u) # Logistic distribution with location u and unit scale, i.e. Logistic(u, 1)
Logistic(u, b) # Logistic distribution with location u ans scale b
params(d) # Get the parameters, i.e. (u, b)
location(d) # Get the location parameter, i.e. u
scale(d) # Get the scale parameter, i.e. b
External links
Distributions.LogitNormal
— TypeLogitNormal(μ,σ)
The logit normal distribution is the distribution of of a random variable whose logit has a Normal
distribution. Or inversely, when applying the logistic function to a Normal random variable then the resulting random variable follows a logit normal distribution.
If $X \sim \operatorname{Normal}(\mu, \sigma)$ then $\operatorname{logistic}(X) \sim \operatorname{LogitNormal}(\mu,\sigma)$.
The probability density function is
\[f(x; \mu, \sigma) = \frac{1}{x \sqrt{2 \pi \sigma^2}} \exp \left( - \frac{(\text{logit}(x) - \mu)^2}{2 \sigma^2} \right), \quad x > 0\]
where the logit-Function is
\[\text{logit}(x) = \ln\left(\frac{x}{1-x}\right) \quad 0 < x < 1\]
LogitNormal() # Logit-normal distribution with zero logit-mean and unit scale
LogitNormal(mu) # Logit-normal distribution with logit-mean mu and unit scale
LogitNormal(mu, sig) # Logit-normal distribution with logit-mean mu and scale sig
params(d) # Get the parameters, i.e. (mu, sig)
median(d) # Get the median, i.e. logistic(mu)
The following properties have no analytical solution but numerical approximations. In order to avoid package dependencies for numerical optimization, they are currently not implemented.
mean(d)
var(d)
std(d)
mode(d)
Similarly, skewness, kurtosis, and entropy are not implemented.
External links
Distributions.LogNormal
— TypeLogNormal(μ,σ)
The log normal distribution is the distribution of the exponential of a Normal
variate: if $X \sim \operatorname{Normal}(\mu, \sigma)$ then $\exp(X) \sim \operatorname{LogNormal}(\mu,\sigma)$. The probability density function is
\[f(x; \mu, \sigma) = \frac{1}{x \sqrt{2 \pi \sigma^2}} \exp \left( - \frac{(\log(x) - \mu)^2}{2 \sigma^2} \right), \quad x > 0\]
LogNormal() # Log-normal distribution with zero log-mean and unit scale
LogNormal(mu) # Log-normal distribution with log-mean mu and unit scale
LogNormal(mu, sig) # Log-normal distribution with log-mean mu and scale sig
params(d) # Get the parameters, i.e. (mu, sig)
meanlogx(d) # Get the mean of log(X), i.e. mu
varlogx(d) # Get the variance of log(X), i.e. sig^2
stdlogx(d) # Get the standard deviation of log(X), i.e. sig
External links
Distributions.NoncentralBeta
— TypeNoncentralBeta(α, β, λ)
Distributions.NoncentralChisq
— TypeNoncentralChisq(ν, λ)
The noncentral chi-squared distribution with ν
degrees of freedom and noncentrality parameter λ
has the probability density function
\[f(x; \nu, \lambda) = \frac{1}{2} e^{-(x + \lambda)/2} \left( \frac{x}{\lambda} \right)^{\nu/4-1/2} I_{\nu/2-1}(\sqrt{\lambda x}), \quad x > 0\]
It is the distribution of the sum of squares of ν
independent Normal
variates with individual means $\mu_i$ and
\[\lambda = \sum_{i=1}^\nu \mu_i^2\]
NoncentralChisq(ν, λ) # Noncentral chi-squared distribution with ν degrees of freedom and noncentrality parameter λ
params(d) # Get the parameters, i.e. (ν, λ)
External links
Distributions.NoncentralF
— TypeNoncentralF(ν1, ν2, λ)
Distributions.NoncentralT
— TypeNoncentralT(ν, λ)
Distributions.Normal
— TypeNormal(μ,σ)
The Normal distribution with mean μ
and standard deviation σ≥0
has probability density function
\[f(x; \mu, \sigma) = \frac{1}{\sqrt{2 \pi \sigma^2}} \exp \left( - \frac{(x - \mu)^2}{2 \sigma^2} \right)\]
Note that if σ == 0
, then the distribution is a point mass concentrated at μ
. Though not technically a continuous distribution, it is allowed so as to account for cases where σ
may have underflowed, and the functions are defined by taking the pointwise limit as $σ → 0$.
Normal() # standard Normal distribution with zero mean and unit variance
Normal(mu) # Normal distribution with mean mu and unit variance
Normal(mu, sig) # Normal distribution with mean mu and variance sig^2
params(d) # Get the parameters, i.e. (mu, sig)
mean(d) # Get the mean, i.e. mu
std(d) # Get the standard deviation, i.e. sig
External links
Distributions.NormalCanon
— TypeNormalCanon(η, λ)
Canonical Form of Normal distribution
Distributions.NormalInverseGaussian
— TypeNormalInverseGaussian(μ,α,β,δ)
The Normal-inverse Gaussian distribution with location μ
, tail heaviness α
, asymmetry parameter β
and scale δ
has probability density function
\[f(x; \mu, \alpha, \beta, \delta) = \frac{\alpha\delta K_1 \left(\alpha\sqrt{\delta^2 + (x - \mu)^2}\right)}{\pi \sqrt{\delta^2 + (x - \mu)^2}} \; e^{\delta \gamma + \beta (x - \mu)}\]
where $K_j$ denotes a modified Bessel function of the third kind.
External links
Distributions.Pareto
— TypePareto(α,θ)
The Pareto distribution with shape α
and scale θ
has probability density function
\[f(x; \alpha, \theta) = \frac{\alpha \theta^\alpha}{x^{\alpha + 1}}, \quad x \ge \theta\]
Pareto() # Pareto distribution with unit shape and unit scale, i.e. Pareto(1, 1)
Pareto(a) # Pareto distribution with shape a and unit scale, i.e. Pareto(a, 1)
Pareto(a, b) # Pareto distribution with shape a and scale b
params(d) # Get the parameters, i.e. (a, b)
shape(d) # Get the shape parameter, i.e. a
scale(d) # Get the scale parameter, i.e. b
External links
Distributions.PGeneralizedGaussian
— TypePGeneralizedGaussian(α, μ, p)
The p-Generalized Gaussian distribution, more commonly known as the exponential power or the generalized normal distribution, with scale α
, location μ
, and shape p
has the probability density function
\[f(x, \mu, \alpha, p) = \frac{p}{2\alpha\Gamma(1/p)} e^{-(\frac{|x-\mu|}{\alpha})^p} \quad x \in (-\infty, +\infty) , \alpha > 0, p > 0\]
The p-Generalized Gaussian (GGD) is a parametric distribution that incorporates the Normal and Laplacian distributions as special cases where p = 1
and p = 2
. As p → ∞
, the distribution approaches the Uniform distribution on [μ-α, μ+α]
.
PGeneralizedGaussian() # GGD with shape 2, scale 1, location 0, (the Normal distribution)
PGeneralizedGaussian(loc,s,sh) # GGD with location loc, scale s, and shape sh
params(d) # Get the parameters, i.e. (loc,s,sh,)
shape(d) # Get the shape parameter, sh
scale(d) # Get the scale parameter, s
location(d) # Get the location parameter, loc
External Links
Distributions.Rayleigh
— TypeRayleigh(σ)
The Rayleigh distribution with scale σ
has probability density function
\[f(x; \sigma) = \frac{x}{\sigma^2} e^{-\frac{x^2}{2 \sigma^2}}, \quad x > 0\]
It is related to the Normal
distribution via the property that if $X, Y \sim \operatorname{Normal}(0,\sigma)$, independently, then $\sqrt{X^2 + Y^2} \sim \operatorname{Rayleigh}(\sigma)$.
Rayleigh() # Rayleigh distribution with unit scale, i.e. Rayleigh(1)
Rayleigh(s) # Rayleigh distribution with scale s
params(d) # Get the parameters, i.e. (s,)
scale(d) # Get the scale parameter, i.e. s
External links
Distributions.Semicircle
— TypeSemicircle(r)
The Wigner semicircle distribution with radius parameter r
has probability density function
\[f(x; r) = \frac{2}{\pi r^2} \sqrt{r^2 - x^2}, \quad x \in [-r, r].\]
Semicircle(r) # Wigner semicircle distribution with radius r
params(d) # Get the radius parameter, i.e. (r,)
External links
Distributions.StudentizedRange
— TypeStudentizedRange(ν, k)
The studentized range distribution has probability density function:
\[f(q; k, \nu) = \frac{\sqrt{2\pi}k(k - 1)\nu^{\nu/2}}{\Gamma{\left(\frac{\nu}{2}\right)}2^{\nu/2 - 1}} \int_{0}^{\infty} {x^{\nu}\phi(\sqrt{\nu}x)} \left[\int_{-\infty}^{\infty} {\phi(u)\phi(u - qx)[\Phi(u) - \Phi(u - qx)]^{k - 2}}du\right]dx\]
where
\[\begin{aligned} \Phi(x) &= \frac{1 + erf(\frac{x}{\sqrt{2}})}{2} &&(\text{Normal Distribution CDF})\\ \phi(x) &= \Phi'(x) &&(\text{Normal Distribution PDF}) \end{aligned}\]
StudentizedRange(ν, k) # Studentized Range Distribution with parameters ν and k
params(d) # Get the parameters, i.e. (ν, k)
External links
Distributions.SymTriangularDist
— TypeSymTriangularDist(μ,σ)
The Symmetric triangular distribution with location μ
and scale σ
has probability density function
\[f(x; \mu, \sigma) = \frac{1}{\sigma} \left( 1 - \left| \frac{x - \mu}{\sigma} \right| \right), \quad \mu - \sigma \le x \le \mu + \sigma\]
SymTriangularDist() # Symmetric triangular distribution with zero location and unit scale
SymTriangularDist(u) # Symmetric triangular distribution with location u and unit scale
SymTriangularDist(u, s) # Symmetric triangular distribution with location u and scale s
params(d) # Get the parameters, i.e. (u, s)
location(d) # Get the location parameter, i.e. u
scale(d) # Get the scale parameter, i.e. s
Distributions.TDist
— TypeTDist(ν)
The Students T distribution with ν
degrees of freedom has probability density function
\[f(x; d) = \frac{1}{\sqrt{d} B(1/2, d/2)} \left( 1 + \frac{x^2}{d} \right)^{-\frac{d + 1}{2}}\]
TDist(d) # t-distribution with d degrees of freedom
params(d) # Get the parameters, i.e. (d,)
dof(d) # Get the degrees of freedom, i.e. d
External links
Distributions.TriangularDist
— TypeTriangularDist(a,b,c)
The triangular distribution with lower limit a
, upper limit b
and mode c
has probability density function
\[f(x; a, b, c)= \begin{cases} 0 & \mathrm{for\ } x < a, \\ \frac{2(x-a)}{(b-a)(c-a)} & \mathrm{for\ } a \le x \leq c, \\[4pt] \frac{2(b-x)}{(b-a)(b-c)} & \mathrm{for\ } c < x \le b, \\[4pt] 0 & \mathrm{for\ } b < x, \end{cases}\]
TriangularDist(a, b) # Triangular distribution with lower limit a, upper limit b, and mode (a+b)/2
TriangularDist(a, b, c) # Triangular distribution with lower limit a, upper limit b, and mode c
params(d) # Get the parameters, i.e. (a, b, c)
minimum(d) # Get the lower bound, i.e. a
maximum(d) # Get the upper bound, i.e. b
mode(d) # Get the mode, i.e. c
External links
Distributions.Triweight
— TypeTriweight(μ, σ)
Distributions.Uniform
— TypeUniform(a,b)
The continuous uniform distribution over an interval $[a, b]$ has probability density function
\[f(x; a, b) = \frac{1}{b - a}, \quad a \le x \le b\]
Uniform() # Uniform distribution over [0, 1]
Uniform(a, b) # Uniform distribution over [a, b]
params(d) # Get the parameters, i.e. (a, b)
minimum(d) # Get the lower bound, i.e. a
maximum(d) # Get the upper bound, i.e. b
location(d) # Get the location parameter, i.e. a
scale(d) # Get the scale parameter, i.e. b - a
External links
Distributions.VonMises
— TypeVonMises(μ, κ)
The von Mises distribution with mean μ
and concentration κ
has probability density function
\[f(x; \mu, \kappa) = \frac{1}{2 \pi I_0(\kappa)} \exp \left( \kappa \cos (x - \mu) \right)\]
VonMises() # von Mises distribution with zero mean and unit concentration
VonMises(κ) # von Mises distribution with zero mean and concentration κ
VonMises(μ, κ) # von Mises distribution with mean μ and concentration κ
External links
Distributions.Weibull
— TypeWeibull(α,θ)
The Weibull distribution with shape α
and scale θ
has probability density function
\[f(x; \alpha, \theta) = \frac{\alpha}{\theta} \left( \frac{x}{\theta} \right)^{\alpha-1} e^{-(x/\theta)^\alpha}, \quad x \ge 0\]
Weibull() # Weibull distribution with unit shape and unit scale, i.e. Weibull(1, 1)
Weibull(a) # Weibull distribution with shape a and unit scale, i.e. Weibull(a, 1)
Weibull(a, b) # Weibull distribution with shape a and scale b
params(d) # Get the parameters, i.e. (a, b)
shape(d) # Get the shape parameter, i.e. a
scale(d) # Get the scale parameter, i.e. b
External links
Discrete Distributions
Distributions.Bernoulli
— TypeBernoulli(p)
A Bernoulli distribution is parameterized by a success rate p
, which takes value 1 with probability p
and 0 with probability 1-p
.
\[P(X = k) = \begin{cases} 1 - p & \quad \text{for } k = 0, \\ p & \quad \text{for } k = 1. \end{cases}\]
Bernoulli() # Bernoulli distribution with p = 0.5
Bernoulli(p) # Bernoulli distribution with success rate p
params(d) # Get the parameters, i.e. (p,)
succprob(d) # Get the success rate, i.e. p
failprob(d) # Get the failure rate, i.e. 1 - p
External links:
Distributions.BetaBinomial
— TypeBetaBinomial(n,α,β)
A Beta-binomial distribution is the compound distribution of the Binomial
distribution where the probability of success p
is distributed according to the Beta
. It has three parameters: n
, the number of trials and two shape parameters α
, β
\[P(X = k) = {n \choose k} B(k + \alpha, n - k + \beta) / B(\alpha, \beta), \quad \text{ for } k = 0,1,2, \ldots, n.\]
BetaBinomial(n, a, b) # BetaBinomial distribution with n trials and shape parameters a, b
params(d) # Get the parameters, i.e. (n, a, b)
ntrials(d) # Get the number of trials, i.e. n
External links:
Distributions.Binomial
— TypeBinomial(n,p)
A Binomial distribution characterizes the number of successes in a sequence of independent trials. It has two parameters: n
, the number of trials, and p
, the probability of success in an individual trial, with the distribution:
\[P(X = k) = {n \choose k}p^k(1-p)^{n-k}, \quad \text{ for } k = 0,1,2, \ldots, n.\]
Binomial() # Binomial distribution with n = 1 and p = 0.5
Binomial(n) # Binomial distribution for n trials with success rate p = 0.5
Binomial(n, p) # Binomial distribution for n trials with success rate p
params(d) # Get the parameters, i.e. (n, p)
ntrials(d) # Get the number of trials, i.e. n
succprob(d) # Get the success rate, i.e. p
failprob(d) # Get the failure rate, i.e. 1 - p
External links:
Distributions.Categorical
— TypeCategorical(p)
A Categorical distribution is parameterized by a probability vector p
(of length K
).
\[P(X = k) = p[k] \quad \text{for } k = 1, 2, \ldots, K.\]
Categorical(p) # Categorical distribution with probability vector p
params(d) # Get the parameters, i.e. (p,)
probs(d) # Get the probability vector, i.e. p
ncategories(d) # Get the number of categories, i.e. K
Here, p
must be a real vector, of which all components are nonnegative and sum to one.
Note: The input vector p
is directly used as a field of the constructed distribution, without being copied.
Categorical
is simply a type alias describing a special case of a DiscreteNonParametric
distribution, so non-specialized methods defined for DiscreteNonParametric
apply to Categorical
as well.
External links:
Distributions.DiscreteUniform
— TypeDiscreteUniform(a,b)
A Discrete uniform distribution is a uniform distribution over a consecutive sequence of integers between a
and b
, inclusive.
\[P(X = k) = 1 / (b - a + 1) \quad \text{for } k = a, a+1, \ldots, b.\]
DiscreteUniform(a, b) # a uniform distribution over {a, a+1, ..., b}
params(d) # Get the parameters, i.e. (a, b)
span(d) # Get the span of the support, i.e. (b - a + 1)
probval(d) # Get the probability value, i.e. 1 / (b - a + 1)
minimum(d) # Return a
maximum(d) # Return b
External links
Distributions.DiscreteNonParametric
— TypeDiscreteNonParametric(xs, ps)
A Discrete nonparametric distribution explicitly defines an arbitrary probability mass function in terms of a list of real support values and their corresponding probabilities
d = DiscreteNonParametric(xs, ps)
params(d) # Get the parameters, i.e. (xs, ps)
support(d) # Get a sorted AbstractVector describing the support (xs) of the distribution
probs(d) # Get a Vector of the probabilities (ps) associated with the support
External links
Distributions.Geometric
— TypeGeometric(p)
A Geometric distribution characterizes the number of failures before the first success in a sequence of independent Bernoulli trials with success rate p
.
\[P(X = k) = p (1 - p)^k, \quad \text{for } k = 0, 1, 2, \ldots.\]
Geometric() # Geometric distribution with success rate 0.5
Geometric(p) # Geometric distribution with success rate p
params(d) # Get the parameters, i.e. (p,)
succprob(d) # Get the success rate, i.e. p
failprob(d) # Get the failure rate, i.e. 1 - p
External links
Distributions.Hypergeometric
— TypeHypergeometric(s, f, n)
A Hypergeometric distribution describes the number of successes in n
draws without replacement from a finite population containing s
successes and f
failures.
\[P(X = k) = {{{s \choose k} {f \choose {n-k}}}\over {s+f \choose n}}, \quad \text{for } k = \max(0, n - f), \ldots, \min(n, s).\]
Hypergeometric(s, f, n) # Hypergeometric distribution for a population with
# s successes and f failures, and a sequence of n trials.
params(d) # Get the parameters, i.e. (s, f, n)
External links
Distributions.NegativeBinomial
— TypeNegativeBinomial(r,p)
A Negative binomial distribution describes the number of failures before the r
th success in a sequence of independent Bernoulli trials. It is parameterized by r
, the number of successes, and p
, the probability of success in an individual trial.
\[P(X = k) = {k + r - 1 \choose k} p^r (1 - p)^k, \quad \text{for } k = 0,1,2,\ldots.\]
The distribution remains well-defined for any positive r
, in which case
\[P(X = k) = \frac{\Gamma(k+r)}{k! \Gamma(r)} p^r (1 - p)^k, \quad \text{for } k = 0,1,2,\ldots.\]
NegativeBinomial() # Negative binomial distribution with r = 1 and p = 0.5
NegativeBinomial(r, p) # Negative binomial distribution with r successes and success rate p
params(d) # Get the parameters, i.e. (r, p)
succprob(d) # Get the success rate, i.e. p
failprob(d) # Get the failure rate, i.e. 1 - p
External links:
Note: The definition of the negative binomial distribution in Wolfram is different from the Wikipedia definition. In Wikipedia, r
is the number of failures and k
is the number of successes.
Distributions.Poisson
— TypePoisson(λ)
A Poisson distribution descibes the number of independent events occurring within a unit time interval, given the average rate of occurrence λ
.
\[P(X = k) = \frac{\lambda^k}{k!} e^{-\lambda}, \quad \text{ for } k = 0,1,2,\ldots.\]
Poisson() # Poisson distribution with rate parameter 1
Poisson(lambda) # Poisson distribution with rate parameter lambda
params(d) # Get the parameters, i.e. (λ,)
mean(d) # Get the mean arrival rate, i.e. λ
External links:
Distributions.PoissonBinomial
— TypePoissonBinomial(p)
A Poisson-binomial distribution describes the number of successes in a sequence of independent trials, wherein each trial has a different success rate. It is parameterized by a vector p
(of length $K$), where $K$ is the total number of trials and p[i]
corresponds to the probability of success of the i
th trial.
\[P(X = k) = \sum\limits_{A\in F_k} \prod\limits_{i\in A} p[i] \prod\limits_{j\in A^c} (1-p[j]), \quad \text{ for } k = 0,1,2,\ldots,K,\]
where $F_k$ is the set of all subsets of $k$ integers that can be selected from $\{1,2,3,...,K\}$.
PoissonBinomial(p) # Poisson Binomial distribution with success rate vector p
params(d) # Get the parameters, i.e. (p,)
succprob(d) # Get the vector of success rates, i.e. p
failprob(d) # Get the vector of failure rates, i.e. 1-p
External links:
Distributions.Skellam
— TypeSkellam(μ1, μ2)
A Skellam distribution describes the difference between two independent Poisson
variables, respectively with rate μ1
and μ2
.
\[P(X = k) = e^{-(\mu_1 + \mu_2)} \left( \frac{\mu_1}{\mu_2} \right)^{k/2} I_k(2 \sqrt{\mu_1 \mu_2}) \quad \text{for integer } k\]
where $I_k$ is the modified Bessel function of the first kind.
Skellam(mu1, mu2) # Skellam distribution for the difference between two Poisson variables,
# respectively with expected values mu1 and mu2.
params(d) # Get the parameters, i.e. (mu1, mu2)
External links:
Vectorized evaluation
Vectorized computation and inplace vectorized computation have been deprecated.
Index
Distributions.Arcsine
Distributions.Bernoulli
Distributions.Beta
Distributions.BetaBinomial
Distributions.BetaPrime
Distributions.Binomial
Distributions.Biweight
Distributions.Categorical
Distributions.Cauchy
Distributions.Chernoff
Distributions.Chi
Distributions.Chisq
Distributions.Cosine
Distributions.DiscreteNonParametric
Distributions.DiscreteUniform
Distributions.Epanechnikov
Distributions.Erlang
Distributions.Exponential
Distributions.FDist
Distributions.Frechet
Distributions.Gamma
Distributions.GeneralizedExtremeValue
Distributions.GeneralizedPareto
Distributions.Geometric
Distributions.Gumbel
Distributions.Hypergeometric
Distributions.InverseGamma
Distributions.InverseGaussian
Distributions.KSDist
Distributions.KSOneSided
Distributions.Kolmogorov
Distributions.Laplace
Distributions.Levy
Distributions.LocationScale
Distributions.LogNormal
Distributions.Logistic
Distributions.LogitNormal
Distributions.NegativeBinomial
Distributions.NoncentralBeta
Distributions.NoncentralChisq
Distributions.NoncentralF
Distributions.NoncentralT
Distributions.Normal
Distributions.NormalCanon
Distributions.NormalInverseGaussian
Distributions.PGeneralizedGaussian
Distributions.Pareto
Distributions.Poisson
Distributions.PoissonBinomial
Distributions.Rayleigh
Distributions.Semicircle
Distributions.Skellam
Distributions.StudentizedRange
Distributions.SymTriangularDist
Distributions.TDist
Distributions.TriangularDist
Distributions.Triweight
Distributions.Uniform
Distributions.VonMises
Distributions.Weibull
Base.extrema
Base.maximum
Base.minimum
Base.rand
Distributions.ccdf
Distributions.cdf
Distributions.cf
Distributions.cquantile
Distributions.failprob
Distributions.insupport
Distributions.invlogccdf
Distributions.invlogcdf
Distributions.isleptokurtic
Distributions.ismesokurtic
Distributions.isplatykurtic
Distributions.location
Distributions.logccdf
Distributions.logcdf
Distributions.logdiffcdf
Distributions.logpdf
Distributions.mgf
Distributions.ncategories
Distributions.ntrials
Distributions.pdf
Distributions.rate
Distributions.scale
Distributions.shape
Distributions.succprob
Random.rand!
Statistics.mean
Statistics.median
Statistics.quantile
Statistics.std
Statistics.var
StatsBase.dof
StatsBase.entropy
StatsBase.entropy
StatsBase.entropy
StatsBase.kurtosis
StatsBase.kurtosis
StatsBase.loglikelihood
StatsBase.mode
StatsBase.modes
StatsBase.params
StatsBase.skewness