Skip to content

Commit fbe32f7

Browse files
dkarraschKristofferC
authored andcommitted
Improve performance of svd and eigen of Diagonals (#43856)
(cherry picked from commit 6bad936)
1 parent 4cfc794 commit fbe32f7

File tree

2 files changed

+52
-10
lines changed

2 files changed

+52
-10
lines changed

stdlib/LinearAlgebra/src/diagonal.jl

Lines changed: 27 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -667,20 +667,38 @@ function eigen(D::Diagonal; permute::Bool=true, scale::Bool=true, sortby::Union{
667667
if any(!isfinite, D.diag)
668668
throw(ArgumentError("matrix contains Infs or NaNs"))
669669
end
670-
Eigen(sorteig!(eigvals(D), eigvecs(D), sortby)...)
670+
Td = Base.promote_op(/, eltype(D), eltype(D))
671+
λ = eigvals(D)
672+
if !isnothing(sortby)
673+
p = sortperm(λ; alg=QuickSort, by=sortby)
674+
λ = λ[p] # make a copy, otherwise this permutes D.diag
675+
evecs = zeros(Td, size(D))
676+
@inbounds for i in eachindex(p)
677+
evecs[p[i],i] = one(Td)
678+
end
679+
else
680+
evecs = Matrix{Td}(I, size(D))
681+
end
682+
Eigen(λ, evecs)
671683
end
672684

673685
#Singular system
674686
svdvals(D::Diagonal{<:Number}) = sort!(abs.(D.diag), rev = true)
675687
svdvals(D::Diagonal) = [svdvals(v) for v in D.diag]
676-
function svd(D::Diagonal{<:Number})
677-
S = abs.(D.diag)
678-
piv = sortperm(S, rev = true)
679-
U = Diagonal(D.diag ./ S)
680-
Up = hcat([U[:,i] for i = 1:length(D.diag)][piv]...)
681-
V = Diagonal(fill!(similar(D.diag), one(eltype(D.diag))))
682-
Vp = hcat([V[:,i] for i = 1:length(D.diag)][piv]...)
683-
return SVD(Up, S[piv], copy(Vp'))
688+
function svd(D::Diagonal{T}) where {T<:Number}
689+
d = D.diag
690+
s = abs.(d)
691+
piv = sortperm(s, rev = true)
692+
S = s[piv]
693+
Td = typeof(oneunit(T)/oneunit(T))
694+
U = zeros(Td, size(D))
695+
Vt = copy(U)
696+
for i in 1:length(d)
697+
j = piv[i]
698+
U[j,i] = d[j] / S[i]
699+
Vt[i,j] = one(Td)
700+
end
701+
return SVD(U, S, Vt)
684702
end
685703

686704
# disambiguation methods: * of Diagonal and Adj/Trans AbsVec

stdlib/LinearAlgebra/test/diagonal.jl

Lines changed: 25 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -5,6 +5,10 @@ module TestDiagonal
55
using Test, LinearAlgebra, SparseArrays, Random
66
using LinearAlgebra: mul!, mul!, rmul!, lmul!, ldiv!, rdiv!, BlasFloat, BlasComplex, SingularException
77

8+
const BASE_TEST_PATH = joinpath(Sys.BINDIR, "..", "share", "julia", "test")
9+
isdefined(Main, :Furlongs) || @eval Main include(joinpath($(BASE_TEST_PATH), "testhelpers", "Furlongs.jl"))
10+
using .Main.Furlongs
11+
812
n=12 #Size of matrix problem to test
913
Random.seed!(1)
1014

@@ -330,8 +334,12 @@ Random.seed!(1)
330334

331335
@testset "Eigensystem" begin
332336
eigD = eigen(D)
333-
@test Diagonal(eigD.values) D
337+
@test Diagonal(eigD.values) == D
334338
@test eigD.vectors == Matrix(I, size(D))
339+
eigsortD = eigen(D, sortby=LinearAlgebra.eigsortby)
340+
@test eigsortD.values !== D.diag
341+
@test eigsortD.values == sort(D.diag, by=LinearAlgebra.eigsortby)
342+
@test Matrix(eigsortD) == D
335343
end
336344

337345
@testset "ldiv" begin
@@ -395,6 +403,22 @@ Random.seed!(1)
395403
@test svd(D).V == V
396404
end
397405

406+
@testset "svd/eigen with Diagonal{Furlong}" begin
407+
Du = Furlong.(D)
408+
@test Du isa Diagonal{<:Furlong{1}}
409+
F = svd(Du)
410+
U, s, V = F
411+
@test map(x -> x.val, Matrix(F)) map(x -> x.val, Du)
412+
@test svdvals(Du) == s
413+
@test U isa AbstractMatrix{<:Furlong{0}}
414+
@test V isa AbstractMatrix{<:Furlong{0}}
415+
@test s isa AbstractVector{<:Furlong{1}}
416+
E = eigen(Du)
417+
vals, vecs = E
418+
@test Matrix(E) == Du
419+
@test vals isa AbstractVector{<:Furlong{1}}
420+
@test vecs isa AbstractMatrix{<:Furlong{0}}
421+
end
398422
end
399423

400424
@testset "svdvals and eigvals (#11120/#11247)" begin

0 commit comments

Comments
 (0)