Skip to content

Probability Distributions

scipy.stats integration

l_rv_generic

Additional methods that are patched into scipy.stats.rv_continuous and scipy.stats.rv_discrete.

l_moment(r, /, *args, trim=0, quad_opts=None, **kwds)

l_moment(
    r: lmt.AnyOrderND,
    /,
    *args: Any,
    trim: lmt.AnyTrim = ...,
    quad_opts: lspt.QuadOptions | None = ...,
    **kwds: Any,
) -> _ArrF8
l_moment(
    r: lmt.AnyOrder,
    /,
    *args: Any,
    trim: lmt.AnyTrim = ...,
    quad_opts: lspt.QuadOptions | None = ...,
    **kwds: Any,
) -> np.float64

Population L-moment(s) \(\lambda^{(s,t)}_r\).

\[ \lambda^{(s, t)}_r = \frac{r+s+t}{r} \frac{B(r,\,r+s+t)}{B(r+s,\,r+t)} \mathbb{E}_X \left[ U^s \left(1 - U\right)^t \,\tilde{P}^{(t, s)}_{r-1}(U) \,X \right] \;, \]

with \(U = F_X(X)\) the rank of \(X\), and \(\tilde{P}^{(a,b)}_n(x)\) the shifted (\(x \mapsto 2x-1\)) Jacobi polynomial.

Examples:

Evaluate the population L-moments of the normally-distributed IQ test:

>>> import lmo
>>> from scipy.stats import norm
>>> norm(100, 15).l_moment([1, 2, 3, 4]).round(6)
array([100.      ,   8.462844,   0.      ,   1.037559])
>>> _[1] * np.sqrt(np.pi)
15.0000004

Discrete distributions are also supported, e.g. the Binomial distribution:

>>> from scipy.stats import binom
>>> binom(10, 0.6).l_moment([1, 2, 3, 4]).round(6)
array([ 6.      ,  0.862238, -0.019729,  0.096461])

Parameters:

Name Type Description Default
r AnyOrder | AnyOrderND

L-moment order(s), non-negative integer or array-like of integers.

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

0
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
**kwds Any

Additional keyword arguments to pass to the distribution.

{}

Raises:

Type Description
TypeError

r is not integer-valued

ValueError

r is empty or negative

Returns:

Name Type Description
lmbda float64 | _ArrF8

The population L-moment(s), a scalar or float array like r.

References
See Also

l_ratio(r, k, /, *args, trim=0, quad_opts=None, **kwds)

l_ratio(
    order: lmt.AnyOrderND,
    order_denom: lmt.AnyOrder | lmt.AnyOrderND,
    /,
    *args: Any,
    trim: lmt.AnyTrim = ...,
    quad_opts: lspt.QuadOptions | None = ...,
    **kwds: Any,
) -> _ArrF8
l_ratio(
    order: lmt.AnyOrder,
    order_denom: lmt.AnyOrder | lmt.AnyOrderND,
    /,
    *args: Any,
    trim: lmt.AnyTrim = ...,
    quad_opts: lspt.QuadOptions | None = ...,
    **kwds: Any,
) -> np.float64

L-moment ratio(’s) \(\tau^{(s,t)}_{r,k}\).

\[ \tau^{(s,t)}_{r,k} = \frac{\lambda^{(s,t)}_r}{\lambda^{(s,t)}_k} \]

Unless explicitly specified, the r-th (\(r>2\)) L-ratio, \(\tau^{(s,t)}_r\) refers to \(\tau^{(s,t)}_{r, 2}\). Another special case is the L-variation, or the L-CV, \(\tau^{(s,t)} = \tau^{(s,t)}_{2, 1}\). This is the L-moment analogue of the coefficient of variation.

Examples:

Evaluate the population L-CV and LL-CV (CV = coefficient of variation) of the standard Rayleigh distribution.

>>> import lmo
>>> from scipy.stats import distributions
>>> X = distributions.rayleigh()
>>> X.std() / X.mean()  # legacy CV
0.5227232
>>> X.l_ratio(2, 1)
0.2928932
>>> X.l_ratio(2, 1, trim=(0, 1))
0.2752551

And similarly, for the (discrete) Poisson distribution with rate parameter set to 2, the L-CF and LL-CV evaluate to:

