From fb0e617686aae88d44bf04a4cbc47aa7b0ecb885 Mon Sep 17 00:00:00 2001 From: Sacha Verweij Date: Thu, 16 Jun 2016 21:05:10 -0700 Subject: [PATCH] Add backslash method specialized for `SparseMatrixCSC`s and a test checking that that specialized method gets called, a simple if inelegant fix for #16548. --- base/sparse/linalg.jl | 26 ++++++++++++++++++++++++++ test/sparsedir/sparse.jl | 3 +++ 2 files changed, 29 insertions(+) diff --git a/base/sparse/linalg.jl b/base/sparse/linalg.jl index 1b24452007c97..05247ccf0343e 100644 --- a/base/sparse/linalg.jl +++ b/base/sparse/linalg.jl @@ -856,6 +856,32 @@ end scale!(A::SparseMatrixCSC, b::Number) = (scale!(A.nzval, b); A) scale!(b::Number, A::SparseMatrixCSC) = (scale!(b, A.nzval); A) +function (\)(A::SparseMatrixCSC, B::AbstractVecOrMat) + m, n = size(A) + if m == n + if istril(A) + if istriu(A) + return Diagonal(A) \ B + else + return LowerTriangular(A) \ B + end + elseif istriu(A) + return UpperTriangular(A) \ B + end + if ishermitian(A) + try + return cholfact(Hermitian(A)) \ B + catch e + isa(e, PosDefException) || rethrow(e) + return ldltfact(Hermitian(A)) \ B + end + end + return lufact(A) \ B + else + return qrfact(A) \ B + end +end + function factorize(A::SparseMatrixCSC) m, n = size(A) if m == n diff --git a/test/sparsedir/sparse.jl b/test/sparsedir/sparse.jl index 5f4ae28caae14..aa3151c2244a1 100644 --- a/test/sparsedir/sparse.jl +++ b/test/sparsedir/sparse.jl @@ -1385,3 +1385,6 @@ end @test !issparse([rand(10,10) rand(10,10)]) @test !issparse([rand(10,10); rand(10,10)]) @test !issparse([rand(10,10) rand(10,10); rand(10,10) rand(10,10)]) + +# Test temporary fix for issue #16548 in PR #16979. Brittle. Expect to remove with `\` revisions. +@test which(\, (SparseMatrixCSC, AbstractVecOrMat)).module == Base.SparseArrays \ No newline at end of file