diff --git a/base/sparse/umfpack.jl b/base/sparse/umfpack.jl index aeaa590dc82db..27cbe194b602a 100644 --- a/base/sparse/umfpack.jl +++ b/base/sparse/umfpack.jl @@ -397,10 +397,9 @@ for (f!, umfpack) in ((:A_ldiv_B!, :UMFPACK_A), end return x end - $f!{T<:UMFVTypes}(lu::UmfpackLU{T}, b::StridedVector{T}) = $f!(b, lu, copy(b)) - $f!{T<:UMFVTypes}(lu::UmfpackLU{T}, b::StridedMatrix{T}) = $f!(b, lu, copy(b)) + $f!{T<:UMFVTypes}(lu::UmfpackLU{T}, b::StridedVecOrMat{T}) = $f!(b, lu, copy(b)) - function $f!{Tb<:Complex}(x::StridedVector{Tb}, lu::UmfpackLU{Float64}, b::StridedVector{Tb}) + function $f!{Tb<:Complex}(x::StridedVecOrMat{Tb}, lu::UmfpackLU{Float64}, b::StridedVecOrMat{Tb}) m, n = size(x, 1), size(x, 2) if n != size(b, 2) throw(DimensionMismatch("in and output arrays must have the same number of columns")) @@ -416,7 +415,7 @@ for (f!, umfpack) in ((:A_ldiv_B!, :UMFPACK_A), end return x end - $f!{Tb<:Complex}(lu::UmfpackLU{Float64}, b::StridedVector{Tb}) = $f!(b, lu, copy(b)) + $f!{Tb<:Complex}(lu::UmfpackLU{Float64}, b::StridedVecOrMat{Tb}) = $f!(b, lu, copy(b)) end end diff --git a/test/sparse/umfpack.jl b/test/sparse/umfpack.jl index dbfe2dfae0387..65f63f6ebdb20 100644 --- a/test/sparse/umfpack.jl +++ b/test/sparse/umfpack.jl @@ -143,3 +143,15 @@ let F = lufact(A) @test F[:p] == [3 ; 4 ; 2 ; 1] end + +# Test that A[c|t]_ldiv_B!{T<:Complex}(X::StridedMatrix{T}, lu::UmfpackLU{Float64}, +# B::StridedMatrix{T}) works as expected. +let N = 10, p = 0.5 + A = N*speye(N) + sprand(N, N, p) + X = zeros(Complex{Float64}, N, N) + B = complex.(rand(N, N), rand(N, N)) + luA, lufA = lufact(A), lufact(Array(A)) + @test A_ldiv_B!(copy(X), luA, B) ≈ A_ldiv_B!(copy(X), lufA, B) + @test At_ldiv_B!(copy(X), luA, B) ≈ At_ldiv_B!(copy(X), lufA, B) + @test Ac_ldiv_B!(copy(X), luA, B) ≈ Ac_ldiv_B!(copy(X), lufA, B) +end