>>> X = distributions.poisson(2)
>>> X.std() / X.mean()
0.7071067
>>> X.l_ratio(2, 1)
0.3857527
>>> X.l_ratio(2, 1, trim=(0, 1))
0.4097538

Note that (untrimmed) L-CV requires a higher (subdivision) limit in the integration routine, otherwise it’ll complain that it didn’t converge (enough) yet. This is because it’s effectively integrating a non-smooth function, which is mathematically iffy, but works fine in this numerical application.

Parameters:

Name Type Description Default
r AnyOrder | AnyOrderND

L-moment ratio order(s), non-negative integer or array-like of integers.

required
k AnyOrder | AnyOrderND

L-moment order of the denominator, e.g. 2.

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

0
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
**kwds Any

Additional keyword arguments to pass to the distribution.

{}
See Also

l_stats(*args, trim=0, moments=4, quad_opts=None, **kwds)

The L-moments (for \(r \le 2\)) and L-ratio’s (for \(r > 2\)).

By default, the first moments = 4 population L-stats are calculated:

  • \(\lambda^{(s,t)}_1\) - L-location
  • \(\lambda^{(s,t)}_2\) - L-scale
  • \(\tau^{(s,t)}_3\) - L-skewness coefficient
  • \(\tau^{(s,t)}_4\) - L-kurtosis coefficient

This method is equivalent to X.l_ratio([1, 2, 3, 4], [0, 0, 2, 2], *, **), for with default moments = 4.

Examples:

Summarize the standard exponential distribution for different trim-orders.

>>> import lmo
>>> from scipy.stats import distributions
>>> X = distributions.expon()
>>> X.l_stats().round(6)
array([1.      , 0.5     , 0.333333, 0.166667])
>>> X.l_stats(trim=(0, 1 / 2)).round(6)
array([0.666667, 0.333333, 0.266667, 0.114286])
>>> X.l_stats(trim=(0, 1)).round(6)
array([0.5     , 0.25    , 0.222222, 0.083333])
Note

This should not be confused with the term L-statistic, which is sometimes used to describe any linear combination of order statistics.

Parameters:

Name Type Description Default
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

0
moments int

The amount of L-moments to return. Defaults to 4.

4
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
**kwds Any

Additional keyword arguments to pass to the distribution.

{}
See Also

l_loc(*args, trim=0, **kwds)

L-location of the distribution, i.e. the 1st L-moment.

Alias for X.l_moment(1, ...).

l_scale(*args, trim=0, **kwds)

L-scale of the distribution, i.e. the 2nd L-moment.

Alias for X.l_moment(2, ...).

l_skew(*args, trim=0, **kwds)

L-skewness coefficient of the distribution; the 3rd L-moment ratio.

Alias for X.l_ratio(3, 2, ...).

l_kurt(*args, trim=0, **kwds)

L-kurtosis coefficient of the distribution; the 4th L-moment ratio.

Alias for X.l_ratio(4, 2, ...).

l_moments_cov(r_max, /, *args, trim=0, quad_opts=None, **kwds)

Variance/covariance matrix of the L-moment estimators.

L-moments that are estimated from \(n\) samples of a distribution with CDF \(F\), converge to the multivariate normal distribution as the sample size \(n \rightarrow \infty\).

\[ \sqrt{n} \left( \vec{l}^{(s, t)} - \vec{\lambda}^{(s, t)} \right) \sim \mathcal{N}( \vec{0}, \mathbf{\Lambda}^{(s, t)} ) \]

Here, \(\vec{l}^{(s, t)} = \left[l^{(s,t)}_r, \dots, l^{(s,t)}_{r_{max}}\right]^T\) is a vector of estimated sample L-moments, and \(\vec{\lambda}^{(s, t)}\) its theoretical (“true”) counterpart.

This function calculates the covariance matrix

\[ \begin{align} \bf{\Lambda}^{(s,t)}_{k, r} &= \mathrm{Cov}[l^{(s, t)}_k, l^{(s, t)}_r] \\ &= c_k c_r \iint\limits_{x < y} \Big[ p_k\big(F(x)\big) \, p_r\big(F(y)\big) + p_r\big(F(x)\big) \, p_k\big(F(y)\big) \Big] w^{(s+1,\, t)}\big(F(x)\big) \, w^{(s,\, t+1)}\big(F(y)\big) \, \mathrm{d}x \, \mathrm{d}y \;, \end{align} \]

