Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add AbstractPDMat constructor #172

Merged
merged 2 commits into from
Feb 27, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 4 additions & 0 deletions src/generics.jl
Original file line number Diff line number Diff line change
@@ -1,5 +1,9 @@
# Generic functions (on top of the type-specific implementations)

## constructors
AbstractPDMat(A::AbstractPDMat) = A
AbstractPDMat(A::AbstractMatrix) = PDMat(A)

## arithmetics

pdadd!(r::Matrix, a::Matrix, b::AbstractPDMat{T}) where {T<:Real} = pdadd!(r, a, b, one(T))
Expand Down
4 changes: 4 additions & 0 deletions src/pdiagmat.jl
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,10 @@ end

PDiagMat(v::AbstractVector{<:Real}) = PDiagMat{eltype(v),typeof(v)}(length(v), v)

AbstractPDMat(A::Diagonal{<:Real}) = PDiagMat(A.diag)
AbstractPDMat(A::Symmetric{<:Real,<:Diagonal{<:Real}}) = PDiagMat(A.data.diag)
AbstractPDMat(A::Hermitian{<:Real,<:Diagonal{<:Real}}) = PDiagMat(A.data.diag)

### Conversion
Base.convert(::Type{PDiagMat{T}}, a::PDiagMat) where {T<:Real} = PDiagMat(convert(AbstractArray{T}, a.diag))
Base.convert(::Type{AbstractArray{T}}, a::PDiagMat) where {T<:Real} = convert(PDiagMat{T}, a)
Expand Down
2 changes: 2 additions & 0 deletions src/pdmat.jl
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,8 @@ end
PDMat(mat::AbstractMatrix) = PDMat(mat, cholesky(mat))
PDMat(fac::Cholesky) = PDMat(AbstractMatrix(fac), fac)

AbstractPDMat(A::Cholesky) = PDMat(A)

### Conversion
Base.convert(::Type{PDMat{T}}, a::PDMat) where {T<:Real} = PDMat(convert(AbstractArray{T}, a.mat))
Base.convert(::Type{AbstractArray{T}}, a::PDMat) where {T<:Real} = convert(PDMat{T}, a)
Expand Down
5 changes: 4 additions & 1 deletion src/pdsparsemat.jl
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,10 @@ function PDSparseMat(mat::AbstractSparseMatrix,chol::CholTypeSparse)
end

PDSparseMat(mat::SparseMatrixCSC) = PDSparseMat(mat, cholesky(mat))
PDSparseMat(fac::CholTypeSparse) = PDSparseMat(sparse(fac) |> x -> x*x', fac)
PDSparseMat(fac::CholTypeSparse) = PDSparseMat(sparse(fac), fac)

AbstractPDMat(A::SparseMatrixCSC) = PDSparseMat(A)
AbstractPDMat(A::CholTypeSparse) = PDSparseMat(A)

### Conversion
Base.convert(::Type{PDSparseMat{T}}, a::PDSparseMat) where {T<:Real} = PDSparseMat(convert(SparseMatrixCSC{T}, a.mat))
Expand Down
38 changes: 38 additions & 0 deletions test/pdmtypes.jl
Original file line number Diff line number Diff line change
Expand Up @@ -121,4 +121,42 @@ using Test
@test Matrix(M) ≈ A
end
end

@testset "AbstractPDMat constructors (#136)" begin
x = rand(10, 10)
A = x' * x + I

M = @inferred AbstractPDMat(A)
@test M isa PDMat
@test Matrix(M) ≈ A

M = @inferred AbstractPDMat(cholesky(A))
@test M isa PDMat
@test Matrix(M) ≈ A

M = @inferred AbstractPDMat(Diagonal(A))
@test M isa PDiagMat
@test Matrix(M) ≈ Diagonal(A)

M = @inferred AbstractPDMat(Symmetric(Diagonal(A)))
@test M isa PDiagMat
@test Matrix(M) ≈ Diagonal(A)

M = @inferred AbstractPDMat(Hermitian(Diagonal(A)))
@test M isa PDiagMat
@test Matrix(M) ≈ Diagonal(A)

M = @inferred AbstractPDMat(sparse(A))
@test M isa PDSparseMat
@test Matrix(M) ≈ A

if VERSION < v"1.6"
# inference fails e.g. on Julia 1.0
M = AbstractPDMat(cholesky(sparse(A)))
else
M = @inferred AbstractPDMat(cholesky(sparse(A)))
end
@test M isa PDSparseMat
@test Matrix(M) ≈ A
end
end
7 changes: 7 additions & 0 deletions test/testutils.jl
Original file line number Diff line number Diff line change
Expand Up @@ -91,6 +91,13 @@ function pdtest_basics(C, Cmat::Matrix, d::Int, verbose::Int)

_pdt(verbose, "ishermitian")
@test ishermitian(C)

_pdt(verbose, "AbstractPDMat")
M = AbstractPDMat(C)
@test M isa AbstractPDMat
if C isa AbstractPDMat
@test M === C
end
end


Expand Down