Skip to content

Removed MulBandedLayout, removed print #33

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 4 commits into from
Jun 6, 2020
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
18 changes: 9 additions & 9 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.3.3"
version = "0.3.4"

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

[compat]
ArrayLayouts = "0.3"
BandedMatrices = "0.15.7"
BlockArrays = "0.12.6"
BlockBandedMatrices = "0.8.4"
FillArrays = "0.8.8"
InfiniteArrays = "0.7"
LazyArrays = "0.16.8"
LazyBandedMatrices = "0.2.8"
ArrayLayouts = "0.3.4"
BandedMatrices = "0.15.11"
BlockArrays = "0.12.7"
BlockBandedMatrices = "0.8.5"
FillArrays = "0.8.10"
InfiniteArrays = "0.7.2"
LazyArrays = "0.16.12"
LazyBandedMatrices = "0.2.11"
MatrixFactorizations = "0.4.1"
SemiseparableMatrices = "0.1"
julia = "1.3"
Expand Down
30 changes: 18 additions & 12 deletions src/InfiniteLinearAlgebra.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,17 +2,17 @@ module InfiniteLinearAlgebra
using BlockArrays, BlockBandedMatrices, BandedMatrices, LazyArrays, LazyBandedMatrices, SemiseparableMatrices,
FillArrays, InfiniteArrays, MatrixFactorizations, ArrayLayouts, LinearAlgebra

import Base: +, -, *, /, \, ^, OneTo, getindex, promote_op, _unsafe_getindex, print_matrix_row, size, axes,
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
import Base.Broadcast: BroadcastStyle, Broadcasted, broadcasted

import ArrayLayouts: colsupport, rowsupport, triangularlayout, MatLdivVec, triangulardata, TriangularLayout, sublayout, _qr,
import ArrayLayouts: colsupport, rowsupport, triangularlayout, MatLdivVec, triangulardata, TriangularLayout, sublayout, _qr,
MatLmulVec, MatLmulMat, AbstractQLayout, materialize!
import BandedMatrices: BandedMatrix, _BandedMatrix, AbstractBandedMatrix, bandeddata, bandwidths, BandedColumns, bandedcolumns,
_default_banded_broadcast
import FillArrays: AbstractFill, getindex_value
import InfiniteArrays: OneToInf, InfUnitRange, Infinity, InfStepRange, AbstractInfUnitRange, InfAxes
import FillArrays: AbstractFill, getindex_value, axes_print_matrix_row
import InfiniteArrays: OneToInf, InfUnitRange, Infinity, InfStepRange, AbstractInfUnitRange, InfAxes
import LinearAlgebra: lmul!, ldiv!, matprod, qr, AbstractTriangular, AbstractQ, adjoint, transpose
import LazyArrays: applybroadcaststyle, CachedArray, CachedMatrix, CachedVector, DenseColumnMajor, FillLayout, ApplyMatrix, check_mul_axes, ApplyStyle, LazyArrayApplyStyle, LazyArrayStyle,
resizedata!, MemoryLayout, mulapplystyle, LmulStyle, RmulStyle,
Expand All @@ -26,7 +26,7 @@ import BlockArrays: AbstractBlockVecOrMat, sizes_from_blocks, _length, BlockedUn

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

import LazyBandedMatrices: MulBandedLayout, BroadcastBandedLayout, _krontrav_axes
import LazyBandedMatrices: ApplyBandedLayout, BroadcastBandedLayout, _krontrav_axes

import BlockBandedMatrices: _BlockSkylineMatrix, _BandedMatrix, _BlockSkylineMatrix, blockstart, blockstride,
BlockSkylineSizes, BlockSkylineMatrix, BlockBandedMatrix, _BlockBandedMatrix, BlockTridiagonal,
Expand All @@ -41,7 +41,7 @@ LazyArrays.@lazymul BandedMatrix{<:Any,<:Any,<:OneToInf}
# BroadcastStyle(::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}) = LazyArrayStyle{2}()