where

\[ c_n = \frac{\Gamma(n) \Gamma(n+s+t+1)}{n \Gamma(n+s) \Gamma(n+t)}\;, \]

the shifted Jacobi polynomial \(p_n(u) = P^{(t, s)}_{n-1}(2u - 1)\), \(P^{(t, s)}_m\), and \(w^{(s,t)}(u) = u^s (1-u)^t\) its weight function.

Notes

This function is not vectorized or parallelized.

For small sample sizes (\(n < 100\)), the covariances of the higher-order L-moments (\(r > 2\)) can be biased. But this bias quickly disappears at roughly \(n > 200\) (depending on the trim- and L-moment orders).

Examples:

>>> import lmo
>>> from scipy.stats import distributions
>>> X = distributions.expon()  # standard exponential distribution
>>> X.l_moments_cov(4).round(6)
array([[1.      , 0.5     , 0.166667, 0.083333],
    [0.5     , 0.333333, 0.166667, 0.083333],
    [0.166667, 0.166667, 0.133333, 0.083333],
    [0.083333, 0.083333, 0.083333, 0.071429]])
>>> X.l_moments_cov(4, trim=(0, 1)).round(6)
array([[0.333333, 0.125   , 0.      , 0.      ],
    [0.125   , 0.075   , 0.016667, 0.      ],
    [0.      , 0.016667, 0.016931, 0.00496 ],
    [0.      , 0.      , 0.00496 , 0.0062  ]])

Parameters:

Name Type Description Default
r_max int

The amount of L-moment orders to consider. If for example r_max = 4, the covariance matrix will be of shape (4, 4), and the columns and rows correspond to the L-moments of order \(r = 1, \dots, r_{max}\).

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float. or floats.

0
quad_opts NQuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
**kwds Any

Additional keyword arguments to pass to the distribution.

{}

Returns:

Name Type Description
cov _ArrF8

Covariance matrix, with shape (r_max, r_max).

Raises:

Type Description
RuntimeError

If the covariance matrix is invalid.

References

l_stats_cov(*args, moments=4, trim=0, quad_opts=None, **kwds)

Similar to l_moments_cov , but for the l_rv_generic.l_stats.

As the sample size \(n \rightarrow \infty\), the L-moment ratio’s are also distributed (multivariate) normally. The L-stats are defined to be L-moments for \(r\le 2\), and L-ratio coefficients otherwise.

The corresponding covariance matrix has been found to be

\[ \bf{T}^{(s, t)}_{k, r} = \begin{cases} \bf{\Lambda}^{(s, t)}_{k, r} & k \le 2 \wedge r \le 2 \\ \frac{ \bf{\Lambda}^{(s, t)}_{k, r} - \tau_r \bf{\Lambda}^{(s, t)}_{k, 2} }{ \lambda^{(s,t)}_{2} } & k \le 2 \wedge r > 2 \\ \frac{ \bf{\Lambda}^{(s, t)}_{k, r} - \tau_k \bf{\Lambda}^{(s, t)}_{2, r} - \tau_r \bf{\Lambda}^{(s, t)}_{k, 2} + \tau_k \tau_r \bf{\Lambda}^{(s, t)}_{2, 2} }{ \Big( \lambda^{(s,t)}_{2} \Big)^2 } & k > 2 \wedge r > 2 \end{cases} \]

where \(\bf{\Lambda}^{(s, t)}\) is the covariance matrix of the L-moments from l_moment_cov_from_cdf, and \(\tau^{(s,t)}_r = \lambda^{(s,t)}_r / \lambda^{(s,t)}_2\) the population L-ratio.

Examples:

Evaluate the LL-stats covariance matrix of the standard exponential distribution, for 0, 1, and 2 degrees of trimming.

>>> import lmo
>>> from scipy.stats import distributions
>>> X = distributions.expon()  # standard exponential distribution
>>> X.l_stats_cov().round(6)
array([[1.      , 0.5     , 0.      , 0.      ],
    [0.5     , 0.333333, 0.111111, 0.055556],
    [0.      , 0.111111, 0.237037, 0.185185],
    [0.      , 0.055556, 0.185185, 0.21164 ]])
