Skip to content

Commit

Permalink
Make concatenations involving combinations of special matrices with s…
Browse files Browse the repository at this point in the history
…pecial matrices, sparse matrices, or dense matrices/vectors yield sparse arrays.
  • Loading branch information
Sacha0 committed Jul 29, 2016
1 parent ba46baf commit b2523ae
Show file tree
Hide file tree
Showing 2 changed files with 47 additions and 5 deletions.
14 changes: 9 additions & 5 deletions base/sparse/sparsematrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3229,19 +3229,23 @@ function hcat(X::SparseMatrixCSC...)
end


# Sparse/dense concatenation
# Sparse/special/dense concatenation

function hcat(Xin::Union{Vector, Matrix, SparseMatrixCSC}...)
# 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, SparseMatrixCSC, SpecialArrays}...)
X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin]
hcat(X...)
end

function vcat(Xin::Union{Vector, Matrix, SparseMatrixCSC}...)
function vcat(Xin::Union{Vector, Matrix, SparseMatrixCSC, SpecialArrays}...)
X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin]
vcat(X...)
end

function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseMatrixCSC}...)
function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseMatrixCSC, SpecialArrays}...)
nbr = length(rows) # number of block rows

tmp_rows = Array{SparseMatrixCSC}(nbr)
Expand All @@ -3253,7 +3257,7 @@ function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseMatrixCS
vcat(tmp_rows...)
end

function cat(catdims, Xin::Union{Vector, Matrix, SparseMatrixCSC}...)
function cat(catdims, Xin::Union{Vector, Matrix, SparseMatrixCSC, SpecialArrays}...)
X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin]
T = promote_eltype(Xin...)
Base.cat_t(catdims, T, X...)
Expand Down
38 changes: 38 additions & 0 deletions test/linalg/special.jl
Original file line number Diff line number Diff line change
Expand Up @@ -128,3 +128,41 @@ for typ in [UpperTriangular,LowerTriangular,Base.LinAlg.UnitUpperTriangular,Base
@test Base.LinAlg.A_mul_Bc(atri,qrb[:Q]) full(atri) * qrb[:Q]'
@test Base.LinAlg.A_mul_Bc!(copy(atri),qrb[:Q]) full(atri) * qrb[:Q]'
end

# Test that concatenations of combinations of special and other matrix types yield sparse arrays
let
N = 4
# Test concatenating pairwise combinations of special matrices
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))
specialmats = (diagmat, bidiagmat, tridiagmat, symtridiagmat)
for specialmata in specialmats, specialmatb in specialmats
@test issparse(hcat(specialmata, specialmatb))
@test issparse(vcat(specialmata, specialmatb))
@test issparse(hvcat((1,1), specialmata, specialmatb))
@test issparse(cat((1,2), specialmata, specialmatb))
end
# Test concatenating pairwise combinations of special matrices with sparse matrices,
# dense matrices, or dense vectors
densevec = ones(N)
densemat = diagm(ones(N))
spmat = spdiagm(ones(N))
for specialmat in specialmats
# --> Tests applicable only to pairs of matrices
for othermat in (spmat, densemat)
@test issparse(vcat(specialmat, othermat))
@test issparse(vcat(othermat, specialmat))
end
# --> Tests applicable also to pairs including vectors
for specialmat in specialmats, othermatorvec in (spmat, densemat, densevec)
@test issparse(hcat(specialmat, othermatorvec))
@test issparse(hcat(othermatorvec, specialmat))
@test issparse(hvcat((2,), specialmat, othermatorvec))
@test issparse(hvcat((2,), othermatorvec, specialmat))
@test issparse(cat((1,2), specialmat, othermatorvec))
@test issparse(cat((1,2), othermatorvec, specialmat))
end
end
end

0 comments on commit b2523ae

Please sign in to comment.