scipy.stats.levy_stable#

scipy.stats.levy_stable = <scipy.stats._levy_stable.levy_stable_gen object>[source]#

A Levy-stable continuous random variable.

As an instance of the rv_continuous class, levy_stable 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.

See also

levy, levy_l, cauchy, norm

Notes

The distribution for levy_stable has characteristic function:

\[\varphi(t, \alpha, \beta, c, \mu) = e^{it\mu -|ct|^{\alpha}(1-i\beta\operatorname{sign}(t)\Phi(\alpha, t))}\]

where two different parameterizations are supported. The first \(S_1\):

\[\begin{split}\Phi = \begin{cases} \tan \left({\frac {\pi \alpha }{2}}\right)&\alpha \neq 1\\ -{\frac {2}{\pi }}\log |t|&\alpha =1 \end{cases}\end{split}\]

The second \(S_0\):

\[\begin{split}\Phi = \begin{cases} -\tan \left({\frac {\pi \alpha }{2}}\right)(|ct|^{1-\alpha}-1) &\alpha \neq 1\\ -{\frac {2}{\pi }}\log |ct|&\alpha =1 \end{cases}\end{split}\]

The probability density function for levy_stable is:

\[f(x) = \frac{1}{2\pi}\int_{-\infty}^\infty \varphi(t)e^{-ixt}\,dt\]

where \(-\infty < t < \infty\). This integral does not have a known closed form.

levy_stable generalizes several distributions. Where possible, they should be used instead. Specifically, when the shape parameters assume the values in the table below, the corresponding equivalent distribution should be used.

alpha

beta

Equivalent

1/2

-1

levy_l

1/2

1

levy

1

0

cauchy

2

any

norm (with scale=sqrt(2))

Evaluation of the pdf uses Nolan’s piecewise integration approach with the Zolotarev \(M\) parameterization by default. There is also the option to use direct numerical integration of the standard parameterization of the characteristic function or to evaluate by taking the FFT of the characteristic function.

The default method can changed by setting the class variable levy_stable.pdf_default_method to one of ‘piecewise’ for Nolan’s approach, ‘dni’ for direct numerical integration, or ‘fft-simpson’ for the FFT based approach. For the sake of backwards compatibility, the methods ‘best’ and ‘zolotarev’ are equivalent to ‘piecewise’ and the method ‘quadrature’ is equivalent to ‘dni’.

The parameterization can be changed by setting the class variable levy_stable.parameterization to either ‘S0’ or ‘S1’. The default is ‘S1’.

To improve performance of piecewise and direct numerical integration one can specify levy_stable.quad_eps (defaults to 1.2e-14). This is used as both the absolute and relative quadrature tolerance for direct numerical integration and as the relative quadrature tolerance for the piecewise method. One can also specify levy_stable.piecewise_x_tol_near_zeta (defaults to 0.005) for how close x is to zeta before it is considered the same as x [NO]. The exact check is abs(x0 - zeta) < piecewise_x_tol_near_zeta*alpha**(1/alpha). One can also specify levy_stable.piecewise_alpha_tol_near_one (defaults to 0.005) for how close alpha is to 1 before being considered equal to 1.

To increase accuracy of FFT calculation one can specify levy_stable.pdf_fft_grid_spacing (defaults to 0.001) and pdf_fft_n_points_two_power (defaults to None which means a value is calculated that sufficiently covers the input range).

Further control over FFT calculation is available by setting pdf_fft_interpolation_degree (defaults to 3) for spline order and pdf_fft_interpolation_level for determining the number of points to use in the Newton-Cotes formula when approximating the characteristic function (considered experimental).

Evaluation of the cdf uses Nolan’s piecewise integration approach with the Zolatarev \(S_0\) parameterization by default. There is also the option to evaluate through integration of an interpolated spline of the pdf calculated by means of the FFT method. The settings affecting FFT calculation are the same as for pdf calculation. The default cdf method can be changed by setting levy_stable.cdf_default_method to either ‘piecewise’ or ‘fft-simpson’. For cdf calculations the Zolatarev method is superior in accuracy, so FFT is disabled by default.

Fitting estimate uses quantile estimation method in [MC]. MLE estimation of parameters in fit method uses this quantile estimate initially. Note that MLE doesn’t always converge if using FFT for pdf calculations; this will be the case if alpha <= 1 where the FFT approach doesn’t give good approximations.

Any non-missing value for the attribute levy_stable.pdf_fft_min_points_threshold will set levy_stable.pdf_default_method to ‘fft-simpson’ if a valid default method is not otherwise set.

