Skip to content

Commit a0ada0f

Browse files
committed
at-deprecate spones(A) fillstored!(copy(A), 1)
1 parent bdae975 commit a0ada0f

File tree

10 files changed

+18
-38
lines changed

10 files changed

+18
-38
lines changed

NEWS.md

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -704,6 +704,9 @@ Deprecated or removed
704704
have been deprecated in favor of `spdiagm(d => x)` and `spdiagm(d[1] => x[1], d[2] => x[2], ...)`
705705
respectively. The new `spdiagm` implementation now always returns a square matrix ([#23757]).
706706

707+
* `spones(A::AbstractSparseArray{T})` has been deprecated in favor of
708+
`LinAlg.fillstored!(copy(A), one(T))` ([#TBD]).
709+
707710
* Constructors for `LibGit2.UserPasswordCredentials` and `LibGit2.SSHCredentials` which take a
708711
`prompt_if_incorrect` argument are deprecated. Instead, prompting behavior is controlled using
709712
the `allow_prompt` keyword in the `LibGit2.CredentialPayload` constructor ([#23690]).

base/deprecated.jl

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1861,6 +1861,12 @@ end
18611861
# PR #25030
18621862
@eval LinAlg @deprecate fillslots! fillstored! false
18631863

1864+
# PR #TBD
1865+
@eval SparseArrays @deprecate spones(A::SparseMatrixCSC{T}) where {T} fillstored!(copy(A), one(T))
1866+
@eval SparseArrays @deprecate spones(A::SparseVector{T}) where {T} fillstored!(copy(A), one(T))
1867+
using .SparseArrays.spones
1868+
export spones
1869+
18641870
function diagm(v::BitVector)
18651871
depwarn(string("diagm(v::BitVector) is deprecated, use diagm(0 => v) or ",
18661872
"BitMatrix(Diagonal(v)) instead"), :diagm)

base/exports.jl

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1258,7 +1258,6 @@ export
12581258
sparse,
12591259
sparsevec,
12601260
spdiagm,
1261-
spones,
12621261
sprand,
12631262
sprandn,
12641263
spzeros,

base/sparse/sparse.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -7,7 +7,7 @@ module SparseArrays
77

88
using Base: ReshapedArray, promote_op, setindex_shape_check, to_shape, tail
99
using Base.Sort: Forward
10-
using Base.LinAlg: AbstractTriangular, PosDefException
10+
using Base.LinAlg: AbstractTriangular, PosDefException, fillstored!
1111

1212
import Base: +, -, *, \, /, &, |, xor, ==
1313
import Base: A_mul_B!, Ac_mul_B, Ac_mul_B!, At_mul_B, At_mul_B!
@@ -31,7 +31,7 @@ import Base: @get!, acos, acosd, acot, acotd, acsch, asech, asin, asind, asinh,
3131

3232
export AbstractSparseArray, AbstractSparseMatrix, AbstractSparseVector,
3333
SparseMatrixCSC, SparseVector, blkdiag, droptol!, dropzeros!, dropzeros,
34-
issparse, nonzeros, nzrange, rowvals, sparse, sparsevec, spdiagm, spones,
34+
issparse, nonzeros, nzrange, rowvals, sparse, sparsevec, spdiagm, fillstored!,
3535
sprand, sprandn, spzeros, nnz, permute
3636

3737
include("abstractsparse.jl")

base/sparse/sparsematrix.jl

Lines changed: 1 addition & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -1431,31 +1431,7 @@ julia> sprandn(rng, 2, 2, 0.75)
14311431
sprandn(r::AbstractRNG, m::Integer, n::Integer, density::AbstractFloat) = sprand(r,m,n,density,randn,Float64)
14321432
sprandn(m::Integer, n::Integer, density::AbstractFloat) = sprandn(GLOBAL_RNG,m,n,density)
14331433

1434-
"""
1435-
spones(S)
1436-
1437-
Create a sparse array with the same structure as that of `S`, but with every nonzero
1438-
element having the value `1.0`.
1439-
1440-
# Examples
1441-
```jldoctest
1442-
julia> A = sparse([1,2,3,4],[2,4,3,1],[5.,4.,3.,2.])
1443-
4×4 SparseMatrixCSC{Float64,Int64} with 4 stored entries:
1444-
[4, 1] = 2.0
1445-
[1, 2] = 5.0
1446-
[3, 3] = 3.0
1447-
[2, 4] = 4.0
1448-
1449-
julia> spones(A)
1450-
4×4 SparseMatrixCSC{Float64,Int64} with 4 stored entries:
1451-
[4, 1] = 1.0
1452-
[1, 2] = 1.0
1453-
[3, 3] = 1.0
1454-
[2, 4] = 1.0
1455-
```
1456-
"""
1457-
spones(S::SparseMatrixCSC{T}) where {T} =
1458-
SparseMatrixCSC(S.m, S.n, copy(S.colptr), copy(S.rowval), ones(T, S.colptr[end]-1))
1434+
LinAlg.fillstored!(S::SparseMatrixCSC, x) = (fill!(S.nzval, x); S)
14591435

14601436
"""
14611437
spzeros([type,]m[,n])

base/sparse/sparsevector.jl

Lines changed: 1 addition & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -99,8 +99,7 @@ spzeros(len::Integer) = spzeros(Float64, len)
9999
spzeros(::Type{T}, len::Integer) where {T} = SparseVector(len, Int[], T[])
100100
spzeros(::Type{Tv}, ::Type{Ti}, len::Integer) where {Tv,Ti<:Integer} = SparseVector(len, Ti[], Tv[])
101101

102-
# Construction of same structure, but with all ones
103-
spones(x::SparseVector{T}) where {T} = SparseVector(x.n, copy(x.nzind), ones(T, length(x.nzval)))
102+
LinAlg.fillstored!(x::SparseVector{T}, y) where {T} = (fill!(x.nzval, y); x)
104103

105104
### Construction from lists of indices and values
106105

doc/src/manual/arrays.md

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -883,7 +883,6 @@ section of the standard library reference.
883883
| Sparse | Dense | Description |
884884
|:-------------------------- |:---------------------- |:--------------------------------------------------------------------------------------------------------------------------------------------------------------------- |
885885
| [`spzeros(m,n)`](@ref) | [`zeros(m,n)`](@ref) | Creates a *m*-by-*n* matrix of zeros. ([`spzeros(m,n)`](@ref) is empty.) |
886-
| [`spones(S)`](@ref) | [`ones(m,n)`](@ref) | Creates a matrix filled with ones. Unlike the dense version, [`spones`](@ref) has the same sparsity pattern as *S*. |
887886
| [`sparse(I, n, n)`](@ref) | [`Matrix(I,n,n)`](@ref)| Creates a *n*-by-*n* identity matrix. |
888887
| [`Array(S)`](@ref) | [`sparse(A)`](@ref) | Interconverts between dense and sparse formats. |
889888
| [`sprand(m,n,d)`](@ref) | [`rand(m,n)`](@ref) | Creates a *m*-by-*n* random matrix (of density *d*) with iid non-zero elements distributed uniformly on the half-open interval ``[0, 1)``. |

doc/src/stdlib/arrays.md

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -198,7 +198,6 @@ Base.SparseArrays.sparsevec
198198
Base.SparseArrays.issparse
199199
Base.SparseArrays.nnz
200200
Base.SparseArrays.spzeros
201-
Base.SparseArrays.spones
202201
Base.SparseArrays.spdiagm
203202
Base.SparseArrays.sprand
204203
Base.SparseArrays.sprandn

test/sparse/sparse.jl

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1737,8 +1737,8 @@ end
17371737
end
17381738
end
17391739

1740-
@testset "spones" begin
1741-
@test spones(sparse(2.0I, 5, 5)) == Matrix(I, 5, 5)
1740+
@testset "fillstored!" begin
1741+
@test SparseArrays.fillstored!(sparse(2.0I, 5, 5), 1) == Matrix(I, 5, 5)
17421742
end
17431743

17441744
@testset "factorization" begin
@@ -2014,7 +2014,7 @@ end
20142014
sprand(5, 5, 1/5)
20152015
end
20162016
A = max.(A, A')
2017-
A = spones(A)
2017+
SparseArrays.fillstored!(A, 1)
20182018
B = A[5:-1:1, 5:-1:1]
20192019
@test issymmetric(B)
20202020
end

test/sparse/sparsevector.jl

Lines changed: 2 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -119,10 +119,9 @@ end
119119
@test exact_equal(sparsevec(d), SparseVector(3, [1, 2, 3], [0.0, 1.0, 2.0]))
120120
end
121121
end
122-
@testset "spones" begin
123-
# copies structure, but replaces nzvals with ones
122+
@testset "fillstored!" begin
124123
x = SparseVector(8, [2, 3, 6], [12.0, 18.0, 25.0])
125-
y = spones(x)
124+
y = SparseArrays.fillstored!(copy(x), 1)
126125
@test (x .!= 0) == (y .!= 0)
127126
@test y == SparseVector(8, [2, 3, 6], [1.0, 1.0, 1.0])
128127
end

0 commit comments

Comments
 (0)