Skip to content

Commit

Permalink
Fix method ambiguity for qr (#931) and lingering ambiguities for `l…
Browse files Browse the repository at this point in the history
…u` (#932)

* fix `qr` method ambiguities (#931) and lingering `lu` ambiguities (#920)

* fix inferrence issues due to using `@invoke` for `lu` keyword arguments

* bump version
  • Loading branch information
thchr committed Jul 2, 2021
1 parent 1a06987 commit 7dc6499
Show file tree
Hide file tree
Showing 5 changed files with 51 additions and 31 deletions.
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "StaticArrays"
uuid = "90137ffa-7385-5640-81b9-e52037218182"
version = "1.2.4"
version = "1.2.5"

[deps]
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
Expand Down
28 changes: 12 additions & 16 deletions src/lu.jl
Original file line number Diff line number Diff line change
Expand Up @@ -30,25 +30,21 @@ function Base.show(io::IO, mime::MIME{Symbol("text/plain")}, F::LU)
end

# LU decomposition
function lu(A::StaticMatrix, pivot::Union{Val{false},Val{true}}=Val(true); check = true)
L, U, p = _lu(A, pivot, check)
LU(L, U, p)
end

# For the square version, return explicit lower and upper triangular matrices.
# We would do this for the rectangular case too, but Base doesn't support that.
function lu(A::StaticMatrix{N,N}, pivot::Union{Val{false},Val{true}}=Val(true);
check = true) where {N}
L, U, p = _lu(A, pivot, check)
LU(LowerTriangular(L), UpperTriangular(U), p)
end
for pv in (:true, :false)
# ... define each `pivot::Val{true/false}` method individually to avoid ambiguties
@eval function lu(A::StaticMatrix, pivot::Val{$pv}; check = true)
L, U, p = _lu(A, pivot, check)
LU(L, U, p)
end

@static if VERSION >= v"1.7-DEV"
# disambiguation
function lu(A::StaticMatrix{N,N}, pivot::Val{true}) where {N}
Base.@invoke lu(A::StaticMatrix{N,N} where N, pivot::Union{Val{false},Val{true}})
# For the square version, return explicit lower and upper triangular matrices.
# We would do this for the rectangular case too, but Base doesn't support that.
@eval function lu(A::StaticMatrix{N,N}, pivot::Val{$pv}; check = true) where {N}
L, U, p = _lu(A, pivot, check)
LU(LowerTriangular(L), UpperTriangular(U), p)
end
end
lu(A::StaticMatrix; check = true) = lu(A, Val(true); check=check)

# location of the first zero on the diagonal, 0 when not found
function _first_zero_on_diagonal(A::StaticMatrix{M,N,T}) where {M,N,T}
Expand Down
33 changes: 19 additions & 14 deletions src/qr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -11,10 +11,26 @@ Base.iterate(S::QR, ::Val{:R}) = (S.R, Val(:p))
Base.iterate(S::QR, ::Val{:p}) = (S.p, Val(:done))
Base.iterate(S::QR, ::Val{:done}) = nothing

for pv in (:true, :false)
@eval begin
@inline function qr(A::StaticMatrix, pivot::Val{$pv})
QRp = _qr(Size(A), A, pivot)
if length(QRp) === 2
# create an identity permutation since that is cheap,
# and much safer since, in the case of isbits types, we can't
# safely leave the field undefined.
p = identity_perm(QRp[2])
return QR(QRp[1], QRp[2], p)
else # length(QRp) === 3
return QR(QRp[1], QRp[2], QRp[3])
end
end
end
end
"""
qr(A::StaticMatrix, pivot=Val(false))
qr(A::StaticMatrix, pivot::Union{Val{true}, Val{false}} = Val(false))
Compute the QR factorization of `A`. The factors can be obtain by iteration:
Compute the QR factorization of `A`. The factors can be obtained by iteration:
```julia
julia> A = @SMatrix rand(3,4);
Expand All @@ -34,18 +50,7 @@ julia> F.Q * F.R ≈ A
true
```
"""
@inline function qr(A::StaticMatrix, pivot::Union{Val{false}, Val{true}} = Val(false))
QRp = _qr(Size(A), A, pivot)
if length(QRp) === 2
# create an identity permutation since that is cheap,
# and much safer since, in the case of isbits types, we can't
# safely leave the field undefined.
p = identity_perm(QRp[2])
return QR(QRp[1], QRp[2], p)
else # length(QRp) === 3
return QR(QRp[1], QRp[2], QRp[3])
end
end
qr(A::StaticMatrix) = qr(A, Val(false))

function identity_perm(R::StaticMatrix{N,M,T}) where {N,M,T}
return similar_type(R, Int, Size((M,)))(ntuple(x -> x, Val{M}()))
Expand Down
11 changes: 11 additions & 0 deletions test/lu.jl
Original file line number Diff line number Diff line change
Expand Up @@ -65,3 +65,14 @@ end
@test_throws SingularException lu(A)
@test !issuccess(lu(A; check = false))
end

@testset "LU method ambiguity" begin
# Issue #920; just test that methods do not throw an ambiguity error when called
for A in ((@SMatrix [1.0 2.0; 3.0 4.0]), (@SMatrix [1.0 2.0 3.0; 4.0 5.0 6.0]))
@test isa(lu(A), StaticArrays.LU)
@test isa(lu(A, Val(true)), StaticArrays.LU)
@test isa(lu(A, Val(false)), StaticArrays.LU)
@test isa(lu(A; check=false), StaticArrays.LU)
@test isa(lu(A; check=true), StaticArrays.LU)
end
end
8 changes: 8 additions & 0 deletions test/qr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -60,3 +60,11 @@ Random.seed!(42)
test_qr(arr)
end
end

@testset "QR method ambiguity" begin
# Issue #931; just test that methods do not throw an ambiguity error when called
A = @SMatrix [1.0 2.0 3.0; 4.0 5.0 6.0]
@test isa(qr(A), StaticArrays.QR)
@test isa(qr(A, Val(true)), StaticArrays.QR)
@test isa(qr(A, Val(false)), StaticArrays.QR)
end

2 comments on commit 7dc6499

@mateuszbaran
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/40103

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v1.2.5 -m "<description of version>" 7dc649948b26fc7e5c7525560d41064f221244a5
git push origin v1.2.5

Please sign in to comment.