>>> X.l_stats_cov(trim=(0, 1)).round(6)
array([[ 0.333333,  0.125   , -0.111111, -0.041667],
    [ 0.125   ,  0.075   ,  0.      , -0.025   ],
    [-0.111111,  0.      ,  0.21164 ,  0.079365],
    [-0.041667, -0.025   ,  0.079365,  0.10754 ]])
>>> X.l_stats_cov(trim=(0, 2)).round(6)
array([[ 0.2     ,  0.066667, -0.114286, -0.02    ],
    [ 0.066667,  0.038095, -0.014286, -0.023333],
    [-0.114286, -0.014286,  0.228571,  0.04    ],
    [-0.02    , -0.023333,  0.04    ,  0.086545]])

Note that with 0 trim the L-location is independent of the L-skewness and L-kurtosis. With 1 trim, the L-scale and L-skewness are independent. And with 2 trim, all L-stats depend on each other.

Parameters:

Name Type Description Default
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
moments int

The amount of L-statistics to consider. Defaults to 4.

4
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float. or floats.

0
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
**kwds Any

Additional keyword arguments to pass to the distribution.

{}
References

l_moment_influence(r, /, *args, trim=0, quad_opts=None, tol=1e-08, **kwds)

Returns the influence function (IF) of an L-moment.

\[ \psi_{\lambda^{(s, t)}_r | F}(x) = c^{(s,t)}_r \, F(x)^s \, \big( 1-{F}(x) \big)^t \, \tilde{P}^{(s,t)}_{r-1} \big( F(x) \big) \, x - \lambda^{(s,t)}_r \;, \]

with \(F\) the CDF, \(\tilde{P}^{(s,t)}_{r-1}\) the shifted Jacobi polynomial, and

\[ c^{(s,t)}_r = \frac{r+s+t}{r} \frac{B(r, \, r+s+t)}{B(r+s, \, r+t)} \;, \]

where \(B\) is the (complete) Beta function.

The proof is trivial, because population L-moments are linear functionals.

Notes

The order parameter r is not vectorized.

Parameters:

Name Type Description Default
r AnyOrder

The L-moment order \(r \in \mathbb{N}^+\)..

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float. or floats.

0
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
tol float

Values that are absolutely smaller than this will be rounded to zero.

1e-08
**kwds Any

Additional keyword arguments to pass to the distribution.

{}

Returns:

Name Type Description
influence_function _Fn1

The (vectorized) influence function, \(\psi_{\lambda^{(s, t)}_r | F}(x)\).

See Also
References

l_ratio_influence(r, k, /, *args, trim=0, quad_opts=None, tol=1e-08, **kwds)

Returns the influence function (IF) of an L-moment ratio.

\[ \psi_{\tau^{(s, t)}_{r,k}|F}(x) = \frac{ \psi_{\lambda^{(s, t)}_r|F}(x) - \tau^{(s, t)}_{r,k} \, \psi_{\lambda^{(s, t)}_k|F}(x) }{ \lambda^{(s,t)}_k } \;, \]

where the L-moment ratio is defined as

\[ \tau^{(s, t)}_{r,k} = \frac{ \lambda^{(s, t)}_r }{ \lambda^{(s, t)}_k } \;. \]

Because IF’s are a special case of the general Gâteuax derivative, the L-ratio IF is derived by applying the chain rule to the L-moment IF.

Parameters:

Name Type Description Default
r AnyOrder

L-moment ratio order, i.e. the order of the numerator L-moment.

required
k AnyOrder

Denominator L-moment order, defaults to 2.

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information)

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float. or floats.

0
quad_opts QuadOptions | None

Optional dict of options to pass to scipy.integrate.quad.

None
tol float

Values that are absolutely smaller than this will be rounded to zero.

1e-08
**kwds Any

Additional keyword arguments to pass to the distribution.

{}

Returns:

Name Type Description
influence_function Callable[[_T_x], _T_x]

The influence function, with vectorized signature () -> ().

See Also
References

l_fit(data, *args, n_extra=0, trim=0, full_output=False, fit_kwargs=None, random_state=None, **kwds)

l_fit(
    data: lnpt.AnyVectorFloat,
    *args: float,
    n_extra: int = ...,
    trim: lmt.AnyTrimInt = ...,
    full_output: Literal[True],
    fit_kwargs: Mapping[str, Any] | None = ...,
    **kwds: Any
) -> tuple[float, ...]
l_fit(
    data: lnpt.AnyVectorFloat,
    *args: float,
    n_extra: int = ...,
    trim: lmt.AnyTrimInt = ...,
    full_output: bool = ...,
    fit_kwargs: Mapping[str, Any] | None = ...,
    **kwds: Any
) -> tuple[float, ...]

