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

QR: handle xtype/dtype returned from LibSuiteSparse that don't match matrix element type #586

Merged
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
9 changes: 9 additions & 0 deletions src/solvers/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -901,6 +901,15 @@

Dense(A::Sparse) = sparse_to_dense(A)

function Dense(ptr::Ptr{cholmod_dense})
if ptr == C_NULL
throw(ArgumentError("dense matrix construction failed for " *

Check warning on line 906 in src/solvers/cholmod.jl

View check run for this annotation

Codecov / codecov/patch

src/solvers/cholmod.jl#L906

Added line #L906 was not covered by tests
"unknown reasons. Please submit a bug report."))
end
s = unsafe_load(ptr)
return Dense{jlxtype(s.xtype, s.dtype)}(ptr)
end

function Base.convert(::Type{Dense{Tnew}}, A::Dense{T}) where {Tnew, T}
GC.@preserve A begin
Ap = unsafe_load(pointer(A))
Expand Down
8 changes: 4 additions & 4 deletions src/solvers/spqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -204,10 +204,10 @@ function LinearAlgebra.qr(A::SparseMatrixCSC{Tv, Ti}; tol=_default_tol(A), order
C_NULL, C_NULL, C_NULL, C_NULL,
R, E, H, HPinv, HTau)

R_ = SparseMatrixCSC(Sparse(R[]))
return QRSparse(SparseMatrixCSC(Sparse(H[])),
vec(Array(CHOLMOD.Dense{Tv}(HTau[]))),
SparseMatrixCSC(min(size(A)...),
R_ = SparseMatrixCSC{Tv, Ti}(Sparse(R[]))
return QRSparse(SparseMatrixCSC{Tv, Ti}(Sparse(H[])),
vec(Array{Tv}(CHOLMOD.Dense(HTau[]))),
SparseMatrixCSC{Tv, Ti}(min(size(A)...),
size(R_, 2),
getcolptr(R_),
rowvals(R_),
Expand Down
6 changes: 6 additions & 0 deletions test/spqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -110,6 +110,12 @@ end
@test (F.Q*F.R)::SparseMatrixCSC == A[F.prow,F.pcol]
end

@testset "Issue #585 for element type: $eltyA" for eltyA in (Float64, Float32)
A = sparse(eltyA[1 0; 0 1])
F = qr(A)
@test eltype(F.Q) == eltype(F.R) == eltyA
end

@testset "select ordering overdetermined" begin
A = sparse([1:n; rand(1:m, nn - n)], [1:n; rand(1:n, nn - n)], randn(nn), m, n)
b = randn(m)
Expand Down
Loading