Skip to content

Remove * overload hacks #26

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 1 commit into from
Apr 16, 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
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.3.0"
version = "0.3.1"

[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.2.1"
BandedMatrices = "0.15.1"
ArrayLayouts = "0.2.4"
BandedMatrices = "0.15.3"
BlockArrays = "0.12"
BlockBandedMatrices = "0.8"
FillArrays = "0.8.6"
InfiniteArrays = "0.7"
LazyArrays = "0.16.4"
LazyBandedMatrices = "0.2.5"
LazyArrays = "0.16.5"
LazyBandedMatrices = "0.2.6"
MatrixFactorizations = "0.3.1"
SemiseparableMatrices = "0.0.1"
julia = "1.3"
Expand Down
32 changes: 1 addition & 31 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, print_matrix_row, size, axes,
AbstractMatrix, AbstractArray, Matrix, Array, Vector, AbstractVector, Slice,
show, getproperty, copy, map
show, getproperty, copy, map, require_one_based_indexing
import Base.Broadcast: BroadcastStyle, Broadcasted, broadcasted

import ArrayLayouts: colsupport, rowsupport, triangularlayout, MatLdivVec, triangulardata, TriangularLayout, sublayout, _qr
Expand Down Expand Up @@ -49,13 +49,6 @@ function ^(A::BandedMatrix{T,<:Any,<:OneToInf}, p::Integer) where T
end
end

if VERSION < v"1.2-"
import Base: has_offset_axes
require_one_based_indexing(A...) = !has_offset_axes(A...) || throw(ArgumentError("offset arrays are not supported but got an array with index other than 1"))
else
import Base: require_one_based_indexing
end

export Vcat, Fill, ql, ql!, ∞, ContinuousSpectrumError, BlockTridiagonal

include("banded/hessenbergq.jl")
Expand All @@ -66,29 +59,6 @@ include("banded/infqltoeplitz.jl")
include("infql.jl")
include("infqr.jl")

###
# temporary work arounds
###

# Fix ∞ BandedMatrix
ApplyStyle(::typeof(*), ::Type{<:Diagonal}, ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}) =
LmulStyle()
ApplyStyle(::typeof(*), ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}, ::Type{<:Diagonal}) =
RmulStyle()
ApplyStyle(::typeof(*), ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}, ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}) =
LazyArrayApplyStyle()
ApplyStyle(::typeof(*), ::Type{<:AbstractArray}, ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}) =
LazyArrayApplyStyle()
ApplyStyle(::typeof(*), ::Type{<:BandedMatrix{<:Any,<:Any,<:OneToInf}}, ::Type{<:AbstractArray}) =
LazyArrayApplyStyle()
ApplyStyle(::typeof(*), ::Type{<:Adjoint{<:Any,<:BandedMatrix{<:Any,<:Any,<:OneToInf}}}, ::Type{<:AbstractArray}) =
LazyArrayApplyStyle()
ApplyStyle(::typeof(*), ::Type{<:Transpose{<:Any,<:BandedMatrix{<:Any,<:Any,<:OneToInf}}}, ::Type{<:AbstractArray}) =
LazyArrayApplyStyle()

ApplyStyle(::typeof(*), ::Type{<:BandedMatrix{<:Any,<:AbstractFill,<:OneToInf}}, ::Type{<:BandedMatrix{<:Any,<:AbstractFill,<:OneToInf}}) =
MulAddStyle()