function ^(A::BandedMatrix{T,<:Any,<:OneToInf}, p::Integer) where T
if p < 0
if p < 0
inv(A)^(-p)
elseif p == 0
Eye{T}(∞)
Expand All @@ -66,18 +66,24 @@ include("infqr.jl")
# block broadcasted
######

const CumsumOneToInf2 = BroadcastArray{Int64,1,typeof(div),Tuple{BroadcastArray{Int64,1,typeof(*),Tuple{InfiniteArrays.OneToInf{Int64},InfiniteArrays.InfUnitRange{Int64}}},Int64}}
const CumsumOneToInf2 = Cumsum{<:Any,1,<:OneToInf}
BlockArrays.sortedunion(a::CumsumOneToInf2, ::CumsumOneToInf2) = a

function BlockArrays.sortedunion(a::Vcat{Int,1,<:Tuple{<:AbstractVector{Int},InfStepRange{Int,Int}}},
b::Vcat{Int,1,<:Tuple{<:AbstractVector{Int},InfStepRange{Int,Int}}})
@assert a == b
a
end


map(::typeof(length), A::BroadcastArray{OneTo{Int},1,Type{OneTo}}) = A.args[1]
map(::typeof(length), A::BroadcastArray{<:Fill,1,Type{Fill}}) = A.args[2]
map(::typeof(length), A::BroadcastArray{<:Zeros,1,Type{Zeros}}) = A.args[1]
map(::typeof(length), A::BroadcastArray{<:Fill,1,Type{Fill}}) = A.args[2]
map(::typeof(length), A::BroadcastArray{<:Zeros,1,Type{Zeros}}) = A.args[1]
map(::typeof(length), A::BroadcastArray{<:Vcat,1,Type{Vcat}}) = broadcast(+,map.(length,A.args)...)
broadcasted(::LazyArrayStyle{1}, ::typeof(length), A::BroadcastArray{OneTo{Int},1,Type{OneTo}}) =
A.args[1]
broadcasted(::LazyArrayStyle{1}, ::typeof(length), A::BroadcastArray{<:Fill,1,Type{Fill}}) =
A.args[2]
A.args[2]

BlockArrays._length(::BlockedUnitRange, ::OneToInf) = ∞
BlockArrays._last(::BlockedUnitRange, ::OneToInf) = ∞
Expand All @@ -86,8 +92,8 @@ BlockArrays._last(::BlockedUnitRange, ::OneToInf) = ∞
# KronTrav
###

_krontrav_axes(A::NTuple{N,OneToInf{Int}}, B::NTuple{N,OneToInf{Int}}) where N =
_krontrav_axes(A::NTuple{N,OneToInf{Int}}, B::NTuple{N,OneToInf{Int}}) where N =
@. blockedrange(OneTo(length(A)))


end # module
2 changes: 1 addition & 1 deletion src/banded/infbanded.jl
Original file line number Diff line number Diff line change
Expand Up @@ -328,7 +328,7 @@ _BandedMatrix(::PertToeplitzLayout, A::AbstractMatrix) =
# end


@inline sub_materialize(::MulBandedLayout, V, ::Tuple{InfAxes,InfAxes}) = V
@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(::BandedColumns, V, ::Tuple{InfAxes,InfAxes}) = BandedMatrix(V)
Expand Down
12 changes: 3 additions & 9 deletions src/blockbanded/infblocktridiagonal.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,16 +34,10 @@ sizes_from_blocks(A::Diagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.diag, 1), si

sizes_from_blocks(A::Tridiagonal, ::NTuple{2,OneToInf{Int}}) = size.(A.d, 1), size.(A.d,2)

print_matrix_row(io::IO,
X::AbstractBlockVecOrMat, A::Vector,
i::Integer, cols::AbstractVector{<:Infinity}, sep::AbstractString) = nothing
axes_print_matrix_row(_, io, X, A, i, ::AbstractVector{<:Infinity}, sep) = nothing
axes_print_matrix_row(::NTuple{2,BlockedUnitRange}, io, X, A, i, ::AbstractVector{<:Infinity}, sep) = nothing


