Skip to content

Commit 6ccc974

Browse files
authored
Add sparse!(I, J, V) and spzeros!(I, J) (#397)
This patch adds convenience wrappers around the "full" calls to `sparse!` and `spzeros` that allocates the intermediate arrays of appropriate size, but reuses the input vectors for the final matrix storage. These methods are useful in the quite common case where the COO representation is constructed with the sole purpose of constructing the sparse matrix and then thrown away.
1 parent 2c8b8b1 commit 6ccc974

File tree

3 files changed

+101
-0
lines changed

3 files changed

+101
-0
lines changed

docs/src/index.md

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -211,12 +211,14 @@ SparseArrays.AbstractSparseMatrix
211211
SparseArrays.SparseVector
212212
SparseArrays.SparseMatrixCSC
213213
SparseArrays.sparse
214+
SparseArrays.sparse!
214215
SparseArrays.sparsevec
215216
Base.similar(::SparseArrays.AbstractSparseMatrixCSC, ::Type)
216217
SparseArrays.issparse
217218
SparseArrays.nnz
218219
SparseArrays.findnz
219220
SparseArrays.spzeros
221+
SparseArrays.spzeros!
220222
SparseArrays.spdiagm
221223
SparseArrays.sparse_hcat
222224
SparseArrays.sparse_vcat

src/sparsematrix.jl

Lines changed: 53 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1270,6 +1270,34 @@ function sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti},
12701270
Vector{Ti}(undef, n+1), Vector{Ti}(), Vector{Tv}())
12711271
end
12721272

1273+
"""
1274+
SparseArrays.sparse!(I, J, V, [m, n, combine]) -> SparseMatrixCSC
1275+
1276+
Variant of `sparse!` that re-uses the input vectors (`I`, `J`, `V`) for the final matrix
1277+
storage. After construction the input vectors will alias the matrix buffers; `S.colptr ===
1278+
I`, `S.rowval === J`, and `S.nzval === V` holds, and they will be `resize!`d as necessary.
1279+
1280+
Note that some work buffers will still be allocated. Specifically, this method is a
1281+
convenience wrapper around `sparse!(I, J, V, m, n, combine, klasttouch, csrrowptr,
1282+
csrcolval, csrnzval, csccolptr, cscrowval, cscnzval)` where this method allocates
1283+
`klasttouch`, `csrrowptr`, `csrcolval`, and `csrnzval` of appropriate size, but reuses `I`,
1284+
`J`, and `V` for `csccolptr`, `cscrowval`, and `cscnzval`.
1285+
1286+
Arguments `m`, `n`, and `combine` defaults to `maximum(I)`, `maximum(J)`, and `+`,
1287+
respectively.
1288+
1289+
!!! compat "Julia 1.10"
1290+
This method requires Julia version 1.10 or later.
1291+
"""
1292+
function sparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti}, V::AbstractVector{Tv},
1293+
m::Integer=dimlub(I), n::Integer=dimlub(J), combine::Function=+) where {Tv, Ti<:Integer}
1294+
klasttouch = Vector{Ti}(undef, n)
1295+
csrrowptr = Vector{Ti}(undef, m + 1)
1296+
csrcolval = Vector{Ti}(undef, length(I))
1297+
csrnzval = Vector{Tv}(undef, length(I))
1298+
sparse!(I, J, V, Int(m), Int(n), combine, klasttouch, csrrowptr, csrcolval, csrnzval, I, J, V)
1299+
end
1300+
12731301
dimlub(I) = isempty(I) ? 0 : Int(maximum(I)) #least upper bound on required sparse matrix dimension
12741302

12751303
sparse(I,J,v::Number) = sparse(I, J, fill(v,length(I)))
@@ -2116,6 +2144,31 @@ function spzeros!(::Type{Tv}, I::AbstractVector{Ti}, J::AbstractVector{Ti}, m::I
21162144
csrrowptr, csrcolval, cscnzval, csccolptr, cscrowval, cscnzval)
21172145
end
21182146

