Skip to content

Commit 1f33cdd

Browse files
authored
LinAlg.fillslots! -> LinAlg.fillstored! (#25030)
* rename LinAlg.fillslots! LinAlg.fillstored! * rewrite fillstored!(::SpecialMatrix, x) without generated function move implementations to more suitable places * implement fillband! and call it from fillstored!(::AbstractTriangular) * review fixes
1 parent e09eca2 commit 1f33cdd

File tree

7 files changed

+92
-87
lines changed

7 files changed

+92
-87
lines changed

NEWS.md

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -471,8 +471,10 @@ Deprecated or removed
471471
`Matrix{Int}(uninitialized, (2, 4))`, and `Array{Float32,3}(11, 13, 17)` is now
472472
`Array{Float32,3}(uninitialized, 11, 13, 17)` ([#24781]).
473473

474+
* `LinAlg.fillslots!` has been renamed `LinAlg.fillstored!` ([#25030]).
475+
474476
* `fill!(A::Diagonal, x)` and `fill!(A::AbstractTriangular, x)` have been deprecated
475-
in favor of `Base.LinAlg.fillslots!(A, x)` ([#24413]).
477+
in favor of `Base.LinAlg.fillstored!(A, x)` ([#24413]).
476478

477479
* `eye` has been deprecated in favor of `I` and `Matrix` constructors. Please see the
478480
deprecation warnings for replacement details ([#24438]).

base/deprecated.jl

Lines changed: 5 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1882,8 +1882,11 @@ end
18821882
# Also un-comment the new definition in base/indices.jl
18831883

18841884
# deprecate odd fill! methods
1885-
@deprecate fill!(D::Diagonal, x) LinAlg.fillslots!(D, x)
1886-
@deprecate fill!(A::Base.LinAlg.AbstractTriangular, x) LinAlg.fillslots!(A, x)
1885+
@deprecate fill!(D::Diagonal, x) LinAlg.fillstored!(D, x)
1886+
@deprecate fill!(A::Base.LinAlg.AbstractTriangular, x) LinAlg.fillstored!(A, x)
1887+
1888+
# PR #25030
1889+
@eval LinAlg @deprecate fillslots! fillstored! false
18871890

18881891
function diagm(v::BitVector)
18891892
depwarn(string("diagm(v::BitVector) is deprecated, use diagm(0 => v) or ",

base/linalg/bidiag.jl

Lines changed: 0 additions & 34 deletions
Original file line numberDiff line numberDiff line change
@@ -636,37 +636,3 @@ function eigvecs(M::Bidiagonal{T}) where T
636636
Q #Actually Triangular
637637
end
638638
eigfact(M::Bidiagonal) = Eigen(eigvals(M), eigvecs(M))
639-
640-
# fill! methods
641-
_valuefields(::Type{<:Diagonal}) = [:diag]
642-
_valuefields(::Type{<:Bidiagonal}) = [:dv, :ev]
643-
_valuefields(::Type{<:Tridiagonal}) = [:dl, :d, :du]
644-
_valuefields(::Type{<:SymTridiagonal}) = [:dv, :ev]
645-
_valuefields(::Type{<:AbstractTriangular}) = [:data]
646-
647-
const SpecialArrays = Union{Diagonal,Bidiagonal,Tridiagonal,SymTridiagonal,AbstractTriangular}
648-
649-
function fillslots!(A::SpecialArrays, x)
650-
xT = convert(eltype(A), x)
651-
if @generated
652-
quote
653-
$([ :(fill!(A.$field, xT)) for field in _valuefields(A) ]...)
654-
end
655-
else
656-
for field in _valuefields(A)
657-
fill!(getfield(A, field), xT)
658-
end
659-
end
660-
return A
661-
end
662-
663-
_small_enough(A::Bidiagonal) = size(A, 1) <= 1
664-
_small_enough(A::Tridiagonal) = size(A, 1) <= 2
665-
_small_enough(A::SymTridiagonal) = size(A, 1) <= 2
666-
667-
function fill!(A::Union{Bidiagonal,Tridiagonal,SymTridiagonal}, x)
668-
xT = convert(eltype(A), x)
669-
(xT == zero(eltype(A)) || _small_enough(A)) && return fillslots!(A, xT)
670-
throw(ArgumentError("array A of type $(typeof(A)) and size $(size(A)) can
671-
not be filled with x=$x, since some of its entries are constrained."))
672-
end

base/linalg/dense.jl

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -242,6 +242,22 @@ function tril!(M::AbstractMatrix, k::Integer)
242242
end
243243
tril(M::Matrix, k::Integer) = tril!(copy(M), k)
244244

245+
"""
246+
fillband!(A::AbstractMatrix, x, l, u)
247+
248+
Fill the band between diagonals `l` and `u` with the value `x`.
249+
"""
250+
function fillband!(A::AbstractMatrix{T}, x, l, u) where T
251+
m, n = size(A)
252+
xT = convert(T, x)
253+
for j in 1:n
254+
for i in max(1,j-u):min(m,j-l)
255+
@inbounds A[i, j] = xT
256+
end
257+
end
258+
return A
259+
end
260+
245261
function diagind(m::Integer, n::Integer, k::Integer=0)
246262
if !(-m <= k <= n)
247263
throw(ArgumentError(string("requested diagonal, $k, must be at least $(-m) and ",

base/linalg/special.jl

Lines changed: 18 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -126,3 +126,21 @@ mul!(A::AbstractTriangular, adjB::Adjoint{<:Any,<:Union{QRCompactWYQ,QRPackedQ}}
126126
(B = adjB.parent; mul!(full!(A), Adjoint(B)))
127127
*(A::AbstractTriangular, adjB::Adjoint{<:Any,<:Union{QRCompactWYQ,QRPackedQ}}) =
128128
(B = adjB.parent; *(copy!(similar(parent(A)), A), Adjoint(B)))
129+
130+
# fill[stored]! methods
131+
fillstored!(A::Diagonal, x) = (fill!(A.diag, x); A)
132+
fillstored!(A::Bidiagonal, x) = (fill!(A.dv, x); fill!(A.ev, x); A)
133+
fillstored!(A::Tridiagonal, x) = (fill!(A.dl, x); fill!(A.d, x); fill!(A.du, x); A)
134+
fillstored!(A::SymTridiagonal, x) = (fill!(A.dv, x); fill!(A.ev, x); A)
135+
136+
_small_enough(A::Bidiagonal) = size(A, 1) <= 1
137+
_small_enough(A::Tridiagonal) = size(A, 1) <= 2
138+
_small_enough(A::SymTridiagonal) = size(A, 1) <= 2
139+
140+
# TODO: Add Diagonal to this method when 0.7 deprecations are removed
141+
function fill!(A::Union{Bidiagonal,Tridiagonal,SymTridiagonal}, x)
142+
xT = convert(eltype(A), x)
143+
(iszero(xT) || _small_enough(A)) && return fillstored!(A, xT)
144+
throw(ArgumentError("array of type $(typeof(A)) and size $(size(A)) can
145+
not be filled with $x, since some of its entries are constrained."))
146+
end

base/linalg/triangular.jl

Lines changed: 5 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -418,6 +418,11 @@ end
418418
scale!(A::Union{UpperTriangular,LowerTriangular}, c::Number) = scale!(A,A,c)
419419
scale!(c::Number, A::Union{UpperTriangular,LowerTriangular}) = scale!(A,c)
420420

421+
fillstored!(A::LowerTriangular, x) = (fillband!(A.data, x, 1-size(A,1), 0); A)
422+
fillstored!(A::UnitLowerTriangular, x) = (fillband!(A.data, x, 1-size(A,1), -1); A)
423+
fillstored!(A::UpperTriangular, x) = (fillband!(A.data, x, 0, size(A,2)-1); A)
424+
fillstored!(A::UnitUpperTriangular, x) = (fillband!(A.data, x, 1, size(A,2)-1); A)
425+
421426
# Binary operations
422427
+(A::UpperTriangular, B::UpperTriangular) = UpperTriangular(A.data + B.data)
423428
+(A::LowerTriangular, B::LowerTriangular) = LowerTriangular(A.data + B.data)

test/linalg/bidiag.jl

Lines changed: 45 additions & 50 deletions
Original file line numberDiff line numberDiff line change
@@ -310,58 +310,53 @@ end
310310
@test promote(C,A) isa Tuple{Tridiagonal, Tridiagonal}
311311
end
312312

313-
import Base.LinAlg: fillslots!, UnitLowerTriangular
314-
@testset "fill! and fillslots!" begin
315-
let #fill!
316-
let # fillslots!
317-
A = Tridiagonal(randn(2), randn(3), randn(2))
318-
@test fillslots!(A, 3) == Tridiagonal([3, 3.], [3, 3, 3.], [3, 3.])
319-
B = Bidiagonal(randn(3), randn(2), :U)
320-
@test fillslots!(B, 2) == Bidiagonal([2.,2,2], [2,2.], :U)
321-
S = SymTridiagonal(randn(3), randn(2))
322-
@test fillslots!(S, 1) == SymTridiagonal([1,1,1.], [1,1.])
323-
Ult = UnitLowerTriangular(randn(3,3))
324-
@test fillslots!(Ult, 3) == UnitLowerTriangular([1 0 0; 3 1 0; 3 3 1])
313+
using Base.LinAlg: fillstored!, UnitLowerTriangular
314+
@testset "fill! and fillstored!" begin
315+
let # fillstored!
316+
A = Tridiagonal(randn(2), randn(3), randn(2))
317+
@test fillstored!(A, 3) == Tridiagonal([3, 3], [3, 3, 3], [3, 3])
318+
B = Bidiagonal(randn(3), randn(2), :U)
319+
@test fillstored!(B, 2) == Bidiagonal([2,2,2], [2,2], :U)
320+
S = SymTridiagonal(randn(3), randn(2))
321+
@test fillstored!(S, 1) == SymTridiagonal([1,1,1], [1,1])
322+
Ult = UnitLowerTriangular(randn(3,3))
323+
@test fillstored!(Ult, 3) == UnitLowerTriangular([1 0 0; 3 1 0; 3 3 1])
324+
end
325+
let # fill!(exotic, 0)
326+
exotic_arrays = Any[Tridiagonal(randn(3), randn(4), randn(3)),
327+
Bidiagonal(randn(3), randn(2), rand([:U,:L])),
328+
SymTridiagonal(randn(3), randn(2)),
329+
sparse(randn(3,4)),
330+
# Diagonal(randn(5)), # Diagonal fill! deprecated, see below
331+
sparse(rand(3)),
332+
# LowerTriangular(randn(3,3)), # AbstractTriangular fill! deprecated, see below
333+
# UpperTriangular(randn(3,3)) # AbstractTriangular fill! deprecated, see below
334+
]
335+
for A in exotic_arrays
336+
@test iszero(fill!(A, 0))
325337
end
326-
let # fill!(exotic, 0)
327-
exotic_arrays = Any[Tridiagonal(randn(3), randn(4), randn(3)),
328-
Bidiagonal(randn(3), randn(2), rand([:U,:L])),
329-
SymTridiagonal(randn(3), randn(2)),
330-
sparse(randn(3,4)),
331-
# Diagonal(randn(5)), # Diagonal fill! deprecated, see below
332-
sparse(rand(3)),
333-
# LowerTriangular(randn(3,3)), # AbstractTriangular fill! deprecated, see below
334-
# UpperTriangular(randn(3,3)) # AbstractTriangular fill! deprecated, see below
335-
]
336-
for A in exotic_arrays
337-
fill!(A, 0)
338-
for a in A
339-
@test a == 0
340-
end
341-
end
342-
# Diagonal and AbstractTriangular fill! were defined as fillslots!,
343-
# not matching the general behavior of fill!, and so have been deprecated.
344-
# In a future dev cycle, these fill! methods should probably be reintroduced
345-
# with behavior matching that of fill! for other structured matrix types.
346-
# In the interm, equivalently test fillslots! below
347-
@test iszero(fillslots!(Diagonal(fill(1, 3)), 0))
348-
@test iszero(fillslots!(LowerTriangular(fill(1, 3, 3)), 0))
349-
@test iszero(fillslots!(UpperTriangular(fill(1, 3, 3)), 0))
338+
# Diagonal and AbstractTriangular fill! were defined as fillstored!,
339+
# not matching the general behavior of fill!, and so have been deprecated.
340+
# In a future dev cycle, these fill! methods should probably be reintroduced
341+
# with behavior matching that of fill! for other structured matrix types.
342+
# In the interm, equivalently test fillstored! below
343+
@test iszero(fillstored!(Diagonal(fill(1, 3)), 0))
344+
@test iszero(fillstored!(LowerTriangular(fill(1, 3, 3)), 0))
345+
@test iszero(fillstored!(UpperTriangular(fill(1, 3, 3)), 0))
346+
end
347+
let # fill!(small, x)
348+
val = randn()
349+
b = Bidiagonal(randn(1,1), :U)
350+
st = SymTridiagonal(randn(1,1))
351+
for x in (b, st)
352+
@test Array(fill!(x, val)) == fill!(Array(x), val)
350353
end
351-
let # fill!(small, x)
352-
val = randn()
353-
b = Bidiagonal(randn(1,1), :U)
354-
st = SymTridiagonal(randn(1,1))
355-
for x in (b, st)
356-
@test Array(fill!(x, val)) == fill!(Array(x), val)
357-
end
358-
b = Bidiagonal(randn(2,2), :U)
359-
st = SymTridiagonal(randn(3), randn(2))
360-
t = Tridiagonal(randn(3,3))
361-
for x in (b, t, st)
362-
@test_throws ArgumentError fill!(x, val)
363-
@test Array(fill!(x, 0)) == fill!(Array(x), 0)
364-
end
354+
b = Bidiagonal(randn(2,2), :U)
355+
st = SymTridiagonal(randn(3), randn(2))
356+
t = Tridiagonal(randn(3,3))
357+
for x in (b, t, st)
358+
@test_throws ArgumentError fill!(x, val)
359+
@test Array(fill!(x, 0)) == fill!(Array(x), 0)
365360
end
366361
end
367362
end

0 commit comments

Comments
 (0)