Skip to content

Commit a01f866

Browse files
committed
Improve Rician distribution
1 parent 71309bc commit a01f866

File tree

4 files changed

+95
-13
lines changed

4 files changed

+95
-13
lines changed

Project.toml

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,13 +1,14 @@
11
name = "Distributions"
22
uuid = "31c24e10-a181-5473-b8eb-7969acd0382f"
3-
version = "0.25.123"
43
authors = ["JuliaStats"]
4+
version = "0.25.124"
55

66
[deps]
77
AliasTables = "66dad0bd-aa9a-41b7-9441-69ab47430ed8"
88
ChainRulesCore = "d360d2e6-b24c-11e9-a2a3-2a2ae2dbcce4"
99
DensityInterface = "b429d917-457f-4dbc-8f4c-0cc954292b1d"
1010
FillArrays = "1a297f60-69ca-5386-bcde-b61e274b549b"
11+
HypergeometricFunctions = "34004b35-14d8-5ef3-9330-4cdb6864b03a"
1112
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
1213
PDMats = "90014a1f-27ba-587c-ab20-58faa44d9150"
1314
Printf = "de0858da-6303-5e67-8744-51eddeeeb8d7"
@@ -42,6 +43,7 @@ Distributed = "<0.0.1, 1"
4243
FillArrays = "0.9, 0.10, 0.11, 0.12, 0.13, 1"
4344
FiniteDifferences = "0.12"
4445
ForwardDiff = "0.10, 1"
46+
HypergeometricFunctions = "0.3.28"
4547
JSON = "0.21, 1"
4648
LinearAlgebra = "<0.0.1, 1"
4749
OffsetArrays = "1"

src/Distributions.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -25,6 +25,7 @@ import StatsBase: kurtosis, skewness, entropy, mode, modes,
2525
import PDMats: dim, PDMat, invquad
2626

2727
using SpecialFunctions
28+
import HypergeometricFunctions
2829
using Base.MathConstants: eulergamma
2930

3031
import AliasTables

src/univariate/continuous/rician.jl

Lines changed: 35 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -63,11 +63,11 @@ partype(d::Rician{T}) where {T<:Real} = T
6363

6464
#### Statistics
6565

66-
# helper
67-
_Lhalf(x) = exp(x/2) * ((1-x) * besseli(zero(x), -x/2) - x * besseli(oneunit(x), -x/2))
66+
# Laguerre function L_{1/2}
67+
_Lhalf(x::Real) = HypergeometricFunctions._₁F₁(-one(x)/2, one(x), x)
6868

69-
mean(d::Rician) = d.σ * sqrthalfπ * _Lhalf(-d.ν^2/(2 * d.σ^2))
70-
var(d::Rician) = 2 * d.σ^2 + d.ν^2 - halfπ * d.σ^2 * _Lhalf(-d.ν^2/(2 * d.σ^2))^2
69+
mean(d::Rician) = d.σ * sqrthalfπ * _Lhalf(-(d.ν / d.σ)^2/2)
70+
var(d::Rician) = 2 * d.σ^2 + d.ν^2 - halfπ * d.σ^2 * _Lhalf(-(d.ν/d.σ)^2/2)^2
7171

7272
function mode(d::Rician)
7373
m = mean(d)
@@ -91,7 +91,21 @@ function _minimize_gss(f, a, b; tol=1e-12)
9191
(b + a) / 2
9292
end
9393

94-
#### PDF/CDF/quantile delegated to NoncentralChisq
94+
#### PDF
95+
96+
function pdf(d::Rician, x::Real)
97+
ν, σ = params(d)
98+
result = x / σ^2 * exp(-((x-ν)/σ)^2/2) * besselix(0, (x * ν)/ σ^2)
99+
return x < 0 || isinf(x) ? zero(result) : result
100+
end
101+
102+
function logpdf(d::Rician, x::Real)
103+
ν, σ = params(d)
104+
result = log(abs(x) / σ^2 * besselix(0, (x * ν)/ σ^2)) - ((x-ν)/σ)^2/2
105+
return x < 0 || isinf(x) ? oftype(result, -Inf) : result
106+
end
107+
108+
#### quantile/CDF delegated to NoncentralChisq
95109

