diff --git a/base/array.jl b/base/array.jl index 8c753b3e3924dc..b9452e62d7e889 100644 --- a/base/array.jl +++ b/base/array.jl @@ -682,15 +682,6 @@ function reverse!(v::AbstractVector, s=1, n=length(v)) end -# concatenations of combinations (homogeneous, heterogeneous) of dense matrices/vectors # -vcat{T}(A::Union{Vector{T},Matrix{T}}...) = typed_vcat(T, A...) -vcat(A::Union{Vector,Matrix}...) = typed_vcat(promote_eltype(A...), A...) -hcat{T}(A::Union{Vector{T},Matrix{T}}...) = typed_hcat(T, A...) -hcat(A::Union{Vector,Matrix}...) = typed_hcat(promote_eltype(A...), A...) -hvcat{T}(rows::Tuple{Vararg{Int}}, xs::Union{Vector{T},Matrix{T}}...) = typed_hvcat(T, rows, xs...) -hvcat(rows::Tuple{Vararg{Int}}, xs::Union{Vector,Matrix}...) = typed_hvcat(promote_eltype(xs...), rows, xs...) -cat{T}(catdims, xs::Union{Vector{T},Matrix{T}}...) = Base.cat_t(catdims, T, xs...) -cat(catdims, xs::Union{Vector,Matrix}...) = Base.cat_t(catdims, promote_eltype(xs...), xs...) # concatenations of homogeneous combinations of vectors, horizontal and vertical function hcat{T}(V::Vector{T}...) height = length(V[1]) diff --git a/base/sparse/sparsematrix.jl b/base/sparse/sparsematrix.jl index d8d3fa47d0f424..c0174fe850a1bb 100644 --- a/base/sparse/sparsematrix.jl +++ b/base/sparse/sparsematrix.jl @@ -265,6 +265,8 @@ function convert{Tv,Ti}(::Type{SparseMatrixCSC{Tv,Ti}}, S::SparseMatrixCSC) end # convert'ing from other matrix types to SparseMatrixCSC (also see sparse()) convert(::Type{SparseMatrixCSC}, M::Matrix) = sparse(M) +convert{Tv}(::Type{SparseMatrixCSC}, M::AbstractMatrix{Tv}) = convert(SparseMatrixCSC{Tv,Int}, M) +convert{Tv}(::Type{SparseMatrixCSC{Tv}}, M::AbstractMatrix{Tv}) = convert(SparseMatrixCSC{Tv,Int}, M) function convert{Tv,Ti}(::Type{SparseMatrixCSC{Tv,Ti}}, M::AbstractMatrix) (I, J, V) = findnz(M) eltypeTiI = convert(Vector{Ti}, I) diff --git a/base/sparse/sparsevector.jl b/base/sparse/sparsevector.jl index cd002a2574bdbf..e94d34a2afb7e8 100644 --- a/base/sparse/sparsevector.jl +++ b/base/sparse/sparsevector.jl @@ -793,27 +793,48 @@ hcat(Xin::Union{Vector, AbstractSparseVector}...) = hcat(map(sparse, Xin)...) vcat(Xin::Union{Vector, AbstractSparseVector}...) = vcat(map(sparse, Xin)...) -### Sparse/special/dense vector/matrix concatenation - -# TODO: These methods should be moved to a more appropriate location, particularly some -# future equivalent of base/linalg/special.jl dedicated to interactions between a broader -# set of matrix types. - -# TODO: A similar definition also exists in base/linalg/bidiag.jl. These definitions should -# be consolidated in a more appropriate location, for example base/linalg/special.jl. -SpecialArrays = Union{Diagonal, Bidiagonal, Tridiagonal, SymTridiagonal} - -function hcat(Xin::Union{Vector, Matrix, SparseVector, SparseMatrixCSC, SpecialArrays}...) +### Concatenation of un/annotated sparse/special/dense vectors/matrices + +# TODO: These methods and definitions should be moved to a more appropriate location, +# particularly some future equivalent of base/linalg/special.jl dedicated to interactions +# between a broader set of matrix types. + +# TODO: A definition similar to the third exists in base/linalg/bidiag.jl. These definitions +# should be consolidated in a more appropriate location, e.g. base/linalg/special.jl. +typealias _SparseArrays Union{SparseVector, SparseMatrixCSC} +typealias _SpecialArrays Union{Diagonal, Bidiagonal, Tridiagonal, SymTridiagonal} +typealias _SparseConcatArrays Union{_SpecialArrays, _SparseArrays} + +typealias _Symmetric_SparseConcatArrays{T,A<:_SparseConcatArrays} Symmetric{T,A} +typealias _Hermitian_SparseConcatArrays{T,A<:_SparseConcatArrays} Hermitian{T,A} +typealias _Triangular_SparseConcatArrays{T,A<:_SparseConcatArrays} Base.LinAlg.AbstractTriangular{T,A} +typealias _Annotated_SparseConcatArrays Union{_Triangular_SparseConcatArrays, _Symmetric_SparseConcatArrays, _Hermitian_SparseConcatArrays} + +typealias _Symmetric_DenseArrays{T,A<:Matrix} Symmetric{T,A} +typealias _Hermitian_DenseArrays{T,A<:Matrix} Hermitian{T,A} +typealias _Triangular_DenseArrays{T,A<:Matrix} Base.LinAlg.AbstractTriangular{T,A} +typealias _Annotated_DenseArrays Union{_Triangular_DenseArrays, _Symmetric_DenseArrays, _Hermitian_DenseArrays} +typealias _Annotated_Typed_DenseArrays{T} Union{_Triangular_DenseArrays{T}, _Symmetric_DenseArrays{T}, _Hermitian_DenseArrays{T}} + +typealias _SparseConcatGroup Union{Vector, Matrix, _SparseConcatArrays, _Annotated_SparseConcatArrays, _Annotated_DenseArrays} +typealias _DenseConcatGroup Union{Vector, Matrix, _Annotated_DenseArrays} +typealias _TypedDenseConcatGroup{T} Union{Vector{T}, Matrix{T}, _Annotated_Typed_DenseArrays{T}} + +# Concatenations involving un/annotated sparse/special matrices/vectors should yield sparse arrays +function cat(catdims, Xin::_SparseConcatGroup...) + X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin] + T = promote_eltype(Xin...) + Base.cat_t(catdims, T, X...) +end +function hcat(Xin::_SparseConcatGroup...) X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin] hcat(X...) end - -function vcat(Xin::Union{Vector, Matrix, SparseVector, SparseMatrixCSC, SpecialArrays}...) +function vcat(Xin::_SparseConcatGroup...) X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin] vcat(X...) end - -function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseVector, SparseMatrixCSC, SpecialArrays}...) +function hvcat(rows::Tuple{Vararg{Int}}, X::_SparseConcatGroup...) nbr = length(rows) # number of block rows tmp_rows = Array{SparseMatrixCSC}(nbr) @@ -825,11 +846,16 @@ function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseVector, vcat(tmp_rows...) end -function cat(catdims, Xin::Union{Vector, Matrix, SparseVector, SparseMatrixCSC, SpecialArrays}...) - X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin] - T = promote_eltype(Xin...) - Base.cat_t(catdims, T, X...) -end +# Concatenations strictly involving un/annotated dense matrices/vectors should yield dense arrays +cat(catdims, xs::_DenseConcatGroup...) = Base.cat_t(catdims, promote_eltype(xs...), xs...) +vcat(A::_DenseConcatGroup...) = Base.typed_vcat(promote_eltype(A...), A...) +hcat(A::_DenseConcatGroup...) = Base.typed_hcat(promote_eltype(A...), A...) +hvcat(rows::Tuple{Vararg{Int}}, xs::_DenseConcatGroup...) = Base.typed_hvcat(promote_eltype(xs...), rows, xs...) +# For performance, specially handle the case where the matrices/vectors have homogeneous eltype +cat{T}(catdims, xs::_TypedDenseConcatGroup{T}...) = Base.cat_t(catdims, T, xs...) +vcat{T}(A::_TypedDenseConcatGroup{T}...) = Base.typed_vcat(T, A...) +hcat{T}(A::_TypedDenseConcatGroup{T}...) = Base.typed_hcat(T, A...) +hvcat{T}(rows::Tuple{Vararg{Int}}, xs::_TypedDenseConcatGroup{T}...) = Base.typed_hvcat(T, rows, xs...) ### math functions diff --git a/test/linalg/special.jl b/test/linalg/special.jl index 05105acd149a86..fc8e4c67d6d84b 100644 --- a/test/linalg/special.jl +++ b/test/linalg/special.jl @@ -166,3 +166,86 @@ let end end end + +# Test that concatenations of annotated sparse/special matrix types with other matrix +# types yield sparse arrays, and that the code which effects that does not make concatenations +# strictly involving un/annotated dense matrices yield sparse arrays +# +# TODO: As with the associated code, these tests should be moved to a more appropriate +# location, particularly some future equivalent of base/linalg/special.jl dedicated to +# intereactions between a broader set of matrix types +let + N = 4 + # The tested annotation types + utriannotations = (UpperTriangular, Base.LinAlg.UnitUpperTriangular) + ltriannotations = (LowerTriangular, Base.LinAlg.UnitLowerTriangular) + annotations = (utriannotations..., ltriannotations..., Symmetric, Hermitian) + # Concatenations involving these types, un/annotated, should yield sparse arrays + spvec = spzeros(N) + spmat = speye(N) + diagmat = Diagonal(ones(N)) + bidiagmat = Bidiagonal(ones(N), ones(N-1), true) + tridiagmat = Tridiagonal(ones(N-1), ones(N), ones(N-1)) + symtridiagmat = SymTridiagonal(ones(N), ones(N-1)) + sparseconcatmats = (spmat, diagmat, bidiagmat, tridiagmat, symtridiagmat) + # Concatenations involving strictly these types, un/annotated, should yield dense arrays + densevec = ones(N) + densemat = ones(N, N) + # Annotated collections + annodmats = [annot(densemat) for annot in annotations] + annospcmats = [annot(spcmat) for annot in annotations, spcmat in sparseconcatmats] + # Test that concatenations of pairwise combinations of annotated sparse/special + # yield sparse matrices + for annospcmata in annospcmats, annospcmatb in annospcmats + @test issparse(vcat(annospcmata, annospcmatb)) + @test issparse(hcat(annospcmata, annospcmatb)) + @test issparse(hvcat((2,), annospcmata, annospcmatb)) + @test issparse(cat((1,2), annospcmata, annospcmatb)) + end + # Test that concatenations of pairwise combinations of annotated sparse/special + # matrices and other matrix/vector types yield sparse matrices + for annospcmat in annospcmats + # --> Tests applicable to pairs including only matrices + for othermat in (densemat, annodmats..., sparseconcatmats...) + @test issparse(vcat(annospcmat, othermat)) + @test issparse(vcat(othermat, annospcmat)) + end + # --> Tests applicable to pairs including other vectors or matrices + for other in (spvec, densevec, densemat, annodmats..., sparseconcatmats...) + @test issparse(hcat(annospcmat, other)) + @test issparse(hcat(other, annospcmat)) + @test issparse(hvcat((2,), annospcmat, other)) + @test issparse(hvcat((2,), other, annospcmat)) + @test issparse(cat((1,2), annospcmat, other)) + @test issparse(cat((1,2), other, annospcmat)) + end + end + # The preceding tests should cover multi-way combinations of those types, but for good + # measure test a few multi-way combinations involving those types + @test issparse(vcat(spmat, densemat, annospcmats[1], annodmats[2])) + @test issparse(vcat(densemat, spmat, annodmats[1], annospcmats[2])) + @test issparse(hcat(spvec, annodmats[3], annospcmats[3], densevec, diagmat)) + @test issparse(hcat(annodmats[4], annospcmats[4], spvec, densevec, diagmat)) + @test issparse(hvcat((5,), diagmat, densevec, spvec, annodmats[1], annospcmats[5])) + @test issparse(hvcat((5,), spvec, annodmats[2], diagmat, densevec, annospcmats[6])) + @test issparse(cat((1,2), annodmats[3], diagmat, annospcmats[7], densevec, spvec)) + @test issparse(cat((1,2), spvec, diagmat, densevec, annospcmats[8], annodmats[4])) + # Test that concatenations strictly involving un/annotated dense matrices/vectors + # yield dense arrays + for densemata in (densemat, annodmats...) + # --> Tests applicable to pairs including only matrices + for densematb in (densemat, annodmats...) + @test !issparse(vcat(densemata, densematb)) + @test !issparse(vcat(densematb, densemata)) + end + # --> Tests applicable to pairs including vectors or matrices + for otherdense in (densevec, densemat, annodmats...) + @test !issparse(hcat(densemata, otherdense)) + @test !issparse(hcat(otherdense, densemata)) + @test !issparse(hvcat((2,), densemata, otherdense)) + @test !issparse(hvcat((2,), otherdense, densemata)) + @test !issparse(cat((1,2), densemata, otherdense)) + @test !issparse(cat((1,2), otherdense, densemata)) + end + end +end