-
Notifications
You must be signed in to change notification settings - Fork 422
/
Copy pathlkjcholesky.jl
267 lines (232 loc) · 8.65 KB
/
lkjcholesky.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
"""
LKJCholesky(d::Int, η::Real, uplo='L')
The `LKJCholesky` distribution of size ``d`` with shape parameter ``\\eta`` is a
distribution over `LinearAlgebra.Cholesky` factorisations of ``d\\times d`` real correlation
matrices (positive-definite matrices with ones on the diagonal).
Variates or samples of the distribution are `LinearAlgebra.Cholesky` objects, as might
be returned by `F = LinearAlgebra.cholesky(R)`, so that `Matrix(F) ≈ R` is a variate or
sample of [`LKJ`](@ref).
Sampling `LKJCholesky` is faster than sampling `LKJ`, and often having the correlation
matrix in factorized form makes subsequent computations cheaper as well.
!!! note
`LinearAlgebra.Cholesky` stores either the upper or lower Cholesky factor, related by
`F.U == F.L'`. Both can be accessed with `F.U` and `F.L`, but if the factor
not stored is requested, then a copy is made. The `uplo` parameter specifies whether the
upper (`'U'`) or lower (`'L'`) Cholesky factor is stored when randomly generating
samples. Set `uplo` to `'U'` if the upper factor is desired to avoid allocating a copy
when calling `F.U`.
See [`LKJ`](@ref) for more details.
External links
* Lewandowski D, Kurowicka D, Joe H.
Generating random correlation matrices based on vines and extended onion method,
Journal of Multivariate Analysis (2009), 100(9): 1989-2001
doi: [10.1016/j.jmva.2009.04.008](https://doi.org/10.1016/j.jmva.2009.04.008)
"""
struct LKJCholesky{T <: Real} <: Distribution{CholeskyVariate,Continuous}
d::Int
η::T
uplo::Char
logc0::T
end
# -----------------------------------------------------------------------------
# Constructors
# -----------------------------------------------------------------------------
function LKJCholesky(d::Int, η::Real, _uplo::Union{Char,Symbol} = 'L'; check_args::Bool=true)
@check_args(
LKJCholesky,
(d, d > 0, "matrix dimension must be positive"),
(η, η > 0, "shape parameter must be positive"),
)
logc0 = lkj_logc0(d, η)
uplo = _char_uplo(_uplo)
T = Base.promote_eltype(η, logc0)
return LKJCholesky(d, T(η), uplo, T(logc0))
end
# adapted from LinearAlgebra.char_uplo
function _char_uplo(uplo::Union{Symbol,Char})
uplo ∈ (:U, 'U') && return 'U'
uplo ∈ (:L, 'L') && return 'L'
throw(ArgumentError("uplo argument must be either 'U' (upper) or 'L' (lower)"))
end
# -----------------------------------------------------------------------------
# REPL display
# -----------------------------------------------------------------------------
Base.show(io::IO, d::LKJCholesky) = show(io, d, (:d, :η, :uplo))
# -----------------------------------------------------------------------------
# Conversion
# -----------------------------------------------------------------------------
function Base.convert(::Type{LKJCholesky{T}}, d::LKJCholesky) where T <: Real
return LKJCholesky{T}(d.d, T(d.η), d.uplo, T(d.logc0))
end
Base.convert(::Type{LKJCholesky{T}}, d::LKJCholesky{T}) where T <: Real = d
function convert(::Type{LKJCholesky{T}}, d::Integer, η::Real, uplo::Char, logc0::Real) where T <: Real
return LKJCholesky{T}(Int(d), T(η), uplo, T(logc0))
end
# -----------------------------------------------------------------------------
# Properties
# -----------------------------------------------------------------------------
Base.eltype(::Type{LKJCholesky{T}}) where {T} = T
function Base.size(d::LKJCholesky)
p = d.d
return (p, p)
end
function insupport(d::LKJCholesky, R::LinearAlgebra.Cholesky)
p = d.d
factors = R.factors
(isreal(factors) && size(factors, 1) == p) || return false
iinds, jinds = axes(factors)
# check that the diagonal of U'*U or L*L' is all ones
@inbounds if R.uplo === 'U'
for (j, jind) in enumerate(jinds)
col_iinds = view(iinds, 1:j)
sum(abs2, view(factors, col_iinds, jind)) ≈ 1 || return false
end
else # R.uplo === 'L'
for (i, iind) in enumerate(iinds)
row_jinds = view(jinds, 1:i)
sum(abs2, view(factors, iind, row_jinds)) ≈ 1 || return false
end
end
return true
end
function StatsBase.mode(d::LKJCholesky)
factors = Matrix{eltype(d)}(LinearAlgebra.I, size(d))
return LinearAlgebra.Cholesky(factors, d.uplo, 0)
end
StatsBase.params(d::LKJCholesky) = (d.d, d.η, d.uplo)
@inline partype(::LKJCholesky{T}) where {T <: Real} = T
# -----------------------------------------------------------------------------
# Evaluation
# -----------------------------------------------------------------------------
function logkernel(d::LKJCholesky, R::LinearAlgebra.Cholesky)
factors = R.factors
p, η = params(d)
c = p + 2(η - 1)
p == 1 && return c * log(first(factors))
# assuming D = diag(factors) with length(D) = p,
# logp = sum(i -> (c - i) * log(D[i]), 2:p)
logp = sum(Iterators.drop(enumerate(diagind(factors)), 1)) do (i, di)
return (c - i) * log(factors[di])
end
return logp
end
function logpdf(d::LKJCholesky, R::LinearAlgebra.Cholesky)
insupport(d, R) || throw(ArgumentError("provided point is not in the support"))
return _logpdf(d, R)
end
_logpdf(d::LKJCholesky, R::LinearAlgebra.Cholesky) = logkernel(d, R) + d.logc0
pdf(d::LKJCholesky, R::LinearAlgebra.Cholesky) = exp(logpdf(d, R))
loglikelihood(d::LKJCholesky, R::LinearAlgebra.Cholesky) = logpdf(d, R)
function loglikelihood(d::LKJCholesky, Rs::AbstractArray{<:LinearAlgebra.Cholesky})
return sum(R -> logpdf(d, R), Rs)
end
# -----------------------------------------------------------------------------
# Sampling
# -----------------------------------------------------------------------------
function Base.rand(rng::AbstractRNG, d::LKJCholesky)
factors = Matrix{eltype(d)}(undef, size(d))
R = LinearAlgebra.Cholesky(factors, d.uplo, 0)
return _lkj_cholesky_onion_sampler!(rng, d, R)
end
function Base.rand(rng::AbstractRNG, d::LKJCholesky, dims::Dims)
p = d.d
uplo = d.uplo
T = eltype(d)
TM = Matrix{T}
Rs = Array{LinearAlgebra.Cholesky{T,TM}}(undef, dims)
for i in eachindex(Rs)
factors = TM(undef, p, p)
Rs[i] = R = LinearAlgebra.Cholesky(factors, uplo, 0)
_lkj_cholesky_onion_sampler!(rng, d, R)
end
return Rs
end
Random.rand!(d::LKJCholesky, R::LinearAlgebra.Cholesky) = Random.rand!(GLOBAL_RNG, d, R)
function Random.rand!(rng::AbstractRNG, d::LKJCholesky, R::LinearAlgebra.Cholesky)
return _lkj_cholesky_onion_sampler!(rng, d, R)
end
function Random.rand!(
rng::AbstractRNG,
d::LKJCholesky,
Rs::AbstractArray{<:LinearAlgebra.Cholesky{T,TM}},
allocate::Bool,
) where {T,TM}
p = d.d
uplo = d.uplo
if allocate
for i in eachindex(Rs)
Rs[i] = _lkj_cholesky_onion_sampler!(
rng,
d,
LinearAlgebra.Cholesky(TM(undef, p, p), uplo, 0),
)
end
else
for i in eachindex(Rs)
_lkj_cholesky_onion_sampler!(rng, d, Rs[i])
end
end
return Rs
end
function Random.rand!(
rng::AbstractRNG,
d::LKJCholesky,
Rs::AbstractArray{<:LinearAlgebra.Cholesky{<:Real}},
)
allocate = any(!isassigned(Rs, i) for i in eachindex(Rs)) || any(R -> size(R, 1) != d.d, Rs)
return Random.rand!(rng, d, Rs, allocate)
end
#
# onion method
#
function _lkj_cholesky_onion_sampler!(
rng::AbstractRNG,
d::LKJCholesky,
R::LinearAlgebra.Cholesky,
)
if R.uplo === 'U'
_lkj_cholesky_onion_tri!(rng, R.factors, d.d, d.η, Val(:U))
else
_lkj_cholesky_onion_tri!(rng, R.factors, d.d, d.η, Val(:L))
end
return R
end
function _lkj_cholesky_onion_tri!(
rng::AbstractRNG,
A::AbstractMatrix,
d::Int,
η::Real,
::Val{uplo},
) where {uplo}
# Section 3.2 in LKJ (2009 JMA)
# reformulated to incrementally construct Cholesky factor as mentioned in Section 5
# equivalent steps in algorithm in reference are marked.
@assert size(A) == (d, d)
A[1, 1] = 1
d > 1 || return R
β = η + (d - 2)//2
# 1. Initialization
w0 = 2 * rand(rng, Beta(β, β)) - 1
@inbounds if uplo === :L
A[2, 1] = w0
else
A[1, 2] = w0
end
@inbounds A[2, 2] = sqrt(1 - w0^2)
# 2. Loop, each iteration k adds row/column k+1
for k in 2:(d - 1)
# (a)
β -= 1//2
# (b)
y = rand(rng, Beta(k//2, β))
# (c)-(e)
# w is directionally uniform vector of length √y
@inbounds w = @views uplo === :L ? A[k + 1, 1:k] : A[1:k, k + 1]
Random.randn!(rng, w)
rmul!(w, sqrt(y) / norm(w))
# normalize so new row/column has unit norm
@inbounds A[k + 1, k + 1] = sqrt(1 - y)
end
# 3.
return A
end