diff --git a/base/linalg/generic.jl b/base/linalg/generic.jl index cd9a04e729fd79..3fbeedf7666bdf 100644 --- a/base/linalg/generic.jl +++ b/base/linalg/generic.jl @@ -237,11 +237,15 @@ function inv{T}(A::AbstractMatrix{T}) A_ldiv_B!(factorize(convert(AbstractMatrix{S}, A)), eye(S, chksquare(A))) end +function \{T}(A::AbstractMatrix{T}, B::AbstractVecOrMat{T}) + size(A,1) == size(B,1) || throw(DimensionMismatch("LHS and RHS should have the same number of rows. LHS has $(size(A,1)) rows, but RHS has $(size(B,1)) rows.")) + factorize(A)\B +end function \{TA,TB}(A::AbstractMatrix{TA}, B::AbstractVecOrMat{TB}) TC = typeof(one(TA)/one(TB)) - size(A,1) == size(B,1) || throw(DimensionMismatch("LHS and RHS should have the same number of rows. LHS has $(size(A,1)) rows, but RHS has $(size(B,1)) rows.")) - \(factorize(TA == TC ? A : convert(AbstractMatrix{TC}, A)), TB == TC ? copy(B) : convert(AbstractArray{TC}, B)) + convert(AbstractMatrix{TC}, A)\convert(AbstractArray{TC}, B) end + \(a::AbstractVector, b::AbstractArray) = reshape(a, length(a), 1) \ b /(A::AbstractVecOrMat, B::AbstractVecOrMat) = (B' \ A')' # \(A::StridedMatrix,x::Number) = inv(A)*x Should be added at some point when the old elementwise version has been deprecated long enough diff --git a/base/sparse.jl b/base/sparse.jl index c7eaee447de11c..51a5258b393ef7 100644 --- a/base/sparse.jl +++ b/base/sparse.jl @@ -19,5 +19,6 @@ include("sparse/csparse.jl") include("sparse/linalg.jl") include("sparse/umfpack.jl") include("sparse/cholmod.jl") +include("sparse/spqr.jl") end # module SparseMatrix diff --git a/base/sparse/cholmod.jl b/base/sparse/cholmod.jl index 12908b4f6beb57..34eaea325833ab 100644 --- a/base/sparse/cholmod.jl +++ b/base/sparse/cholmod.jl @@ -13,7 +13,7 @@ export Factor, Sparse -using Base.SparseMatrix: AbstractSparseMatrix, SparseMatrixCSC, increment, increment!, indtype, decrement, decrement! +using Base.SparseMatrix: AbstractSparseMatrix, SparseMatrixCSC, increment, indtype ######### # Setup # @@ -716,7 +716,7 @@ function Sparse(filename::ASCIIString) end ## convertion back to base Julia types -function convert{T<:VTypes}(::Type{Matrix}, D::Dense{T}) +function convert{T}(::Type{Matrix{T}}, D::Dense{T}) s = unsafe_load(D.p) a = Array(T, s.nrow, s.ncol) if s.d == s.nrow @@ -730,6 +730,14 @@ function convert{T<:VTypes}(::Type{Matrix}, D::Dense{T}) end a end +convert{T}(::Type{Matrix}, D::Dense{T}) = convert(Matrix{T}, D) +function convert{T}(::Type{Vector{T}}, D::Dense{T}) + if size(D, 2) > 1 + throw(DimensionMismatch("input must be a vector but had $(size(D, 2)) columnds")) + end + reshape(convert(Matrix, D), size(D, 1)) +end +convert{T}(::Type{Vector}, D::Dense{T}) = convert(Vector{T}, D) function convert{Tv,Ti}(::Type{SparseMatrixCSC{Tv,Ti}}, A::Sparse{Tv,Ti}) s = unsafe_load(A.p) diff --git a/base/sparse/linalg.jl b/base/sparse/linalg.jl index 3b5249c3136994..bda9591b123ebc 100644 --- a/base/sparse/linalg.jl +++ b/base/sparse/linalg.jl @@ -686,8 +686,10 @@ function factorize(A::SparseMatrixCSC) end end return lufact(A) + elseif m > n + return qrfact(A) else - throw(ArgumentError("sparse least squares problems by QR are not handled yet")) + throw(ArgumentError("underdetermined systemed are not implemented yet")) end end diff --git a/test/runtests.jl b/test/runtests.jl index f499145432a14e..1f393f7148a817 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -3,7 +3,7 @@ testnames = [ "linalg", "core", "keywordargs", "numbers", "strings", "dates", "dict", "hashing", "remote", "iobuffer", "staged", "arrayops", "subarray", "reduce", "reducedim", "random", "intfuncs", - "simdloop", "blas", "fft", "dsp", "sparse", "bitarray", "copy", "math", + "simdloop", "blas", "fft", "dsp", "sparsetests", "bitarray", "copy", "math", "fastmath", "functional", "bigint", "sorting", "statistics", "spawn", "backtrace", "priorityqueue", "arpack", "file", "version", "resolve", "pollfd", "mpfr", "broadcast", "complex", "socket", @@ -25,12 +25,6 @@ push!(testnames, "parallel") tests = (ARGS==["all"] || isempty(ARGS)) ? testnames : ARGS -if "sparse" in tests - # specifically selected case - filter!(x -> x != "sparse", tests) - prepend!(tests, ["sparse/sparse", "sparse/cholmod", "sparse/umfpack"]) -end - if "linalg" in tests # specifically selected case filter!(x -> x != "linalg", tests)