feat(distributions): add lograyleigh and logrice distributions

Add Log-Rayleigh and Log-Rice continuous distributions as
scipy rv_continuous subclasses with PDF, CDF, SF, PPF, ISF,
moments, entropy, and RVS methods.

Log-Rice reduces to Log-Rayleigh when nu=0. Both are derived
via the change-of-variable Y = ln X on their respective parent
distributions. Includes unit tests verifying numerical
correctness and the change-of-variable identity.
This commit is contained in:
2026-04-27 11:11:32 -03:00
parent 9a3f5959cd
commit e780bb956e
2 changed files with 347 additions and 57 deletions

View File

@@ -1,4 +1,4 @@
from scipy.stats import rv_continuous
from scipy.stats import rv_continuous, ncx2
from scipy.special import kve, gammaln
from scipy.integrate import quad
import scipy.special as sc
@@ -9,7 +9,7 @@ import numpy as np
class k_gen(rv_continuous):
"""Generalized K distribution for radar clutter modeling.
The probability density function is::
The probability density function is:
f(x; mu, alpha, beta) = 2 / (Gamma(alpha) * Gamma(beta))
* (alpha*beta/mu)^((alpha+beta)/2)
@@ -103,22 +103,16 @@ k_dist = k_gen(a=0.0, name="k_distribution", shapes="mu, alpha, beta")
class lognakagami_gen(rv_continuous):
r"""A Log-Nakagami continuous random variable.
"""Log-Nakagami continuous random variable.
If X follows a Nakagami distribution with shape parameter m >= 0.5
and spread Omega > 0, then Y = ln(X) follows the Log-Nakagami distribution.
%(before_notes)s
Notes
-----
The probability density function for `lognakagami` is:
The probability density function is:
f(y; m, Omega) = 2*m^m / (Gamma(m)*Omega^m)
* exp(2*m*y) * exp(-(m/Omega)*exp(2*y))
for y in (-inf, +inf). Derived via the change of variable Y = ln X,
whose Jacobian e^y combines with x^{2m-1} to give the e^{2my} term.
for y in (-inf, +inf), m >= 0.5, Omega > 0. Derived via the change of
variable Y = ln X on a Nakagami(m, Omega) variate; the Jacobian e^y
combines with x^(2m-1) to give the exp(2*m*y) term.
The mean and variance are:
@@ -130,10 +124,6 @@ class lognakagami_gen(rv_continuous):
unchanged; entropy is also Omega-independent.
`lognakagami` takes ``m`` and ``Omega`` as shape parameters.
See Also
--------
nakagami : The Nakagami distribution (before the log transform).
"""
def _shape_info(self):
@@ -145,7 +135,6 @@ class lognakagami_gen(rv_continuous):
def _argcheck(self, m, Omega):
return (m >= 0.5) & (Omega > 0)
# ---------- PDF / log-PDF ----------
def _pdf(self, y, m, Omega):
return np.exp(self._logpdf(y, m, Omega))
@@ -159,7 +148,6 @@ class lognakagami_gen(rv_continuous):
+ 2.0 * m * y
- (m / Omega) * np.exp(2.0 * y))
# ---------- CDF / SF / log-CDF / log-SF ----------
def _cdf(self, y, m, Omega):
# P(Y <= y) = P(X <= e^y) = gammainc(m, (m/Omega) * e^{2y})
@@ -174,7 +162,6 @@ class lognakagami_gen(rv_continuous):
def _logsf(self, y, m, Omega):
return np.log(sc.gammaincc(m, (m / Omega) * np.exp(2.0 * y)))
# ---------- PPF / ISF ----------
def _ppf(self, q, m, Omega):
# q = gammainc(m, (m/Omega)*e^{2y})
@@ -185,8 +172,6 @@ class lognakagami_gen(rv_continuous):
def _isf(self, q, m, Omega):
return 0.5 * np.log(Omega * sc.gammainccinv(m, q) / m)
# ---------- Moments & statistics ----------
def _stats(self, m, Omega):
# E[Y] = 0.5*(psi(m) - ln(m) + ln(Omega))
mu = 0.5 * (sc.digamma(m) - np.log(m) + np.log(Omega))
@@ -200,9 +185,7 @@ class lognakagami_gen(rv_continuous):
# H = -E[ln f(Y)] = -ln2 + ln Γ(m) - m*psi(m) + m
# Omega cancels exactly (translation invariance of differential entropy)
return -np.log(2.0) + sc.gammaln(m) - m * sc.digamma(m) + m
# ---------- Random variate generation ----------
def _rvs(self, m, Omega, size=None, random_state=None):
# X ~ Nakagami(m, Omega) => X^2 ~ Gamma(m, Omega/m)
g = random_state.gamma(m, Omega / m, size=size)
@@ -212,34 +195,27 @@ class lognakagami_gen(rv_continuous):
lognakagami = lognakagami_gen(name='lognakagami', shapes="m, Omega")
class loggamma_gen(rv_continuous):
r"""A Log-Gamma continuous random variable.
"""Log-Gamma continuous random variable.
If X follows a Gamma distribution with shape a
and unit scale (beta = 1), then Y = ln X follows
the Log-Gamma distribution.
The probability density function is:
%(before_notes)s
f(y; a) = 1 / Gamma(a) * exp(a*y - exp(y))
Notes
-----
The probability density function for `loggamma_dist` is:
for y in (-inf, +inf), a > 0. Derived via the change of variable
Y = ln X on a Gamma(a, 1) variate; the Jacobian e^y combines with
x^(a-1) to give the exp(a*y) term.
f(y, a) = 1 / Gamma(a) * exp(a*y - exp(y))
The mean and variance are:
for y in (-inf, +inf), a > 0.
E[Y] = psi(a)
Var[Y] = psi_1(a)
This is the standardised form (scale beta = 1).
where psi and psi_1 are the digamma and trigamma functions.
The non-unit scale beta is recovered through scipy's ``loc`` parameter:
loc = ln(beta), since ln(beta * X) = ln(beta) + ln(X).
`loggamma_dist` takes ``a`` as a shape parameter.
The non-unit scale beta is recovered through scipy's
``loc`` parameter: loc = ln(beta), since
ln(beta * X) = ln(beta) + ln(X).
See Also
--------
gamma : The Gamma distribution (before the log transform).
"""
def _shape_info(self):
@@ -248,8 +224,6 @@ class loggamma_gen(rv_continuous):
def _argcheck(self, a):
return a > 0
# ---------- PDF / log-PDF ----------
def _pdf(self, y, a):
return np.exp(self._logpdf(y, a))
@@ -257,8 +231,6 @@ class loggamma_gen(rv_continuous):
# ln f = a*y - e^y - ln Γ(a)
return a * y - np.exp(y) - sc.gammaln(a)
# ---------- CDF / SF / log-CDF / log-SF ----------
def _cdf(self, y, a):
# P(Y <= y) = P(X <= e^y) = gammainc(a, e^y)
return sc.gammainc(a, np.exp(y))
@@ -272,8 +244,6 @@ class loggamma_gen(rv_continuous):
def _logsf(self, y, a):
return np.log(sc.gammaincc(a, np.exp(y)))
# ---------- PPF / ISF ----------
def _ppf(self, q, a):
# Invert CDF: q = gammainc(a, e^y)
# => e^y = gammaincinv(a, q)
@@ -283,8 +253,6 @@ class loggamma_gen(rv_continuous):
def _isf(self, q, a):
return np.log(sc.gammainccinv(a, q))
# ---------- Moments & statistics ----------
def _stats(self, a):
# Y = ln(X), X ~ Gamma(a, 1)
# E[Y] = psi(a) (digamma)
@@ -303,8 +271,6 @@ class loggamma_gen(rv_continuous):
# E[e^Y] = E[X] = a (for Gamma(a, 1))
return sc.gammaln(a) - a * sc.digamma(a) + a
# ---------- Random variate generation ----------
def _rvs(self, a, size=None, random_state=None):
# Generate Gamma(a, 1) variates, then take log.
# For small a, use Gamma(a) ~ Gamma(a+1) * U^(1/a)
@@ -314,4 +280,154 @@ class loggamma_gen(rv_continuous):
+ np.log(random_state.uniform(size=size)) / a)
loggamma_dist = loggamma_gen(name='loggamma_dist')
loggamma_dist = loggamma_gen(name='loggamma_dist')
class lograyleigh_gen(rv_continuous):
"""Log-Rayleigh continuous random variable.
The probability density function is:
f(y; sigma) = (e^{2y} / sigma^2) * exp(-e^{2y} / (2*sigma^2))
for y in (-inf, +inf), sigma > 0. Derived via the change of variable
Y = ln X on a Rayleigh(sigma) variate; the Jacobian e^y combines with
x to give the e^{2y} term.
The mean and variance are:
E[Y] = 0.5 * (ln(2*sigma^2) + psi(1))
Var[Y] = pi^2 / 24
where psi is the digamma function. The scale sigma shifts the mean by
ln(sigma) but leaves all higher central moments unchanged; entropy is
also sigma-independent.
`lograyleigh` takes ``sigma`` as a shape parameter.
"""
def _shape_info(self):
return [_ShapeInfo("sigma", False, (0, np.inf), (False, False))]
def _argcheck(self, sigma):
return sigma > 0
def _pdf(self, y, sigma):
return np.exp(self._logpdf(y, sigma))
def _logpdf(self, y, sigma):
# ln f = 2y - 2*ln(sigma) - e^{2y} / (2*sigma^2)
return 2.0 * y - 2.0 * np.log(sigma) - np.exp(2.0 * y) / (2.0 * sigma**2)
def _cdf(self, y, sigma):
return 1.0 - np.exp(-np.exp(2.0 * y) / (2.0 * sigma**2))
def _sf(self, y, sigma):
return np.exp(-np.exp(2.0 * y) / (2.0 * sigma**2))
def _logcdf(self, y, sigma):
return np.log1p(-np.exp(-np.exp(2.0 * y) / (2.0 * sigma**2)))
def _logsf(self, y, sigma):
return -np.exp(2.0 * y) / (2.0 * sigma**2)
def _ppf(self, q, sigma):
# q = 1 - exp(-e^{2y}/(2*sigma^2))
# => e^{2y} = -2*sigma^2 * ln(1-q)
# => y = 0.5 * ln(-2*sigma^2 * ln(1-q))
return 0.5 * np.log(-2.0 * sigma**2 * np.log1p(-q))
def _isf(self, q, sigma):
# q = exp(-e^{2y}/(2*sigma^2))
# => y = 0.5 * ln(-2*sigma^2 * ln(q))
return 0.5 * np.log(-2.0 * sigma**2 * np.log(q))
def _stats(self, sigma):
# Y = 0.5*ln(2) + 0.5*ln(W) + ln(sigma), W ~ Exp(1) = Gamma(1, 1)
# E[Y] = 0.5*(ln(2*sigma^2) + psi(1))
# Var[Y] = 0.25*psi_1(1) = pi^2/24
# Standardised skewness/kurtosis are sigma-independent (pure shift)
mu = 0.5 * (np.log(2.0 * sigma**2) + sc.digamma(1))
mu2 = 0.25 * sc.polygamma(1, 1)
g1 = sc.polygamma(2, 1) / sc.polygamma(1, 1)**1.5
g2 = sc.polygamma(3, 1) / sc.polygamma(1, 1)**2
return mu, mu2, g1, g2
def _entropy(self, sigma):
# H(Y) = 1 - ln(2) - psi(1) (sigma-independent; derived via H(X) - E[ln X])
return 1.0 - np.log(2.0) - sc.digamma(1)
def _rvs(self, sigma, size=None, random_state=None):
# X ~ Rayleigh(sigma) = sigma*sqrt(-2*ln(U)), U ~ Uniform(0,1)
u = random_state.uniform(size=size)
return np.log(sigma) + 0.5 * np.log(-2.0 * np.log(u))
lograyleigh = lograyleigh_gen(name='lograyleigh', shapes="sigma")
class logrice_gen(rv_continuous):
"""Log-Rice continuous random variable.
The probability density function is:
f(y; nu, sigma) = (e^{2y} / sigma^2)
* exp(-(e^{2y} + nu^2) / (2*sigma^2))
* I_0(e^y * nu / sigma^2)
for y in (-inf, +inf), nu >= 0, sigma > 0. Derived via the change of
variable Y = ln X on a Rice(nu, sigma) variate; the Jacobian e^y
combines with x to give the e^{2y} factor. I_0 is the modified Bessel
function of the first kind of order zero.
The CDF involves the Marcum Q-function of order 1:
F(y; nu, sigma) = 1 - Q_1(nu/sigma, e^y/sigma)
When nu = 0, the distribution reduces to Log-Rayleigh.
`logrice` takes ``nu`` and ``sigma`` as shape parameters.
"""
def _shape_info(self):
return [
_ShapeInfo("nu", False, (0, np.inf), (True, False)),
_ShapeInfo("sigma", False, (0, np.inf), (False, False)),
]
def _argcheck(self, nu, sigma):
return (nu >= 0) & (sigma > 0)
def _pdf(self, y, nu, sigma):
return np.exp(self._logpdf(y, nu, sigma))
def _logpdf(self, y, nu, sigma):
# ln f = 2y - 2*ln(sigma) - (e^{2y} + nu^2)/(2*sigma^2) + ln(I_0(z))
# Use ive(0, z)*exp(z) = I_0(z) for stability: ln(I_0(z)) = ln(ive(0,z)) + z
z = np.exp(y) * nu / sigma**2
return (2.0 * y
- 2.0 * np.log(sigma)
- (np.exp(2.0 * y) + nu**2) / (2.0 * sigma**2)
+ np.log(sc.ive(0, z)) + z)
def _cdf(self, y, nu, sigma):
# P(Y <= y) = 1 - Q_1(nu/sigma, e^y/sigma) = ncx2.cdf(e^{2y}/sigma^2, 2, nu^2/sigma^2)
return ncx2.cdf(np.exp(2.0 * y) / sigma**2, 2, (nu / sigma)**2)
def _sf(self, y, nu, sigma):
return ncx2.sf(np.exp(2.0 * y) / sigma**2, 2, (nu / sigma)**2)
def _logcdf(self, y, nu, sigma):
return ncx2.logcdf(np.exp(2.0 * y) / sigma**2, 2, (nu / sigma)**2)
def _logsf(self, y, nu, sigma):
return ncx2.logsf(np.exp(2.0 * y) / sigma**2, 2, (nu / sigma)**2)
def _rvs(self, nu, sigma, size=None, random_state=None):
# X ~ Rice(nu, sigma): norm of 2D Gaussian with mean (nu, 0) and std sigma
z1 = random_state.normal(nu, sigma, size=size)
z2 = random_state.normal(0.0, sigma, size=size)
return np.log(np.hypot(z1, z2))
logrice = logrice_gen(name='logrice', shapes="nu, sigma")