Skip to content

Commit 589e9e4

Browse files
committed
support user-specified permutation in CHOLMOD factorizations, change β optional argument to a shift keyword arg, some code consolidation; change & to new Ref syntax in cholmod
1 parent f443d48 commit 589e9e4

File tree

4 files changed

+84
-103
lines changed

4 files changed

+84
-103
lines changed

base/deprecated.jl

+3
Original file line numberDiff line numberDiff line change
@@ -439,6 +439,9 @@ end
439439
@deprecate dlsym_e Libdl.dlsym_e
440440
@deprecate find_library Libdl.find_library
441441

442+
@deprecate cholfact(A::AbstractMatrix, β::Number) cholfact(A, shift=β)
443+
@deprecate ldltfact(A::AbstractMatrix, β::Number) ldltfact(A, shift=β)
444+
442445
# 0.4 discontinued functions
443446

444447
@noinline function subtypetree(x::DataType, level=-1)

base/sparse/cholmod.jl

+55-89
Original file line numberDiff line numberDiff line change
@@ -33,6 +33,8 @@ const common_final_ll = (1:4) + cholmod_com_offsets[7]
3333
const common_print = (1:4) + cholmod_com_offsets[13]
3434
const common_itype = (1:4) + cholmod_com_offsets[18]
3535
const common_dtype = (1:4) + cholmod_com_offsets[19]
36+
const common_nmethods = (1:4) + cholmod_com_offsets[15]
37+
const common_postorder = (1:4) + cholmod_com_offsets[17]
3638