Return estimates of shape (if applicable), location, and scale parameters from data. The default estimation method is Method of L-moments (L-MM), but the Generalized Method of L-Moments (L-GMM) is also available (see the n_extra parameter).

See ‘lmo.inference.fit’ for details.

Examples:

Fitting standard normal samples Using scipy’s default MLE (Maximum Likelihood Estimation) method:

>>> import lmo
>>> import numpy as np
>>> from scipy.stats import norm
>>> rng = np.random.default_rng(12345)
>>> x = rng.standard_normal(200)
>>> norm.fit(x)
(0.0033254, 0.95554)

Better results can be obtained different by using Lmo’s L-MM (Method of L-moment):

>>> norm.l_fit(x, random_state=rng)
(0.0033145, 0.96179)
>>> norm.l_fit(x, trim=1, random_state=rng)
(0.0196385, 0.96861)

To use more L-moments than the number of parameters, two in this case, n_extra can be used. This will use the L-GMM (Generalized Method of L-Moments), which results in slightly better estimates:

>>> norm.l_fit(x, n_extra=1, random_state=rng)
(0.0039747, .96233)
>>> norm.l_fit(x, trim=1, n_extra=1, random_state=rng)
(-0.00127874, 0.968547)

Parameters:

Name Type Description Default
data AnyVectorFloat

1-D array-like data to use in estimating the distribution parameters.

required
*args float

Starting value(s) for any shape-characterizing arguments ( those not provided will be determined by a call to fit(data)).

()
trim AnyTrimInt

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

0
n_extra int

The amount of extra L-moment conditions to use than the amount of parameters. If 0 (default), L-MM will be used. If >0, \(k\)-step L-GMM will be used.

0
full_output bool

If set to True, a LGMMResult instance will be returned, instead of only a tuple with parameters.

False
fit_kwargs Mapping[str, Any] | None

Additional keyword arguments to be passed to ‘lmo.inference.fit’ or ‘scipy.optimize.minimize’.

None
random_state int | Generator | None

Integer or numpy.random.Generator instance, used for Monte-Carlo simulation when n_extra > 0. If None (default), the random_state of this distribution will be used.

None
**kwds Any

Special keyword arguments are recognized as holding certain parameters fixed:

- `f0...fn`: hold respective shape parameters fixed.
Alternatively, shape parameters to fix can be specified by
name. For example, if `self.shapes == "a, b"`, `fa` and
`fix_a` are equivalent to `f0`, and `fb` and `fix_b` are
equivalent to `f1`.
- `floc`: hold location parameter fixed to specified value.
- `fscale`: hold scale parameter fixed to specified value.
{}

Returns:

Name Type Description
result tuple[float, ...] | GMMResult

Named tuple with estimates for any shape parameters (if applicable), followed by those for location and scale. For most random variables, shape statistics will be returned, but there are exceptions (e.g. norm). If full_output=True, an instance of LGMMResult will be returned instead.

See Also
References
Todo
  • Support integral parameters.

l_fit_loc_scale(data, *args, trim=0, **kwds)

Estimate loc and scale parameters from data using the first two L-moments.

Notes

The implementation mimics that of fit_loc_scale()

Parameters:

Name Type Description Default
data AnyArrayFloat

Data to fit.

required
*args Any

The shape parameter(s) for the distribution (see docstring of the instance object for more information).

()
trim AnyTrim

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

0
**kwds Any

Additional keyword arguments to pass to the distribution.

{}

Returns:

Name Type Description
loc_hat float

Estimated location parameter for the data.

scale_hat float

Estimated scale parameter for the data.

Parametric

lmo.distributions.kumaraswamy = kumaraswamy_gen(name='kumaraswamy', a=0.0, b=1.0)

A Kumaraswamy random variable, similar to scipy.stats.beta.

The probability density function for kumaraswamy is:

\[ f(x, a, b) = a x^{a - 1} b \left(1 - x^a\right)^{b - 1} \]

for \( 0 < x < 1,\ a > 0,\ b > 0 \).

