Skip to content

Array layouts v0.6 #70

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

Merged
merged 12 commits into from
Feb 18, 2021
10 changes: 5 additions & 5 deletions Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "InfiniteLinearAlgebra"
uuid = "cde9dba0-b1de-11e9-2c62-0bab9446c55c"
version = "0.5.0"
version = "0.5.1"

[deps]
ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a"
Expand All @@ -16,16 +16,16 @@ MatrixFactorizations = "a3b82374-2e81-5b9e-98ce-41277c0e4c87"
SemiseparableMatrices = "f8ebbe35-cbfb-4060-bf7f-b10e4670cf57"

[compat]
ArrayLayouts = "0.5.1, 0.6"
ArrayLayouts = "0.6"
BandedMatrices = "0.16"
BlockArrays = "0.14.5"
BlockBandedMatrices = "0.10"
FillArrays = "0.11"
InfiniteArrays = "0.10"
LazyArrays = "0.20.2"
LazyBandedMatrices = "0.4.5"
MatrixFactorizations = "0.7.1, 0.8"
SemiseparableMatrices = "0.2.2"
LazyBandedMatrices = "0.5"
MatrixFactorizations = "0.8"
SemiseparableMatrices = "0.2.3"
julia = "1.5"

[extras]
Expand Down
4 changes: 2 additions & 2 deletions src/InfiniteLinearAlgebra.jl
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,7 @@ using BlockArrays, BlockBandedMatrices, BandedMatrices, LazyArrays, LazyBandedMa

import Base: +, -, *, /, \, ^, OneTo, getindex, promote_op, _unsafe_getindex, size, axes,
AbstractMatrix, AbstractArray, Matrix, Array, Vector, AbstractVector, Slice,
show, getproperty, copy, map, require_one_based_indexing, similar, inv
show, getproperty, copy, copyto!, map, require_one_based_indexing, similar, inv
import Base.Broadcast: BroadcastStyle, Broadcasted, broadcasted

import ArrayLayouts: colsupport, rowsupport, triangularlayout, MatLdivVec, triangulardata, TriangularLayout, TridiagonalLayout,
Expand All @@ -27,7 +27,7 @@ import BlockArrays: AbstractBlockVecOrMat, sizes_from_blocks, _length, BlockedUn

import BandedMatrices: BandedMatrix, bandwidths, AbstractBandedLayout, _banded_qr!, _banded_qr, _BandedMatrix

import LazyBandedMatrices: ApplyBandedLayout, BroadcastBandedLayout, _krontrav_axes, _block_interlace_axes
import LazyBandedMatrices: ApplyBandedLayout, BroadcastBandedLayout, _krontrav_axes, _block_interlace_axes, LazyBandedLayout,AbstractLazyBandedLayout

import BlockBandedMatrices: _BlockSkylineMatrix, _BandedMatrix, _BlockSkylineMatrix, blockstart, blockstride,
BlockSkylineSizes, BlockSkylineMatrix, BlockBandedMatrix, _BlockBandedMatrix, BlockTridiagonal,
Expand Down
79 changes: 55 additions & 24 deletions src/banded/infbanded.jl
Original file line number Diff line number Diff line change
Expand Up @@ -144,47 +144,47 @@ end