Warning

For pdf calculations FFT calculation is considered experimental.

For cdf calculations FFT calculation is considered experimental. Use Zolatarev’s method instead (default).

The probability density above is defined in the “standardized” form. To shift and/or scale the distribution use the loc and scale parameters. Generally levy_stable.pdf(x, alpha, beta, loc, scale) is identically equivalent to levy_stable.pdf(y, alpha, beta) / scale with y = (x - loc) / scale, except in the S1 parameterization if alpha == 1. In that case levy_stable.pdf(x, alpha, beta, loc, scale) is identically equivalent to levy_stable.pdf(y, alpha, beta) / scale with y = (x - loc - 2 * beta * scale * np.log(scale) / np.pi) / scale. See [NO2] Definition 1.8 for more information. Note that shifting the location of a distribution does not make it a “noncentral” distribution.

References

[MC]

McCulloch, J., 1986. Simple consistent estimators of stable distribution parameters. Communications in Statistics - Simulation and Computation 15, 11091136.

[WZ]

Wang, Li and Zhang, Ji-Hong, 2008. Simpson’s rule based FFT method to compute densities of stable distribution.

[NO]

Nolan, J., 1997. Numerical Calculation of Stable Densities and distributions Functions.

[NO2]

Nolan, J., 2018. Stable Distributions: Models for Heavy Tailed Data.

[HO]

Hopcraft, K. I., Jakeman, E., Tanner, R. M. J., 1999. Lévy random walks with fluctuating step number and multiscale behavior.

Examples

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

Calculate the first four moments:

>>> alpha, beta = 1.8, -0.5
>>> mean, var, skew, kurt = levy_stable.stats(alpha, beta, moments='mvsk')

Display the probability density function (pdf):

>>> x = np.linspace(levy_stable.ppf(0.01, alpha, beta),
...                 levy_stable.ppf(0.99, alpha, beta), 100)
>>> ax.plot(x, levy_stable.pdf(x, alpha, beta),
...        'r-', lw=5, alpha=0.6, label='levy_stable 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 = levy_stable(alpha, beta)
>>> ax.plot(x, rv.pdf(x), 'k-', lw=2, label='frozen pdf')

Check accuracy of cdf and ppf:

>>> vals = levy_stable.ppf([0.001, 0.5, 0.999], alpha, beta)
>>> np.allclose([0.001, 0.5, 0.999], levy_stable.cdf(vals, alpha, beta))
True

Generate random numbers:

>>> r = levy_stable.rvs(alpha, beta, size=1000)

And compare the histogram:

>>> ax.hist(r, density=True, bins='auto', histtype='stepfilled', alpha=0.2)
>>> ax.set_xlim([x[0], x[-1]])
>>> ax.legend(loc='best', frameon=False)
>>> plt.show()
../../_images/scipy-stats-levy_stable-1.png

Methods

rvs(alpha, beta, loc=0, scale=1, size=1, random_state=None)

Random variates.

pdf(x, alpha, beta, loc=0, scale=1)

Probability density function.

logpdf(x, alpha, beta, loc=0, scale=1)

Log of the probability density function.

cdf(x, alpha, beta, loc=0, scale=1)

Cumulative distribution function.

logcdf(x, alpha, beta, loc=0, scale=1)

Log of the cumulative distribution function.

sf(x, alpha, beta, loc=0, scale=1)

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

logsf(x, alpha, beta, loc=0, scale=1)

Log of the survival function.

ppf(q, alpha, beta, loc=0, scale=1)

Percent point function (inverse of cdf — percentiles).

isf(q, alpha, beta, loc=0, scale=1)

Inverse survival function (inverse of sf).

moment(order, alpha, beta, loc=0, scale=1)

Non-central moment of the specified order.

stats(alpha, beta, loc=0, scale=1, moments=’mv’)

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

entropy(alpha, beta, loc=0, scale=1)

(Differential) entropy of the RV.

fit(data)

Parameter estimates for generic data. See scipy.stats.rv_continuous.fit for detailed documentation of the keyword arguments.

expect(func, args=(alpha, beta), 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(alpha, beta, loc=0, scale=1)

Median of the distribution.

mean(alpha, beta, loc=0, scale=1)

Mean of the distribution.

var(alpha, beta, loc=0, scale=1)

Variance of the distribution.

std(alpha, beta, loc=0, scale=1)

Standard deviation of the distribution.

interval(confidence, alpha, beta, loc=0, scale=1)

Confidence interval with equal areas around the median.