kumaraswamy takes \( a \) and \( b \) as shape parameters.

See Also

lmo.distributions.wakeby = wakeby_gen(name='wakeby', a=0.0)

A Wakeby random variable, a generalization of scipy.stats.genpareto.

wakeby takes \( b \), \( d \) and \( f \) as shape parameters.

For a detailed description of the Wakeby distribution, refer to Distributions - Wakeby.

lmo.distributions.genlambda = genlambda_gen(name='genlambda')

A generalized Tukey-Lambda random variable.

genlambda takes b, d and f as shape parameters. b and d can be any float, and f requires -1 <= f <= 1.

If f == 0 and b == d, genlambda is equivalent to scipy.stats.tukeylambda, with b (or d) as shape parameter.

For a detailed description of the GLD, refer to Distributions - GLD.

Nonparametric

lmo.distributions.l_poly(lmbda, /, trim=0, *, seed=None)

Polynomial quantile distribution with (only) the given L-moments.

Parameters:

Name Type Description Default
lmbda AnyVectorFloat

1-d array-like of L-moments \( \tlmoment{s,t}{r} \) for \( r = 1, 2, \ldots, R \). At least 2 L-moments are required. All remaining L-moments with \( r > R \) are considered zero.

required
trim AnyTrim

The trim-length(s) of L-moments lmbda.

0
seed Seed | None

Random number generator.

None

lmo.distributions.l_poly.fit(data, /, moments=None, trim=0)

Fit distribution using the (trimmed) L-moment estimates of the given data.

Parameters:

Name Type Description Default
data _AnyRealND

1-d array-like with sample observations.

required
moments int | None

How many sample L-moments to use, 2 <= moments < len(data). Defaults to \(\sqrt[3]{n}\), where \(n\) is len(data).

None
trim AnyTrim

The left and right trim-lengths \((s, t)\) to use. Defaults to \((0, 0)\).

0

Returns:

Type Description
Self

A fitted l_poly instance.

Raises:

Type Description
TypeError

Invalid data shape.

ValueError

Not enough moments.

ValueError

If the L-moments of the data do not result in strictly monotinically increasing quantile function (PPF).

This generally means that either the left, the right, or both trim-orders are too small.

lmo.distributions.l_poly.rvs(size=None, random_state=None)

rvs(
    size: Literal[1] | None = None, random_state: lnpt.Seed | None = None
) -> np.float64
rvs(
    size: int | tuple[int, ...], random_state: lnpt.Seed | None = None
) -> npt.NDArray[np.float64]

Draw random variates from the relevant distribution.

Parameters:

Name Type Description Default
size int | tuple[int, ...] | None

Defining number of random variates, defaults to 1.

None
random_state Seed | None

Seed or numpy.random.Generator instance. Defaults to l_poly.random_state.

None

Returns:

Type Description
float64 | NDArray[float64]

A scalar or array with shape like size.

lmo.distributions.l_poly.ppf(p)

Percent point function \( Q(p) \) (inverse of CDF, a.k.a. the quantile function) at \( p \) of the given distribution.

Parameters:

Name Type Description Default
p _T_x

Scalar or array-like of lower tail probability values in \( [0, 1] \).

required
See Also

lmo.distributions.l_poly.isf(q)

Inverse survival function \( \bar{Q}(q) = Q(1 - q) \) (inverse of sf) at \( q \).

Parameters:

Name Type Description Default
q _T_x

Scalar or array-like of upper tail probability values in \( [0, 1] \).

required

lmo.distributions.l_poly.qdf(p)

Quantile density function \( q \equiv \frac{\dd{Q}}{\dd{p}} \) ( derivative of the PPF) at \( p \) of the given distribution.

Parameters:

Name Type Description Default
p _T_x

Scalar or array-like of lower tail probability values in \( [0, 1] \).

required
See Also

lmo.distributions.l_poly.cdf(x)

Cumulative distribution function \( F(x) = \mathrm{P}(X \le x) \) at \( x \) of the given distribution.

Note

Because the CDF of l_poly is not analytically expressible, it is evaluated numerically using a root-finding algorithm.

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.logcdf(x)

Logarithm of the cumulative distribution function (CDF) at \( x \), i.e. \( \ln F(x) \).

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.sf(x)

