Line | Exclusive | Inclusive | Code |
---|---|---|---|
1 | # This file is a part of Julia. License is MIT: https://julialang.org/license | ||
2 | |||
3 | ## Matrix factorizations and decompositions | ||
4 | """ | ||
5 | LinearAlgebra.Factorization | ||
6 | |||
7 | Abstract type for [matrix factorizations](https://en.wikipedia.org/wiki/Matrix_decomposition) | ||
8 | a.k.a. matrix decompositions. | ||
9 | See [online documentation](@ref man-linalg-factorizations) for a list of available | ||
10 | matrix factorizations. | ||
11 | """ | ||
12 | abstract type Factorization{T} end | ||
13 | |||
14 | """ | ||
15 | AdjointFactorization | ||
16 | |||
17 | Lazy wrapper type for the adjoint of the underlying `Factorization` object. Usually, the | ||
18 | `AdjointFactorization` constructor should not be called directly, use | ||
19 | [`adjoint(:: Factorization)`](@ref) instead. | ||
20 | """ | ||
21 | struct AdjointFactorization{T,S<:Factorization} <: Factorization{T} | ||
22 | parent::S | ||
23 | end | ||
24 | AdjointFactorization(F::Factorization) = | ||
25 | AdjointFactorization{Base.promote_op(adjoint,eltype(F)),typeof(F)}(F) | ||
26 | |||
27 | """ | ||
28 | TransposeFactorization | ||
29 | |||
30 | Lazy wrapper type for the transpose of the underlying `Factorization` object. Usually, the | ||
31 | `TransposeFactorization` constructor should not be called directly, use | ||
32 | [`transpose(:: Factorization)`](@ref) instead. | ||
33 | """ | ||
34 | struct TransposeFactorization{T,S<:Factorization} <: Factorization{T} | ||
35 | parent::S | ||
36 | end | ||
37 | TransposeFactorization(F::Factorization) = | ||
38 | TransposeFactorization{Base.promote_op(adjoint,eltype(F)),typeof(F)}(F) | ||
39 | |||
40 | eltype(::Type{<:Factorization{T}}) where {T} = T | ||
41 | size(F::AdjointFactorization) = reverse(size(parent(F))) | ||
42 | size(F::TransposeFactorization) = reverse(size(parent(F))) | ||
43 | size(F::Union{AdjointFactorization,TransposeFactorization}, d::Integer) = d in (1, 2) ? size(F)[d] : 1 | ||
44 | parent(F::Union{AdjointFactorization,TransposeFactorization}) = F.parent | ||
45 | |||
46 | """ | ||
47 | adjoint(F::Factorization) | ||
48 | |||
49 | Lazy adjoint of the factorization `F`. By default, returns an | ||
50 | [`AdjointFactorization`](@ref) wrapper. | ||
51 | """ | ||
52 | adjoint(F::Factorization) = AdjointFactorization(F) | ||
53 | """ | ||
54 | transpose(F::Factorization) | ||
55 | |||
56 | Lazy transpose of the factorization `F`. By default, returns a [`TransposeFactorization`](@ref), | ||
57 | except for `Factorization`s with real `eltype`, in which case returns an [`AdjointFactorization`](@ref). | ||
58 | """ | ||
59 | transpose(F::Factorization) = TransposeFactorization(F) | ||
60 | transpose(F::Factorization{<:Real}) = AdjointFactorization(F) | ||
61 | adjoint(F::AdjointFactorization) = F.parent | ||
62 | transpose(F::TransposeFactorization) = F.parent | ||
63 | transpose(F::AdjointFactorization{<:Real}) = F.parent | ||
64 | conj(A::TransposeFactorization) = adjoint(A.parent) | ||
65 | conj(A::AdjointFactorization) = transpose(A.parent) | ||
66 | |||
67 | checkpositivedefinite(info) = info == 0 || throw(PosDefException(info)) | ||
68 | checknonsingular(info, ::RowMaximum) = info == 0 || throw(SingularException(info)) | ||
69 | checknonsingular(info, ::RowNonZero) = info == 0 || throw(SingularException(info)) | ||
70 | checknonsingular(info, ::NoPivot) = info == 0 || throw(ZeroPivotException(info)) | ||
71 | checknonsingular(info) = checknonsingular(info, RowMaximum()) | ||
72 | |||
73 | """ | ||
74 | issuccess(F::Factorization) | ||
75 | |||
76 | Test that a factorization of a matrix succeeded. | ||
77 | |||
78 | !!! compat "Julia 1.6" | ||
79 | `issuccess(::CholeskyPivoted)` requires Julia 1.6 or later. | ||
80 | |||
81 | ```jldoctest | ||
82 | julia> F = cholesky([1 0; 0 1]); | ||
83 | |||
84 | julia> issuccess(F) | ||
85 | true | ||
86 | |||
87 | julia> F = lu([1 0; 0 0]; check = false); | ||
88 | |||
89 | julia> issuccess(F) | ||
90 | false | ||
91 | ``` | ||
92 | """ | ||
93 | issuccess(F::Factorization) | ||
94 | |||
95 | function logdet(F::Factorization) | ||
96 | d, s = logabsdet(F) | ||
97 | return d + log(s) | ||
98 | end | ||
99 | |||
100 | function det(F::Factorization) | ||
101 | d, s = logabsdet(F) | ||
102 | return exp(d)*s | ||
103 | end | ||
104 | |||
105 | convert(::Type{T}, f::T) where {T<:Factorization} = f | ||
106 | convert(::Type{T}, f::Factorization) where {T<:Factorization} = T(f)::T | ||
107 | |||
108 | convert(::Type{T}, f::Factorization) where {T<:AbstractArray} = T(f)::T | ||
109 | |||
110 | ### General promotion rules | ||
111 | Factorization{T}(F::Factorization{T}) where {T} = F | ||
112 | # This no longer looks odd since the return _is_ a Factorization! | ||
113 | Factorization{T}(A::AdjointFactorization) where {T} = | ||
114 | adjoint(Factorization{T}(parent(A))) | ||
115 | Factorization{T}(A::TransposeFactorization) where {T} = | ||
116 | transpose(Factorization{T}(parent(A))) | ||
117 | inv(F::Factorization{T}) where {T} = (n = size(F, 1); ldiv!(F, Matrix{T}(I, n, n))) | ||
118 | |||
119 | Base.hash(F::Factorization, h::UInt) = mapreduce(f -> hash(getfield(F, f)), hash, 1:nfields(F); init=h) | ||
120 | Base.:(==)( F::T, G::T) where {T<:Factorization} = all(f -> getfield(F, f) == getfield(G, f), 1:nfields(F)) | ||
121 | Base.isequal(F::T, G::T) where {T<:Factorization} = all(f -> isequal(getfield(F, f), getfield(G, f)), 1:nfields(F))::Bool | ||
122 | |||
123 | function Base.show(io::IO, x::AdjointFactorization) | ||
124 | print(io, "adjoint of ") | ||
125 | show(io, parent(x)) | ||
126 | end | ||
127 | function Base.show(io::IO, x::TransposeFactorization) | ||
128 | print(io, "transpose of ") | ||
129 | show(io, parent(x)) | ||
130 | end | ||
131 | function Base.show(io::IO, ::MIME"text/plain", x::AdjointFactorization) | ||
132 | print(io, "adjoint of ") | ||
133 | show(io, MIME"text/plain"(), parent(x)) | ||
134 | end | ||
135 | function Base.show(io::IO, ::MIME"text/plain", x::TransposeFactorization) | ||
136 | print(io, "transpose of ") | ||
137 | show(io, MIME"text/plain"(), parent(x)) | ||
138 | end | ||
139 | |||
140 | function (\)(F::Factorization, B::AbstractVecOrMat) | ||
141 | require_one_based_indexing(B) | ||
142 | TFB = typeof(oneunit(eltype(F)) \ oneunit(eltype(B))) | ||
143 | ldiv!(F, copy_similar(B, TFB)) | ||
144 | end | ||
145 | (\)(F::TransposeFactorization, B::AbstractVecOrMat) = conj!(adjoint(F.parent) \ conj.(B)) | ||
146 | # With a real lhs and complex rhs with the same precision, we can reinterpret | ||
147 | # the complex rhs as a real rhs with twice the number of columns or rows | ||
148 | function (\)(F::Factorization{T}, B::VecOrMat{Complex{T}}) where {T<:BlasReal} | ||
149 | require_one_based_indexing(B) | ||
150 | c2r = reshape(copy(transpose(reinterpret(T, reshape(B, (1, length(B)))))), size(B, 1), 2*size(B, 2)) | ||
151 | x = ldiv!(F, c2r) | ||
152 | return reshape(copy(reinterpret(Complex{T}, copy(transpose(reshape(x, div(length(x), 2), 2))))), _ret_size(F, B)) | ||
153 | end | ||
154 | # don't do the reinterpretation for [Adjoint/Transpose]Factorization | ||
155 | (\)(F::TransposeFactorization{T}, B::VecOrMat{Complex{T}}) where {T<:BlasReal} = | ||
156 | conj!(adjoint(parent(F)) \ conj.(B)) | ||
157 | (\)(F::AdjointFactorization{T}, B::VecOrMat{Complex{T}}) where {T<:BlasReal} = | ||
158 | @invoke \(F::typeof(F), B::VecOrMat) | ||
159 | |||
160 | function ldiv!(Y::AbstractVector, A::Factorization, B::AbstractVector) | ||
161 | require_one_based_indexing(Y, B) | ||
162 | m, n = size(A) | ||
163 | if m > n | ||
164 | Bc = copy(B) | ||
165 | ldiv!(A, Bc) | ||
166 | return copyto!(Y, 1, Bc, 1, n) | ||
167 | else | ||
168 | 9 (3 %) |
9 (100 %)
samples spent calling
ldiv!
return ldiv!(A, copyto!(Y, B))
|
|
169 | end | ||
170 | end | ||
171 | function ldiv!(Y::AbstractMatrix, A::Factorization, B::AbstractMatrix) | ||
172 | require_one_based_indexing(Y, B) | ||
173 | m, n = size(A) | ||
174 | if m > n | ||
175 | Bc = copy(B) | ||
176 | ldiv!(A, Bc) | ||
177 | return copyto!(Y, view(Bc, 1:n, :)) | ||
178 | else | ||
179 | copyto!(view(Y, 1:m, :), view(B, 1:m, :)) | ||
180 | return ldiv!(A, Y) | ||
181 | end | ||
182 | end | ||
183 | |||
184 | function (/)(B::AbstractMatrix, F::Factorization) | ||
185 | require_one_based_indexing(B) | ||
186 | TFB = typeof(oneunit(eltype(B)) / oneunit(eltype(F))) | ||
187 | rdiv!(copy_similar(B, TFB), F) | ||
188 | end | ||
189 | # reinterpretation trick for complex lhs and real factorization | ||
190 | function (/)(B::Union{Matrix{Complex{T}},AdjOrTrans{Complex{T},Vector{Complex{T}}}}, F::Factorization{T}) where {T<:BlasReal} | ||
191 | require_one_based_indexing(B) | ||
192 | x = rdiv!(copy(reinterpret(T, B)), F) | ||
193 | return copy(reinterpret(Complex{T}, x)) | ||
194 | end | ||
195 | # don't do the reinterpretation for [Adjoint/Transpose]Factorization | ||
196 | (/)(B::Union{Matrix{Complex{T}},AdjOrTrans{Complex{T},Vector{Complex{T}}}}, F::TransposeFactorization{T}) where {T<:BlasReal} = | ||
197 | @invoke /(B::AbstractMatrix, F::Factorization) | ||
198 | (/)(B::Matrix{Complex{T}}, F::AdjointFactorization{T}) where {T<:BlasReal} = | ||
199 | @invoke /(B::AbstractMatrix, F::Factorization) | ||
200 | (/)(B::Adjoint{Complex{T},Vector{Complex{T}}}, F::AdjointFactorization{T}) where {T<:BlasReal} = | ||
201 | (F' \ B')' | ||
202 | (/)(B::Transpose{Complex{T},Vector{Complex{T}}}, F::TransposeFactorization{T}) where {T<:BlasReal} = | ||
203 | transpose(transpose(F) \ transpose(B)) | ||
204 | |||
205 | rdiv!(B::AbstractMatrix, A::TransposeFactorization) = transpose(ldiv!(A.parent, transpose(B))) | ||
206 | rdiv!(B::AbstractMatrix, A::AdjointFactorization) = adjoint(ldiv!(A.parent, adjoint(B))) |