Line | Exclusive | Inclusive | Code |
---|---|---|---|
1 | # This file is a part of Julia. License is MIT: https://julialang.org/license | ||
2 | |||
3 | ########################## | ||
4 | # Cholesky Factorization # | ||
5 | ########################## | ||
6 | |||
7 | # The dispatch structure in the cholesky, and cholesky! methods is a bit | ||
8 | # complicated and some explanation is therefore provided in the following | ||
9 | # | ||
10 | # In the methods below, LAPACK is called when possible, i.e. StridedMatrices with Float32, | ||
11 | # Float64, ComplexF32, and ComplexF64 element types. For other element or | ||
12 | # matrix types, the unblocked Julia implementation in _chol! is used. For cholesky | ||
13 | # and cholesky! pivoting is supported through a RowMaximum() argument. A type argument is | ||
14 | # necessary for type stability since the output of cholesky and cholesky! is either | ||
15 | # Cholesky or CholeskyPivoted. The latter is only | ||
16 | # supported for the four LAPACK element types. For other types, e.g. BigFloats RowMaximum() will | ||
17 | # give an error. It is required that the input is Hermitian (including real symmetric) either | ||
18 | # through the Hermitian and Symmetric views or exact symmetric or Hermitian elements which | ||
19 | # is checked for and an error is thrown if the check fails. | ||
20 | |||
21 | # The internal structure is as follows | ||
22 | # - _chol! returns the factor and info without checking positive definiteness | ||
23 | # - cholesky/cholesky! returns Cholesky without checking positive definiteness | ||
24 | |||
25 | # FixMe? The dispatch below seems overly complicated. One simplification could be to | ||
26 | # merge the two Cholesky types into one. It would remove the need for Val completely but | ||
27 | # the cost would be extra unnecessary/unused fields for the unpivoted Cholesky and runtime | ||
28 | # checks of those fields before calls to LAPACK to check which version of the Cholesky | ||
29 | # factorization the type represents. | ||
30 | """ | ||
31 | Cholesky <: Factorization | ||
32 | |||
33 | Matrix factorization type of the Cholesky factorization of a dense symmetric/Hermitian | ||
34 | positive definite matrix `A`. This is the return type of [`cholesky`](@ref), | ||
35 | the corresponding matrix factorization function. | ||
36 | |||
37 | The triangular Cholesky factor can be obtained from the factorization `F::Cholesky` | ||
38 | via `F.L` and `F.U`, where `A ≈ F.U' * F.U ≈ F.L * F.L'`. | ||
39 | |||
40 | The following functions are available for `Cholesky` objects: [`size`](@ref), [`\\`](@ref), | ||
41 | [`inv`](@ref), [`det`](@ref), [`logdet`](@ref) and [`isposdef`](@ref). | ||
42 | |||
43 | Iterating the decomposition produces the components `L` and `U`. | ||
44 | |||
45 | # Examples | ||
46 | ```jldoctest | ||
47 | julia> A = [4. 12. -16.; 12. 37. -43.; -16. -43. 98.] | ||
48 | 3×3 Matrix{Float64}: | ||
49 | 4.0 12.0 -16.0 | ||
50 | 12.0 37.0 -43.0 | ||
51 | -16.0 -43.0 98.0 | ||
52 | |||
53 | julia> C = cholesky(A) | ||
54 | Cholesky{Float64, Matrix{Float64}} | ||
55 | U factor: | ||
56 | 3×3 UpperTriangular{Float64, Matrix{Float64}}: | ||
57 | 2.0 6.0 -8.0 | ||
58 | ⋅ 1.0 5.0 | ||
59 | ⋅ ⋅ 3.0 | ||
60 | |||
61 | julia> C.U | ||
62 | 3×3 UpperTriangular{Float64, Matrix{Float64}}: | ||
63 | 2.0 6.0 -8.0 | ||
64 | ⋅ 1.0 5.0 | ||
65 | ⋅ ⋅ 3.0 | ||
66 | |||
67 | julia> C.L | ||
68 | 3×3 LowerTriangular{Float64, Matrix{Float64}}: | ||
69 | 2.0 ⋅ ⋅ | ||
70 | 6.0 1.0 ⋅ | ||
71 | -8.0 5.0 3.0 | ||
72 | |||
73 | julia> C.L * C.U == A | ||
74 | true | ||
75 | |||
76 | julia> l, u = C; # destructuring via iteration | ||
77 | |||
78 | julia> l == C.L && u == C.U | ||
79 | true | ||
80 | ``` | ||
81 | """ | ||
82 | struct Cholesky{T,S<:AbstractMatrix} <: Factorization{T} | ||
83 | factors::S | ||
84 | uplo::Char | ||
85 | info::BlasInt | ||
86 | |||
87 | function Cholesky{T,S}(factors, uplo, info) where {T,S<:AbstractMatrix} | ||
88 | require_one_based_indexing(factors) | ||
89 | new(factors, uplo, info) | ||
90 | end | ||
91 | end | ||
92 | Cholesky(A::AbstractMatrix{T}, uplo::Symbol, info::Integer) where {T} = | ||
93 | Cholesky{T,typeof(A)}(A, char_uplo(uplo), info) | ||
94 | Cholesky(A::AbstractMatrix{T}, uplo::AbstractChar, info::Integer) where {T} = | ||
95 | Cholesky{T,typeof(A)}(A, uplo, info) | ||
96 | Cholesky(U::UpperTriangular{T}) where {T} = Cholesky{T,typeof(U.data)}(U.data, 'U', 0) | ||
97 | Cholesky(L::LowerTriangular{T}) where {T} = Cholesky{T,typeof(L.data)}(L.data, 'L', 0) | ||
98 | |||
99 | # iteration for destructuring into components | ||
100 | Base.iterate(C::Cholesky) = (C.L, Val(:U)) | ||
101 | Base.iterate(C::Cholesky, ::Val{:U}) = (C.U, Val(:done)) | ||
102 | Base.iterate(C::Cholesky, ::Val{:done}) = nothing | ||
103 | |||
104 | |||
105 | """ | ||
106 | CholeskyPivoted | ||
107 | |||
108 | Matrix factorization type of the pivoted Cholesky factorization of a dense symmetric/Hermitian | ||
109 | positive semi-definite matrix `A`. This is the return type of [`cholesky(_, ::RowMaximum)`](@ref), | ||
110 | the corresponding matrix factorization function. | ||
111 | |||
112 | The triangular Cholesky factor can be obtained from the factorization `F::CholeskyPivoted` | ||
113 | via `F.L` and `F.U`, and the permutation via `F.p`, where `A[F.p, F.p] ≈ Ur' * Ur ≈ Lr * Lr'` | ||
114 | with `Ur = F.U[1:F.rank, :]` and `Lr = F.L[:, 1:F.rank]`, or alternatively | ||
115 | `A ≈ Up' * Up ≈ Lp * Lp'` with `Up = F.U[1:F.rank, invperm(F.p)]` and | ||
116 | `Lp = F.L[invperm(F.p), 1:F.rank]`. | ||
117 | |||
118 | The following functions are available for `CholeskyPivoted` objects: | ||
119 | [`size`](@ref), [`\\`](@ref), [`inv`](@ref), [`det`](@ref), and [`rank`](@ref). | ||
120 | |||
121 | Iterating the decomposition produces the components `L` and `U`. | ||
122 | |||
123 | # Examples | ||
124 | ```jldoctest | ||
125 | julia> X = [1.0, 2.0, 3.0, 4.0]; | ||
126 | |||
127 | julia> A = X * X'; | ||
128 | |||
129 | julia> C = cholesky(A, RowMaximum(), check = false) | ||
130 | CholeskyPivoted{Float64, Matrix{Float64}, Vector{Int64}} | ||
131 | U factor with rank 1: | ||
132 | 4×4 UpperTriangular{Float64, Matrix{Float64}}: | ||
133 | 4.0 2.0 3.0 1.0 | ||
134 | ⋅ 0.0 6.0 2.0 | ||
135 | ⋅ ⋅ 9.0 3.0 | ||
136 | ⋅ ⋅ ⋅ 1.0 | ||
137 | permutation: | ||
138 | 4-element Vector{Int64}: | ||
139 | 4 | ||
140 | 2 | ||
141 | 3 | ||
142 | 1 | ||
143 | |||
144 | julia> C.U[1:C.rank, :]' * C.U[1:C.rank, :] ≈ A[C.p, C.p] | ||
145 | true | ||
146 | |||
147 | julia> l, u = C; # destructuring via iteration | ||
148 | |||
149 | julia> l == C.L && u == C.U | ||
150 | true | ||
151 | ``` | ||
152 | """ | ||
153 | struct CholeskyPivoted{T,S<:AbstractMatrix,P<:AbstractVector{<:Integer}} <: Factorization{T} | ||
154 | factors::S | ||
155 | uplo::Char | ||
156 | piv::P | ||
157 | rank::BlasInt | ||
158 | tol::Real | ||
159 | info::BlasInt | ||
160 | |||
161 | function CholeskyPivoted{T,S,P}(factors, uplo, piv, rank, tol, info) where {T,S<:AbstractMatrix,P<:AbstractVector} | ||
162 | require_one_based_indexing(factors) | ||
163 | new{T,S,P}(factors, uplo, piv, rank, tol, info) | ||
164 | end | ||
165 | end | ||
166 | CholeskyPivoted(A::AbstractMatrix{T}, uplo::AbstractChar, piv::AbstractVector{<:Integer}, | ||
167 | rank::Integer, tol::Real, info::Integer) where T = | ||
168 | CholeskyPivoted{T,typeof(A),typeof(piv)}(A, uplo, piv, rank, tol, info) | ||
169 | # backwards-compatible constructors (remove with Julia 2.0) | ||
170 | @deprecate(CholeskyPivoted{T,S}(factors, uplo, piv, rank, tol, info) where {T,S<:AbstractMatrix}, | ||
171 | CholeskyPivoted{T,S,typeof(piv)}(factors, uplo, piv, rank, tol, info), false) | ||
172 | |||
173 | |||
174 | # iteration for destructuring into components | ||
175 | Base.iterate(C::CholeskyPivoted) = (C.L, Val(:U)) | ||
176 | Base.iterate(C::CholeskyPivoted, ::Val{:U}) = (C.U, Val(:done)) | ||
177 | Base.iterate(C::CholeskyPivoted, ::Val{:done}) = nothing | ||
178 | |||
179 | |||
180 | # make a copy that allow inplace Cholesky factorization | ||
181 | choltype(A) = promote_type(typeof(sqrt(oneunit(eltype(A)))), Float32) | ||
182 | 88 (31 %) |
88 (100 %)
samples spent calling
eigencopy_oftype
cholcopy(A::AbstractMatrix) = eigencopy_oftype(A, choltype(A))
|
|
183 | |||
184 | # _chol!. Internal methods for calling unpivoted Cholesky | ||
185 | ## BLAS/LAPACK element types | ||
186 | function _chol!(A::StridedMatrix{<:BlasFloat}, ::Type{UpperTriangular}) | ||
187 | C, info = LAPACK.potrf!('U', A) | ||
188 | return UpperTriangular(C), info | ||
189 | end | ||
190 | function _chol!(A::StridedMatrix{<:BlasFloat}, ::Type{LowerTriangular}) | ||
191 | C, info = LAPACK.potrf!('L', A) | ||
192 | return LowerTriangular(C), info | ||
193 | end | ||
194 | |||
195 | ## Non BLAS/LAPACK element types (generic) | ||
196 | function _chol!(A::AbstractMatrix, ::Type{UpperTriangular}) | ||
197 | require_one_based_indexing(A) | ||
198 | n = checksquare(A) | ||
199 | realdiag = eltype(A) <: Complex | ||
200 | @inbounds begin | ||
201 | for k = 1:n | ||
202 | Akk = realdiag ? real(A[k,k]) : A[k,k] | ||
203 | for i = 1:k - 1 | ||
204 | Akk -= realdiag ? abs2(A[i,k]) : A[i,k]'A[i,k] | ||
205 | end | ||
206 | A[k,k] = Akk | ||
207 | Akk, info = _chol!(Akk, UpperTriangular) | ||
208 | if info != 0 | ||
209 | return UpperTriangular(A), convert(BlasInt, k) | ||
210 | end | ||
211 | A[k,k] = Akk | ||
212 | AkkInv = inv(copy(Akk')) | ||
213 | for j = k + 1:n | ||
214 | for i = 1:k - 1 | ||
215 | A[k,j] -= A[i,k]'A[i,j] | ||
216 | end | ||
217 | A[k,j] = AkkInv*A[k,j] | ||
218 | end | ||
219 | end | ||
220 | end | ||
221 | return UpperTriangular(A), convert(BlasInt, 0) | ||
222 | end | ||
223 | function _chol!(A::AbstractMatrix, ::Type{LowerTriangular}) | ||
224 | require_one_based_indexing(A) | ||
225 | n = checksquare(A) | ||
226 | realdiag = eltype(A) <: Complex | ||
227 | @inbounds begin | ||
228 | for k = 1:n | ||
229 | Akk = realdiag ? real(A[k,k]) : A[k,k] | ||
230 | for i = 1:k - 1 | ||
231 | Akk -= realdiag ? abs2(A[k,i]) : A[k,i]*A[k,i]' | ||
232 | end | ||
233 | A[k,k] = Akk | ||
234 | Akk, info = _chol!(Akk, LowerTriangular) | ||
235 | if info != 0 | ||
236 | return LowerTriangular(A), convert(BlasInt, k) | ||
237 | end | ||
238 | A[k,k] = Akk | ||
239 | AkkInv = inv(Akk) | ||
240 | for j = 1:k - 1 | ||
241 | @simd for i = k + 1:n | ||
242 | A[i,k] -= A[i,j]*A[k,j]' | ||
243 | end | ||
244 | end | ||
245 | for i = k + 1:n | ||
246 | A[i,k] *= AkkInv' | ||
247 | end | ||
248 | end | ||
249 | end | ||
250 | return LowerTriangular(A), convert(BlasInt, 0) | ||
251 | end | ||
252 | |||
253 | ## Numbers | ||
254 | function _chol!(x::Number, _) | ||
255 | rx = real(x) | ||
256 | iszero(rx) && return (rx, convert(BlasInt, 1)) | ||
257 | rxr = sqrt(abs(rx)) | ||
258 | rval = convert(promote_type(typeof(x), typeof(rxr)), rxr) | ||
259 | return (rval, convert(BlasInt, rx != abs(x))) | ||
260 | end | ||
261 | |||
262 | ## for StridedMatrices, check that matrix is symmetric/Hermitian | ||
263 | |||
264 | # cholesky!. Destructive methods for computing Cholesky factorization of real symmetric | ||
265 | # or Hermitian matrix | ||
266 | ## No pivoting (default) | ||
267 | function cholesky!(A::RealHermSymComplexHerm, ::NoPivot = NoPivot(); check::Bool = true) | ||
268 | C, info = _chol!(A.data, A.uplo == 'U' ? UpperTriangular : LowerTriangular) | ||
269 | check && checkpositivedefinite(info) | ||
270 | return Cholesky(C.data, A.uplo, info) | ||
271 | end | ||
272 | |||
273 | ### for AbstractMatrix, check that matrix is symmetric/Hermitian | ||
274 | """ | ||
275 | cholesky!(A::AbstractMatrix, NoPivot(); check = true) -> Cholesky | ||
276 | |||
277 | The same as [`cholesky`](@ref), but saves space by overwriting the input `A`, | ||
278 | instead of creating a copy. An [`InexactError`](@ref) exception is thrown if | ||
279 | the factorization produces a number not representable by the element type of | ||
280 | `A`, e.g. for integer types. | ||
281 | |||
282 | # Examples | ||
283 | ```jldoctest | ||
284 | julia> A = [1 2; 2 50] | ||
285 | 2×2 Matrix{Int64}: | ||
286 | 1 2 | ||
287 | 2 50 | ||
288 | |||
289 | julia> cholesky!(A) | ||
290 | ERROR: InexactError: Int64(6.782329983125268) | ||
291 | Stacktrace: | ||
292 | [...] | ||
293 | ``` | ||
294 | """ | ||
295 | function cholesky!(A::AbstractMatrix, ::NoPivot = NoPivot(); check::Bool = true) | ||
296 | checksquare(A) | ||
297 | if !ishermitian(A) # return with info = -1 if not Hermitian | ||
298 | check && checkpositivedefinite(-1) | ||
299 | return Cholesky(A, 'U', convert(BlasInt, -1)) | ||
300 | else | ||
301 | return cholesky!(Hermitian(A), NoPivot(); check = check) | ||
302 | end | ||
303 | end | ||
304 | @deprecate cholesky!(A::StridedMatrix, ::Val{false}; check::Bool = true) cholesky!(A, NoPivot(); check) false | ||
305 | @deprecate cholesky!(A::RealHermSymComplexHerm, ::Val{false}; check::Bool = true) cholesky!(A, NoPivot(); check) false | ||
306 | |||
307 | ## With pivoting | ||
308 | ### BLAS/LAPACK element types | ||
309 | function cholesky!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}, | ||
310 | ::RowMaximum; tol = 0.0, check::Bool = true) | ||
311 | AA, piv, rank, info = LAPACK.pstrf!(A.uplo, A.data, tol) | ||
312 | C = CholeskyPivoted{eltype(AA),typeof(AA),typeof(piv)}(AA, A.uplo, piv, rank, tol, info) | ||
313 | check && chkfullrank(C) | ||
314 | return C | ||
315 | end | ||
316 | @deprecate cholesky!(A::RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix}, ::Val{true}; kwargs...) cholesky!(A, RowMaximum(); kwargs...) false | ||
317 | |||
318 | ### Non BLAS/LAPACK element types (generic). Since generic fallback for pivoted Cholesky | ||
319 | ### is not implemented yet we throw an error | ||
320 | cholesky!(A::RealHermSymComplexHerm{<:Real}, ::RowMaximum; tol = 0.0, check::Bool = true) = | ||
321 | throw(ArgumentError("generic pivoted Cholesky factorization is not implemented yet")) | ||
322 | @deprecate cholesky!(A::RealHermSymComplexHerm{<:Real}, ::Val{true}; kwargs...) cholesky!(A, RowMaximum(); kwargs...) false | ||
323 | |||
324 | ### for AbstractMatrix, check that matrix is symmetric/Hermitian | ||
325 | """ | ||
326 | cholesky!(A::AbstractMatrix, RowMaximum(); tol = 0.0, check = true) -> CholeskyPivoted | ||
327 | |||
328 | The same as [`cholesky`](@ref), but saves space by overwriting the input `A`, | ||
329 | instead of creating a copy. An [`InexactError`](@ref) exception is thrown if the | ||
330 | factorization produces a number not representable by the element type of `A`, | ||
331 | e.g. for integer types. | ||
332 | """ | ||
333 | function cholesky!(A::AbstractMatrix, ::RowMaximum; tol = 0.0, check::Bool = true) | ||
334 | checksquare(A) | ||
335 | if !ishermitian(A) | ||
336 | C = CholeskyPivoted(A, 'U', Vector{BlasInt}(),convert(BlasInt, 1), | ||
337 | tol, convert(BlasInt, -1)) | ||
338 | check && chkfullrank(C) | ||
339 | return C | ||
340 | else | ||
341 | return cholesky!(Hermitian(A), RowMaximum(); tol = tol, check = check) | ||
342 | end | ||
343 | end | ||
344 | @deprecate cholesky!(A::StridedMatrix, ::Val{true}; kwargs...) cholesky!(A, RowMaximum(); kwargs...) false | ||
345 | |||
346 | # cholesky. Non-destructive methods for computing Cholesky factorization of real symmetric | ||
347 | # or Hermitian matrix | ||
348 | ## No pivoting (default) | ||
349 | """ | ||
350 | cholesky(A, NoPivot(); check = true) -> Cholesky | ||
351 | |||
352 | Compute the Cholesky factorization of a dense symmetric positive definite matrix `A` | ||
353 | and return a [`Cholesky`](@ref) factorization. The matrix `A` can either be a [`Symmetric`](@ref) or [`Hermitian`](@ref) | ||
354 | [`AbstractMatrix`](@ref) or a *perfectly* symmetric or Hermitian `AbstractMatrix`. | ||
355 | |||
356 | The triangular Cholesky factor can be obtained from the factorization `F` via `F.L` and `F.U`, | ||
357 | where `A ≈ F.U' * F.U ≈ F.L * F.L'`. | ||
358 | |||
359 | The following functions are available for `Cholesky` objects: [`size`](@ref), [`\\`](@ref), | ||
360 | [`inv`](@ref), [`det`](@ref), [`logdet`](@ref) and [`isposdef`](@ref). | ||
361 | |||
362 | If you have a matrix `A` that is slightly non-Hermitian due to roundoff errors in its construction, | ||
363 | wrap it in `Hermitian(A)` before passing it to `cholesky` in order to treat it as perfectly Hermitian. | ||
364 | |||
365 | When `check = true`, an error is thrown if the decomposition fails. | ||
366 | When `check = false`, responsibility for checking the decomposition's | ||
367 | validity (via [`issuccess`](@ref)) lies with the user. | ||
368 | |||
369 | # Examples | ||
370 | ```jldoctest | ||
371 | julia> A = [4. 12. -16.; 12. 37. -43.; -16. -43. 98.] | ||
372 | 3×3 Matrix{Float64}: | ||
373 | 4.0 12.0 -16.0 | ||
374 | 12.0 37.0 -43.0 | ||
375 | -16.0 -43.0 98.0 | ||
376 | |||
377 | julia> C = cholesky(A) | ||
378 | Cholesky{Float64, Matrix{Float64}} | ||
379 | U factor: | ||
380 | 3×3 UpperTriangular{Float64, Matrix{Float64}}: | ||
381 | 2.0 6.0 -8.0 | ||
382 | ⋅ 1.0 5.0 | ||
383 | ⋅ ⋅ 3.0 | ||
384 | |||
385 | julia> C.U | ||
386 | 3×3 UpperTriangular{Float64, Matrix{Float64}}: | ||
387 | 2.0 6.0 -8.0 | ||
388 | ⋅ 1.0 5.0 | ||
389 | ⋅ ⋅ 3.0 | ||
390 | |||
391 | julia> C.L | ||
392 | 3×3 LowerTriangular{Float64, Matrix{Float64}}: | ||
393 | 2.0 ⋅ ⋅ | ||
394 | 6.0 1.0 ⋅ | ||
395 | -8.0 5.0 3.0 | ||
396 | |||
397 | julia> C.L * C.U == A | ||
398 | true | ||
399 | ``` | ||
400 | """ | ||
401 | 176 (61 %) |
176 (61 %)
samples spent in cholesky
cholesky(A::AbstractMatrix, ::NoPivot=NoPivot(); check::Bool = true) =
88 (50 %) (incl.) when called from cholesky_instance line 453 88 (50 %) (incl.) when called from cholesky line 401 |
|
402 | cholesky!(cholcopy(A); check) | ||
403 | @deprecate cholesky(A::Union{StridedMatrix,RealHermSymComplexHerm{<:Real,<:StridedMatrix}}, ::Val{false}; check::Bool = true) cholesky(A, NoPivot(); check) false | ||
404 | |||
405 | function cholesky(A::AbstractMatrix{Float16}, ::NoPivot=NoPivot(); check::Bool = true) | ||
406 | X = cholesky!(cholcopy(A); check = check) | ||
407 | return Cholesky{Float16}(X) | ||
408 | end | ||
409 | @deprecate cholesky(A::Union{StridedMatrix{Float16},RealHermSymComplexHerm{Float16,<:StridedMatrix}}, ::Val{false}; check::Bool = true) cholesky(A, NoPivot(); check) false | ||
410 | |||
411 | ## With pivoting | ||
412 | """ | ||
413 | cholesky(A, RowMaximum(); tol = 0.0, check = true) -> CholeskyPivoted | ||
414 | |||
415 | Compute the pivoted Cholesky factorization of a dense symmetric positive semi-definite matrix `A` | ||
416 | and return a [`CholeskyPivoted`](@ref) factorization. The matrix `A` can either be a [`Symmetric`](@ref) | ||
417 | or [`Hermitian`](@ref) [`AbstractMatrix`](@ref) or a *perfectly* symmetric or Hermitian `AbstractMatrix`. | ||
418 | |||
419 | The triangular Cholesky factor can be obtained from the factorization `F` via `F.L` and `F.U`, | ||
420 | and the permutation via `F.p`, where `A[F.p, F.p] ≈ Ur' * Ur ≈ Lr * Lr'` with `Ur = F.U[1:F.rank, :]` | ||
421 | and `Lr = F.L[:, 1:F.rank]`, or alternatively `A ≈ Up' * Up ≈ Lp * Lp'` with | ||
422 | `Up = F.U[1:F.rank, invperm(F.p)]` and `Lp = F.L[invperm(F.p), 1:F.rank]`. | ||
423 | |||
424 | The following functions are available for `CholeskyPivoted` objects: | ||
425 | [`size`](@ref), [`\\`](@ref), [`inv`](@ref), [`det`](@ref), and [`rank`](@ref). | ||
426 | |||
427 | The argument `tol` determines the tolerance for determining the rank. | ||
428 | For negative values, the tolerance is the machine precision. | ||
429 | |||
430 | If you have a matrix `A` that is slightly non-Hermitian due to roundoff errors in its construction, | ||
431 | wrap it in `Hermitian(A)` before passing it to `cholesky` in order to treat it as perfectly Hermitian. | ||
432 | |||
433 | When `check = true`, an error is thrown if the decomposition fails. | ||
434 | When `check = false`, responsibility for checking the decomposition's | ||
435 | validity (via [`issuccess`](@ref)) lies with the user. | ||
436 | |||
437 | # Examples | ||
438 | ```jldoctest | ||
439 | julia> X = [1.0, 2.0, 3.0, 4.0]; | ||
440 | |||
441 | julia> A = X * X'; | ||
442 | |||
443 | julia> C = cholesky(A, RowMaximum(), check = false) | ||
444 | CholeskyPivoted{Float64, Matrix{Float64}, Vector{Int64}} | ||
445 | U factor with rank 1: | ||
446 | 4×4 UpperTriangular{Float64, Matrix{Float64}}: | ||
447 | 4.0 2.0 3.0 1.0 | ||
448 | ⋅ 0.0 6.0 2.0 | ||
449 | ⋅ ⋅ 9.0 3.0 | ||
450 | ⋅ ⋅ ⋅ 1.0 | ||
451 | permutation: | ||
452 | 4-element Vector{Int64}: | ||
453 | 4 | ||
454 | 2 | ||
455 | 3 | ||
456 | 1 | ||
457 | |||
458 | julia> C.U[1:C.rank, :]' * C.U[1:C.rank, :] ≈ A[C.p, C.p] | ||
459 | true | ||
460 | |||
461 | julia> l, u = C; # destructuring via iteration | ||
462 | |||
463 | julia> l == C.L && u == C.U | ||
464 | true | ||
465 | ``` | ||
466 | """ | ||
467 | cholesky(A::AbstractMatrix, ::RowMaximum; tol = 0.0, check::Bool = true) = | ||
468 | cholesky!(cholcopy(A), RowMaximum(); tol, check) | ||
469 | @deprecate cholesky(A::Union{StridedMatrix,RealHermSymComplexHerm{<:Real,<:StridedMatrix}}, ::Val{true}; tol = 0.0, check::Bool = true) cholesky(A, RowMaximum(); tol, check) false | ||
470 | |||
471 | function cholesky(A::AbstractMatrix{Float16}, ::RowMaximum; tol = 0.0, check::Bool = true) | ||
472 | X = cholesky!(cholcopy(A), RowMaximum(); tol, check) | ||
473 | return CholeskyPivoted{Float16}(X) | ||
474 | end | ||
475 | |||
476 | ## Number | ||
477 | function cholesky(x::Number, uplo::Symbol=:U) | ||
478 | C, info = _chol!(x, uplo) | ||
479 | xf = fill(C, 1, 1) | ||
480 | Cholesky(xf, uplo, info) | ||
481 | end | ||
482 | |||
483 | |||
484 | function Cholesky{T}(C::Cholesky) where T | ||
485 | Cnew = convert(AbstractMatrix{T}, C.factors) | ||
486 | Cholesky{T, typeof(Cnew)}(Cnew, C.uplo, C.info) | ||
487 | end | ||
488 | Factorization{T}(C::Cholesky{T}) where {T} = C | ||
489 | Factorization{T}(C::Cholesky) where {T} = Cholesky{T}(C) | ||
490 | CholeskyPivoted{T}(C::CholeskyPivoted{T}) where {T} = C | ||
491 | CholeskyPivoted{T}(C::CholeskyPivoted) where {T} = | ||
492 | CholeskyPivoted(AbstractMatrix{T}(C.factors),C.uplo,C.piv,C.rank,C.tol,C.info) | ||
493 | Factorization{T}(C::CholeskyPivoted{T}) where {T} = C | ||
494 | Factorization{T}(C::CholeskyPivoted) where {T} = CholeskyPivoted{T}(C) | ||
495 | |||
496 | AbstractMatrix(C::Cholesky) = C.uplo == 'U' ? C.U'C.U : C.L*C.L' | ||
497 | AbstractArray(C::Cholesky) = AbstractMatrix(C) | ||
498 | Matrix(C::Cholesky) = Array(AbstractArray(C)) | ||
499 | Array(C::Cholesky) = Matrix(C) | ||
500 | |||
501 | function AbstractMatrix(F::CholeskyPivoted) | ||
502 | ip = invperm(F.p) | ||
503 | U = F.U[1:F.rank,ip] | ||
504 | U'U | ||
505 | end | ||
506 | AbstractArray(F::CholeskyPivoted) = AbstractMatrix(F) | ||
507 | Matrix(F::CholeskyPivoted) = Array(AbstractArray(F)) | ||
508 | Array(F::CholeskyPivoted) = Matrix(F) | ||
509 | |||
510 | copy(C::Cholesky) = Cholesky(copy(C.factors), C.uplo, C.info) | ||
511 | copy(C::CholeskyPivoted) = CholeskyPivoted(copy(C.factors), C.uplo, C.piv, C.rank, C.tol, C.info) | ||
512 | |||
513 | size(C::Union{Cholesky, CholeskyPivoted}) = size(C.factors) | ||
514 | size(C::Union{Cholesky, CholeskyPivoted}, d::Integer) = size(C.factors, d) | ||
515 | |||
516 | function getproperty(C::Cholesky, d::Symbol) | ||
517 | Cfactors = getfield(C, :factors) | ||
518 | Cuplo = getfield(C, :uplo) | ||
519 | if d === :U | ||
520 | return UpperTriangular(Cuplo === char_uplo(d) ? Cfactors : copy(Cfactors')) | ||
521 | elseif d === :L | ||
522 | return LowerTriangular(Cuplo === char_uplo(d) ? Cfactors : copy(Cfactors')) | ||
523 | elseif d === :UL | ||
524 | return (Cuplo === 'U' ? UpperTriangular(Cfactors) : LowerTriangular(Cfactors)) | ||
525 | else | ||
526 | return getfield(C, d) | ||
527 | end | ||
528 | end | ||
529 | Base.propertynames(F::Cholesky, private::Bool=false) = | ||
530 | (:U, :L, :UL, (private ? fieldnames(typeof(F)) : ())...) | ||
531 | |||
532 | function getproperty(C::CholeskyPivoted{T}, d::Symbol) where {T} | ||
533 | Cfactors = getfield(C, :factors) | ||
534 | Cuplo = getfield(C, :uplo) | ||
535 | if d === :U | ||
536 | return UpperTriangular(sym_uplo(Cuplo) == d ? Cfactors : copy(Cfactors')) | ||
537 | elseif d === :L | ||
538 | return LowerTriangular(sym_uplo(Cuplo) == d ? Cfactors : copy(Cfactors')) | ||
539 | elseif d === :p | ||
540 | return getfield(C, :piv) | ||
541 | elseif d === :P | ||
542 | n = size(C, 1) | ||
543 | P = zeros(T, n, n) | ||
544 | for i = 1:n | ||
545 | P[getfield(C, :piv)[i], i] = one(T) | ||
546 | end | ||
547 | return P | ||
548 | else | ||
549 | return getfield(C, d) | ||
550 | end | ||
551 | end | ||
552 | Base.propertynames(F::CholeskyPivoted, private::Bool=false) = | ||
553 | (:U, :L, :p, :P, (private ? fieldnames(typeof(F)) : ())...) | ||
554 | |||
555 | issuccess(C::Union{Cholesky,CholeskyPivoted}) = C.info == 0 | ||
556 | |||
557 | adjoint(C::Union{Cholesky,CholeskyPivoted}) = C | ||
558 | |||
559 | function show(io::IO, mime::MIME{Symbol("text/plain")}, C::Cholesky) | ||
560 | if issuccess(C) | ||
561 | summary(io, C); println(io) | ||
562 | println(io, "$(C.uplo) factor:") | ||
563 | show(io, mime, C.UL) | ||
564 | else | ||
565 | print(io, "Failed factorization of type $(typeof(C))") | ||
566 | end | ||
567 | end | ||
568 | |||
569 | function show(io::IO, mime::MIME{Symbol("text/plain")}, C::CholeskyPivoted) | ||
570 | summary(io, C); println(io) | ||
571 | println(io, "$(C.uplo) factor with rank $(rank(C)):") | ||
572 | show(io, mime, C.uplo == 'U' ? C.U : C.L) | ||
573 | println(io, "\npermutation:") | ||
574 | show(io, mime, C.p) | ||
575 | end | ||
576 | |||
577 | ldiv!(C::Cholesky{T,<:StridedMatrix}, B::StridedVecOrMat{T}) where {T<:BlasFloat} = | ||
578 | LAPACK.potrs!(C.uplo, C.factors, B) | ||
579 | |||
580 | function ldiv!(C::Cholesky, B::AbstractVecOrMat) | ||
581 | if C.uplo == 'L' | ||
582 | return ldiv!(adjoint(LowerTriangular(C.factors)), ldiv!(LowerTriangular(C.factors), B)) | ||
583 | else | ||
584 | return ldiv!(UpperTriangular(C.factors), ldiv!(adjoint(UpperTriangular(C.factors)), B)) | ||
585 | end | ||
586 | end | ||
587 | |||
588 | function ldiv!(C::CholeskyPivoted{T,<:StridedMatrix}, B::StridedVector{T}) where T<:BlasFloat | ||
589 | invpermute!(LAPACK.potrs!(C.uplo, C.factors, permute!(B, C.piv)), C.piv) | ||
590 | end | ||
591 | function ldiv!(C::CholeskyPivoted{T,<:StridedMatrix}, B::StridedMatrix{T}) where T<:BlasFloat | ||
592 | n = size(C, 1) | ||
593 | for i=1:size(B, 2) | ||
594 | permute!(view(B, 1:n, i), C.piv) | ||
595 | end | ||
596 | LAPACK.potrs!(C.uplo, C.factors, B) | ||
597 | for i=1:size(B, 2) | ||
598 | invpermute!(view(B, 1:n, i), C.piv) | ||
599 | end | ||
600 | B | ||
601 | end | ||
602 | |||
603 | function ldiv!(C::CholeskyPivoted, B::AbstractVector) | ||
604 | if C.uplo == 'L' | ||
605 | ldiv!(adjoint(LowerTriangular(C.factors)), | ||
606 | ldiv!(LowerTriangular(C.factors), permute!(B, C.piv))) | ||
607 | else | ||
608 | ldiv!(UpperTriangular(C.factors), | ||
609 | ldiv!(adjoint(UpperTriangular(C.factors)), permute!(B, C.piv))) | ||
610 | end | ||
611 | invpermute!(B, C.piv) | ||
612 | end | ||
613 | |||
614 | function ldiv!(C::CholeskyPivoted, B::AbstractMatrix) | ||
615 | n = size(C, 1) | ||
616 | for i in 1:size(B, 2) | ||
617 | permute!(view(B, 1:n, i), C.piv) | ||
618 | end | ||
619 | if C.uplo == 'L' | ||
620 | ldiv!(adjoint(LowerTriangular(C.factors)), | ||
621 | ldiv!(LowerTriangular(C.factors), B)) | ||
622 | else | ||
623 | ldiv!(UpperTriangular(C.factors), | ||
624 | ldiv!(adjoint(UpperTriangular(C.factors)), B)) | ||
625 | end | ||
626 | for i in 1:size(B, 2) | ||
627 | invpermute!(view(B, 1:n, i), C.piv) | ||
628 | end | ||
629 | B | ||
630 | end | ||
631 | |||
632 | function rdiv!(B::AbstractMatrix, C::Cholesky) | ||
633 | if C.uplo == 'L' | ||
634 | return rdiv!(rdiv!(B, adjoint(LowerTriangular(C.factors))), LowerTriangular(C.factors)) | ||
635 | else | ||
636 | return rdiv!(rdiv!(B, UpperTriangular(C.factors)), adjoint(UpperTriangular(C.factors))) | ||
637 | end | ||
638 | end | ||
639 | |||
640 | function LinearAlgebra.rdiv!(B::AbstractMatrix, C::CholeskyPivoted) | ||
641 | n = size(C, 2) | ||
642 | for i in 1:size(B, 1) | ||
643 | permute!(view(B, i, 1:n), C.piv) | ||
644 | end | ||
645 | if C.uplo == 'L' | ||
646 | rdiv!(rdiv!(B, adjoint(LowerTriangular(C.factors))), | ||
647 | LowerTriangular(C.factors)) | ||
648 | else | ||
649 | rdiv!(rdiv!(B, UpperTriangular(C.factors)), | ||
650 | adjoint(UpperTriangular(C.factors))) | ||
651 | end | ||
652 | for i in 1:size(B, 1) | ||
653 | invpermute!(view(B, i, 1:n), C.piv) | ||
654 | end | ||
655 | B | ||
656 | end | ||
657 | |||
658 | isposdef(C::Union{Cholesky,CholeskyPivoted}) = C.info == 0 | ||
659 | |||
660 | function det(C::Cholesky) | ||
661 | dd = one(real(eltype(C))) | ||
662 | @inbounds for i in 1:size(C.factors,1) | ||
663 | dd *= real(C.factors[i,i])^2 | ||
664 | end | ||
665 | return dd | ||
666 | end | ||
667 | |||
668 | function logdet(C::Cholesky) | ||
669 | dd = zero(real(eltype(C))) | ||
670 | @inbounds for i in 1:size(C.factors,1) | ||
671 | dd += log(real(C.factors[i,i])) | ||
672 | end | ||
673 | dd + dd # instead of 2.0dd which can change the type | ||
674 | end | ||
675 | |||
676 | function det(C::CholeskyPivoted) | ||
677 | if C.rank < size(C.factors, 1) | ||
678 | return zero(real(eltype(C))) | ||
679 | else | ||
680 | dd = one(real(eltype(C))) | ||
681 | for i in 1:size(C.factors,1) | ||
682 | dd *= real(C.factors[i,i])^2 | ||
683 | end | ||
684 | return dd | ||
685 | end | ||
686 | end | ||
687 | |||
688 | function logdet(C::CholeskyPivoted) | ||
689 | if C.rank < size(C.factors, 1) | ||
690 | return real(eltype(C))(-Inf) | ||
691 | else | ||
692 | dd = zero(real(eltype(C))) | ||
693 | for i in 1:size(C.factors,1) | ||
694 | dd += log(real(C.factors[i,i])) | ||
695 | end | ||
696 | return dd + dd # instead of 2.0dd which can change the type | ||
697 | end | ||
698 | end | ||
699 | |||
700 | logabsdet(C::Union{Cholesky, CholeskyPivoted}) = logdet(C), one(eltype(C)) # since C is p.s.d. | ||
701 | |||
702 | inv!(C::Cholesky{<:BlasFloat,<:StridedMatrix}) = | ||
703 | copytri!(LAPACK.potri!(C.uplo, C.factors), C.uplo, true) | ||
704 | |||
705 | inv(C::Cholesky{<:BlasFloat,<:StridedMatrix}) = inv!(copy(C)) | ||
706 | |||
707 | function inv(C::CholeskyPivoted{<:BlasFloat,<:StridedMatrix}) | ||
708 | ipiv = invperm(C.piv) | ||
709 | copytri!(LAPACK.potri!(C.uplo, copy(C.factors)), C.uplo, true)[ipiv, ipiv] | ||
710 | end | ||
711 | |||
712 | function chkfullrank(C::CholeskyPivoted) | ||
713 | if C.rank < size(C.factors, 1) | ||
714 | throw(RankDeficientException(C.info)) | ||
715 | end | ||
716 | end | ||
717 | |||
718 | rank(C::CholeskyPivoted) = C.rank | ||
719 | |||
720 | """ | ||
721 | lowrankupdate!(C::Cholesky, v::AbstractVector) -> CC::Cholesky | ||
722 | |||
723 | Update a Cholesky factorization `C` with the vector `v`. If `A = C.U'C.U` then | ||
724 | `CC = cholesky(C.U'C.U + v*v')` but the computation of `CC` only uses `O(n^2)` | ||
725 | operations. The input factorization `C` is updated in place such that on exit `C == CC`. | ||
726 | The vector `v` is destroyed during the computation. | ||
727 | """ | ||
728 | function lowrankupdate!(C::Cholesky, v::AbstractVector) | ||
729 | A = C.factors | ||
730 | n = length(v) | ||
731 | if size(C, 1) != n | ||
732 | throw(DimensionMismatch("updating vector must fit size of factorization")) | ||
733 | end | ||
734 | if C.uplo == 'U' | ||
735 | conj!(v) | ||
736 | end | ||
737 | |||
738 | for i = 1:n | ||
739 | |||
740 | # Compute Givens rotation | ||
741 | c, s, r = givensAlgorithm(A[i,i], v[i]) | ||
742 | |||
743 | # Store new diagonal element | ||
744 | A[i,i] = r | ||
745 | |||
746 | # Update remaining elements in row/column | ||
747 | if C.uplo == 'U' | ||
748 | for j = i + 1:n | ||
749 | Aij = A[i,j] | ||
750 | vj = v[j] | ||
751 | A[i,j] = c*Aij + s*vj | ||
752 | v[j] = -s'*Aij + c*vj | ||
753 | end | ||
754 | else | ||
755 | for j = i + 1:n | ||
756 | Aji = A[j,i] | ||
757 | vj = v[j] | ||
758 | A[j,i] = c*Aji + s*vj | ||
759 | v[j] = -s'*Aji + c*vj | ||
760 | end | ||
761 | end | ||
762 | end | ||
763 | return C | ||
764 | end | ||
765 | |||
766 | """ | ||
767 | lowrankdowndate!(C::Cholesky, v::AbstractVector) -> CC::Cholesky | ||
768 | |||
769 | Downdate a Cholesky factorization `C` with the vector `v`. If `A = C.U'C.U` then | ||
770 | `CC = cholesky(C.U'C.U - v*v')` but the computation of `CC` only uses `O(n^2)` | ||
771 | operations. The input factorization `C` is updated in place such that on exit `C == CC`. | ||
772 | The vector `v` is destroyed during the computation. | ||
773 | """ | ||
774 | function lowrankdowndate!(C::Cholesky, v::AbstractVector) | ||
775 | A = C.factors | ||
776 | n = length(v) | ||
777 | if size(C, 1) != n | ||
778 | throw(DimensionMismatch("updating vector must fit size of factorization")) | ||
779 | end | ||
780 | if C.uplo == 'U' | ||
781 | conj!(v) | ||
782 | end | ||
783 | |||
784 | for i = 1:n | ||
785 | |||
786 | Aii = A[i,i] | ||
787 | |||
788 | # Compute Givens rotation | ||
789 | s = conj(v[i]/Aii) | ||
790 | s2 = abs2(s) | ||
791 | if s2 > 1 | ||
792 | throw(LinearAlgebra.PosDefException(i)) | ||
793 | end | ||
794 | c = sqrt(1 - abs2(s)) | ||
795 | |||
796 | # Store new diagonal element | ||
797 | A[i,i] = c*Aii | ||
798 | |||
799 | # Update remaining elements in row/column | ||
800 | if C.uplo == 'U' | ||
801 | for j = i + 1:n | ||
802 | vj = v[j] | ||
803 | Aij = (A[i,j] - s*vj)/c | ||
804 | A[i,j] = Aij | ||
805 | v[j] = -s'*Aij + c*vj | ||
806 | end | ||
807 | else | ||
808 | for j = i + 1:n | ||
809 | vj = v[j] | ||
810 | Aji = (A[j,i] - s*vj)/c | ||
811 | A[j,i] = Aji | ||
812 | v[j] = -s'*Aji + c*vj | ||
813 | end | ||
814 | end | ||
815 | end | ||
816 | return C | ||
817 | end | ||
818 | |||
819 | """ | ||
820 | lowrankupdate(C::Cholesky, v::AbstractVector) -> CC::Cholesky | ||
821 | |||
822 | Update a Cholesky factorization `C` with the vector `v`. If `A = C.U'C.U` | ||
823 | then `CC = cholesky(C.U'C.U + v*v')` but the computation of `CC` only uses | ||
824 | `O(n^2)` operations. | ||
825 | """ | ||
826 | lowrankupdate(C::Cholesky, v::AbstractVector) = lowrankupdate!(copy(C), copy(v)) | ||
827 | |||
828 | """ | ||
829 | lowrankdowndate(C::Cholesky, v::AbstractVector) -> CC::Cholesky | ||
830 | |||
831 | Downdate a Cholesky factorization `C` with the vector `v`. If `A = C.U'C.U` | ||
832 | then `CC = cholesky(C.U'C.U - v*v')` but the computation of `CC` only uses | ||
833 | `O(n^2)` operations. | ||
834 | """ | ||
835 | lowrankdowndate(C::Cholesky, v::AbstractVector) = lowrankdowndate!(copy(C), copy(v)) |