Skip to content

Commit

Permalink
Merge pull request #17690 from Sacha0/sparsegencat
Browse files Browse the repository at this point in the history
Fix #17680 (sparse general concatenation) and #17686 (concatenation of heterogeneous Matrix-Vector combinations)
  • Loading branch information
JeffBezanson authored Jul 29, 2016
2 parents 6b41538 + 9a67c16 commit ba46baf
Show file tree
Hide file tree
Showing 6 changed files with 75 additions and 39 deletions.
2 changes: 1 addition & 1 deletion base/abstractarray.jl
Original file line number Diff line number Diff line change
Expand Up @@ -998,7 +998,7 @@ function cat_t(catdims, typeC::Type, X...)
end
end

C = similar(isa(X[1],AbstractArray) ? full(X[1]) : [X[1]], typeC, tuple(dimsC...))
C = similar(isa(X[1],AbstractArray) ? X[1] : [X[1]], typeC, tuple(dimsC...))
if length(catdims)>1
fill!(C,0)
end
Expand Down
48 changes: 20 additions & 28 deletions base/array.jl
Original file line number Diff line number Diff line change
Expand Up @@ -643,6 +643,26 @@ function reverse!(v::AbstractVector, s=1, n=length(v))
return 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])
for j = 2:length(V)
if length(V[j]) != height
throw(DimensionMismatch("vectors must have same lengths"))
end
end
return [ V[j][i]::T for i=1:length(V[1]), j=1:length(V) ]
end
function vcat{T}(arrays::Vector{T}...)
n = 0
for a in arrays
Expand Down Expand Up @@ -670,34 +690,6 @@ function vcat{T}(arrays::Vector{T}...)
return arr
end

function hcat{T}(V::Vector{T}...)
height = length(V[1])
for j = 2:length(V)
if length(V[j]) != height
throw(DimensionMismatch("vectors must have same lengths"))
end
end
return [ V[j][i]::T for i=1:length(V[1]), j=1:length(V) ]
end

hcat(A::Matrix...) = typed_hcat(promote_eltype(A...), A...)
hcat{T}(A::Matrix{T}...) = typed_hcat(T, A...)

vcat(A::Matrix...) = typed_vcat(promote_eltype(A...), A...)
vcat{T}(A::Matrix{T}...) = typed_vcat(T, A...)

hcat(A::Union{Matrix, Vector}...) = typed_hcat(promote_eltype(A...), A...)
hcat{T}(A::Union{Matrix{T}, Vector{T}}...) = typed_hcat(T, A...)


vcat(A::Union{Matrix, Vector}...) = typed_vcat(promote_eltype(A...), A...)
vcat{T}(A::Union{Matrix{T}, Vector{T}}...) = typed_vcat(T, A...)

hvcat(rows::Tuple{Vararg{Int}}, xs::Vector...) = typed_hvcat(promote_eltype(xs...), rows, xs...)
hvcat{T}(rows::Tuple{Vararg{Int}}, xs::Vector{T}...) = typed_hvcat(T, rows, xs...)

hvcat(rows::Tuple{Vararg{Int}}, xs::Matrix...) = typed_hvcat(promote_eltype(xs...), rows, xs...)
hvcat{T}(rows::Tuple{Vararg{Int}}, xs::Matrix{T}...) = typed_hvcat(T, rows, xs...)

## find ##

Expand Down
4 changes: 2 additions & 2 deletions base/sparse/sparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -21,10 +21,10 @@ import Base: @get!, acos, acosd, acot, acotd, acsch, asech, asin, asind, asinh,
tand, tanh, trace, transpose!, tril!, triu!, trunc, vecnorm, abs, abs2,
broadcast, ceil, complex, cond, conj, convert, copy, copy!, ctranspose, diagm,
exp, expm1, factorize, find, findmax, findmin, findnz, float, full, getindex,
hcat, hvcat, imag, indmax, ishermitian, kron, length, log, log1p, max, min,
vcat, hcat, hvcat, cat, imag, indmax, ishermitian, kron, length, log, log1p, max, min,
maximum, minimum, norm, one, promote_eltype, real, reinterpret, reshape, rot180,
rotl90, rotr90, round, scale!, setindex!, similar, size, transpose, tril,
triu, vcat, vec, permute!
triu, vec, permute!