2147+
"""
2148+
SparseArrays.spzeros!(::Type{Tv}, I, J, [m, n]) -> SparseMatrixCSC{Tv}
2149+
2150+
Variant of `spzeros!` that re-uses the input vectors `I` and `J` for the final matrix
2151+
storage. After construction the input vectors will alias the matrix buffers; `S.colptr ===
2152+
I` and `S.rowval === J` holds, and they will be `resize!`d as necessary.
2153+
2154+
Note that some work buffers will still be allocated. Specifically, this method is a
2155+
convenience wrapper around `spzeros!(Tv, I, J, m, n, klasttouch, csrrowptr, csrcolval,
2156+
csccolptr, cscrowval)` where this method allocates `klasttouch`, `csrrowptr`, and
2157+
`csrcolval` of appropriate size, but reuses `I` and `J` for `csccolptr` and `cscrowval`.
2158+
2159+
Arguments `m` and `n` defaults to `maximum(I)` and `maximum(J)`.
2160+
2161+
!!! compat "Julia 1.10"
2162+
This method requires Julia version 1.10 or later.
2163+
"""
2164+
function spzeros!(::Type{Tv}, I::AbstractVector{Ti}, J::AbstractVector{Ti},
2165+
m::Integer=dimlub(I), n::Integer=dimlub(J)) where {Tv, Ti <: Integer}
2166+
klasttouch = Vector{Ti}(undef, n)
2167+
csrrowptr = Vector{Ti}(undef, m + 1)
2168+
csrcolval = Vector{Ti}(undef, length(I))
2169+
return spzeros!(Tv, I, J, Int(m), Int(n), klasttouch, csrrowptr, csrcolval, I, J)
2170+
end
2171+
21192172
import Base._one
21202173
function Base._one(unit::T, S::AbstractSparseMatrixCSC) where T
21212174
size(S, 1) == size(S, 2) || throw(DimensionMismatch("multiplicative identity only defined for square matrices"))

test/sparsematrix_constructors_indexing.jl

Lines changed: 46 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1751,6 +1751,52 @@ end
17511751
@test getcolptr(S!) === I
17521752
@test getrowval(S!) === J
17531753
@test nonzeros(S!) === V
1754+
1755+
# Test reuse of I, J, V for the matrix buffers in
1756+
# sparse!(I, J, V), sparse!(I, J, V, m, n), sparse!(I, J, V, m, n, combine),
1757+
# spzeros!(T, I, J), and spzeros!(T, I, J, m, n).
1758+
I, J, V = allocate_arrays(m, n)
1759+
S = sparse(I, J, V)
1760+
S! = sparse!(I, J, V)
1761+
@test S == S!
1762+
@test same_structure(S, S!)
1763+
@test getcolptr(S!) === I
1764+
@test getrowval(S!) === J
1765+
@test nonzeros(S!) === V
1766+
I, J, V = allocate_arrays(m, n)
1767+
S = sparse(I, J, V, 2m, 2n)
1768+
S! = sparse!(I, J, V, 2m, 2n)
1769+
@test S == S!
1770+
@test same_structure(S, S!)
1771+
@test getcolptr(S!) === I
1772+
@test getrowval(S!) === J
1773+
@test nonzeros(S!) === V
1774+
I, J, V = allocate_arrays(m, n)
1775+
S = sparse(I, J, V, 2m, 2n, *)
1776+
S! = sparse!(I, J, V, 2m, 2n, *)
1777+
@test S == S!
1778+
@test same_structure(S, S!)
1779+
@test getcolptr(S!) === I
1780+
@test getrowval(S!) === J
1781+
@test nonzeros(S!) === V
1782+
for T in (Float32, Float64)
1783+
I, J, = allocate_arrays(m, n)
1784+
S = spzeros(T, I, J)
1785+
S! = spzeros!(T, I, J)
1786+
@test S == S!
1787+
@test same_structure(S, S!)
1788+
@test eltype(S) == eltype(S!) == T
1789+
@test getcolptr(S!) === I
1790+
@test getrowval(S!) === J
1791+
I, J, = allocate_arrays(m, n)
1792+
S = spzeros(T, I, J, 2m, 2n)
1793+
S! = spzeros!(T, I, J, 2m, 2n)
1794+
@test S == S!
1795+
@test same_structure(S, S!)
1796+
@test eltype(S) == eltype(S!) == T
1797+
@test getcolptr(S!) === I
1798+
@test getrowval(S!) === J
1799+
end
17541800
end
17551801
end
17561802

0 commit comments

Comments
 (0)