Survival function \(S(x) = \mathrm{P}(X > x) = 1 - \mathrm{P}(X \le x) = 1 - F(x) \) (the complement of the CDF).

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.logsf(x)

Logarithm of the survical function (SF) at \( x \), i.e. \( \ln \left( S(x) \right) \).

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.pdf(x)

Probability density function \( f \equiv \frac{\dd{F}}{\dd{x}} \) (derivative of the CDF) at \( x \).

By applying the inverse function rule, the PDF can also defined using the QDF as \( f(x) = 1 / q\big(F(x)\big) \).

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.logpdf(x)

Logarithm of the PDF.

lmo.distributions.l_poly.hf(x)

Hazard function \( h(x) = f(x) / S(x) \) at \( x \), with \( f \) and \( S \) the PDF and SF, respectively.

Parameters:

Name Type Description Default
x _T_x

Scalar or array-like of quantiles.

required

lmo.distributions.l_poly.median()

Median (50th percentile) of the distribution. Alias for ppf(1 / 2).

See Also

lmo.distributions.l_poly.mean()

The mean \( \mu = \E[X] \) of random varianble \( X \) of the relevant distribution.

See Also

lmo.distributions.l_poly.var()

The variance \( \Var[X] = \E\bigl[(X - \E[X])^2\bigr] = \E\bigl[X^2\bigr] - \E[X]^2 = \sigma^2 \) (2nd central product moment) of random varianble \( X \) of the relevant distribution.

See Also

lmo.distributions.l_poly.std()

The standard deviation \( \Std[X] = \sqrt{\Var[X]} = \sigma \) of random varianble \( X \) of the relevant distribution.

See Also

lmo.distributions.l_poly.entropy()

Differential entropy \( \mathrm{H}[X] \) of random varianble \( X \) of the relevant distribution.

It is defined as

\[ \mathrm{H}[X] = \E \bigl[ -\ln f(X) \bigr] = -\int_{Q(0)}^{Q(1)} \ln f(x) \dd x = \int_0^1 \ln q(p) \dd p , \]