print_matrix_row(io::IO,
X::Union{AbstractTriangular{<:Any,<:AbstractBlockMatrix},
Symmetric{<:Any,<:AbstractBlockMatrix},
Hermitian{<:Any,<:AbstractBlockMatrix}}, A::Vector,
i::Integer, cols::AbstractVector{<:Infinity}, sep::AbstractString) = nothing

function BlockSkylineSizes(A::BlockTriPertToeplitz, (l,u)::NTuple{2,Int})
N = max(length(A.blocks.du.args[1])+1,length(A.blocks.d.args[1]),length(A.blocks.dl.args[1]))
block_sizes = Vector{Int}(undef, N) # assume square
Expand Down
9 changes: 2 additions & 7 deletions src/infqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -87,7 +87,7 @@ end
partialqr!(F::AdaptiveQRData{<:Any,<:BlockSkylineMatrix}, n::Int) =
partialqr!(F, findblock(axes(F.data,2), n))

struct AdaptiveQRFactors{T,DM<:AbstractMatrix{T},M<:AbstractMatrix{T}} <: AbstractMatrix{T}
struct AdaptiveQRFactors{T,DM<:AbstractMatrix{T},M<:AbstractMatrix{T}} <: LayoutMatrix{T}
data::AdaptiveQRData{T,DM,M}
end

Expand Down Expand Up @@ -128,12 +128,7 @@ rowsupport(F::QRPackedQ{<:Any,<:AdaptiveQRFactors}, j) = rowsupport(F.factors, j
blockcolsupport(F::QRPackedQ{<:Any,<:AdaptiveQRFactors}, j) = blockcolsupport(F.factors, j)


Base.replace_in_print_matrix(A::AdaptiveQRFactors, i::Integer, j::Integer, s::AbstractString) =
i in colsupport(A,j) ? s : Base.replace_with_centered_mark(s)
Base.replace_in_print_matrix(A::UpperTriangular{<:Any,<:AdaptiveQRFactors}, i::Integer, j::Integer, s::AbstractString) =
i in colsupport(A,j) ? s : Base.replace_with_centered_mark(s)

struct AdaptiveQRTau{T,DM<:AbstractMatrix{T},M<:AbstractMatrix{T}} <: AbstractVector{T}
struct AdaptiveQRTau{T,DM<:AbstractMatrix{T},M<:AbstractMatrix{T}} <: LayoutVector{T}
data::AdaptiveQRData{T,DM,M}
end

Expand Down
3 changes: 2 additions & 1 deletion test/test_infqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -221,6 +221,7 @@ import SemiseparableMatrices: AlmostBandedLayout, VcatAlmostBandedLayout
@test (F.Q*[1;zeros(∞)])[1:6] ≈ [-0.9701425001453321,0,0.24253562503633297,0,0,0]

u = F \ [1; zeros(∞)]
@test blockisequal(axes(A,2),axes(u,1))
@test (A*u)[1:10] ≈ [1; zeros(9)]

x = 0.1
Expand All @@ -232,14 +233,14 @@ import SemiseparableMatrices: AlmostBandedLayout, VcatAlmostBandedLayout

@test dot(u[getindex.(Block.(1:50),1)], sin.((1:50) .* θ)/sin(θ)) ≈ 1/(x-2)


L = A+B;
@test MemoryLayout(L) isa BroadcastBandedBlockBandedLayout{typeof(+)}
V = view(L,Block.(1:400),Block.(1:400))
@time u = L \ [1;zeros(∞)]

x,y = 0.1,0.2
θ,φ = acos(x),acos(y)
@test u[Block.(1:50)] isa PseudoBlockArray
@test (sin.((1:50) .* φ)/sin(φ))' * InvDiagTrav(u[Block.(1:50)]) * sin.((1:50) .* θ)/sin(θ) ≈ 1/(x+y-4)
@test (L*u)[1:10] ≈ [1; zeros(9)]

Expand Down