for op in (:-, :+)
@eval begin
function $op(A::SymTriPertToeplitz{T}, λ::UniformScaling) where T
function $op(A::SymTriPertToeplitz{T}, λ::UniformScaling) where T
TV = promote_type(T,eltype(λ))
dv = Vcat(convert.(AbstractVector{TV}, A.dv.args)...)
ev = Vcat(convert.(AbstractVector{TV}, A.ev.args)...)
SymTridiagonal(broadcast($op, dv, Ref(λ.λ)), ev)
end
function $op(λ::UniformScaling, A::SymTriPertToeplitz{V}) where V
TV = promote_type(eltype(λ),V)
SymTridiagonal(convert(AbstractVector{TV}, broadcast($op, Ref(λ.λ), A.dv)),
SymTridiagonal(convert(AbstractVector{TV}, broadcast($op, Ref(λ.λ), A.dv)),
convert(AbstractVector{TV}, broadcast($op, A.ev)))
end
function $op(A::SymTridiagonal{T,<:AbstractFill}, λ::UniformScaling) where T
function $op(A::SymTridiagonal{T,<:AbstractFill}, λ::UniformScaling) where T
TV = promote_type(T,eltype(λ))
SymTridiagonal(convert(AbstractVector{TV}, broadcast($op, A.dv, Ref(λ.λ))),
convert(AbstractVector{TV}, A.ev))
end

function $op(A::TriPertToeplitz{T}, λ::UniformScaling) where T
function $op(A::TriPertToeplitz{T}, λ::UniformScaling) where T
TV = promote_type(T,eltype(λ))
Tridiagonal(Vcat(convert.(AbstractVector{TV}, A.dl.args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.d, λ.λ).args)...),
Tridiagonal(Vcat(convert.(AbstractVector{TV}, A.dl.args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.d, λ.λ).args)...),
Vcat(convert.(AbstractVector{TV}, A.du.args)...))
end
function $op(λ::UniformScaling, A::TriPertToeplitz{V}) where V
TV = promote_type(eltype(λ),V)
Tridiagonal(Vcat(convert.(AbstractVector{TV}, broadcast($op, A.dl.args))...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, λ.λ, A.d).args)...),
Tridiagonal(Vcat(convert.(AbstractVector{TV}, broadcast($op, A.dl.args))...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, λ.λ, A.d).args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.du.args))...))
end
function $op(adjA::AdjTriPertToeplitz{T}, λ::UniformScaling) where T
function $op(adjA::AdjTriPertToeplitz{T}, λ::UniformScaling) where T
A = parent(adjA)
TV = promote_type(T,eltype(λ))
Tridiagonal(Vcat(convert.(AbstractVector{TV}, A.du.args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.d, λ.λ).args)...),
Tridiagonal(Vcat(convert.(AbstractVector{TV}, A.du.args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.d, λ.λ).args)...),
Vcat(convert.(AbstractVector{TV}, A.dl.args)...))
end
function $op(λ::UniformScaling, adjA::AdjTriPertToeplitz{V}) where V
A = parent(adjA)
TV = promote_type(eltype(λ),V)
Tridiagonal(Vcat(convert.(AbstractVector{TV}, broadcast($op, A.du.args))...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, λ.λ, A.d).args)...),
Tridiagonal(Vcat(convert.(AbstractVector{TV}, broadcast($op, A.du.args))...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, λ.λ, A.d).args)...),
Vcat(convert.(AbstractVector{TV}, broadcast($op, A.dl.args))...))
end

Expand Down Expand Up @@ -231,7 +231,7 @@ end

####
# Conversions to BandedMatrix
####
####

function BandedMatrix(A::PertToeplitz{T}, (l,u)::Tuple{Int,Int}) where T
@assert A.u == u # Not implemented
Expand Down Expand Up @@ -323,7 +323,7 @@ _constrows(A::PertConstRowMatrix) = _constrows(A.args[2])
_constrows(A::SubArray) = _constrows(parent(A))[parentindices(A)[1]]

ConstRowMatrix(A::AbstractMatrix{T}) where T = ApplyMatrix(*, A[:,1], Ones{T}(1,size(A,2)))
PertConstRowMatrix(A::AbstractMatrix{T}) where T =
PertConstRowMatrix(A::AbstractMatrix{T}) where T =
Hcat(_pertdata(A), ApplyMatrix(*, _constrows(A), Ones{T}(1,size(A,2))))

struct ConstRows <: MemoryLayout end
Expand All @@ -342,40 +342,40 @@ for Typ in (:ConstRows, :PertConstRows)
end
end

const TridiagonalToeplitzLayout = Union{SymTridiagonalLayout{FillLayout},TridiagonalLayout{FillLayout}}
struct TridiagonalToeplitzLayout <: AbstractLazyBandedLayout end
const BandedToeplitzLayout = BandedColumns{ConstRows}
const PertToeplitzLayout = BandedColumns{PertConstRows}
const PertTriangularToeplitzLayout{UPLO,UNIT} = TriangularLayout{UPLO,UNIT,BandedColumns{PertConstRows}}


_BandedMatrix(::BandedToeplitzLayout, A::AbstractMatrix) =
_BandedMatrix(::BandedToeplitzLayout, A::AbstractMatrix) =
_BandedMatrix(ConstRowMatrix(bandeddata(A)), size(A,1), bandwidths(A)...)
_BandedMatrix(::PertToeplitzLayout, A::AbstractMatrix) =
_BandedMatrix(PertConstRowMatrix(bandeddata(A)), size(A,1), bandwidths(A)...)
_BandedMatrix(::PertToeplitzLayout, A::AbstractMatrix) =
_BandedMatrix(PertConstRowMatrix(bandeddata(A)), size(A,1), bandwidths(A)...)