3739
## macro to generate the name of the C function according to the integer type
3840
macro cholmod_name(nm,typ) string("cholmod_", eval(typ) == SuiteSparse_long ? "l_" : "", nm) end
@@ -283,7 +285,7 @@ function allocate_dense(nrow::Integer, ncol::Integer, d::Integer, ::Type{Complex
283285
d
284286
end
285287

286-
free_dense!{T}(p::Ptr{C_Dense{T}}) = ccall((:cholmod_l_free_dense, :libcholmod), Cint, (Ptr{Ptr{C_Dense{T}}}, Ptr{Void}), &p, common(Cint))
288+
free_dense!{T}(p::Ptr{C_Dense{T}}) = ccall((:cholmod_l_free_dense, :libcholmod), Cint, (Ref{Ptr{C_Dense{T}}}, Ptr{Void}), p, common(Cint))
287289

288290
function zeros{T<:VTypes}(m::Integer, n::Integer, ::Type{T})
289291
d = Dense(ccall((:cholmod_l_zeros, :libcholmod), Ptr{C_Dense{T}},
@@ -568,18 +570,18 @@ for Ti in IndexTypes
568570
A
569571
end
570572

571-
function sdmult!{Tv<:VTypes}(A::Sparse{Tv,$Ti}, transpose::Bool, α::Tv, β::Tv, X::Dense{Tv}, Y::Dense{Tv})
573+
function sdmult!{Tv<:VTypes}(A::Sparse{Tv,$Ti}, transpose::Bool, α::Number, β::Number, X::Dense{Tv}, Y::Dense{Tv})
572574
m, n = size(A)
573575
nc = transpose ? m : n
574576
nr = transpose ? n : m
575577
if nc != size(X, 1)
576578
throw(DimensionMismatch("incompatible dimensions, $nc and $(size(X,1))"))
577579
end
578580
@isok ccall((@cholmod_name("sdmult", $Ti),:libcholmod), Cint,
579-
(Ptr{C_Sparse{Tv,$Ti}}, Cint, Ptr{Cdouble}, Ptr{Cdouble},
580-
Ptr{C_Dense{Tv}}, Ptr{C_Dense{Tv}}, Ptr{UInt8}),
581-
A.p, transpose, &α, &β,
582-
X.p, Y.p, common($Ti))
581+
(Ptr{C_Sparse{Tv,$Ti}}, Cint,
582+
Ref{Complex128}, Ref{Complex128},
583+
Ptr{C_Dense{Tv}}, Ptr{C_Dense{Tv}}, Ptr{UInt8}),
584+
A.p, transpose, α, β, X.p, Y.p, common($Ti))
583585
Y
584586
end
585587

@@ -605,7 +607,7 @@ for Ti in IndexTypes
605607
end
606608

607609
# cholmod_cholesky.h
608-
# For analyze, factorize and factorize_p, the Common argument must be
610+
# For analyze, analyze_p, and factorize_p!, the Common argument must be
609611
# supplied in order to control if the factorization is LLt or LDLt
610612
function analyze{Tv<:VTypes}(A::Sparse{Tv,$Ti}, cmmn::Vector{UInt8})
611613
f = Factor(ccall((@cholmod_name("analyze", $Ti),:libcholmod), Ptr{C_Factor{Tv,$Ti}},
@@ -614,18 +616,22 @@ for Ti in IndexTypes
614616
finalizer(f, free!)
615617
f
616618
end
617-
function factorize!{Tv<:VTypes}(A::Sparse{Tv,$Ti}, F::Factor{Tv,$Ti}, cmmn::Vector{UInt8})
618-
@isok ccall((@cholmod_name("factorize", $Ti),:libcholmod), Cint,
619-
(Ptr{C_Sparse{Tv,$Ti}}, Ptr{C_Factor{Tv,$Ti}}, Ptr{UInt8}),
620-
A.p, F.p, cmmn)
621-
F
619+
function analyze_p{Tv<:VTypes}(A::Sparse{Tv,$Ti}, perm::Vector{$Ti},
620+
cmmn::Vector{UInt8})
621+
length(perm) != size(A,1) && throw(BoundsError())
622+
f = Factor(ccall((@cholmod_name("analyze_p", $Ti),:libcholmod), Ptr{C_Factor{Tv,$Ti}},
623+
(Ptr{C_Sparse{Tv,$Ti}}, Ptr{$Ti}, Ptr{$Ti}, Csize_t, Ptr{UInt8}),
624+
A.p, perm, C_NULL, 0, cmmn))
625+
finalizer(f, free!)
626+
f
622627
end
623-
function factorize_p!{Tv<:VTypes}(A::Sparse{Tv,$Ti}, β::Tv, fset::Vector{$Ti}, F::Factor{Tv,$Ti}, cmmn::Vector{UInt8})
628+
function factorize_p!{Tv<:VTypes}(A::Sparse{Tv,$Ti}, β::Real, F::Factor{Tv,$Ti}, cmmn::Vector{UInt8})
629+
# note that β is passed as a complex number (double beta[2]),
630+
# but the CHOLMOD manual says that only beta[0] (real part) is used
624631
@isok ccall((@cholmod_name("factorize_p", $Ti),:libcholmod), Cint,
625-
(Ptr{C_Sparse{Tv,$Ti}}, Ptr{Cdouble}, Ptr{$Ti}, Csize_t,
632+
(Ptr{C_Sparse{Tv,$Ti}}, Ref{Complex128}, Ptr{$Ti}, Csize_t,
626633
Ptr{C_Factor{Tv,$Ti}}, Ptr{UInt8}),
627-
A.p, &β, fset, length(fset),
628-
F.p, cmmn)
634+
A.p, β, C_NULL, 0, F.p, cmmn)
629635
F
630636
end
631637

@@ -984,42 +990,45 @@ Ac_mul_B(A::Sparse, B::VecOrMat) = Ac_mul_B(A, Dense(B))
984990

985991

986992
## Factorization methods
987-
function cholfact(A::Sparse)
993+
994+
function fact_{Tv<:VTypes,Ti<:ITypes,Ti2<:Integer}(A::Sparse{Tv,Ti}, cm::Array{UInt8};
995+
shift::Real=0.0,
996+
perm::AbstractVector{Ti2}=Int[],
997+
postorder::Bool=true,
998+
userperm_only::Bool=true)
988999
sA = unsafe_load(A.p)
9891000
sA.stype == 0 && throw(ArgumentError("sparse matrix is not symmetric/Hermitian"))
9901001

991-
cm = common(indtype(A))
1002+
if !postorder
1003+
cm[common_postorder] = reinterpret(UInt8, [zero(Cint)])
1004+
end
9921005

993-
# Hack! makes it a llt
994-
cm[common_final_ll] = reinterpret(UInt8, [one(Cint)])
995-
F = analyze(A, cm)
996-
factorize!(A, F, cm)
997-
s = unsafe_load(F.p)
998-
s.minor < size(A, 1) && throw(Base.LinAlg.PosDefException(s.minor))
1006+
if isempty(perm)
1007+
F = analyze(A, cm)
1008+
else # user permutation provided
1009+
if userperm_only # use perm even if it is worse than AMD
1010+
cm[common_nmethods] = reinterpret(UInt8, [one(Cint)])
1011+
end
1012+
F = analyze_p(A, Ti[p-1 for p in perm], cm)
1013+
end
1014+
1015+
factorize_p!(A, shift, F, cm)
9991016
return F
10001017
end
10011018

1002-
function cholfact{Tv<:VTypes,Ti<:ITypes}(A::Sparse{Tv,Ti}, β::Tv)
1003-
sA = unsafe_load(A.p)
1004-
sA.stype == 0 && throw(ArgumentError("sparse matrix is not symmetric/Hermitian"))
1005-
1006-
cm = common(Ti)
1019+
function cholfact(A::Sparse; kws...)
1020+
cm = common(indtype(A))
10071021

10081022
# Hack! makes it a llt
10091023
cm[common_final_ll] = reinterpret(UInt8, [one(Cint)])
10101024

1011-
F = analyze(A, cm)
1012-
factorize_p!(A, β, Ti[0:sA.ncol - 1;], F, cm)
1013-
1025+
F = fact_(A, cm; kws...)
10141026
s = unsafe_load(F.p)
10151027
s.minor < size(A, 1) && throw(Base.LinAlg.PosDefException(s.minor))
10161028
return F
10171029
end
10181030

1019-
function ldltfact(A::Sparse)
1020-
sA = unsafe_load(A.p)
1021-
sA.stype == 0 && throw(ArgumentError("sparse matrix is not symmetric/Hermitian"))
1022-
1031+
function ldltfact(A::Sparse; kws...)
10231032
cm = common(indtype(A))
10241033

10251034
# Hack! makes it a ldlt
@@ -1028,72 +1037,29 @@ function ldltfact(A::Sparse)
10281037
# Hack! really make sure it's a ldlt by avoiding supernodal factorisation
10291038
cm[common_supernodal] = reinterpret(UInt8, [zero(Cint)])
10301039

1031-
F = analyze(A, cm)
1032-
factorize!(A, F, cm)
1033-
1034-
# Check if decomposition failed
1035-
s = unsafe_load(F.p)
1036-
s.minor < size(A, 1) && throw(Base.LinAlg.ArgumentError("matrix has one or more zero pivots"))
1037-
1038-
return F
1039-
end
1040-
1041-
function ldltfact{Tv<:VTypes,Ti<:ITypes}(A::Sparse{Tv,Ti}, β::Tv)
1042-
sA = unsafe_load(A.p)
1043-
sA.stype == 0 && throw(ArgumentError("sparse matrix is not symmetric/Hermitian"))
1044-
1045-
cm = common(Ti)
1046-
1047-
# Hack! makes it a ldlt
1048-
cm[common_final_ll] = reinterpret(UInt8, [zero(Cint)])
1049-
1050-
# Hack! really make sure it's a ldlt by avoiding supernodal factorisation
1051-
cm[common_supernodal] = reinterpret(UInt8, [zero(Cint)])
1052-
1053-
F = analyze(A, cm)
1054-
factorize_p!(A, β, Ti[0:sA.ncol - 1;], F, cm)
1055-
1056-
# Check if decomposition failed
1040+
F = fact_(A, cm; kws...)
10571041
s = unsafe_load(F.p)
10581042
s.minor < size(A, 1) && throw(Base.LinAlg.ArgumentError("matrix has one or more zero pivots"))
1059-
10601043
return F
10611044
end
10621045

1063-
cholfact{T<:VTypes}(A::Sparse{T}, β::Number) = cholfact(A, convert(T, β))
1064-
cholfact(A::SparseMatrixCSC) = cholfact(Sparse(A))
1065-
cholfact(A::SparseMatrixCSC, β::Number) = cholfact(Sparse(A), β)
1066-
cholfact{Ti}(A::Symmetric{Float64,SparseMatrixCSC{Float64,Ti}}) = cholfact(Sparse(A))
1067-
cholfact{Ti}(A::Symmetric{Float64,SparseMatrixCSC{Float64,Ti}}, β::Number) = cholfact(Sparse(A), β)
1068-
cholfact{Ti}(A::Hermitian{Complex{Float64},SparseMatrixCSC{Complex{Float64},Ti}}) = cholfact(Sparse(A))
1069-
cholfact{Ti}(A::Hermitian{Complex{Float64},SparseMatrixCSC{Complex{Float64},Ti}}, β::Number) = cholfact(Sparse(A), β)
1070-
1071-
ldltfact{T<:VTypes}(A::Sparse{T}, β::Number) = ldltfact(A, convert(T, β))
1072-
ldltfact(A::SparseMatrixCSC) = ldltfact(Sparse(A))
1073-
ldltfact(A::SparseMatrixCSC, β::Number) = ldltfact(Sparse(A), β)
1074-
ldltfact{Ti}(A::Symmetric{Float64,SparseMatrixCSC{Float64,Ti}}) = ldltfact(Sparse(A))
1075-
ldltfact{Ti}(A::Symmetric{Float64,SparseMatrixCSC{Float64,Ti}}, β::Number) = ldltfact(Sparse(A), β)
1076-
ldltfact{Ti}(A::Hermitian{Complex{Float64},SparseMatrixCSC{Complex{Float64},Ti}}) = ldltfact(Sparse(A))
1077-
ldltfact{Ti}(A::Hermitian{Complex{Float64},SparseMatrixCSC{Complex{Float64},Ti}}, β::Number) = ldltfact(Sparse(A), β)
1078-
1079-
function update!{Tv<:VTypes,Ti<:ITypes}(F::Factor{Tv,Ti}, A::Sparse{Tv,Ti})
1080-
cm = common(Ti)
1081-
s = unsafe_load(F.p)
1082-
if Bool(s.is_ll)
1083-
cm[common_final_ll] = reinterpret(UInt8, [one(Cint)]) # Hack! makes it a llt
1046+
for f in (:cholfact, :ldltfact)
1047+
@eval begin
1048+
$f(A::SparseMatrixCSC; kws...) = $f(Sparse(A); kws...)
1049+
$f{Ti}(A::Symmetric{Float64,SparseMatrixCSC{Float64,Ti}}; kws...) = $f(Sparse(A); kws...)
1050+
$f{Ti}(A::Hermitian{Complex{Float64},SparseMatrixCSC{Complex{Float64},Ti}}; kws...) = $f(Sparse(A); kws...)
10841051
end
1085-
factorize!(A, F, cm)
10861052
end
1087-
function update!{Tv<:VTypes,Ti<:ITypes}(F::Factor{Tv,Ti}, A::Sparse{Tv,Ti}, β::Tv)
1053+
1054+
function update!{Tv<:VTypes,Ti<:ITypes}(F::Factor{Tv,Ti}, A::Sparse{Tv,Ti}; shift::Real=0.0)
10881055
cm = common(Ti)
10891056
s = unsafe_load(F.p)
10901057
if Bool(s.is_ll)
10911058
cm[common_final_ll] = reinterpret(UInt8, [one(Cint)]) # Hack! makes it a llt
10921059
end
1093-
factorize_p!(A, β, Ti[0:size(F, 1) - 1;], F, cm)
1060+
factorize_p!(A, shift, F, cm)
10941061
end
1095-
update!{T<:VTypes}(F::Factor{T}, A::SparseMatrixCSC{T}) = update!(F, Sparse(A))
1096-
update!{T<:VTypes}(F::Factor{T}, A::SparseMatrixCSC{T}, β::Number) = update!(F, Sparse(A), convert(T, β))
1062+
update!{T<:VTypes}(F::Factor{T}, A::SparseMatrixCSC{T}; kws...) = update!(F, Sparse(A); kws...)
10971063

10981064
## Solvers
10991065

doc/stdlib/linalg.rst

+12-2
Original file line numberDiff line numberDiff line change
@@ -95,10 +95,15 @@ Linear algebra functions in Julia are largely implemented by calling functions f
9595

9696
Compute the Cholesky factorization of a dense symmetric positive (semi)definite matrix ``A`` and return either a ``Cholesky`` if ``pivot==Val{false}`` or ``CholeskyPivoted`` if ``pivot==Val{true}``. ``LU`` may be ``:L`` for using the lower part or ``:U`` for the upper part. The default is to use ``:U``. The triangular matrix can be obtained from the factorization ``F`` with: ``F[:L]`` and ``F[:U]``. The following functions are available for ``Cholesky`` objects: ``size``, ``\``, ``inv``, ``det``. For ``CholeskyPivoted`` there is also defined a ``rank``. If ``pivot==Val{false}`` a ``PosDefException`` exception is thrown in case the matrix is not positive definite. The argument ``tol`` determines the tolerance for determining the rank. For negative values, the tolerance is the machine precision.
9797

98-
.. function:: cholfact(A) -> CHOLMOD.Factor
98+
.. function:: cholfact(A; shift=0, perm=Int[]) -> CHOLMOD.Factor
9999

100100
Compute the Cholesky factorization of a sparse positive definite matrix ``A``. A fill-reducing permutation is used. The main application of this type is to solve systems of equations with ``\``, but also the methods ``diag``, ``det``, ``logdet`` are defined. The function calls the C library CHOLMOD and many other functions from the library are wrapped but not exported.
101101

102+
Setting optional ``shift`` keyword argument computes the factorization
103+
of ``A+shift*I`` instead of ``A``. If the ``perm`` argument is nonempty,
104+
it should be a permutation of `1:size(A,1)` giving the ordering to use
105+
(instead of CHOLMOD's default AMD ordering).
106+
102107
.. function:: cholfact!(A [,LU=:U [,pivot=Val{false}]][;tol=-1.0]) -> Cholesky
103108

104109
``cholfact!`` is the same as :func:`cholfact`, but saves space by overwriting the input ``A``, instead of creating a copy. ``cholfact!`` can also reuse the symbolic factorization from a different matrix ``F`` with the same structure when used as: ``cholfact!(F::CholmodFactor, A)``.
@@ -107,10 +112,15 @@ Linear algebra functions in Julia are largely implemented by calling functions f
107112

108113
Compute a factorization of a positive definite matrix ``A`` such that ``A=L*Diagonal(d)*L'`` where ``L`` is a unit lower triangular matrix and ``d`` is a vector with non-negative elements.
109114

110-
.. function:: ldltfact(A) -> CHOLMOD.Factor
115+
.. function:: ldltfact(A; shift=0, perm=Int[]) -> CHOLMOD.Factor
111116

112117
Compute the LDLt factorization of a sparse symmetric or Hermitian matrix ``A``. A fill-reducing permutation is used. The main application of this type is to solve systems of equations with ``\``, but also the methods ``diag``, ``det``, ``logdet`` are defined. The function calls the C library CHOLMOD and many other functions from the library are wrapped but not exported.
113118

119+
Setting optional ``shift`` keyword argument computes the factorization
120+
of ``A+shift*I`` instead of ``A``. If the ``perm`` argument is nonempty,
121+
it should be a permutation of `1:size(A,1)` giving the ordering to use
122+
(instead of CHOLMOD's default AMD ordering).
123+
114124
.. function:: qr(A [,pivot=Val{false}][;thin=true]) -> Q, R, [p]
115125

116126
Compute the (pivoted) QR factorization of ``A`` such that either ``A = Q*R`` or ``A[:,p] = Q*R``. Also see ``qrfact``. The default is to compute a thin factorization. Note that ``R`` is not extended with zeros when the full ``Q`` is requested.

test/sparsedir/cholmod.jl

+14-12
Original file line numberDiff line numberDiff line change
@@ -98,13 +98,15 @@ chma = ldltfact(A) # LDL' form
9898
@test unsafe_load(chma.p).is_ll == 0 # check that it is in fact an LDLt
9999
x = chma\B
100100
@test_approx_eq x ones(size(x))
101+
@test nnz(ldltfact(A, perm=1:size(A,1))) > nnz(chma)
101102

102103
chma = cholfact(A) # LL' form
103104
@test CHOLMOD.isvalid(chma)
104105
@test unsafe_load(chma.p).is_ll == 1 # check that it is in fact an LLt
105106
x = chma\B
106107
@test_approx_eq x ones(size(x))
107108
@test nnz(chma) == 489
109+
@test nnz(cholfact(A, perm=1:size(A,1))) > nnz(chma)
108110

109111
#lp_afiro example
110112
afiro = CHOLMOD.Sparse(27, 51,
@@ -360,13 +362,13 @@ for elty in (Float64, Complex{Float64})
360362
# Factor
361363
@test_throws ArgumentError cholfact(A1)
362364
@test_throws Base.LinAlg.PosDefException cholfact(A1 + A1' - 2eigmax(full(A1 + A1'))I)
363-
@test_throws Base.LinAlg.PosDefException cholfact(A1 + A1', -2eigmax(full(A1 + A1')))
365+
@test_throws Base.LinAlg.PosDefException cholfact(A1 + A1', shift=-2eigmax(full(A1 + A1')))
364366
@test_throws ArgumentError ldltfact(A1 + A1' - 2real(A1[1,1])I)
365-
@test_throws ArgumentError ldltfact(A1 + A1', -2real(A1[1,1]))
367+
@test_throws ArgumentError ldltfact(A1 + A1', shift=-2real(A1[1,1]))
366368
@test_throws ArgumentError cholfact(A1)
367-
@test_throws ArgumentError cholfact(A1, 1.0)
369+
@test_throws ArgumentError cholfact(A1, shift=1.0)
368370
@test_throws ArgumentError ldltfact(A1)
369-
@test_throws ArgumentError ldltfact(A1, 1.0)
371+
@test_throws ArgumentError ldltfact(A1, shift=1.0)
370372
F = cholfact(A1pd)
371373
tmp = IOBuffer()
372374
show(tmp, F)
@@ -391,16 +393,16 @@ for elty in (Float64, Complex{Float64})
391393
if elty <: Real
392394
@test CHOLMOD.issym(Sparse(A1pd, 0))
393395
@test CHOLMOD.Sparse(cholfact(Symmetric(A1pd, :L))) == CHOLMOD.Sparse(cholfact(A1pd))
394-
@test CHOLMOD.Sparse(cholfact(Symmetric(A1pd, :L), 2)) == CHOLMOD.Sparse(cholfact(A1pd, 2))
396+
@test CHOLMOD.Sparse(cholfact(Symmetric(A1pd, :L), shift=2)) == CHOLMOD.Sparse(cholfact(A1pd, shift=2))
395397
@test CHOLMOD.Sparse(ldltfact(Symmetric(A1pd, :L))) == CHOLMOD.Sparse(ldltfact(A1pd))
396-
@test CHOLMOD.Sparse(ldltfact(Symmetric(A1pd, :L), 2)) == CHOLMOD.Sparse(ldltfact(A1pd, 2))
398+
@test CHOLMOD.Sparse(ldltfact(Symmetric(A1pd, :L), shift=2)) == CHOLMOD.Sparse(ldltfact(A1pd, shift=2))
397399
else
398400
@test !CHOLMOD.issym(Sparse(A1pd, 0))
399401
@test CHOLMOD.ishermitian(Sparse(A1pd, 0))
400402
@test CHOLMOD.Sparse(cholfact(Hermitian(A1pd, :L))) == CHOLMOD.Sparse(cholfact(A1pd))
401-
@test CHOLMOD.Sparse(cholfact(Hermitian(A1pd, :L), 2)) == CHOLMOD.Sparse(cholfact(A1pd, 2))
403+
@test CHOLMOD.Sparse(cholfact(Hermitian(A1pd, :L), shift=2)) == CHOLMOD.Sparse(cholfact(A1pd, shift=2))
402404
@test CHOLMOD.Sparse(ldltfact(Hermitian(A1pd, :L))) == CHOLMOD.Sparse(ldltfact(A1pd))
403-
@test CHOLMOD.Sparse(ldltfact(Hermitian(A1pd, :L), 2)) == CHOLMOD.Sparse(ldltfact(A1pd, 2))
405+
@test CHOLMOD.Sparse(ldltfact(Hermitian(A1pd, :L), shift=2)) == CHOLMOD.Sparse(ldltfact(A1pd, shift=2))
404406
end
405407

406408

@@ -414,17 +416,17 @@ for elty in (Float64, Complex{Float64})
414416
@test size(F, 3) == 1
415417
@test_throws ArgumentError size(F, 0)
416418

417-
F = cholfact(A1pdSparse, 2)
419+
F = cholfact(A1pdSparse, shift=2)
418420
@test isa(CHOLMOD.Sparse(F), CHOLMOD.Sparse{elty, CHOLMOD.SuiteSparse_long})
419-
@test_approx_eq CHOLMOD.Sparse(CHOLMOD.update!(copy(F), A1pd, 2.0)) CHOLMOD.Sparse(F) # surprisingly, this can cause small ulp size changes so we cannot test exact equality
421+
@test_approx_eq CHOLMOD.Sparse(CHOLMOD.update!(copy(F), A1pd, shift=2.0)) CHOLMOD.Sparse(F) # surprisingly, this can cause small ulp size changes so we cannot test exact equality
420422

421423
F = ldltfact(A1pd)
422424
@test isa(CHOLMOD.Sparse(F), CHOLMOD.Sparse{elty, CHOLMOD.SuiteSparse_long})
423425
@test_approx_eq CHOLMOD.Sparse(CHOLMOD.update!(copy(F), A1pd)) CHOLMOD.Sparse(F) # surprisingly, this can cause small ulp size changes so we cannot test exact equality
424426

425-
F = ldltfact(A1pdSparse, 2)
427+
F = ldltfact(A1pdSparse, shift=2)
426428
@test isa(CHOLMOD.Sparse(F), CHOLMOD.Sparse{elty, CHOLMOD.SuiteSparse_long})
427-
@test_approx_eq CHOLMOD.Sparse(CHOLMOD.update!(copy(F), A1pd, 2.0)) CHOLMOD.Sparse(F) # surprisingly, this can cause small ulp size changes so we cannot test exact equality
429+
@test_approx_eq CHOLMOD.Sparse(CHOLMOD.update!(copy(F), A1pd, shift=2.0)) CHOLMOD.Sparse(F) # surprisingly, this can cause small ulp size changes so we cannot test exact equality
428430

429431
@test isa(CHOLMOD.factor_to_sparse!(F), CHOLMOD.Sparse)
430432
@test_throws CHOLMOD.CHOLMODException CHOLMOD.factor_to_sparse!(F)

0 commit comments

Comments
 (0)