Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Check isbanded in conversions #289

Merged
merged 2 commits into from
Feb 7, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
38 changes: 31 additions & 7 deletions src/banded/BandedMatrix.jl
Original file line number Diff line number Diff line change
Expand Up @@ -103,20 +103,38 @@ for MAT in (:AbstractBandedMatrix, :AbstractMatrix, :AbstractArray)
end
end

function copyto_bandeddata!(B::BandedMatrix, M::AbstractMatrix)
copyto_bandeddata!(B, MemoryLayout(B), M, MemoryLayout(M))
end
function copyto_bandeddata!(B, ::BandedColumns{DenseColumnMajor},
M, ::Union{BandedColumns{DenseColumnMajor}, DiagonalLayout{DenseColumnMajor}})
copyto!(B.data, bandeddata(M))
return B
end
copyto_bandeddata!(B, @nospecialize(_), M, @nospecialize(_)) = copyto!(B, M)

@inline function _convert_common(T, Container, M)
B = BandedMatrix{T, Container}(undef, size(M), bandwidths(M))
copyto_bandeddata!(B, M)
end

function convert(::Type{BandedMatrix{<:,C,OneTo{Int}}}, M::AbstractMatrix) where {C}
Container = typeof(convert(C, similar(M, 0, 0)))
T = eltype(Container)
copyto!(BandedMatrix{T, Container}(undef, size(M), bandwidths(M)), M)
_convert_common(T, Container, M)
end

function convert(::Type{BandedMatrix{T,C,OneTo{Int}}}, M::AbstractMatrix) where {T, C}
@inline function _convert_common_container(T, C, M)
Container = typeof(convert(C, similar(M, T, 0, 0)))
copyto!(BandedMatrix{T, Container}(undef, size(M), bandwidths(M)), M)
_convert_common(T, Container, M)
end

function convert(::Type{BandedMatrix{T,C,OneTo{Int}}}, M::AbstractMatrix) where {T, C}
_convert_common_container(T, C, M)
end

function convert(::Type{BandedMatrix{T,C,OneTo{Int}}}, M::AbstractMatrix) where {T, C<:AbstractMatrix{T}}
Container = typeof(convert(C, similar(M, T, 0, 0)))
copyto!(BandedMatrix{T, Container}(undef, size(M), bandwidths(M)), M)
_convert_common_container(T, C, M)
end

convert(::Type{BandedMatrix{<:, C}}, M::AbstractMatrix) where {C} = convert(BandedMatrix{<:,C,OneTo{Int}}, M)
Expand Down Expand Up @@ -198,9 +216,9 @@ function BandedMatrix{T}(E::Eye, bnds::NTuple{2,Integer}) where T
end

BandedMatrix{T,C}(A::AbstractMatrix) where {T, C<:AbstractMatrix{T}} =
copyto!(BandedMatrix{T, C}(undef, size(A), bandwidths(A)), A)
copyto_bandeddata!(BandedMatrix{T, C}(undef, size(A), bandwidths(A)), A)
BandedMatrix{T}(A::AbstractMatrix) where T =
copyto!(BandedMatrix{T}(undef, size(A), bandwidths(A)), A)
copyto_bandeddata!(BandedMatrix{T}(undef, size(A), bandwidths(A)), A)


_BandedMatrix(_, A::AbstractMatrix) = BandedMatrix(A, bandwidths(A))
Expand All @@ -210,6 +228,12 @@ BandedMatrix(A::AbstractMatrix) = _BandedMatrix(MemoryLayout(A), A)
## specialised
# use bandeddata if possible
_BandedMatrix(::BandedColumns, A::AbstractMatrix) = _BandedMatrix(copy(bandeddata(A)), axes(A,1), bandwidths(A)...)
function _BandedMatrix(::DiagonalLayout, A::AbstractMatrix{T}) where T
m,n = size(A)
dat = Matrix{T}(undef, 1, n)
copyto!(vec(dat), diagonaldata(A))
_BandedMatrix(dat, m, 0, 0)
end
function _BandedMatrix(::BidiagonalLayout, A::AbstractMatrix{T}) where T
m,n = size(A)
dat = Matrix{T}(undef, 2, n)
Expand Down
4 changes: 3 additions & 1 deletion src/interfaceimpl.jl
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,9 @@ isbanded(::Diagonal) = true
bandwidths(::Diagonal) = (0,0)
inbands_getindex(D::Diagonal, k::Integer, j::Integer) = D.diag[k]
inbands_setindex!(D::Diagonal, v, k::Integer, j::Integer) = (D.diag[k] = v)
bandeddata(D::Diagonal) = permutedims(D.diag)
function bandeddata(D::Union{Diagonal, Transpose{<:Any, <:Diagonal}, Adjoint{<:Any,<:Diagonal}})
permutedims(diagonaldata(D))
end

# treat subinds as banded
sublayout(::DiagonalLayout{L}, inds::Type) where L = sublayout(bandedcolumns(L()), inds)
Expand Down
13 changes: 13 additions & 0 deletions test/test_banded.jl
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,19 @@ Base.similar(::MyMatrix, ::Type{T}, m::Int, n::Int) where T = MyMatrix{T}(undef,
@test BandedMatrix{Int64, typeof(matrix)}(matrix, (1, 1)) isa BandedMatrix{Int64, typeof(matrix)}
@test BandedMatrix{Int64, typeof(matrix)}(matrix, (1, 1)).data isa typeof(matrix)
@test_throws UndefRefError BandedMatrix{Vector{Float64}}(undef, (5,5), (1,1))[1,1]

@testset "Construction from Diagonal" begin
D = Diagonal(1:4)
B1 = BandedMatrix(D)
@test B1 isa BandedMatrix{Int}
@test B1 == D
B2 = BandedMatrix{Float64}(D)
@test B2 isa BandedMatrix{Float64}
@test B2 == D
B3 = BandedMatrix{Float32,Matrix{Float32}}(D)
@test B3 isa BandedMatrix{Float32,Matrix{Float32}}
@test B3 == D
end
end

@testset "BandedMatrix arithmetic" begin
Expand Down
4 changes: 3 additions & 1 deletion test/test_interface.jl
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
using BandedMatrices, LinearAlgebra, ArrayLayouts, FillArrays, Test, Base64
import BandedMatrices: banded_mul!, isbanded, AbstractBandedLayout, BandedStyle,
rowsupport, colsupport, _BandedMatrix, BandedColumns
rowsupport, colsupport, _BandedMatrix, BandedColumns, bandeddata
import ArrayLayouts: OnesLayout, UnknownLayout

struct PseudoBandedMatrix{T} <: AbstractMatrix{T}
Expand Down Expand Up @@ -70,6 +70,8 @@ LinearAlgebra.fill!(A::PseudoBandedMatrix, v) = fill!(A.data,v)
@test A[1,2] == 0
@test BandedMatrices.@inbands(A[1,2]) == 2

@test bandeddata(A) == bandeddata(Adjoint(A)) == bandeddata(Transpose(A)) == diag(A)'

@test MemoryLayout(view(A, 1:3,2:4)) isa BandedColumns{DenseColumnMajor}
@test MemoryLayout(view(A, [1,2,3],2:4)) isa UnknownLayout

Expand Down