# for Lay in (:BandedToeplitzLayout, :PertToeplitzLayout)
# @eval begin
# @eval begin
# sublayout(::$Lay, ::Type{<:Tuple{AbstractInfUnitRange{Int},AbstractInfUnitRange{Int}}}) = $Lay()
# sublayout(::$Lay, ::Type{<:Tuple{Slice,AbstractInfUnitRange{Int}}}) = $Lay()
# sublayout(::$Lay, ::Type{<:Tuple{AbstractInfUnitRange{Int},Slice}}) = $Lay()
# sublayout(::$Lay, ::Type{<:Tuple{Slice,Slice}}) = $Lay()

# sub_materialize(::$Lay, V) = BandedMatrix(V)
# sub_materialize(::$Lay, V) = BandedMatrix(V)
# end
# end


@inline sub_materialize(::ApplyBandedLayout{typeof(*)}, V, ::Tuple{InfAxes,InfAxes}) = V
@inline sub_materialize(::BroadcastBandedLayout, V, ::Tuple{InfAxes,InfAxes}) = V
@inline sub_materialize(::AbstractBandedLayout, V, ::Tuple{InfAxes,InfAxes}) = BandedMatrix(V)
@inline sub_materialize(::AbstractBandedLayout, V, ::Tuple{InfAxes,InfAxes}) = V
@inline sub_materialize(::BandedColumns, V, ::Tuple{InfAxes,InfAxes}) = BandedMatrix(V)


##
# UniformScaling
##