#######
# block broadcasted
######
Expand Down
11 changes: 9 additions & 2 deletions src/banded/infbanded.jl
Original file line number Diff line number Diff line change
Expand Up @@ -354,7 +354,7 @@ _default_banded_broadcast(bc::Broadcasted, ::Tuple{<:OneToInf,<:Any}) = copy(Bro
# BandedFill * BandedFill
###

copy(M::MulAdd{BandedMatrices.BandedColumns{FillLayout},BandedMatrices.BandedColumns{FillLayout},ZerosLayout}) =
copy(M::MulAdd{BandedColumns{FillLayout},BandedColumns{FillLayout},ZerosLayout}) =
_bandedfill_mul(M, axes(M.A), axes(M.B))

_bandedfill_mul(M, _, _) = copyto!(similar(M), M)
Expand All @@ -368,4 +368,11 @@ function _bandedfill_mul(M::MulAdd, ::Tuple{InfAxes,InfAxes}, ::Tuple{InfAxes,In
ret = _BandedMatrix(Hcat(Array{typeof(λ)}(undef, l+u+1,u), [1:m-1; Fill(m,l+u-2m+3); m-1:-1:1]*Fill(λ,1,∞)), ∞, l, u)
mul!(view(ret, 1:l+u,1:u), view(A,1:l+u,1:u+Bl), view(B,1:u+Bl,1:u))
ret
end
end

mulapplystyle(::BandedColumns{FillLayout}, ::PertToeplitzLayout) = LazyArrayApplyStyle()
mulapplystyle(::PertToeplitzLayout, ::BandedColumns{FillLayout}) = LazyArrayApplyStyle()
mulapplystyle(::BandedColumns{FillLayout}, ::BandedToeplitzLayout) = LazyArrayApplyStyle()
mulapplystyle(::BandedToeplitzLayout, ::BandedColumns{FillLayout}) = LazyArrayApplyStyle()
mulapplystyle(::QLayout, ::BandedToeplitzLayout) = LazyArrayApplyStyle()
mulapplystyle(::QLayout, ::PertToeplitzLayout) = LazyArrayApplyStyle()
11 changes: 9 additions & 2 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,7 @@ import InfiniteLinearAlgebra: qltail, toeptail, tailiterate , tailiterate!, tail
InfToeplitz, PertToeplitz, TriToeplitz, InfBandedMatrix,
rightasymptotics, QLHessenberg, ConstRows, PertConstRows, BandedToeplitzLayout, PertToeplitzLayout
import Base: BroadcastStyle
import BlockArrays: _BlockArray
import BlockArrays: _BlockArray
import BlockBandedMatrices: isblockbanded, _BlockBandedMatrix
import MatrixFactorizations: QLPackedQ
import BandedMatrices: bandeddata, _BandedMatrix, BandedStyle
Expand All @@ -16,6 +16,12 @@ import InfiniteArrays: OneToInf
@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 ApplyStyle(*, typeof(A), typeof(x)) isa LazyArrays.MulAddStyle
@test A*x isa Vcat
@test (A*x)[1:10] == A[1:10,1:10]*x[1:10]
end

@testset "∞-block arrays" begin
Expand Down Expand Up @@ -91,7 +97,7 @@ end
A = _BandedMatrix(Fill(1,4,∞),∞,1,2)
@test A*A isa BandedMatrix
@test (A^2)[1:10,1:10] == (A*A)[1:10,1:10] == (A[1:100,1:100]^2)[1:10,1:10]
@test (A^3)[1:10,1:10] == (A*A*A)[1:10,1:10] == (A[1:100,1:100]^3)[1:10,1:10]
@test (A^3)[1:10,1:10] == (A*A*A)[1:10,1:10] == ((A*A)*A)[1:10,1:10] == (A*(A*A))[1:10,1:10] == (A[1:100,1:100]^3)[1:10,1:10]
end

@testset "BlockTridiagonal" begin
Expand Down Expand Up @@ -152,6 +158,7 @@ end
C = Diagonal( 2 ./ (1:2:∞))
@test A*(B*C) isa MulMatrix
@test bandwidths(A*(B*C)) == (-1,1)
@test bandwidths((A*B)*C) == (-1,1)
end

@testset "Triangle OP recurrences" begin
Expand Down
1 change: 1 addition & 0 deletions test/test_infql.jl
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@ import BandedMatrices: _BandedMatrix
L = _BandedMatrix(Hcat([e; X[2,2]; X[2,1]], X[2,end:-1:1] * Ones{Float64}(1,∞)), ∞, 2, 0)
Qn,Ln = ql(A[1:1000,1:1000])
@test Q[1:10,1:10] ≈ Qn[1:10,1:10]
@test Q'A isa MulMatrix
@test Array((Q'A)[1:10,1:10]) ≈ [(Q'A)[k,j] for k=1:10,j=1:10]
@test (Q'A)[1:10,1:10] ≈ Ln[1:10,1:10] ≈ L[1:10,1:10]

Expand Down