Skip to content

Commit 25de324

Browse files
authored
Fix infQR for upper triangular (#9)
* Fix infQR for upper triangular * Update test_infqr.jl * v0.1.1
1 parent 7e4896c commit 25de324

File tree

4 files changed

+19
-8
lines changed

4 files changed

+19
-8
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "InfiniteLinearAlgebra"
22
uuid = "cde9dba0-b1de-11e9-2c62-0bab9446c55c"
3-
version = "0.1"
3+
version = "0.1.1"
44

55
[deps]
66
ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a"

src/InfiniteLinearAlgebra.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -16,7 +16,7 @@ import LinearAlgebra: lmul!, ldiv!, matprod, qr, AbstractTriangular, AbstractQ,
1616
import LazyArrays: applybroadcaststyle, CachedArray, CachedMatrix, CachedVector, DenseColumnMajor, FillLayout, ApplyMatrix, check_mul_axes, ApplyStyle, LazyArrayApplyStyle, LazyArrayStyle,
1717
CachedMatrix, CachedArray, resizedata!, MemoryLayout, mulapplystyle, LmulStyle, RmulStyle,
1818
factorize, sub_materialize, LazyLayout, LazyArrayStyle,
19-
@lazymul, applylayout, ApplyLayout, PaddedLayout, materialize!
19+
@lazymul, applylayout, ApplyLayout, PaddedLayout, materialize!, zero!
2020
import MatrixFactorizations: ql, ql!, QLPackedQ, getL, getR, reflector!, reflectorApply!, QL, QR, QRPackedQ
2121

2222
import BlockArrays: BlockSizes, cumulsizes, _find_block, AbstractBlockVecOrMat, sizes_from_blocks

src/infqr.jl

Lines changed: 10 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -18,13 +18,17 @@ function partialqr!(F::AdaptiveQRData{<:Any,<:BandedMatrix}, n::Int)
1818
resizedata!(F.data,n+l,n+u+l);
1919
resize!(F.τ,n);
2020
= F.ncols
21-
factors = view(F.data.data,ñ+1:n+l,ñ+1:n);
2221
τ = view(F.τ,ñ+1:n);
23-
_banded_qr!(factors, τ);
24-
# multiply remaining columns
25-
= max(ñ+1,n-l-u+1) # first column interacting with extra data
26-
Q = QRPackedQ(view(F.data.data,n̄:n+l,n̄:n), view(F.τ,n̄:n))
27-
lmul!(Q',view(F.data.data,n̄:n+l,n+1:n+u+l))
22+
if l 0
23+
zero!(τ)
24+
else
25+
factors = view(F.data.data,ñ+1:n+l,ñ+1:n);
26+
_banded_qr!(factors, τ);
27+
# multiply remaining columns
28+
= max(ñ+1,n-l-u+1) # first column interacting with extra data
29+
Q = QRPackedQ(view(F.data.data,n̄:n+l,n̄:n), view(F.τ,n̄:n))
30+
lmul!(Q',view(F.data.data,n̄:n+l,n+1:n+u+l))
31+
end
2832
F.ncols = n
2933
end
3034
F

test/test_infqr.jl

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -125,4 +125,11 @@ import InfiniteLinearAlgebra: partialqr!, AdaptiveQRData, AdaptiveLayout
125125
@time x = qr(AB) \ b;
126126
@test x[1:300] AB[1:300,1:300] \ b[1:300]
127127
end
128+
129+
@testset "triangular infqr" begin
130+
A = BandedMatrix(0 => 1:∞, 2 => Ones(∞))
131+
F = qr(A)
132+
@test F.Q[1:10,1:10] == Eye(10)
133+
@test F.R[1:10,1:10] == A[1:10,1:10]
134+
end
128135
end

0 commit comments

Comments
 (0)