Skip to content

Commit

Permalink
Merge pull request #10834 from stevengj/sparsefact-nnz
Browse files Browse the repository at this point in the history
nnz methods for sparse factorizations
  • Loading branch information
andreasnoack committed Apr 15, 2015
2 parents 4234487 + 1600f6a commit f443d48
Show file tree
Hide file tree
Showing 4 changed files with 12 additions and 3 deletions.
6 changes: 4 additions & 2 deletions base/sparse/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,7 @@ import Base.LinAlg: (\), A_mul_Bc, A_mul_Bt, Ac_ldiv_B, Ac_mul_B, At_ldiv_B, At_
cholfact, cholfact!, det, diag, ishermitian, isposdef,
issym, ldltfact, logdet

import Base.SparseMatrix: sparse
import Base.SparseMatrix: sparse, nnz

export
Dense,
Expand Down Expand Up @@ -884,13 +884,15 @@ eltype{T<:VTypes}(A::Dense{T}) = T
eltype{T<:VTypes}(A::Factor{T}) = T
eltype{T<:VTypes}(A::Sparse{T}) = T

nnz(F::Factor) = nnz(Sparse(F))

function show(io::IO, F::Factor)
s = unsafe_load(F.p)
println(io, typeof(F))
@printf(io, "type: %12s\n", Bool(s.is_ll) ? "LLt" : "LDLt")
@printf(io, "method: %10s\n", Bool(s.is_super) ? "supernodal" : "simplicial")
@printf(io, "maxnnz: %10d\n", Int(s.nzmax))
@printf(io, "nnz: %13d\n", nnz(Sparse(F)))
@printf(io, "nnz: %13d\n", nnz(F))
end

isvalid(A::Dense) = check_dense(A)
Expand Down
7 changes: 6 additions & 1 deletion base/sparse/umfpack.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,7 @@ import Base: (\), Ac_ldiv_B, At_ldiv_B, findnz, getindex, show, size
import Base.LinAlg: A_ldiv_B!, Ac_ldiv_B!, At_ldiv_B!, Factorization, det, lufact

importall Base.SparseMatrix
import Base.SparseMatrix: increment, increment!, decrement, decrement!
import Base.SparseMatrix: increment, increment!, decrement, decrement!, nnz

include("umfpack_h.jl")
type MatrixIllConditionedException <: Exception
Expand Down Expand Up @@ -323,6 +323,11 @@ for itype in UmfpackIndexTypes
end
end

function nnz(lu::UmfpackLU)
lnz, unz, = umf_lunz(lu)
return Int(lnz + unz)
end

### Solve with Factorization
A_ldiv_B!{T<:UMFVTypes}(lu::UmfpackLU{T}, b::Vector{T}) = solve(lu, b, UMFPACK_A)
A_ldiv_B!{T<:UMFVTypes}(lu::UmfpackLU{T}, b::Matrix{T}) = solve(lu, b, UMFPACK_A)
Expand Down
1 change: 1 addition & 0 deletions test/sparsedir/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -104,6 +104,7 @@ chma = cholfact(A) # LL' form
@test unsafe_load(chma.p).is_ll == 1 # check that it is in fact an LLt
x = chma\B
@test_approx_eq x ones(size(x))
@test nnz(chma) == 489

#lp_afiro example
afiro = CHOLMOD.Sparse(27, 51,
Expand Down
1 change: 1 addition & 0 deletions test/sparsedir/umfpack.jl
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@ for Tv in (Float64, Complex128)
for Ti in Base.SparseMatrix.UMFPACK.UMFITypes.types
A = convert(SparseMatrixCSC{Tv,Ti}, A0)
lua = lufact(A)
@test nnz(lua) == 18
L,U,p,q,Rs = lua[:(:)]
@test_approx_eq scale(Rs,A)[p,q] L*U

Expand Down

0 comments on commit f443d48

Please sign in to comment.