96110
function quantile(d::Rician, x::Real)
97111
ν, σ = params(d)
@@ -103,16 +117,14 @@ function cquantile(d::Rician, x::Real)
103117
return sqrt(cquantile(NoncentralChisq(2, (ν / σ)^2), x)) * σ
104118
end
105119

106-
function pdf(d::Rician, x::Real)
120+
function invlogcdf(d::Rician, lx::Real)
107121
ν, σ = params(d)
108-
result = 2 * x / σ^2 * pdf(NoncentralChisq(2, (ν / σ)^2), (x / σ)^2)
109-
return x < 0 || isinf(x) ? zero(result) : result
122+
return sqrt(invlogcdf(NoncentralChisq(2, (ν / σ)^2), lx)) * σ
110123
end
111124

112-
function logpdf(d::Rician, x::Real)
125+
function invlogccdf(d::Rician, lx::Real)
113126
ν, σ = params(d)
114-
result = log(2 * abs(x) / σ^2) + logpdf(NoncentralChisq(2, (ν / σ)^2), (x / σ)^2)
115-
return x < 0 || isinf(x) ? oftype(result, -Inf) : result
127+
return sqrt(invlogccdf(NoncentralChisq(2, (ν / σ)^2), lx)) * σ
116128
end
117129

118130
function cdf(d::Rician, x::Real)
@@ -127,6 +139,18 @@ function logcdf(d::Rician, x::Real)
127139
return x < 0 ? oftype(result, -Inf) : result
128140
end
129141

142+
function ccdf(d::Rician, x::Real)
143+
ν, σ = params(d)
144+
result = ccdf(NoncentralChisq(2, (ν / σ)^2), (x / σ)^2)
145+
return x < 0 ? zero(result) : result
146+
end
147+
148+
function logccdf(d::Rician, x::Real)
149+
ν, σ = params(d)
150+
result = logccdf(NoncentralChisq(2, (ν / σ)^2), (x / σ)^2)
151+
return x < 0 ? oftype(result, -Inf) : result
152+
end
153+
130154
#### Sampling
131155

132156
function rand(rng::AbstractRNG, d::Rician)

test/univariate/continuous/rician.jl

Lines changed: 56 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,5 @@
1+
using StatsFuns
2+
13
@testset "Rician" begin
24

35
d1 = Rician(0.0, 10.0)
@@ -9,6 +11,24 @@
911
@test eltype(d1) === Float64
1012
@test rand(d1) isa Float64
1113

14+
# Reference: WolframAlpha
15+
@test @inferred(mean(d1))::Float64 5 * sqrt(2 * π)
16+
@test @inferred(std(d1))::Float64 sqrt(200 - 50 * π)
17+
@test @inferred(var(d1))::Float64 200 - 50 * π
18+
@test @inferred(mode(d1))::Float64 10
19+
ps = 0:0.1:1
20+
@test all(splat(isapprox), zip(quantile.(d1, ps), @. sqrt(-200 * log1p(-ps))))
21+
@test all(splat(isapprox), zip(cquantile.(d1, ps), @. sqrt(-200 * log(ps))))
22+
@test all(splat(isapprox), zip(invlogcdf.(d1, log.(ps)), @. sqrt(-200 * log1p(-ps))))
23+
@test all(splat(isapprox), zip(invlogccdf.(d1, log.(ps)), @. sqrt(-200 * log(ps))))
24+
xs = 0:0.5:30 # 99th percentile is approx 30
25+
@test all(splat(isapprox), zip(pdf.(d1, xs), map(x -> x / 100 * exp(-x^2/200), xs)))
26+
@test all(splat(isapprox), zip(logpdf.(d1, xs), map(x -> log(x / 100) - x^2/200, xs)))
27+
@test all(splat(isapprox), zip(cdf.(d1, xs), @. 1 - exp(-xs^2/200)))
28+
@test all(splat(isapprox), zip(logcdf.(d1, xs), @. log1mexp(-xs^2/200)))
29+
@test all(splat(isapprox), zip(ccdf.(d1, xs), @. exp(-xs^2/200)))
30+
@test all(splat(isapprox), zip(logccdf.(d1, xs), @. -xs^2/200))
31+
1232
d2 = Rayleigh(10.0)
1333
@test mean(d1) mean(d2)
1434
@test var(d1) var(d2)
@@ -23,7 +43,22 @@
2343
x = Base.Fix1(quantile, d1).([0.25, 0.45, 0.60, 0.80, 0.90])
2444
@test all(Base.Fix1(cdf, d1).(x) .≈ [0.25, 0.45, 0.60, 0.80, 0.90])
2545