# for op in (:+, :-), Typ in (:(BandedMatrix{<:Any,<:Any,OneToInf{Int}}),
# for op in (:+, :-), Typ in (:(BandedMatrix{<:Any,<:Any,OneToInf{Int}}),
# :(Adjoint{<:Any,<:BandedMatrix{<:Any,<:Any,OneToInf{Int}}}),
# :(Transpose{<:Any,<:BandedMatrix{<:Any,<:Any,OneToInf{Int}}}))
# @eval begin
Expand All @@ -395,7 +395,7 @@ _default_banded_broadcast(bc::Broadcasted, ::Tuple{<:OneToInf,<:Any}) = copy(Bro
# Banded * Banded
###

BandedMatrix{T}(::UndefInitializer, axes::Tuple{OneToInf{Int},OneTo{Int}}, lu::NTuple{2,Integer}) where T =
BandedMatrix{T}(::UndefInitializer, axes::Tuple{OneToInf{Int},OneTo{Int}}, lu::NTuple{2,Integer}) where T =
BandedMatrix{T}(undef, map(length,axes), lu)

similar(M::MulAdd{<:AbstractBandedLayout,<:AbstractBandedLayout}, ::Type{T}, axes::Tuple{OneTo{Int},OneToInf{Int}}) where T =
Expand Down Expand Up @@ -471,3 +471,34 @@ function _bidiag_forwardsub!(M::Ldiv{<:Any,<:PaddedLayout})
b_in
end

###
# Inf-Toeplitz layout
# this could possibly be avoided via an InfFillLayout
###

const InfFill = AbstractFill{<:Any,1,<:Tuple{OneToInf}}

for Typ in (:(LinearAlgebra.Tridiagonal{<:Any,<:InfFill}),
:(LinearAlgebra.SymTridiagonal{<:Any,<:InfFill}),
:(LazyBandedMatrices.Tridiagonal{<:Any,<:InfFill,<:InfFill,<:InfFill}),
:(LazyBandedMatrices.SymTridiagonal{<:Any,<:InfFill,<:InfFill}))
@eval begin
MemoryLayout(::Type{<:$Typ}) = TridiagonalToeplitzLayout()
Base.BroadcastStyle(::Type{<:$Typ}) = LazyArrayStyle{2}()
end
end

struct BidiagonalToeplitzLayout <: AbstractLazyBandedLayout end

for Typ in (:(LinearAlgebra.Bidiagonal{<:Any,<:InfFill}),
:(LazyBandedMatrices.Bidiagonal{<:Any,<:InfFill,<:InfFill}))
@eval begin
MemoryLayout(::Type{<:$Typ}) = BidiagonalToeplitzLayout()
Base.BroadcastStyle(::Type{<:$Typ}) = LazyArrayStyle{2}()
end
end

# fall back for Ldiv
triangularlayout(::Type{<:TriangularLayout{UPLO,'N'}}, ::TridiagonalToeplitzLayout) where UPLO = BidiagonalToeplitzLayout()
materialize!(L::MatLdivVec{BidiagonalToeplitzLayout,Lay}) where Lay = materialize!(Ldiv{BidiagonalLayout{FillLayout,FillLayout},Lay}(L.A, L.B))
copyto!(dest::AbstractArray, L::Ldiv{BidiagonalToeplitzLayout,Lay}) where Lay = copyto!(dest, Ldiv{BidiagonalLayout{FillLayout,FillLayout},Lay}(L.A, L.B))
4 changes: 3 additions & 1 deletion src/blockbanded/infblocktridiagonal.jl
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
const BlockTriPertToeplitz{T} = BlockMatrix{T,Tridiagonal{Matrix{T},Vcat{Matrix{T},1,Tuple{Vector{Matrix{T}},Fill{Matrix{T},1,Tuple{OneToInf{Int}}}}}},
NTuple{2,BlockedUnitRange{Vcat{Int,1,Tuple{Vector{Int},InfStepRange{Int,Int}}}}}}

const BlockTridiagonalToeplitzLayout = BlockLayout{TridiagonalLayout{FillLayout},DenseColumnMajor}
const BlockTridiagonalToeplitzLayout = BlockLayout{TridiagonalToeplitzLayout,DenseColumnMajor}

function BlockTridiagonal(adjA::Adjoint{T,BlockTriPertToeplitz{T}}) where T
A = parent(adjA)
Expand Down Expand Up @@ -33,7 +33,9 @@ end

sizes_from_blocks(A::Diagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.diag, 1), size.(A.diag,2)
sizes_from_blocks(A::Tridiagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.d, 1), size.(A.d,2)
sizes_from_blocks(A::LazyBandedMatrices.Tridiagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.d, 1), size.(A.d,2)
sizes_from_blocks(A::Bidiagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.dv, 1), size.(A.dv,2)
sizes_from_blocks(A::LazyBandedMatrices.Bidiagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.dv, 1), size.(A.dv,2)

axes_print_matrix_row(_, io, X, A, i, ::AbstractVector{<:PosInfinity}, sep) = nothing
axes_print_matrix_row(::NTuple{2,BlockedUnitRange}, io, X, A, i, ::AbstractVector{<:PosInfinity}, sep) = nothing
Expand Down
108 changes: 80 additions & 28 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,42 +2,45 @@ using InfiniteLinearAlgebra, BlockBandedMatrices, BlockArrays, BandedMatrices, I
MatrixFactorizations, ArrayLayouts, LinearAlgebra, Random, LazyBandedMatrices, StaticArrays
import InfiniteLinearAlgebra: qltail, toeptail, tailiterate , tailiterate!, tail_de, ql_X!,
InfToeplitz, PertToeplitz, TriToeplitz, InfBandedMatrix, InfBandCartesianIndices,
rightasymptotics, QLHessenberg, ConstRows, PertConstRows, BandedToeplitzLayout, PertToeplitzLayout
rightasymptotics, QLHessenberg, ConstRows, PertConstRows,
BandedToeplitzLayout, PertToeplitzLayout, TridiagonalToeplitzLayout, BidiagonalToeplitzLayout
import Base: BroadcastStyle
import BlockArrays: _BlockArray
import BlockBandedMatrices: isblockbanded, _BlockBandedMatrix
import MatrixFactorizations: QLPackedQ
import BandedMatrices: bandeddata, _BandedMatrix, BandedStyle
import LazyArrays: colsupport, ApplyStyle, MemoryLayout, ApplyLayout, LazyArrayStyle, arguments
import InfiniteArrays: OneToInf, oneto
import LazyBandedMatrices: BroadcastBandedBlockBandedLayout, BroadcastBandedLayout
import LazyBandedMatrices: BroadcastBandedBlockBandedLayout, BroadcastBandedLayout, LazyBandedLayout


@testset "∞-banded" begin
D = Diagonal(Fill(2,∞))
@testset "Diagonal and BandedMatrix" begin
D = Diagonal(Fill(2,∞))

B = D[1:∞,2:∞]
@test B isa BandedMatrix
@test B[1:10,1:10] == diagm(-1 => Fill(2,9))
@test B[1:∞,2:∞] isa BandedMatrix
B = D[1:∞,2:∞]
@test B isa BandedMatrix
@test B[1:10,1:10] == diagm(-1 => Fill(2,9))
@test B[1:∞,2:∞] isa BandedMatrix

A = BandedMatrix(0 => 1:∞, 1=> Fill(2.0,∞), -1 => Fill(3.0,∞))
x = [1; 2; zeros(∞)]
@test A*x isa Vcat
@test (A*x)[1:10] == A[1:10,1:10]*x[1:10]
A = BandedMatrix(0 => 1:∞, 1=> Fill(2.0,∞), -1 => Fill(3.0,∞))
x = [1; 2; zeros(∞)]
@test A*x isa Vcat
@test (A*x)[1:10] == A[1:10,1:10]*x[1:10]

@test InfBandCartesianIndices(0)[1:5] == CartesianIndex.(1:5,1:5)
@test InfBandCartesianIndices(1)[1:5] == CartesianIndex.(1:5,2:6)
@test InfBandCartesianIndices(-1)[1:5] == CartesianIndex.(2:6,1:5)
@test InfBandCartesianIndices(0)[1:5] == CartesianIndex.(1:5,1:5)
@test InfBandCartesianIndices(1)[1:5] == CartesianIndex.(1:5,2:6)
@test InfBandCartesianIndices(-1)[1:5] == CartesianIndex.(2:6,1:5)

@test D[band(0)] ≡ Fill(2,∞)
@test D[band(1)] ≡ Fill(0,∞)
@test A[band(0)][2:10] == 2:10
@test A[band(0)][2:10] == 2:10
@test D[band(0)] ≡ Fill(2,∞)
@test D[band(1)] ≡ Fill(0,∞)

@test B*A*x isa Vcat
@test (B*A*x)[1:10] == [0; 10; 14; 12; zeros(6)]
@test B*A*x isa Vcat
@test (B*A*x)[1:10] == [0; 10; 14; 12; zeros(6)]

@test _BandedMatrix((1:∞)', ∞, -1, 1) isa BandedMatrix
@test _BandedMatrix((1:∞)', ∞, -1, 1) isa BandedMatrix
end

@testset "∞-Toeplitz" begin
A = BandedMatrix(1 => Fill(2im,∞), 2 => Fill(-1,∞), 3 => Fill(2,∞), -2 => Fill(-4,∞), -3 => Fill(-2im,∞))
Expand All @@ -51,19 +54,56 @@ import LazyBandedMatrices: BroadcastBandedBlockBandedLayout, BroadcastBandedLayo
@test A[:,3:end] isa InfToeplitz

@test (A + 2I)[1:10,1:10] == (2I + A)[1:10,1:10] == A[1:10,1:10] + 2I

@test MemoryLayout(Tridiagonal(Fill(1,∞), Fill(2,∞), Fill(3,∞))) isa TridiagonalToeplitzLayout
@test MemoryLayout(Bidiagonal(Fill(1,∞), Fill(2,∞), :U)) isa BidiagonalToeplitzLayout
@test MemoryLayout(SymTridiagonal(Fill(1,∞), Fill(2,∞))) isa TridiagonalToeplitzLayout
@test MemoryLayout(LazyBandedMatrices.Tridiagonal(Fill(1,∞), Zeros(∞), Fill(3,∞))) isa TridiagonalToeplitzLayout
@test MemoryLayout(LazyBandedMatrices.Bidiagonal(Fill(1,∞), Zeros(∞), :U)) isa BidiagonalToeplitzLayout
@test MemoryLayout(LazyBandedMatrices.SymTridiagonal(Fill(1,∞), Zeros(∞))) isa TridiagonalToeplitzLayout

T = LazyBandedMatrices.Tridiagonal(Fill(1,∞), Zeros(∞), Fill(3,∞))
@test T[2:∞,3:∞] isa SubArray
@test exp.(T) isa BroadcastMatrix
@test exp.(T)[2:∞,3:∞] isa SubArray

B = LazyBandedMatrices.Bidiagonal(Fill(1,∞), Zeros(∞), :U)
@test B[2:∞,3:∞] isa SubArray
@test exp.(B) isa BroadcastMatrix
@test exp.(B)[2:∞,3:∞] isa SubArray

@testset "algebra" begin
T = Tridiagonal(Fill(1,∞), Fill(2,∞), Fill(3,∞))
@test T isa InfiniteLinearAlgebra.TriToeplitz
@test (T + 2I)[1:10,1:10] == (2I + T)[1:10,1:10] == T[1:10,1:10] + 2I
end
end

@testset "Pert-Toeplitz" begin
A = BandedMatrix(-2 => Vcat(Float64[], Fill(1/4,∞)), 0 => Vcat([1.0+im,2,3],Fill(0,∞)), 1 => Vcat(Float64[], Fill(1,∞)))
@test A isa PertToeplitz
@test MemoryLayout(typeof(A)) == PertToeplitzLayout()
V = view(A,2:∞,2:∞)
@test MemoryLayout(typeof(V)) == PertToeplitzLayout()
@test BandedMatrix(V) isa PertToeplitz
@test A[2:∞,2:∞] isa PertToeplitz
@testset "Inf Pert" begin
A = BandedMatrix(-2 => Vcat(Float64[], Fill(1/4,∞)), 0 => Vcat([1.0+im,2,3],Fill(0,∞)), 1 => Vcat(Float64[], Fill(1,∞)))
@test A isa PertToeplitz
@test MemoryLayout(typeof(A)) == PertToeplitzLayout()
V = view(A,2:∞,2:∞)
@test MemoryLayout(typeof(V)) == PertToeplitzLayout()
@test BandedMatrix(V) isa PertToeplitz
@test A[2:∞,2:∞] isa PertToeplitz

@test (A + 2I)[1:10,1:10] == (2I + A)[1:10,1:10] == A[1:10,1:10] + 2I
end

@test (A + 2I)[1:10,1:10] == (2I + A)[1:10,1:10] == A[1:10,1:10] + 2I
@testset "TriPert" begin
A = SymTridiagonal(Vcat([1,2.], Fill(2.,∞)), Vcat([3.,4.], Fill.(0.5,∞)))
@test A isa InfiniteLinearAlgebra.SymTriPertToeplitz
@test (A + 2I)[1:10,1:10] == (2I + A)[1:10,1:10] == A[1:10,1:10] + 2I

A = Tridiagonal(Vcat([3.,4.], Fill.(0.5,∞)), Vcat([1,2.], Fill(2.,∞)), Vcat([3.,4.], Fill.(0.5,∞)))
@test A isa InfiniteLinearAlgebra.TriPertToeplitz
@test Adjoint(A) isa InfiniteLinearAlgebra.AdjTriPertToeplitz
@test (A + 2I)[1:10,1:10] == (2I + A)[1:10,1:10] == A[1:10,1:10] + 2I
@test (Adjoint(A) + 2I)[1:10,1:10] == (2I + Adjoint(A))[1:10,1:10] == Adjoint(A)[1:10,1:10] + 2I
end


@testset "InfBanded" begin
A = _BandedMatrix(Fill(2,4,∞),ℵ₀,2,1)
Expand Down Expand Up @@ -296,6 +336,18 @@ end
@test (A + im*I)[1:100,1:100] == A[1:100,1:100]+im*I
@test (im*I + A)[1:100,1:100] == im*I+A[1:100,1:100]
@test (im*I - A)[1:100,1:100] == im*I-A[1:100,1:100]

T = mortar(LazyBandedMatrices.Tridiagonal(Fill([1 2; 3 4], ∞), Fill([1 2; 3 4], ∞), Fill([1 2; 3 4], ∞)))
#TODO: copy BlockBidiagonal code from BlockBandedMatrices to LazyBandedMatrices
@test T[Block(2,2)] == [1 2; 3 4]
@test_broken T[Block(1,3)] == Zeros(2,2)
end

@testset "BlockBidiagonal" begin
B = mortar(LazyBandedMatrices.Bidiagonal(Fill([1 2; 3 4], ∞), Fill([1 2; 3 4], ∞), :U))
#TODO: copy BlockBidiagonal code from BlockBandedMatrices to LazyBandedMatrices
@test B[Block(2,3)] == [1 2; 3 4]
@test_broken B[Block(1,3)] == Zeros(2,2)
end

@testset "Fill" begin
Expand Down
2 changes: 1 addition & 1 deletion test/test_inful.jl
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
using InfiniteLinearAlgebra, BlockArrays, ArrayLayouts, Test
import InfiniteLinearAlgebra: BlockTridiagonalToeplitzLayout, ul, adaptiveqr

@testset "Block Toeplitz UL" begin
@testset "∞-UL" begin
@testset "Toeplitz" begin
Δ = SymTridiagonal(Fill(-2,∞), Fill(1,∞))
h = 0.01
Expand Down