diff --git a/stdlib/LinearAlgebra/src/tridiag.jl b/stdlib/LinearAlgebra/src/tridiag.jl index 497e265b280bc..bdedcd1b5e023 100644 --- a/stdlib/LinearAlgebra/src/tridiag.jl +++ b/stdlib/LinearAlgebra/src/tridiag.jl @@ -174,7 +174,7 @@ Base.copy(S::Adjoint{<:Any,<:SymTridiagonal}) = SymTridiagonal(map(x -> copy.(ad ishermitian(S::SymTridiagonal) = isreal(S.dv) && isreal(_evview(S)) issymmetric(S::SymTridiagonal) = true -tr(S::SymTridiagonal) = sum(S.dv) +tr(S::SymTridiagonal) = sum(symmetric, S.dv) @noinline function throw_diag_outofboundserror(n, sz) sz1, sz2 = sz diff --git a/stdlib/LinearAlgebra/test/tridiag.jl b/stdlib/LinearAlgebra/test/tridiag.jl index 487f808076707..00414fcdc8b56 100644 --- a/stdlib/LinearAlgebra/test/tridiag.jl +++ b/stdlib/LinearAlgebra/test/tridiag.jl @@ -468,7 +468,7 @@ end end @testset "SymTridiagonal/Tridiagonal block matrix" begin - M = [1 2; 2 4] + M = [1 2; 3 4] n = 5 A = SymTridiagonal(fill(M, n), fill(M, n-1)) @test @inferred A[1,1] == Symmetric(M) @@ -482,6 +482,9 @@ end @test_throws ArgumentError diag(A, n+1) @test_throws ArgumentError diag(A, -n-1) + @test tr(A) == sum(diag(A)) + @test issymmetric(tr(A)) + A = Tridiagonal(fill(M, n-1), fill(M, n), fill(M, n-1)) @test @inferred A[1,1] == M @test @inferred A[1,2] == M