Univariate Distributions

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

StatsBase.paramsMethod.
params(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$.

source
succprob(d::UnivariateDistribution)

Get the probability of success.

source
failprob(d::UnivariateDistribution)

Get the probability of failure.

source
scale(d::UnivariateDistribution)

Get the scale parameter.

source
location(d::UnivariateDistribution)

Get the location parameter.

source
shape(d::UnivariateDistribution)

Get the shape parameter.

source
Distributions.rateMethod.
rate(d::UnivariateDistribution)

Get the rate parameter.

source
ncategories(d::UnivariateDistribution)

Get the number of categories.

source
ntrials(d::UnivariateDistribution)

Get the number of trials.

source
StatsBase.dofMethod.
dof(d::UnivariateDistribution)

Get the degrees of freedom.

source

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.

Computation of statistics

Base.maximumMethod.
maximum(d::Distribution)

Return the maximum of the support of d.

source
Base.minimumMethod.
minimum(d::Distribution)

Return the minimum of the support of d.

source
Base.meanMethod.
mean(d::UnivariateDistribution)

Compute the expectation.

source
Base.varMethod.
var(d::UnivariateDistribution)

Compute the variance. (A generic std is provided as std(d) = sqrt(var(d)))

source
Base.stdMethod.
std(d::UnivariateDistribution)

Return the standard deviation of distribution d, i.e. sqrt(var(d)).

source
Base.medianMethod.
median(d::UnivariateDistribution)

Return the median value of distribution d.

source
StatsBase.modesMethod.
modes(d::UnivariateDistribution)

Get all modes (if this makes sense).

source
StatsBase.modeMethod.
mode(d::UnivariateDistribution)

Returns the first mode.

source
StatsBase.skewnessMethod.
skewness(d::UnivariateDistribution)

Compute the skewness.

source
StatsBase.kurtosisMethod.
kurtosis(d::UnivariateDistribution)

Compute the excessive kurtosis.

source
isplatykurtic(d)

Return whether d is platykurtic (i.e kurtosis(d) > 0).

source
isleptokurtic(d)

Return whether d is leptokurtic (i.e kurtosis(d) < 0).

source
ismesokurtic(d)

Return whether d is mesokurtic (i.e kurtosis(d) == 0).

source
StatsBase.entropyMethod.
entropy(d::UnivariateDistribution)

Compute the entropy value of distribution d.

source
StatsBase.entropyMethod.
entropy(d::UnivariateDistribution, b::Real)

Compute the entropy value of distribution d, w.r.t. a given base.

source
Distributions.mgfMethod.
mgf(d::UnivariateDistribution, t)

Evaluate the moment generating function of distribution d.

source
Distributions.cfMethod.
cf(d::UnivariateDistribution, t)

Evaluate the characteristic function of distribution d.

source

Probability Evaluation

insupport(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.

source
Distributions.pdfMethod.
pdf(d::UnivariateDistribution, x::Real)

Evaluate the probability density (mass) at x.

Note: The package implements the following generic methods to evaluate pdf values in batch.

  • pdf!(dst::AbstractArray, d::Distribution, x::AbstractArray)

  • pdf(d::UnivariateDistribution, x::AbstractArray)

If there exists more efficient routine to evaluate pdf in batch (faster than repeatedly calling the scalar version of pdf), then one can also provide a specialized method of pdf!. The vectorized version of pdf simply delegats to pdf!.

source
logpdf(d::UnivariateDistribution, x::Real)

Evaluate the logarithm of probability density (mass) at x. Whereas there is a fallback implemented logpdf(d, x) = log(pdf(d, x)). Relying on this fallback is not recommended in general, as it is prone to overflow or underflow. Again, the package provides vectorized version of logpdf! and logpdf. One may override logpdf! to provide more efficient vectorized evaluation. Furthermore, the generic loglikelihood function delegates to _loglikelihood, which repeatedly calls logpdf. If there is a better way to compute log-likelihood, one should override _loglikelihood.

source
loglikelihood(d::UnivariateDistribution, X::AbstractArray)

The log-likelihood of distribution d w.r.t. all samples contained in array x.

source
Distributions.cdfMethod.
cdf(d::UnivariateDistribution, x::Real)

Evaluate the cumulative probability at x. The package provides generic functions to compute ccdf, logcdf, and logccdf in both scalar and vectorized forms. One may override these generic fallbacks if the specialized versions provide better numeric stability or higher efficiency.

source
logcdf(d::UnivariateDistribution, x::Real)

The logarithm of the cumulative function value(s) evaluated at x, i.e. log(cdf(x)).

source
Distributions.ccdfMethod.
ccdf(d::UnivariateDistribution, x::Real)

The complementary cumulative function evaluated at x, i.e. 1 - cdf(d, x).

source
logccdf(d::UnivariateDistribution, x::Real)

The logarithm of the complementary cumulative function values evaluated at x, i.e. log(ccdf(x)).

source
Base.quantileMethod.
quantile(d::UnivariateDistribution, q::Real)

Evaluate the inverse cumulative distribution function at q. The package provides generic functions to compute cquantile, invlogcdf, and invlogccdf in both scalar and vectorized forms. One may override these generic fallbacks if the specialized versions provide better numeric stability or higher efficiency. A generic median is provided, as median(d) = quantile(d, 0.5). However, one should implement a specialized version of median if it can be computed faster than $quantile$.

source
cquantile(d::UnivariateDistribution, q::Real)

The complementary quantile value, i.e. quantile(d, 1-q).

source
invlogcdf(d::UnivariateDistribution, lp::Real)

The inverse function of logcdf.

source
invlogcdf(d::UnivariateDistribution, lp::Real)

The inverse function of logcdf.

source

Vectorized evaluation

Vectorized computation and inplace vectorized computation are supported for the following functions:

For example, when x is an array, then r = pdf(d, x) returns an array r of the same size, such that r[i] = pdf(d, x[i]). One can also use pdf! to write results to pre-allocated storage, as pdf!(r, d, x).

Sampling (Random number generation)

Base.Random.randMethod.
rand(d::UnivariateDistribution)

Generate a scalar sample from d. The general fallback is quantile(d, rand()).

rand(d::UnivariateDistribution, n::Int) -> Vector

Generates a vector of n random scalar samples from d. The general fallback is to pick random samples from sampler(d).

source
Base.Random.rand!Method.
rand!(d::UnivariateDistribution, A::AbstractArray)

Populates the array A with scalar samples from d. The general fallback is to pick random samples from sampler(d).

source

Continuous Distributions

Arcsine(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

source
Beta(α,β)

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

source
BetaPrime(α,β)

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

source
Biweight(μ, σ)
source
Cauchy(μ, σ)

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

source
Chi(ν)

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

source
Chisq(ν)

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

source
Cosine(μ, σ)

A raised Cosine distribution.

External link:

source
Epanechnikov(μ, σ)
source
Erlang(α,θ)

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 b

External links

source
Exponential(θ)

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

source
FDist(ν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

source
Frechet(α,θ)

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

source
Gamma(α,θ)

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

source
GeneralizedExtremeValue(μ, σ, ξ)

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

source
GeneralizedPareto(μ, σ, ξ)

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

source
Gumbel(μ, θ)

The Gumbel 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

source
InverseGamma(α, θ)

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(a)       # Inverse Gamma distribution with shape a and unit scale, i.e. InverseGamma(a, 1)
InverseGamma(a, b)    # Inverse Gamma 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

source
InverseGaussian(μ,λ)

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

source
Kolmogorov()

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.

source
KSDist(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.

source
KSOneSided(n)

Distribution of the one-sided Kolmogorov-Smirnov test statistic:

\[D^+_n = \sup_x (\hat{F}_n(x) -F(x))\]
source
Laplace(μ,θ)

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(u)      # Laplace distribution with location u and unit scale, i.e. Laplace(u, 1)
Laplace(u, b)   # Laplace 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

source
Levy(μ, σ)

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

source
Logistic(μ,θ)

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

source
LogNormal(μ,σ)

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

source
NoncentralBeta(α, β, λ)
source
NoncentralChisq(ν, λ)

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

source
NoncentralF(ν1, ν2, λ)
source
NoncentralT(ν, λ)
source
Normal(μ,σ)

The Normal distribution with mean μ and standard deviation σ 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)\]
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

source
NormalCanon(η, λ)

Canonical Form of Normal distribution

source
NormalInverseGaussian(μ,α,β,δ)

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

source
Pareto(α,θ)

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

source
Rayleigh(σ)

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

source
SymTriangularDist(μ,σ)

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
source
TDist(ν)

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

Student's T distribution on Wikipedia

source
TriangularDist(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

source
Triweight(μ, σ)
source
Uniform(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

source
VonMises(μ, κ)

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

source
Weibull(α,θ)

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

source

Discrete Distributions

Bernoulli(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:

source
BetaBinomial(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:

source
Binomial(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:

source
Categorical(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. External links:

source
DiscreteUniform(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

source
Geometric(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

source
Hypergeometric(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

source
NegativeBinomial(r,p)

A Negative binomial distribution describes the number of failures before the rth 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:

source
Poisson(λ)

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:

source
PoissonBinomial(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 ith 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:

source
Skellam(μ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:

source