From e4a42704d69383676a7a504f15b402a381038a75 Mon Sep 17 00:00:00 2001 From: Jishnu Bhattacharya Date: Sun, 8 Sep 2024 18:12:38 +0530 Subject: [PATCH] Avoid stack overflow in triangular eigvecs (#55497) MIME-Version: 1.0 Content-Type: text/plain; charset=UTF-8 Content-Transfer-Encoding: 8bit This fixes a stack overflow in ```julia julia> using LinearAlgebra, StaticArrays julia> U = UpperTriangular(SMatrix{2,2}(1:4)) 2×2 UpperTriangular{Int64, SMatrix{2, 2, Int64, 4}} with indices SOneTo(2)×SOneTo(2): 1 3 ⋅ 4 julia> eigvecs(U) Warning: detected a stack overflow; program state may be corrupted, so further execution might be unreliable. ERROR: StackOverflowError: Stacktrace: [1] eigvecs(A::UpperTriangular{Float32, SMatrix{2, 2, Float32, 4}}) (repeats 79984 times) @ LinearAlgebra ~/.julia/juliaup/julia-nightly/share/julia/stdlib/v1.12/LinearAlgebra/src/triangular.jl:2749 ``` After this, ```julia julia> eigvecs(U) 2×2 Matrix{Float32}: 1.0 1.0 0.0 1.0 ``` --- stdlib/LinearAlgebra/src/triangular.jl | 8 ++++++++ stdlib/LinearAlgebra/test/triangular.jl | 16 ++++++++++++++++ 2 files changed, 24 insertions(+) diff --git a/stdlib/LinearAlgebra/src/triangular.jl b/stdlib/LinearAlgebra/src/triangular.jl index 71473e0dc11743..df0d0d4fd0d8b4 100644 --- a/stdlib/LinearAlgebra/src/triangular.jl +++ b/stdlib/LinearAlgebra/src/triangular.jl @@ -2782,6 +2782,14 @@ end # Generic eigensystems eigvals(A::AbstractTriangular) = diag(A) +# fallback for unknown types +function eigvecs(A::AbstractTriangular{<:BlasFloat}) + if istriu(A) + eigvecs(UpperTriangular(Matrix(A))) + else # istril(A) + eigvecs(LowerTriangular(Matrix(A))) + end +end function eigvecs(A::AbstractTriangular{T}) where T TT = promote_type(T, Float32) if TT <: BlasFloat diff --git a/stdlib/LinearAlgebra/test/triangular.jl b/stdlib/LinearAlgebra/test/triangular.jl index 8748d11bd1da4d..a09d0092e9f394 100644 --- a/stdlib/LinearAlgebra/test/triangular.jl +++ b/stdlib/LinearAlgebra/test/triangular.jl @@ -1198,6 +1198,22 @@ end end end +@testset "eigvecs for AbstractTriangular" begin + S = SizedArrays.SizedArray{(3,3)}(reshape(1:9,3,3)) + for T in (UpperTriangular, UnitUpperTriangular, + LowerTriangular, UnitLowerTriangular) + U = T(S) + V = eigvecs(U) + λ = eigvals(U) + @test U * V ≈ V * Diagonal(λ) + + MU = MyTriangular(U) + V = eigvecs(U) + λ = eigvals(U) + @test MU * V ≈ V * Diagonal(λ) + end +end + @testset "(l/r)mul! and (l/r)div! for generic triangular" begin @testset for T in (UpperTriangular, LowerTriangular, UnitUpperTriangular, UnitLowerTriangular) M = MyTriangular(T(rand(4,4)))