Skip to content

Commit

Permalink
Add qrfact(SparseMatrixCSC) by wrapping SPQR. Make \(SparseMatrixCSC)…
Browse files Browse the repository at this point in the history
… work for least squares problems. Some tests.
  • Loading branch information
andreasnoack committed Feb 12, 2015
1 parent 3451d8f commit 74c5b11
Show file tree
Hide file tree
Showing 5 changed files with 21 additions and 12 deletions.
8 changes: 6 additions & 2 deletions base/linalg/generic.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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
Expand Down
1 change: 1 addition & 0 deletions base/sparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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
12 changes: 10 additions & 2 deletions base/sparse/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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 #
Expand Down Expand Up @@ -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
Expand All @@ -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"))

This comment has been minimized.

Copy link
@jiahao

jiahao Feb 12, 2015

Member

*columns

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)
Expand Down
4 changes: 3 additions & 1 deletion base/sparse/linalg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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"))

This comment has been minimized.

Copy link
@jiahao

jiahao Feb 12, 2015

Member

*systems

end
end

Expand Down
8 changes: 1 addition & 7 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -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",
Expand All @@ -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)
Expand Down

0 comments on commit 74c5b11

Please sign in to comment.