with \( f(x) \) the PDF, \( Q(p) \) the PPF, and \( q(p) = Q'(p) \) the QDF.

See Also

lmo.distributions.l_poly.support()

The support \( (Q(0), Q(1)) \) of the distribution, where \( Q(p) \) is the PPF.

lmo.distributions.l_poly.interval(confidence)

interval(confidence: _AnyReal0D) -> tuple[np.float64, np.float64]
interval(
    confidence: _AnyRealND,
) -> tuple[npt.NDArray[np.float64], npt.NDArray[np.float64]]

Confidence interval with equal areas around the median.

For confidence level \( \alpha \in [0, 1] \), this function evaluates

\[ \left[ Q\left( \frac{1 - \alpha}{2} \right) , Q\left( \frac{1 + \alpha}{2} \right) \right], \]

where \( Q(p) \) is the PPF.

Parameters:

Name Type Description Default
confidence _AnyReal0D | _AnyRealND

Scalar or array-like. The Probability that a random varianble will be drawn from the returned range.

Each confidence value should be between 0 and 1.

required

lmo.distributions.l_poly.moment(n)

Non-central product moment \( \E[X^n] \) of \( X \) of specified order \( n \).

Note

The product moment is evaluated using numerical integration (scipy.integrate.quad), which cannot check whether the product-moment actually exists for the distribution, in which case an invalid result will be returned.

Parameters:

Name Type Description Default
n int | integer[Any]

Order \( n \ge 0 \) of the moment.

required
See Also
Todo
  • For n>=2, attempt tot infer from _l_moments if the 2nd moment condition holds, using diagnostics.l_moment_bounds.

lmo.distributions.l_poly.stats(moments='mv')

stats() -> _Tuple2[float]
stats(moments: _Stats0) -> tuple[]
stats(moments: _Stats1) -> _Tuple1[float]
stats(moments: _Stats2) -> _Tuple2[float]
stats(moments: _Stats3) -> _Tuple3[float]
stats(moments: _Stats4) -> _Tuple4[float]

Some product-moment statistics of the given distribution.

Parameters:

Name Type Description Default
moments _Stats

Composed of letters mvsk defining which product-moment statistic to compute:

'm':
Mean \( \mu = \E[X] \)
'v':
Variance \( \sigma^2 = \E[(X - \mu)^2] \)
's':
Skewness \( \E[(X - \mu)^3] / \sigma^3 \)
'k':
Ex. Kurtosis \( \E[(X - \mu)^4] / \sigma^4 - 3 \)
'mv'

lmo.distributions.l_poly.expect(g)

Calculate expected value of a function with respect to the distribution by numerical integration.

The expected value of a function \( g(x) \) with respect to a random variable \( X \) is defined as

\[ \E\left[ g(X) \right] = \int_{Q(0)}^{Q(1)} g(x) f(x) \dd x = \int_0^1 g\big(Q(u)\big) \dd u , \]

with \( f(x) \) the PDF, and \( Q \) the PPF.

Parameters:

Name Type Description Default
g (float) -> float

Continuous and deterministic function \( g: \reals \mapsto \reals \).

required

lmo.distributions.l_poly.l_moment(r, /, trim=None)

l_moment(r: lmt.AnyOrder, /, trim: lmt.AnyTrim | None = None) -> np.float64
l_moment(
    r: lmt.AnyOrderND, /, trim: lmt.AnyTrim | None = None
) -> npt.NDArray[np.float64]

Evaluate the population L-moment(s) \(\lambda^{(s,t)}_r\).

Parameters:

Name Type Description Default
r AnyOrder | AnyOrderND

L-moment order(s), non-negative integer or array-like of integers.

required
trim AnyTrim | None

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

None

lmo.distributions.l_poly.l_ratio(r, k, /, trim=None)

l_ratio(
    r: lmt.AnyOrder, k: lmt.AnyOrder, /, trim: lmt.AnyTrim | None = None
) -> np.float64
l_ratio(
    r: lmt.AnyOrderND,
    k: lmt.AnyOrder | lmt.AnyOrderND,
    /,
    trim: lmt.AnyTrim | None = None,
) -> npt.NDArray[np.float64]
l_ratio(
    r: lmt.AnyOrder | lmt.AnyOrderND,
    k: lmt.AnyOrderND,
    /,
    trim: lmt.AnyTrim | None = None,
) -> npt.NDArray[np.float64]

Evaluate the population L-moment ratio(’s) \(\tau^{(s,t)}_{r,k}\).

Parameters:

Name Type Description Default
r AnyOrder | AnyOrderND

L-moment order(s), non-negative integer or array-like of integers.

required
k AnyOrder | AnyOrderND

L-moment order of the denominator, e.g. 2.

required
trim AnyTrim | None

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

None

lmo.distributions.l_poly.l_stats(trim=None, moments=4)

Evaluate the L-moments (for \(r \le 2\)) and L-ratio’s (for \(r > 2\)).

Parameters:

Name Type Description Default
trim AnyTrim | None

Left- and right- trim. Can be scalar or 2-tuple of non-negative int or float.

None
moments int

The amount of L-moments to return. Defaults to 4.

4

lmo.distributions.l_poly.l_loc(trim=None)

L-location of the distribution, i.e. the 1st L-moment.

Alias for l_poly.l_moment(1, ...).

See Also

lmo.distributions.l_poly.l_scale(trim=None)

L-scale of the distribution, i.e. the 2nd L-moment.

Alias for l_poly.l_moment(2, ...).

See Also

lmo.distributions.l_poly.l_skew(trim=None)

L-skewness coefficient of the distribution; the 3rd L-moment ratio.

Alias for l_poly.l_ratio(3, 2, ...).

See Also

lmo.distributions.l_poly.l_kurtosis(trim=None)

L-kurtosis coefficient of the distribution; the 4th L-moment ratio.

Alias for l_poly.l_ratio(4, 2, ...).

See Also

lmo.distributions.l_poly.nnlf(theta, x)

Negative loglikelihood function.

This is calculated as -sum(log pdf(x, *theta), axis=0), where theta are the vector of L-moments, and optionally the trim.

Notes

This is mostly for compatibility rv_generic, and is impractically slow (due to the numerical inversion of the ppf).

Parameters:

Name Type Description Default
theta _LPolyParams

Tuple of size 1 or 2, with the L-moments vector, and optionally the trim (defaults to 0).

required
x NDArray[float64]

Array-like with observations of shape (n) or (n, *ks).

required

Returns:

Type Description
float | NDArray[float64]

Scalar or array of shape (*ks) with negative loglikelihoods.