Skip to content

Commit

Permalink
deprecate diagm(A::SparseMatrixCSC) in favor of diagm(sparsevec(A))
Browse files Browse the repository at this point in the history
  • Loading branch information
fredrikekre committed Aug 18, 2017
1 parent feaa2f6 commit c8e505c
Show file tree
Hide file tree
Showing 5 changed files with 37 additions and 50 deletions.
3 changes: 3 additions & 0 deletions base/deprecated.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1687,6 +1687,9 @@ export hex2num
# issue #5148, PR #23259
# warning for `const` on locals should be changed to an error in julia-syntax.scm

# PR 23341
@deprecate diagm(A::SparseMatrixCSC) diagm(sparsevec(A))

# END 0.7 deprecations

# BEGIN 1.0 deprecations
Expand Down
43 changes: 0 additions & 43 deletions base/sparse/sparsematrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3413,49 +3413,6 @@ function trace(A::SparseMatrixCSC{Tv}) where Tv
return s
end

function diagm(v::SparseMatrixCSC{Tv,Ti}) where {Tv,Ti}
if size(v,1) != 1 && size(v,2) != 1
throw(DimensionMismatch("input should be nx1 or 1xn"))
end

n = length(v)
numnz = nnz(v)
colptr = Vector{Ti}(n+1)
rowval = Vector{Ti}(numnz)
nzval = Vector{Tv}(numnz)

if size(v,1) == 1
copy!(colptr, 1, v.colptr, 1, n+1)
ptr = 1
for col = 1:n
if colptr[col] != colptr[col+1]
rowval[ptr] = col
nzval[ptr] = v.nzval[ptr]
ptr += 1
end
end
else
copy!(rowval, 1, v.rowval, 1, numnz)
copy!(nzval, 1, v.nzval, 1, numnz)
colptr[1] = 1
ptr = 1
col = 1
while col <= n && ptr <= numnz
while rowval[ptr] > col
colptr[col+1] = colptr[col]
col += 1
end
colptr[col+1] = colptr[col] + 1
ptr += 1
col += 1
end
if col <= n
colptr[(col+1):(n+1)] = colptr[col]
end
end

return SparseMatrixCSC(n, n, colptr, rowval, nzval)
end

# Sort all the indices in each column of a CSC sparse matrix
# sortSparseMatrixCSC!(A, sortindices = :sortcols) # Sort each column with sort()
Expand Down
27 changes: 27 additions & 0 deletions base/sparse/sparsevector.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2004,3 +2004,30 @@ function fill!(A::Union{SparseVector, SparseMatrixCSC}, x)
end
return A
end

function diagm(v::SparseVector{Tv,Ti}) where {Tv,Ti}
n = length(v)
numnz = nnz(v)
colptr = Vector{Ti}(n+1)
rowval = Vector{Ti}(numnz)
nzval = Vector{Tv}(numnz)

copy!(rowval, 1, v.nzind, 1, numnz)
copy!(nzval, 1, v.nzval, 1, numnz)
colptr[1] = 1
ptr = 1
col = 1
while col <= n && ptr <= numnz
while rowval[ptr] > col
colptr[col+1] = colptr[col]
col += 1
end
colptr[col+1] = colptr[col] + 1
ptr += 1
col += 1
end
if col <= n
colptr[(col+1):(n+1)] = colptr[col]
end
return SparseMatrixCSC(n, n, colptr, rowval, nzval)
end
7 changes: 0 additions & 7 deletions test/sparse/sparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1316,13 +1316,6 @@ end
@test trace(speye(5)) == 5
end

@testset "diagm on a matrix" begin
@test_throws DimensionMismatch diagm(sparse(ones(5,2)))
@test_throws DimensionMismatch diagm(sparse(ones(2,5)))
@test diagm(sparse(ones(1,5))) == speye(5)
@test diagm(sparse(ones(5,1))) == speye(5)
end

@testset "diag" begin
for T in (Float64, Complex128)
S1 = sprand(T, 5, 5, 0.5)
Expand Down
7 changes: 7 additions & 0 deletions test/sparse/sparsevector.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1150,3 +1150,10 @@ end
@testset "spzeros with index type" begin
@test typeof(spzeros(Float32, Int16, 3)) == SparseVector{Float32,Int16}
end

@testset "diagm" begin
v = sprand(10, 0.4)
@test diagm(v)::SparseMatrixCSC == diagm(Vector(v))
@test diagm(sparse(ones(5)))::SparseMatrixCSC == speye(5)
@test diagm(sparse(zeros(5)))::SparseMatrixCSC == spzeros(5,5)
end

0 comments on commit c8e505c

Please sign in to comment.