Skip to content

Commit 84df00b

Browse files
Improve chol_lower and chol_upper (#149)
* Improve `chol_lower` and `chol_upper` * Increase allowed overhead for Ubuntu x64 on Julia 1.0 Co-authored-by: Andreas Noack <[email protected]>
1 parent 8727724 commit 84df00b

File tree

3 files changed

+13
-6
lines changed

3 files changed

+13
-6
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "PDMats"
22
uuid = "90014a1f-27ba-587c-ab20-58faa44d9150"
3-
version = "0.11.7"
3+
version = "0.11.8"
44

55
[deps]
66
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"

src/chol.jl

Lines changed: 10 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,14 @@
11
# Accessing a.L directly might involve an extra copy();
2-
# instead, always use the stored Cholesky factor:
3-
chol_lower(a::Cholesky) = a.uplo === 'L' ? a.L : a.U'
4-
chol_upper(a::Cholesky) = a.uplo === 'U' ? a.U : a.L'
2+
# instead, always use the stored Cholesky factor
3+
# Using `a.factors` instead of `a.L` or `a.U` avoids one
4+
# additional `LowerTriangular` or `UpperTriangular` wrapper and
5+
# leads to better performance
6+
function chol_lower(a::Cholesky)
7+
return a.uplo === 'L' ? LowerTriangular(a.factors) : LowerTriangular(a.factors')
8+
end
9+
function chol_upper(a::Cholesky)
10+
return a.uplo === 'U' ? UpperTriangular(a.factors) : UpperTriangular(a.factors')
11+
end
512

613
# For a dense Matrix, the following allows us to avoid the Adjoint wrapper:
714
chol_lower(a::Matrix) = cholesky(Symmetric(a, :L)).L

test/chol.jl

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -13,8 +13,8 @@ using PDMats: chol_lower, chol_upper
1313
for uplo in (:L, :U)
1414
ch = cholesky(Symmetric(C, uplo))
1515
chol_lower(ch)
16-
@test (@allocated chol_lower(ch)) < 50 # allow small overhead for wrapper types
16+
@test (@allocated chol_lower(ch)) < 33 # allow small overhead for wrapper types
1717
chol_upper(ch)
18-
@test (@allocated chol_upper(ch)) < 50 # allow small overhead for wrapper types
18+
@test (@allocated chol_upper(ch)) < 33 # allow small overhead for wrapper types
1919
end
2020
end

0 commit comments

Comments
 (0)