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

Avoid reading inactive rowvals when computing the rank of a QRSparse #29925

Merged
merged 1 commit into from
Nov 5, 2018
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: 2 additions & 2 deletions stdlib/SuiteSparse/src/spqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -23,7 +23,7 @@ const ORDERING_BESTAMD = Int32(9) # try COLAMD and AMD; pick best#
# tried. If there is a high fill-in with AMD then try METIS(A'A) and take
# the best of AMD and METIS. METIS is not tried if it isn't installed.

using SparseArrays: SparseMatrixCSC
using SparseArrays: SparseMatrixCSC, nnz
using ..SuiteSparse.CHOLMOD
using ..SuiteSparse.CHOLMOD: change_stype!, free!

Expand Down Expand Up @@ -342,7 +342,7 @@ end
_ret_size(F::QRSparse, b::AbstractVector) = (size(F, 2),)
_ret_size(F::QRSparse, B::AbstractMatrix) = (size(F, 2), size(B, 2))

LinearAlgebra.rank(F::QRSparse) = maximum(F.R.rowval)
LinearAlgebra.rank(F::QRSparse) = reduce(max, view(F.R.rowval, 1:nnz(F.R)), init = eltype(F.R.rowval)(0))
Copy link
Member

Choose a reason for hiding this comment

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

Use SparseArrays.nzvalview instead?

Copy link
Member Author

Choose a reason for hiding this comment

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

I'm using the rowvals here, not the nzvals. We could add rowvalview but I think we should wait until there are more use cases.


function (\)(F::QRSparse{T}, B::VecOrMat{Complex{T}}) where T<:LinearAlgebra.BlasReal
# |z1|z3| reinterpret |x1|x2|x3|x4| transpose |x1|y1| reshape |x1|y1|x3|y3|
Expand Down
6 changes: 6 additions & 0 deletions stdlib/SuiteSparse/test/spqr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -93,4 +93,10 @@ end
@test propertynames(F) == (:R, :Q, :prow, :pcol)
@test propertynames(F, true) == (:R, :Q, :prow, :pcol, :factors, :τ, :cpiv, :rpivinv)
end

@testset "rank" begin
@test rank(qr(sprandn(10, 5, 1.0)*sprandn(5, 10, 1.0))) == 5
@test all(iszero, (rank(qr(spzeros(10, i))) for i in 1:10))
end

end