import Base.Broadcast: broadcast_shape

Expand Down
6 changes: 6 additions & 0 deletions base/sparse/sparsematrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3253,6 +3253,12 @@ function hvcat(rows::Tuple{Vararg{Int}}, X::Union{Vector, Matrix, SparseMatrixCS
vcat(tmp_rows...)
end

function cat(catdims, Xin::Union{Vector, Matrix, SparseMatrixCSC}...)
X = SparseMatrixCSC[issparse(x) ? x : sparse(x) for x in Xin]
T = promote_eltype(Xin...)
Base.cat_t(catdims, T, X...)
end

"""
blkdiag(A...)
Expand Down
22 changes: 22 additions & 0 deletions test/arrayops.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1703,3 +1703,25 @@ for op in (:.+, :.*, :.÷, :.%, :.<<, :.>>, :.-, :./, :.\, :.//, :.^)
end

end

# Test that concatenations of dense matrices/vectors yield dense matrices/vectors
let
N = 4
densevec = ones(N)
densemat = diagm(ones(N))
# Test that concatenations of homogeneous pairs of either dense matrices or dense vectors
# (i.e., Matrix-Matrix concatenations, and Vector-Vector concatenations) yield dense arrays
for densearray in (densevec, densemat)
@test isa(vcat(densearray, densearray), Array)
@test isa(hcat(densearray, densearray), Array)
@test isa(hvcat((2,), densearray, densearray), Array)
@test isa(cat((1,2), densearray, densearray), Array)
end
# Test that concatenations of heterogeneous Matrix-Vector pairs yield dense matrices
@test isa(hcat(densemat, densevec), Array)
@test isa(hcat(densevec, densemat), Array)
@test isa(hvcat((2,), densemat, densevec), Array)
@test isa(hvcat((2,), densevec, densemat), Array)
@test isa(cat((1,2), densemat, densevec), Array)
@test isa(cat((1,2), densevec, densemat), Array)
end
32 changes: 24 additions & 8 deletions test/sparsedir/sparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1554,14 +1554,30 @@ end
@inferred sprand(1, 1, 1.0, rand, Float64)
@inferred sprand(1, 1, 1.0, x->round(Int,rand(x)*100))

# dense sparse concatenation -> sparse return type
@test issparse([sprand(10,10,.1) rand(10,10)])
@test issparse([sprand(10,10,.1); rand(10,10)])
@test issparse([sprand(10,10,.1) rand(10,10); rand(10,10) rand(10,10)])
# ---
@test !issparse([rand(10,10) rand(10,10)])
@test !issparse([rand(10,10); rand(10,10)])
@test !issparse([rand(10,10) rand(10,10); rand(10,10) rand(10,10)])
# Test that concatenations of combinations of sparse matrices with sparse matrices or dense
# matrices/vectors yield sparse arrays
let
N = 4
densevec = ones(N)
densemat = diagm(ones(N))
spmat = spdiagm(ones(N))
# Test that concatenations of pairs of sparse matrices yield sparse arrays
@test issparse(vcat(spmat, spmat))
@test issparse(hcat(spmat, spmat))
@test issparse(hvcat((2,), spmat, spmat))
@test issparse(cat((1,2), spmat, spmat))
# Test that concatenations of a sparse matrice with a dense matrix/vector yield sparse arrays
@test issparse(vcat(spmat, densemat))
@test issparse(vcat(densemat, spmat))
for densearg in (densevec, densemat)
@test issparse(hcat(spmat, densearg))
@test issparse(hcat(densearg, spmat))
@test issparse(hvcat((2,), spmat, densearg))
@test issparse(hvcat((2,), densearg, spmat))
@test issparse(cat((1,2), spmat, densearg))
@test issparse(cat((1,2), densearg, spmat))
end
end

# issue #14816
let m = 5
Expand Down

2 comments on commit ba46baf

@nanosoldier
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Executing the daily benchmark build, I will reply here when finished:

@nanosoldier runbenchmarks(ALL, isdaily = true)

@nanosoldier
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Your benchmark job has completed - possible performance regressions were detected. A full report can be found here. cc @jrevels

Please sign in to comment.