SciPy

scipy.stats.exponweib

scipy.stats.exponweib(*args, **kwds) = <scipy.stats._continuous_distns.exponweib_gen object>[source]

An exponentiated Weibull continuous random variable.

As an instance of the rv_continuous class, exponweib object inherits from it a collection of generic methods (see below for the full list), and completes them with details specific for this particular distribution.

Notes

The probability density function for exponweib is:

\[f(x, a, c) = a c [1-\exp(-x^c)]^{a-1} \exp(-x^c) x^{c-1}\]

and its cumulative distribution function is:

\[F(x, a, c) = [1-\exp(-x^c)]^a\]

for \(x > 0\), \(a > 0\), \(c > 0\).

exponweib takes \(a\) and \(c\) as shape parameters:

  • \(a\) is the exponentiation parameter, with the special case \(a=1\) corresponding to the (non-exponentiated) Weibull distribution weibull_min.

  • \(c\) is the shape parameter of the non-exponentiated Weibull law.

The probability density above is defined in the “standardized” form. To shift and/or scale the distribution use the loc and scale parameters. Specifically, exponweib.pdf(x, a, c, loc, scale) is identically equivalent to exponweib.pdf(y, a, c) / scale with y = (x - loc) / scale.

References

https://en.wikipedia.org/wiki/Exponentiated_Weibull_distribution

Examples

>>> from scipy.stats import exponweib
>>> import matplotlib.pyplot as plt
>>> fig, ax = plt.subplots(1, 1)

Calculate a few first moments:

>>> a, c = 2.89, 1.95
>>> mean, var, skew, kurt = exponweib.stats(a, c, moments='mvsk')

Display the probability density function (pdf):

>>> x = np.linspace(exponweib.ppf(0.01, a, c),
...                 exponweib.ppf(0.99, a, c), 100)
>>> ax.plot(x, exponweib.pdf(x, a, c),
...        'r-', lw=5, alpha=0.6, label='exponweib pdf')

Alternatively, the distribution object can be called (as a function) to fix the shape, location and scale parameters. This returns a “frozen” RV object holding the given parameters fixed.

Freeze the distribution and display the frozen pdf:

>>> rv = exponweib(a, c)
>>> ax.plot(x, rv.pdf(x), 'k-', lw=2, label='frozen pdf')

Check accuracy of cdf and ppf:

>>> vals = exponweib.ppf([0.001, 0.5, 0.999], a, c)
>>> np.allclose([0.001, 0.5, 0.999], exponweib.cdf(vals, a, c))
True

Generate random numbers:

>>> r = exponweib.rvs(a, c, size=1000)

And compare the histogram:

>>> ax.hist(r, density=True, histtype='stepfilled', alpha=0.2)
>>> ax.legend(loc='best', frameon=False)
>>> plt.show()
../_images/scipy-stats-exponweib-1.png

Methods

rvs(a, c, loc=0, scale=1, size=1, random_state=None)

Random variates.

pdf(x, a, c, loc=0, scale=1)

Probability density function.

logpdf(x, a, c, loc=0, scale=1)

Log of the probability density function.

cdf(x, a, c, loc=0, scale=1)

Cumulative distribution function.

logcdf(x, a, c, loc=0, scale=1)

Log of the cumulative distribution function.

sf(x, a, c, loc=0, scale=1)

Survival function (also defined as 1 - cdf, but sf is sometimes more accurate).

logsf(x, a, c, loc=0, scale=1)

Log of the survival function.

ppf(q, a, c, loc=0, scale=1)

Percent point function (inverse of cdf — percentiles).

isf(q, a, c, loc=0, scale=1)

Inverse survival function (inverse of sf).

moment(n, a, c, loc=0, scale=1)

Non-central moment of order n

stats(a, c, loc=0, scale=1, moments=’mv’)

Mean(‘m’), variance(‘v’), skew(‘s’), and/or kurtosis(‘k’).

entropy(a, c, loc=0, scale=1)

(Differential) entropy of the RV.

fit(data, a, c, loc=0, scale=1)

Parameter estimates for generic data.

expect(func, args=(a, c), loc=0, scale=1, lb=None, ub=None, conditional=False, **kwds)

Expected value of a function (of one argument) with respect to the distribution.

median(a, c, loc=0, scale=1)

Median of the distribution.

mean(a, c, loc=0, scale=1)

Mean of the distribution.

var(a, c, loc=0, scale=1)

Variance of the distribution.

std(a, c, loc=0, scale=1)

Standard deviation of the distribution.

interval(alpha, a, c, loc=0, scale=1)

Endpoints of the range that contains alpha percent of the distribution

Previous topic

scipy.stats.exponnorm

Next topic

scipy.stats.exponpow