Skip to content

Commit b85c6d6

Browse files
authored
Merge branch 'master' into st/remove_chol_copy
2 parents 35e8d2e + 443aa0f commit b85c6d6

File tree

3 files changed

+33
-7
lines changed

3 files changed

+33
-7
lines changed

.ci/Manifest.toml

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -52,7 +52,7 @@ version = "1.11.0"
5252

5353
[[deps.JuliaSyntaxHighlighting]]
5454
deps = ["StyledStrings"]
55-
uuid = "dc6e5ff7-fb65-4e79-a425-ec3bc9c03011"
55+
uuid = "ac6e5ff7-fb65-4e79-a425-ec3bc9c03011"
5656
version = "1.12.0"
5757

5858
[[deps.LazyArtifacts]]
@@ -93,7 +93,7 @@ version = "1.11.0"
9393
deps = ["Libdl", "OpenBLAS_jll", "libblastrampoline_jll"]
9494
path = ".."
9595
uuid = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
96-
version = "1.11.0"
96+
version = "1.12.0"
9797

9898
[[deps.Logging]]
9999
uuid = "56ddb016-857b-54e1-b83d-db4d58db5568"

src/dense.jl

Lines changed: 12 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -1587,24 +1587,25 @@ factorize(A::Transpose) = transpose(factorize(parent(A)))
15871587
factorize(a::Number) = a # same as how factorize behaves on Diagonal types
15881588

15891589
function getstructure(A::StridedMatrix)
1590+
require_one_based_indexing(A)
15901591
m, n = size(A)
15911592
if m == 1 return A[1] end
15921593
utri = true
15931594
utri1 = true
15941595
herm = true
15951596
sym = true
1596-
for j = 1:n-1, i = j+1:m
1597-
if utri1
1597+
for j = 1:n, i = j:m
1598+
if (j < n) && (i > j) && utri1 # indices are off-diagonal
15981599
if A[i,j] != 0
15991600
utri1 = i == j + 1
16001601
utri = false
16011602
end
16021603
end
16031604
if sym
1604-
sym &= A[i,j] == A[j,i]
1605+
sym &= A[i,j] == transpose(A[j,i])
16051606
end
16061607
if herm
1607-
herm &= A[i,j] == conj(A[j,i])
1608+
herm &= A[i,j] == adjoint(A[j,i])
16081609
end
16091610
if !(utri1|herm|sym) break end
16101611
end
@@ -1617,10 +1618,12 @@ function getstructure(A::StridedMatrix)
16171618
if ltri1
16181619
for i = 1:n-1
16191620
if A[i,i+1] != 0
1620-
ltri &= false
1621+
ltri = false
16211622
break
16221623
end
16231624
end
1625+
else
1626+
ltri = false
16241627
end
16251628
return (utri, utri1, ltri, ltri1, sym, herm)
16261629
end
@@ -1779,6 +1782,10 @@ Condition number of the matrix `M`, computed using the operator `p`-norm. Valid
17791782
"""
17801783
function cond(A::AbstractMatrix, p::Real=2)
17811784
if p == 2
1785+
if isempty(A)
1786+
checksquare(A)
1787+
return zero(real(eigtype(eltype(A))))
1788+
end
17821789
v = svdvals(A)
17831790
maxv = maximum(v)
17841791
return iszero(maxv) ? oftype(real(maxv), Inf) : maxv / minimum(v)

test/dense.jl

Lines changed: 19 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -56,6 +56,18 @@ Random.seed!(1234323)
5656
@test cond(Mars, 2) 6.181867355918493
5757
@test cond(Mars, Inf) 7.1
5858
end
59+
@testset "Empty matrices" begin
60+
for p in (1,2,Inf)
61+
# zero for square (i.e. 0×0) matrices
62+
@test cond(zeros(Int, 0, 0), p) === 0.0
63+
@test cond(zeros(0, 0), p) === 0.0
64+
@test cond(zeros(ComplexF64, 0, 0), p) === 0.0
65+
# error for non-square matrices
66+
for size in ((10,0), (0,10))
67+
@test_throws DimensionMismatch cond(zeros(size...), p)
68+
end
69+
end
70+
end
5971
end
6072

6173
areal = randn(n,n)/2
@@ -1347,4 +1359,11 @@ end
13471359
end
13481360
end
13491361

1362+
@testset "structure of dense matrices" begin
1363+
# A is neither triangular nor symmetric/Hermitian
1364+
A = [1 im 2; -im 0 3; 2 3 im]
1365+
@test factorize(A) isa LU{ComplexF64, Matrix{ComplexF64}, Vector{Int}}
1366+
@test !any(LinearAlgebra.getstructure(A))
1367+
end
1368+
13501369
end # module TestDense

0 commit comments

Comments
 (0)