26-
x = rand(Rician(5.0, 5.0), 100000)
46+
# Reference: WolframAlpha
47+
@test mean(d1) 15.4857246055
48+
@test std(d1) 7.75837182934
49+
@test var(d1) 60.1923334423
50+
@test median(d1) 14.7547909179
51+
@test all(xy -> isapprox(xy...; rtol=1e-5), zip(quantile.(d1, [0.1, 0.25, 0.5, 0.75, 0.9]), [5.86808, 9.62923, 14.7548, 20.5189, 26.0195]))
52+
53+
d1 = Rician(5.0, 5.0)
54+
# Reference: WolframAlpha
55+
@test mean(d1) 7.7428623027557
56+
@test std(d1) 3.8791859146687
57+
@test var(d1) 15.0480833605606
58+
@test median(d1) 7.37739545894
59+
@test all(xy -> isapprox(xy...; rtol=1e-5), zip(quantile.(d1, [0.1, 0.25, 0.5, 0.75, 0.9]), [2.93404, 4.81462, 7.3774, 10.2595, 13.0097]))
60+
61+
x = rand(d1, 100000)
2762
d1 = fit(Rician, x)
2863
@test d1 isa Rician{Float64}
2964
@test params(d1)[1] 5.0 atol=0.2
@@ -37,6 +72,11 @@
3772
@test partype(d1) === Float32
3873
@test eltype(d1) === Float64
3974
@test rand(d1) isa Float64
75+
d1_f64 = Rician(10.0, 10.0)
76+
@test @inferred(mean(d1))::Float32 Float32(mean(d1_f64))
77+
@test @inferred(std(d1))::Float32 Float32(std(d1_f64))
78+
@test @inferred(var(d1))::Float32 Float32(var(d1_f64))
79+
@test @inferred(median(d1))::Float32 Float32(median(d1_f64))
4080

4181
d1 = Rician()
4282
@test d1 isa Rician{Float64}
@@ -94,4 +134,19 @@
94134
@inferred logcdf(d1, 1//2)
95135
@inferred logcdf(d1, Inf)
96136

137+
# issue #2035
138+
# reference values computed with WolframAlpha
139+
d = Rician(1.0, 0.12)
140+
@test mean(d) 1.00722650704
141+
@test std(d) 0.119560710568
142+
@test var(d) 0.0142947635114
143+
@test median(d) 1.00719144719
144+
@test pdf(d, 0.5) 0.000400758853946
145+
@test pdf(d, 1) 3.33055235263
146+
@test pdf(d, 1.5) 0.000692437774661
147+
@test pdf(d, 2) 3.91711741719e-15
148+
@test logpdf(d, 0.5) -7.82215067328
149+
@test logpdf(d, 1) 1.2031381619
150+
@test logpdf(d, 1.5) -7.27529218003
151+
@test logpdf(d, 2) -33.1734203644
97152
end

0 commit comments

Comments
 (0)