Skip to content

Commit 2fe4190

Browse files
authored
Reduce code duplication in eigen(::AbstractMatrix) methods (#52450)
This also fixes the following: ```julia julia> D = Diagonal(Float16[1,2,4]) 3×3 Diagonal{Float16, Vector{Float16}}: 1.0 ⋅ ⋅ ⋅ 2.0 ⋅ ⋅ ⋅ 4.0 julia> eigen(Array(D)) # returns Float32 Eigen{Float32, Float32, Matrix{Float32}, Vector{Float32}} values: 3-element Vector{Float32}: 1.0 2.0 4.0 vectors: 3×3 Matrix{Float32}: 1.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 1.0 julia> eigen(D) # returns Float16 as expected Eigen{Float16, Float16, Matrix{Float16}, Vector{Float16}} values: 3-element Vector{Float16}: 1.0 2.0 4.0 vectors: 3×3 Matrix{Float16}: 1.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 1.0 ``` After this, ```julia julia> eigen(Array(D)) Eigen{Float16, Float16, Matrix{Float16}, Vector{Float16}} values: 3-element Vector{Float16}: 1.0 2.0 4.0 vectors: 3×3 Matrix{Float16}: 1.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 1.0 ```
1 parent 5e9cd58 commit 2fe4190

File tree

2 files changed

+24
-8
lines changed

2 files changed

+24
-8
lines changed

stdlib/LinearAlgebra/src/eigen.jl

Lines changed: 9 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -236,22 +236,23 @@ true
236236
```
237237
"""
238238
function eigen(A::AbstractMatrix{T}; permute::Bool=true, scale::Bool=true, sortby::Union{Function,Nothing}=eigsortby) where T
239-
isdiag(A) && return eigen(Diagonal{eigtype(T)}(diag(A)); sortby)
240-
ishermitian(A) && return eigen!(eigencopy_oftype(Hermitian(A), eigtype(T)); sortby)
241-
AA = eigencopy_oftype(A, eigtype(T))
242-
return eigen!(AA; permute, scale, sortby)
239+
_eigen(A; permute, scale, sortby)
243240
end
244241
function eigen(A::AbstractMatrix{T}; permute::Bool=true, scale::Bool=true, sortby::Union{Function,Nothing}=eigsortby) where {T <: Union{Float16,Complex{Float16}}}
242+
E = _eigen(A; permute, scale, sortby)
243+
values = convert(AbstractVector{isreal(E.values) ? Float16 : Complex{Float16}}, E.values)
244+
vectors = convert(AbstractMatrix{isreal(E.vectors) ? Float16 : Complex{Float16}}, E.vectors)
245+
return Eigen(values, vectors)
246+
end
247+
function _eigen(A::AbstractMatrix{T}; permute=true, scale=true, sortby=eigsortby) where {T}
245248
isdiag(A) && return eigen(Diagonal{eigtype(T)}(diag(A)); sortby)
246-
E = if ishermitian(A)
249+
if ishermitian(A)
247250
eigen!(eigencopy_oftype(Hermitian(A), eigtype(T)); sortby)
248251
else
249252
eigen!(eigencopy_oftype(A, eigtype(T)); permute, scale, sortby)
250253
end
251-
values = convert(AbstractVector{isreal(E.values) ? Float16 : Complex{Float16}}, E.values)
252-
vectors = convert(AbstractMatrix{isreal(E.vectors) ? Float16 : Complex{Float16}}, E.vectors)
253-
return Eigen(values, vectors)
254254
end
255+
255256
eigen(x::Number) = Eigen([x], fill(one(x), 1, 1))
256257

257258
"""

stdlib/LinearAlgebra/test/eigen.jl

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -241,6 +241,21 @@ end
241241
@test F.vectors isa Matrix{ComplexF16}
242242
@test F.values F32.values
243243
@test F.vectors F32.vectors
244+
245+
for T in (Float16, ComplexF16)
246+
D = Diagonal(T[1,2,4])
247+
A = Array(D)
248+
B = eigen(A)
249+
@test B isa Eigen{Float16, Float16, Matrix{Float16}, Vector{Float16}}
250+
@test B.values isa Vector{Float16}
251+
@test B.vectors isa Matrix{Float16}
252+
end
253+
D = Diagonal(ComplexF16[im,2,4])
254+
A = Array(D)
255+
B = eigen(A)
256+
@test B isa Eigen{Float16, ComplexF16, Matrix{Float16}, Vector{ComplexF16}}
257+
@test B.values isa Vector{ComplexF16}
258+
@test B.vectors isa Matrix{Float16}
244259
end
245260

246261
@testset "complex eigen inference (#52289)" begin

0 commit comments

Comments
 (0)