Skip to content

Commit fd5a2de

Browse files
authored
Faster QR (#75)
* Faster QR * Update Project.toml * Update Project.toml
1 parent fa51908 commit fd5a2de

File tree

5 files changed

+62
-9
lines changed

5 files changed

+62
-9
lines changed

src/AbstractBlockBandedMatrix.jl

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -79,9 +79,10 @@ blockbanded_blockrowstart(A, i::BlockRange) = blockbanded_blockrowstart(A, minim
7979
blockbanded_blockcolstop(A, i::BlockRange) = blockbanded_blockcolstop(A, maximum(i))
8080
blockbanded_blockrowstop(A, i::BlockRange) = blockbanded_blockrowstop(A, maximum(i))
8181

82+
const AllBlockBandedLayout{UPLO,UNIT} = Union{AbstractBlockBandedLayout,TriangularLayout{UPLO,UNIT,<:AbstractBlockBandedLayout}}
8283

83-
@inline blockcolsupport(::AbstractBlockBandedLayout, A, i) = blockbanded_blockcolstart(A,i):blockbanded_blockcolstop(A,i)
84-
@inline blockrowsupport(::AbstractBlockBandedLayout, A, i) = blockbanded_blockrowstart(A,i):blockbanded_blockrowstop(A,i)
84+
@inline blockcolsupport(::AllBlockBandedLayout, A, i) = blockbanded_blockcolstart(A,i):blockbanded_blockcolstop(A,i)
85+
@inline blockrowsupport(::AllBlockBandedLayout, A, i) = blockbanded_blockrowstart(A,i):blockbanded_blockrowstop(A,i)
8586

8687

8788
# length of i-the column/row

src/BlockBandedMatrices.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -32,7 +32,7 @@ import BlockArrays: blocksize, blockcheckbounds, BlockedUnitRange, blockisequal,
3232
Block, BlockSlice, getblock, unblock, setblock!, block, blockindex,
3333
_blocklengths2blocklasts, BlockIndexRange, sizes_from_blocks, BlockSlice1,
3434
blockcolsupport, blockrowsupport, blockcolstart, blockcolstop, blockrowstart, blockrowstop,
35-
AbstractBlockLayout, BlockLayout, blocks
35+
AbstractBlockLayout, BlockLayout, blocks, hasmatchingblocks
3636

3737
import BandedMatrices: isbanded, bandwidths, bandwidth, banded_getindex, colrange,
3838
inbands_setindex!, inbands_getindex, banded_setindex!,

src/blockskylineqr.jl

Lines changed: 26 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -109,15 +109,38 @@ function materialize!(Mul::MatLmulVec{<:AdjQLPackedQLayout{<:AbstractBlockBanded
109109
Mul.B
110110
end
111111

112+
function materialize!(Mul::MatLmulMat{<:AdjQRPackedQLayout{<:AbstractBlockBandedLayout}})
113+
adjQ,Bin = Mul.A,Mul.B
114+
Q = parent(adjQ)
115+
A = Q.factors
116+
l,u = blockbandwidths(A)
117+
N,M = blocksize(A)
118+
# impose block structure
119+
ax1 = (axes(A,1),)
120+
τ = PseudoBlockArray(Q.τ, ax1)
121+
B = PseudoBlockArray(Bin, (axes(A,1),axes(Bin,2)))
122+
for K = 1:min(N,M), J = 1:blocksize(Bin,2)
123+
KR = Block.(K:min(K+l,N))
124+
V = view(A,KR,Block(K))
125+
t = view(τ,Block(K))
126+
apply_qr!(V, t, view(B,KR,Block(J)))
127+
end
128+
Bin
129+
end
130+
112131
# avoid LinearALgebra Strided obsession
113132

114133
for Typ in (:StridedVector, :StridedMatrix, :AbstractVector, :AbstractMatrix, :LayoutMatrix)
115134
@eval function ldiv!(A::QR{<:Any,<:BlockSkylineMatrix}, B::$Typ)
116135
lmul!(adjoint(A.Q), B)
117136
M,N = blocksize(A.factors)
118-
MN = min(M,N)
119-
V = view(A.factors,Block.(1:MN), Block.(1:MN))
120-
materialize!(Ldiv(UpperTriangular(V), view(B,1:size(V,1),:)))
137+
if M == N
138+
materialize!(Ldiv(UpperTriangular(A.factors), view(B,1:size(A.factors,1),:)))
139+
else
140+
MN = min(M,N)
141+
V = view(A.factors,Block.(1:MN), Block.(1:MN))
142+
materialize!(Ldiv(UpperTriangular(V), view(B,1:size(V,1),:)))
143+
end
121144
B
122145
end
123146
end

src/triblockbanded.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -10,7 +10,7 @@ isbandedblockbanded(A::AbstractTriangular) =
1010
function blockbandwidths(A::Union{UpperTriangular,UnitUpperTriangular})
1111
P = parent(A)
1212
if hasmatchingblocks(P)
13-
(min(0,blockbandwidths(P,1)), blockbandwidth(P,2))
13+
(min(0,blockbandwidth(P,1)), blockbandwidth(P,2))
1414
else
1515
blockbandwidths(P)
1616
end

test/test_blockskylineqr.jl

Lines changed: 31 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,5 @@
1-
using BlockBandedMatrices, LinearAlgebra, MatrixFactorizations
1+
using BlockBandedMatrices, BlockArrays, LinearAlgebra, MatrixFactorizations, Test
2+
import BlockBandedMatrices: blockcolsupport
23

34
@testset "BlockBandedMatrix QR/QL" begin
45
@testset "Square QR" begin
@@ -118,7 +119,7 @@ using BlockBandedMatrices, LinearAlgebra, MatrixFactorizations
118119
N = 5
119120
A = BlockBandedMatrix{Float64}(undef, 1:N,1:N+1, (2,1))
120121
A.data .= randn.()
121-
122+
122123
@test_throws ArgumentError ql(A)
123124
end
124125

@@ -148,6 +149,34 @@ using BlockBandedMatrices, LinearAlgebra, MatrixFactorizations
148149
end
149150
end
150151
end
152+
153+
@testset "Fast QR \\" begin
154+
A = [1. 2; 3 4]; A = A + A'
155+
B = [5. 6; 7 8]
156+
N = 10_000;
157+
T = BlockBandedMatrix(mortar(Tridiagonal(fill(Matrix(B'),N-1), fill(zeros(2,2),N), fill(B,N-1))))
158+
z = 1+2im
159+
@time F = qr(T - z*I);
160+
@test blockbandwidths(UpperTriangular(F.factors)) == (0,2)
161+
@test blockcolsupport(UpperTriangular(F.factors),Block(4)) == Block.(2:4)
162+
V = view(F.factors,Block.(1:N), Block.(1:N))
163+
@test blocksize(V) == (N,N)
164+
@test @allocated(blocksize(V)) 40
165+
@test blockbandwidths(V) == (1,2)
166+
@test @allocated(blockbandwidths(V)) 40
167+
R = UpperTriangular(V)
168+
@test axes(R) == axes(V) == axes(F.factors)
169+
@test blocksize(R) == (N,N)
170+
@test @allocated(blocksize(R))  40
171+
@test blockbandwidths(R) == (0,2)
172+
@test @allocated(blockbandwidths(R)) 40
173+
@test blockcolsupport(R,Block(4)) == Block.(2:4)
174+
@test @allocated(blockcolsupport(R,Block(4))) 40
175+
176+
b = [1; zeros(size(T,1)-1)]
177+
B = [Matrix(I,2,2); zeros(size(T,1)-2,2)]
178+
@test ((T - z*I)\b)[1] (F\b)[1] (F \ B)[1,1] ((T - z*I)\B)[1,1] -0.1309123477325813 + 0.28471699370329884im
179+
end
151180
end
152181

153182

0 commit comments

Comments
 (0)