From d8e10c7edead31916ab2d970b338438033022119 Mon Sep 17 00:00:00 2001 From: Daniel VandenHeuvel <95613936+DanielVandH@users.noreply.github.com> Date: Mon, 19 Aug 2024 11:07:57 +0100 Subject: [PATCH] Fix layout of `cholesky(Symmetric(view(BandedMatrix)))` (#451) * Fix 450 * Change tests for LTS * Change tests for LTS * Should also test that U is actually correct * Fix check of LTS --- Project.toml | 2 +- src/symbanded/BandedCholesky.jl | 2 ++ test/test_symbanded.jl | 17 +++++++++++++++++ 3 files changed, 20 insertions(+), 1 deletion(-) diff --git a/Project.toml b/Project.toml index 388aa786..a1238e8b 100644 --- a/Project.toml +++ b/Project.toml @@ -1,6 +1,6 @@ name = "BandedMatrices" uuid = "aae01518-5342-5314-be14-df237901396f" -version = "1.7.2" +version = "1.7.3" [deps] ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a" diff --git a/src/symbanded/BandedCholesky.jl b/src/symbanded/BandedCholesky.jl index 5f25e811..0e8549a1 100644 --- a/src/symbanded/BandedCholesky.jl +++ b/src/symbanded/BandedCholesky.jl @@ -98,5 +98,7 @@ ArrayLayouts._cholesky(::Union{SymmetricLayout{<:AbstractBandedLayout},Hermitian cholcopy(A::RealHermSymComplexHerm{<:Any,<:BandedMatrix}) = copyto!(similar(A, LinearAlgebra.choltype(A)), A) +cholcopy(A::RealHermSymComplexHerm{<:Any,<:SubArray{<:Number, 2, <:BandedMatrix}}) = + copyto!(similar(A, LinearAlgebra.choltype(A)), A) # cholesky(A::Hermitian{T,<:BandedMatrix{T}}, # ::Val{false}=Val(false); check::Bool = true) where T = cholesky!(cholcopy(A); check = check) diff --git a/test/test_symbanded.jl b/test/test_symbanded.jl index c1ac0fb1..2cf802a3 100644 --- a/test/test_symbanded.jl +++ b/test/test_symbanded.jl @@ -402,6 +402,23 @@ end @test isposdef(B) == isposdef(Matrix(B)) == false B = BandedMatrix(1=>Float64[1:4;], 0=>Float64[2:2:10;], -1=>Float64[1:4;]) @test isposdef(B) == isposdef(Matrix(B)) == true + + @testset "Issue #450: Cholesky of a BandedMatrix's view" begin + L = brand(10, 10, 2, 0) + B = BandedMatrix(Symmetric(L * L')); + Bv = Symmetric(view(B, :, :)) + chol = cholesky(Bv) + if VERSION >= v"1.7" + @test MemoryLayout(chol.U) isa TriangularLayout{'U', 'N', typeof(MemoryLayout(B))} + @test MemoryLayout(chol.L) isa TriangularLayout{'L', 'N', typeof(MemoryLayout(B'))} + @test isbanded(chol.L) + end + @test isbanded(chol.U) + cc = LinearAlgebra.cholcopy(Bv) + @test cc isa Symmetric{Float64, typeof(B)} + @test chol.U ≈ cholesky(Matrix(B)).U + @test cc == B + end end end # module