This repository has been archived by the owner on Dec 16, 2024. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 6
/
Copy pathcholmod.jl
1746 lines (1501 loc) · 62 KB
/
cholmod.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# This file is a part of Julia. License is MIT: https://julialang.org/license
# Theoretically CHOLMOD supports both Int32 and Int64 indices on 64-bit.
# However experience suggests that using both in the same session causes memory
# leaks, so we restrict indices to be `SuiteSparse_long`.
# Ref: https://github.com/JuliaLang/julia/issues/12664
# Additionally, only Float64/ComplexF64 are supported in practice.
# Ref: https://github.com/JuliaLang/julia/issues/25986
module CHOLMOD
import Base: (*), convert, copy, eltype, getindex, getproperty, show, size,
IndexStyle, IndexLinear, IndexCartesian, adjoint, axes
using Base: require_one_based_indexing
using LinearAlgebra
import LinearAlgebra: (\),
cholesky, cholesky!, det, diag, ishermitian, isposdef,
issuccess, issymmetric, ldlt, ldlt!, logdet
using SparseArrays
using SparseArrays: getcolptr
import Libdl
export
Dense,
Factor,
Sparse
import SparseArrays: AbstractSparseMatrix, SparseMatrixCSC, indtype, sparse, spzeros, nnz
import ..increment, ..increment!, ..decrement, ..decrement!
using ..LibSuiteSparse
import ..LibSuiteSparse: SuiteSparse_long, TRUE, FALSE
# # itype defines the types of integer used:
# CHOLMOD_INT, # all integer arrays are int
# CHOLMOD_INTLONG, # most are int, some are SuiteSparse_long
# CHOLMOD_LONG, # all integer arrays are SuiteSparse_long
# # dtype defines what the numerical type is (double or float):
# CHOLMOD_DOUBLE, # all numerical values are double
# CHOLMOD_SINGLE, # all numerical values are float
# # xtype defines the kind of numerical values used:
# CHOLMOD_PATTERN, # pattern only, no numerical values
# CHOLMOD_REAL, # a real matrix
# CHOLMOD_COMPLEX, # a complex matrix (ANSI C99 compatible)
# CHOLMOD_ZOMPLEX, # a complex matrix (MATLAB compatible)
# # Scaling modes, selected by the scale input parameter:
# CHOLMOD_SCALAR, # A = s*A
# CHOLMOD_ROW, # A = diag(s)*A
# CHOLMOD_COL, # A = A*diag(s)
# CHOLMOD_SYM, # A = diag(s)*A*diag(s)
# # Types of systems to solve
# CHOLMOD_A, # solve Ax=b
# CHOLMOD_LDLt, # solve LDL'x=b
# CHOLMOD_LD, # solve LDx=b
# CHOLMOD_DLt, # solve DL'x=b
# CHOLMOD_L, # solve Lx=b
# CHOLMOD_Lt, # solve L'x=b
# CHOLMOD_D, # solve Dx=b
# CHOLMOD_P, # permute x=Px
# CHOLMOD_Pt, # permute x=P'x
# # Symmetry types
# CHOLMOD_MM_RECTANGULAR,
# CHOLMOD_MM_UNSYMMETRIC,
# CHOLMOD_MM_SYMMETRIC,
# CHOLMOD_MM_HERMITIAN,
# CHOLMOD_MM_SKEW_SYMMETRIC,
# CHOLMOD_MM_SYMMETRIC_POSDIAG,
# CHOLMOD_MM_HERMITIAN_POSDIAG
dtyp(::Type{Float32}) = CHOLMOD_SINGLE
dtyp(::Type{Float64}) = CHOLMOD_DOUBLE
dtyp(::Type{ComplexF32}) = CHOLMOD_SINGLE
dtyp(::Type{ComplexF64}) = CHOLMOD_DOUBLE
xtyp(::Type{Float32}) = CHOLMOD_REAL
xtyp(::Type{Float64}) = CHOLMOD_REAL
xtyp(::Type{ComplexF32}) = CHOLMOD_COMPLEX
xtyp(::Type{ComplexF64}) = CHOLMOD_COMPLEX
# check the size of SuiteSparse_long
if sizeof(SuiteSparse_long) == 4
const IndexTypes = (:Int32,)
const ITypes = Union{Int32}
else
const IndexTypes = (:Int32, :Int64)
const ITypes = Union{Int32, Int64}
end
ityp(::Type{SuiteSparse_long}) = CHOLMOD_LONG
const VTypes = Union{ComplexF64, Float64}
const VRealTypes = Union{Float64}
# overload field access methods
function Base.getproperty(x::cholmod_sparse, f::Symbol)
f === :p && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :i && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :nz && return Ptr{SuiteSparse_long}(getfield(x, f))
return getfield(x, f)
end
function Base.getproperty(x::cholmod_factor, f::Symbol)
f === :Perm && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :ColCount && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :IPerm && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :p && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :i && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :nz && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :next && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :prev && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :super && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :pi && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :px && return Ptr{SuiteSparse_long}(getfield(x, f))
f === :s && return Ptr{SuiteSparse_long}(getfield(x, f))
return getfield(x, f)
end
# exception
struct CHOLMODException <: Exception
msg::String
end
function error_handler(status::Cint, file::Cstring, line::Cint, message::Cstring)::Cvoid
status < 0 && throw(CHOLMODException(unsafe_string(message)))
nothing
end
const CHOLMOD_MIN_VERSION = v"2.1.1"
# Set a `common` field, execute some code and then safely reset the field to
# its initial value
macro cholmod_param(kwarg, code)
@assert kwarg.head == :(=)
param = kwarg.args[1]
value = kwarg.args[2]
common_param = # Read `common.param`
Expr(:., :(COMMONS[Threads.threadid()][]), QuoteNode(param))
return quote
default_value = $common_param
try
$common_param = $(esc(value))
$(esc(code))
finally
$common_param = default_value
end
end
end
const COMMONS = Vector{Ref{cholmod_common}}()
const BUILD_VERSION = VersionNumber(CHOLMOD_MAIN_VERSION, CHOLMOD_SUB_VERSION, CHOLMOD_SUBSUB_VERSION)
function __init__()
try
### Check if the linked library is compatible with the Julia code
if Libdl.dlsym_e(Libdl.dlopen("libcholmod"), :cholmod_version) != C_NULL
current_version_array = Vector{Cint}(undef, 3)
cholmod_version(current_version_array)
current_version = VersionNumber(current_version_array...)
else # CHOLMOD < 2.1.1 does not include cholmod_version()
current_version = v"0.0.0"
end
if current_version < CHOLMOD_MIN_VERSION
@warn """
CHOLMOD version incompatibility
Julia was compiled with CHOLMOD version $BUILD_VERSION. It is
currently linked with a version older than
$(CHOLMOD_MIN_VERSION). This might cause Julia to
terminate when working with sparse matrix factorizations,
e.g. solving systems of equations with \\.
It is recommended that you use Julia with a recent version
of CHOLMOD, or download the generic binaries
from www.julialang.org, which ship with the correct
versions of all dependencies.
"""
elseif BUILD_VERSION.major != current_version.major
@warn """
CHOLMOD version incompatibility
Julia was compiled with CHOLMOD version $BUILD_VERSION. It is
currently linked with version $current_version.
This might cause Julia to terminate when working with
sparse matrix factorizations, e.g. solving systems of
equations with \\.
It is recommended that you use Julia with the same major
version of CHOLMOD as the one used during the build, or
download the generic binaries from www.julialang.org,
which ship with the correct versions of all dependencies.
"""
end
intsize = sizeof(SuiteSparse_long)
if intsize != 4length(IndexTypes)
@error """
CHOLMOD integer size incompatibility
Julia was compiled with a version of CHOLMOD that
supported $(32length(IndexTypes)) bit integers. It is
currently linked with version that supports $(8intsize)
integers. This might cause Julia to terminate when
working with sparse matrix factorizations, e.g. solving
systems of equations with \\.
This problem can be fixed by modifying the Julia build
configuration or by downloading the OS X or generic
Linux binary from www.julialang.org, which include
the correct versions of all dependencies.
"""
end
### Initiate CHOLMOD
### common controls the type of factorization and keeps pointers
### to temporary memory. We need to manage a copy for each thread.
nt = Threads.nthreads()
resize!(COMMONS, nt)
errorhandler = @cfunction(error_handler, Cvoid, (Cint, Cstring, Cint, Cstring))
for i in 1:nt
COMMONS[i] = Ref(cholmod_common())
result = cholmod_l_start(COMMONS[i])
@assert result == TRUE "failed to run `cholmod_l_start`!"
COMMONS[i][].print = 0 # no printing from CHOLMOD by default
COMMONS[i][].error_handler = errorhandler
end
# Register gc tracked allocator if CHOLMOD is new enough
if current_version >= v"3.0.0"
cnfg = cglobal((:SuiteSparse_config, :libsuitesparseconfig), Ptr{Cvoid})
unsafe_store!(cnfg, cglobal(:jl_malloc, Ptr{Cvoid}), 1)
unsafe_store!(cnfg, cglobal(:jl_calloc, Ptr{Cvoid}), 2)
unsafe_store!(cnfg, cglobal(:jl_realloc, Ptr{Cvoid}), 3)
unsafe_store!(cnfg, cglobal(:jl_free, Ptr{Cvoid}), 4)
end
catch ex
@error "Error during initialization of module CHOLMOD" exception=ex,catch_backtrace()
end
end
####################
# Type definitions #
####################
# The three core data types for CHOLMOD: Dense, Sparse and Factor.
# CHOLMOD manages the memory, so the Julia versions only wrap a
# pointer to a struct. Therefore finalizers should be registered each
# time a pointer is returned from CHOLMOD.
mutable struct Dense{Tv<:VTypes} <: DenseMatrix{Tv}
ptr::Ptr{cholmod_dense}
function Dense{Tv}(ptr::Ptr{cholmod_dense}) where Tv<:VTypes
if ptr == C_NULL
throw(ArgumentError("dense matrix construction failed for " *
"unknown reasons. Please submit a bug report."))
end
s = unsafe_load(ptr)
if s.xtype != xtyp(Tv)
free!(ptr)
throw(CHOLMODException("xtype=$(s.xtype) not supported"))
elseif s.dtype != dtyp(Tv)
free!(ptr)
throw(CHOLMODException("dtype=$(s.dtype) not supported"))
end
obj = new(ptr)
finalizer(free!, obj)
return obj
end
end
mutable struct Sparse{Tv<:VTypes} <: AbstractSparseMatrix{Tv,SuiteSparse_long}
ptr::Ptr{cholmod_sparse}
function Sparse{Tv}(ptr::Ptr{cholmod_sparse}) where Tv<:VTypes
if ptr == C_NULL
throw(ArgumentError("sparse matrix construction failed for " *
"unknown reasons. Please submit a bug report."))
end
s = unsafe_load(ptr)
if s.itype != ityp(SuiteSparse_long)
free!(ptr)
throw(CHOLMODException("itype=$(s.itype) not supported"))
elseif s.xtype != xtyp(Tv)
free!(ptr)
throw(CHOLMODException("xtype=$(s.xtype) not supported"))
elseif s.dtype != dtyp(Tv)
free!(ptr)
throw(CHOLMODException("dtype=$(s.dtype) not supported"))
end
A = new(ptr)
finalizer(free!, A)
return A
end
end
# Useful when reading in files, but not type stable
function Sparse(p::Ptr{cholmod_sparse})
if p == C_NULL
throw(ArgumentError("sparse matrix construction failed for " *
"unknown reasons. Please submit a bug report."))
end
s = unsafe_load(p)
Tv = s.xtype == CHOLMOD_REAL ? Float64 : ComplexF64
Sparse{Tv}(p)
end
mutable struct Factor{Tv<:VTypes} <: Factorization{Tv}
ptr::Ptr{cholmod_factor}
function Factor{Tv}(ptr::Ptr{cholmod_factor}, register_finalizer = true) where Tv
if ptr == C_NULL
throw(ArgumentError("factorization construction failed for " *
"unknown reasons. Please submit a bug report."))
end
s = unsafe_load(ptr)
if s.itype != ityp(SuiteSparse_long)
free!(ptr)
throw(CHOLMODException("itype=$(s.itype) not supported"))
elseif s.xtype != xtyp(Tv) && s.xtype != CHOLMOD_PATTERN
free!(ptr)
throw(CHOLMODException("xtype=$(s.xtype) not supported"))
elseif s.dtype != dtyp(Tv)
free!(ptr)
throw(CHOLMODException("dtype=$(s.dtype) not supported"))
end
F = new(ptr)
if register_finalizer
finalizer(free!, F)
end
return F
end
end
Base.adjoint(F::Factor) = Adjoint(F)
Base.transpose(F::Factor) = Transpose(F)
const SuiteSparseStruct = Union{cholmod_dense, cholmod_sparse, cholmod_factor}
# All pointer loads should be checked to make sure that SuiteSparse is not called with
# a C_NULL pointer which could cause a segfault. Pointers are set to null
# when serialized so this can happen when multiple processes are in use.
function Base.unsafe_convert(::Type{Ptr{T}}, x::Union{Dense,Sparse,Factor}) where T<:SuiteSparseStruct
xp = getfield(x, :ptr)
if xp == C_NULL
throw(ArgumentError("pointer to the $T object is null. This can " *
"happen if the object has been serialized."))
else
return xp
end
end
Base.pointer(x::Dense{Tv}) where {Tv} = Base.unsafe_convert(Ptr{cholmod_dense}, x)
Base.pointer(x::Sparse{Tv}) where {Tv} = Base.unsafe_convert(Ptr{cholmod_sparse}, x)
Base.pointer(x::Factor{Tv}) where {Tv} = Base.unsafe_convert(Ptr{cholmod_factor}, x)
# FactorComponent, for encoding particular factors from a factorization
mutable struct FactorComponent{Tv,S} <: AbstractMatrix{Tv}
F::Factor{Tv}
function FactorComponent{Tv,S}(F::Factor{Tv}) where {Tv,S}
s = unsafe_load(pointer(F))
if s.is_ll != 0
if !(S === :L || S === :U || S === :PtL || S === :UP)
throw(CHOLMODException(string(S, " not supported for sparse ",
"LLt matrices; try :L, :U, :PtL, or :UP")))
end
elseif !(S === :L || S === :U || S === :PtL || S === :UP ||
S === :D || S === :LD || S === :DU || S === :PtLD || S === :DUP)
throw(CHOLMODException(string(S, " not supported for sparse LDLt ",
"matrices; try :L, :U, :PtL, :UP, :D, :LD, :DU, :PtLD, or :DUP")))
end
new(F)
end
end
function FactorComponent(F::Factor{Tv}, sym::Symbol) where Tv
FactorComponent{Tv,sym}(F)
end
Factor(FC::FactorComponent) = FC.F
#################
# Thin wrappers #
#################
# Dense wrappers
function allocate_dense(m::Integer, n::Integer, d::Integer, ::Type{Tv}) where {Tv<:VTypes}
Dense{Tv}(cholmod_l_allocate_dense(m, n, d, xtyp(Tv), COMMONS[Threads.threadid()]))
end
function free!(p::Ptr{cholmod_dense})
cholmod_l_free_dense(Ref(p), COMMONS[Threads.threadid()]) == TRUE
end
function zeros(m::Integer, n::Integer, ::Type{Tv}) where Tv<:VTypes
Dense{Tv}(cholmod_l_zeros(m, n, xtyp(Tv), COMMONS[Threads.threadid()]))
end
zeros(m::Integer, n::Integer) = zeros(m, n, Float64)
function ones(m::Integer, n::Integer, ::Type{Tv}) where Tv<:VTypes
Dense{Tv}(cholmod_l_ones(m, n, xtyp(Tv), COMMONS[Threads.threadid()]))
end
ones(m::Integer, n::Integer) = ones(m, n, Float64)
function eye(m::Integer, n::Integer, ::Type{Tv}) where Tv<:VTypes
Dense{Tv}(cholmod_l_eye(m, n, xtyp(Tv), COMMONS[Threads.threadid()]))
end
eye(m::Integer, n::Integer) = eye(m, n, Float64)
eye(n::Integer) = eye(n, n, Float64)
function copy(A::Dense{Tv}) where Tv<:VTypes
Dense{Tv}(cholmod_l_copy_dense(A, COMMONS[Threads.threadid()]))
end
function sort!(S::Sparse{Tv}) where Tv<:VTypes
cholmod_l_sort(S, COMMONS[Threads.threadid()])
return S
end
function norm_dense(D::Dense{Tv}, p::Integer) where Tv<:VTypes
s = unsafe_load(pointer(D))
if p == 2
if s.ncol > 1
throw(ArgumentError("2 norm only supported when matrix has one column"))
end
elseif p != 0 && p != 1
throw(ArgumentError("second argument must be either 0 (Inf norm), 1, or 2"))
end
cholmod_l_norm_dense(D, p, COMMONS[Threads.threadid()])
end
function check_dense(A::Dense{Tv}) where Tv<:VTypes
cholmod_l_check_dense(pointer(A), COMMONS[Threads.threadid()]) != 0
end
# Non-Dense wrappers
function allocate_sparse(nrow::Integer, ncol::Integer, nzmax::Integer,
sorted::Bool, packed::Bool, stype::Integer, ::Type{Tv}) where {Tv<:VTypes}
Sparse{Tv}(cholmod_l_allocate_sparse(nrow, ncol, nzmax, sorted, packed, stype,
xtyp(Tv), COMMONS[Threads.threadid()]))
end
function free!(ptr::Ptr{cholmod_sparse})
cholmod_l_free_sparse(Ref(ptr), COMMONS[Threads.threadid()]) == TRUE
end
function free!(ptr::Ptr{cholmod_factor})
# Warning! Important that finalizer doesn't modify the global Common struct.
cholmod_l_free_factor(Ref(ptr), COMMONS[Threads.threadid()]) == TRUE
end
function aat(A::Sparse{Tv}, fset::Vector{SuiteSparse_long}, mode::Integer) where Tv<:VRealTypes
Sparse{Tv}(cholmod_l_aat(A, fset, length(fset), mode, COMMONS[Threads.threadid()]))
end
function sparse_to_dense(A::Sparse{Tv}) where Tv<:VTypes
Dense{Tv}(cholmod_l_sparse_to_dense(A, COMMONS[Threads.threadid()]))
end
function dense_to_sparse(D::Dense{Tv}, ::Type{SuiteSparse_long}) where Tv<:VTypes
Sparse{Tv}(cholmod_l_dense_to_sparse(D, true, COMMONS[Threads.threadid()]))
end
function factor_to_sparse!(F::Factor{Tv}) where Tv<:VTypes
ss = unsafe_load(pointer(F))
ss.xtype == CHOLMOD_PATTERN && throw(CHOLMODException("only numeric factors are supported"))
Sparse{Tv}(cholmod_l_factor_to_sparse(F, COMMONS[Threads.threadid()]))
end
function change_factor!(F::Factor{Tv}, to_ll::Bool, to_super::Bool, to_packed::Bool,
to_monotonic::Bool) where Tv<:VTypes
cholmod_l_change_factor(xtyp(Tv), to_ll, to_super, to_packed, to_monotonic, F, COMMONS[Threads.threadid()]) == TRUE
end
function check_sparse(A::Sparse{Tv}) where Tv<:VTypes
cholmod_l_check_sparse(A, COMMONS[Threads.threadid()]) != 0
end
function check_factor(F::Factor{Tv}) where Tv<:VTypes
cholmod_l_check_factor(F, COMMONS[Threads.threadid()]) != 0
end
nnz(A::Sparse{<:VTypes}) = cholmod_l_nnz(A, COMMONS[Threads.threadid()])
function speye(m::Integer, n::Integer, ::Type{Tv}) where Tv<:VTypes
Sparse{Tv}(cholmod_l_speye(m, n, xtyp(Tv), COMMONS[Threads.threadid()]))
end
function spzeros(m::Integer, n::Integer, nzmax::Integer, ::Type{Tv}) where Tv<:VTypes
Sparse{Tv}(cholmod_l_spzeros(m, n, nzmax, xtyp(Tv), COMMONS[Threads.threadid()]))
end
function transpose_(A::Sparse{Tv}, values::Integer) where Tv<:VTypes
Sparse{Tv}(cholmod_l_transpose(A, values, COMMONS[Threads.threadid()]))
end
function copy(F::Factor{Tv}) where Tv<:VTypes
Factor{Tv}(cholmod_l_copy_factor(F, COMMONS[Threads.threadid()]))
end
function copy(A::Sparse{Tv}) where Tv<:VTypes
Sparse{Tv}(cholmod_l_copy_sparse(A, COMMONS[Threads.threadid()]))
end
function copy(A::Sparse{Tv}, stype::Integer, mode::Integer) where Tv<:VRealTypes
Sparse{Tv}(cholmod_l_copy(A, stype, mode, COMMONS[Threads.threadid()]))
end
function print_sparse(A::Sparse{Tv}, name::String) where Tv<:VTypes
isascii(name) || error("non-ASCII name: $name")
@cholmod_param print = 3 begin
cholmod_l_print_sparse(A, name, COMMONS[Threads.threadid()])
end
nothing
end
function print_factor(F::Factor{Tv}, name::String) where Tv<:VTypes
@cholmod_param print = 3 begin
cholmod_l_print_factor(F, name, COMMONS[Threads.threadid()])
end
nothing
end
function ssmult(A::Sparse{Tv}, B::Sparse{Tv}, stype::Integer,
values::Bool, sorted::Bool) where Tv<:VRealTypes
lA = unsafe_load(pointer(A))
lB = unsafe_load(pointer(B))
if lA.ncol != lB.nrow
throw(DimensionMismatch("inner matrix dimensions do not fit"))
end
Sparse{Tv}(cholmod_l_ssmult(A, B, stype, values, sorted, COMMONS[Threads.threadid()]))
end
function norm_sparse(A::Sparse{Tv}, norm::Integer) where Tv<:VTypes
if norm != 0 && norm != 1
throw(ArgumentError("norm argument must be either 0 or 1"))
end
cholmod_l_norm_sparse(A, norm, COMMONS[Threads.threadid()])
end
function horzcat(A::Sparse{Tv}, B::Sparse{Tv}, values::Bool) where Tv<:VRealTypes
Sparse{Tv}(cholmod_l_horzcat(A, B, values, COMMONS[Threads.threadid()]))
end
function scale!(S::Dense{Tv}, scale::Integer, A::Sparse{Tv}) where Tv<:VRealTypes
sS = unsafe_load(pointer(S))
sA = unsafe_load(pointer(A))
if sS.ncol != 1 && sS.nrow != 1
throw(DimensionMismatch("first argument must be a vector"))
end
if scale == CHOLMOD_SCALAR && sS.nrow != 1
throw(DimensionMismatch("scaling argument must have length one"))
elseif scale == CHOLMOD_ROW && sS.nrow*sS.ncol != sA.nrow
throw(DimensionMismatch("scaling vector has length $(sS.nrow*sS.ncol), " *
"but matrix has $(sA.nrow) rows."))
elseif scale == CHOLMOD_COL && sS.nrow*sS.ncol != sA.ncol
throw(DimensionMismatch("scaling vector has length $(sS.nrow*sS.ncol), " *
"but matrix has $(sA.ncol) columns"))
elseif scale == CHOLMOD_SYM
if sA.nrow != sA.ncol
throw(DimensionMismatch("matrix must be square"))
elseif sS.nrow*sS.ncol != sA.nrow
throw(DimensionMismatch("scaling vector has length $(sS.nrow*sS.ncol), " *
"but matrix has $(sA.ncol) columns and rows"))
end
end
sA = unsafe_load(pointer(A))
cholmod_l_scale(S, scale, A, COMMONS[Threads.threadid()])
A
end
function sdmult!(A::Sparse{Tv}, transpose::Bool,
α::Number, β::Number, X::Dense{Tv}, Y::Dense{Tv}) where Tv<:VTypes
m, n = size(A)
nc = transpose ? m : n
nr = transpose ? n : m
if nc != size(X, 1)
throw(DimensionMismatch("incompatible dimensions, $nc and $(size(X,1))"))
end
cholmod_l_sdmult(A, transpose, Ref(α), Ref(β), X, Y, COMMONS[Threads.threadid()])
Y
end
function vertcat(A::Sparse{Tv}, B::Sparse{Tv}, values::Bool) where Tv<:VRealTypes
Sparse{Tv}(cholmod_l_vertcat(A, B, values, COMMONS[Threads.threadid()]))
end
function symmetry(A::Sparse{Tv}, option::Integer) where Tv<:VTypes
xmatched = Ref{SuiteSparse_long}()
pmatched = Ref{SuiteSparse_long}()
nzoffdiag = Ref{SuiteSparse_long}()
nzdiag = Ref{SuiteSparse_long}()
rv = cholmod_l_symmetry(A, option, xmatched, pmatched,
nzoffdiag, nzdiag, COMMONS[Threads.threadid()])
rv, xmatched[], pmatched[], nzoffdiag[], nzdiag[]
end
# For analyze, analyze_p, and factorize_p!, the Common argument must be
# supplied in order to control if the factorization is LLt or LDLt
function analyze(A::Sparse{Tv}) where Tv<:VTypes
Factor{Tv}(cholmod_l_analyze(A, COMMONS[Threads.threadid()]))
end
function analyze_p(A::Sparse{Tv}, perm::Vector{SuiteSparse_long}) where Tv<:VTypes
length(perm) != size(A,1) && throw(BoundsError())
Factor{Tv}(cholmod_l_analyze_p(A, perm, C_NULL, 0, COMMONS[Threads.threadid()]))
end
function factorize!(A::Sparse{Tv}, F::Factor{Tv}) where Tv<:VTypes
cholmod_l_factorize(A, F, COMMONS[Threads.threadid()])
F
end
function factorize_p!(A::Sparse{Tv}, β::Real, F::Factor{Tv}) where Tv<:VTypes
# note that β is passed as a complex number (double beta[2]),
# but the CHOLMOD manual says that only beta[0] (real part) is used
cholmod_l_factorize_p(A, Ref{Cdouble}(β), C_NULL, 0, F, COMMONS[Threads.threadid()])
F
end
function solve(sys::Integer, F::Factor{Tv}, B::Dense{Tv}) where Tv<:VTypes
if size(F,1) != size(B,1)
throw(DimensionMismatch("LHS and RHS should have the same number of rows. " *
"LHS has $(size(F,1)) rows, but RHS has $(size(B,1)) rows."))
end
if !issuccess(F)
s = unsafe_load(pointer(F))
if s.is_ll == 1
throw(LinearAlgebra.PosDefException(s.minor))
else
throw(LinearAlgebra.ZeroPivotException(s.minor))
end
end
Dense{Tv}(cholmod_l_solve(sys, F, B, COMMONS[Threads.threadid()]))
end
function spsolve(sys::Integer, F::Factor{Tv}, B::Sparse{Tv}) where Tv<:VTypes
if size(F,1) != size(B,1)
throw(DimensionMismatch("LHS and RHS should have the same number of rows. " *
"LHS has $(size(F,1)) rows, but RHS has $(size(B,1)) rows."))
end
Sparse{Tv}(cholmod_l_spsolve(sys, F, B, COMMONS[Threads.threadid()]))
end
# Autodetects the types
function read_sparse(file::Libc.FILE, ::Type{SuiteSparse_long})
Sparse(cholmod_l_read_sparse(file.ptr, COMMONS[Threads.threadid()]))
end
function read_sparse(file::IO, T)
cfile = Libc.FILE(file)
try return read_sparse(cfile, T)
finally close(cfile)
end
end
function get_perm(F::Factor)
s = unsafe_load(pointer(F))
p = unsafe_wrap(Array, s.Perm, s.n, own = false)
p .+ 1
end
get_perm(FC::FactorComponent) = get_perm(Factor(FC))
#########################
# High level interfaces #
#########################
# Conversion/construction
function Dense{T}(A::StridedVecOrMat) where T<:VTypes
d = allocate_dense(size(A, 1), size(A, 2), stride(A, 2), T)
GC.@preserve d begin
s = unsafe_load(pointer(d))
for (i, c) in enumerate(eachindex(A))
unsafe_store!(Ptr{T}(s.x), A[c], i)
end
end
d
end
function Dense{T}(A::Union{Adjoint{<:Any, <:StridedVecOrMat}, Transpose{<:Any, <:StridedVecOrMat}}) where T<:VTypes
d = allocate_dense(size(A, 1), size(A, 2), size(A, 1), T)
GC.@preserve d begin
s = unsafe_load(pointer(d))
for (i, c) in enumerate(eachindex(A))
unsafe_store!(Ptr{T}(s.x), A[c], i)
end
end
d
end
function Dense(A::Union{StridedVecOrMat, Adjoint{<:Any, <:StridedVecOrMat}, Transpose{<:Any, <:StridedVecOrMat}})
T = promote_type(eltype(A), Float64)
return Dense{T}(A)
end
Dense(A::Sparse) = sparse_to_dense(A)
# This constructior assumes zero based colptr and rowval
function Sparse(m::Integer, n::Integer,
colptr0::Vector{SuiteSparse_long}, rowval0::Vector{SuiteSparse_long},
nzval::Vector{Tv}, stype) where Tv<:VTypes
# checks
## length of input
if length(colptr0) <= n
throw(ArgumentError("length of colptr0 must be at least n + 1 = $(n + 1) but was $(length(colptr0))"))
end
if colptr0[n + 1] > length(rowval0)
throw(ArgumentError("length of rowval0 is $(length(rowval0)) but value of colptr0 requires length to be at least $(colptr0[n + 1])"))
end
if colptr0[n + 1] > length(nzval)
throw(ArgumentError("length of nzval is $(length(nzval)) but value of colptr0 requires length to be at least $(colptr0[n + 1])"))
end
## columns are sorted
iss = true
for i = 2:length(colptr0)
if !issorted(view(rowval0, colptr0[i - 1] + 1:colptr0[i]))
iss = false
break
end
end
o = allocate_sparse(m, n, colptr0[n + 1], iss, true, stype, Tv)
s = unsafe_load(pointer(o))
unsafe_copyto!(s.p, pointer(colptr0), n + 1)
unsafe_copyto!(s.i, pointer(rowval0), colptr0[n + 1])
unsafe_copyto!(Ptr{Tv}(s.x), pointer(nzval) , colptr0[n + 1])
check_sparse(o)
return o
end
function Sparse(m::Integer, n::Integer,
colptr0::Vector{SuiteSparse_long},
rowval0::Vector{SuiteSparse_long},
nzval::Vector{<:VTypes})
o = Sparse(m, n, colptr0, rowval0, nzval, 0)
# sort indices
sort!(o)
# check if array is symmetric and change stype if it is
if ishermitian(o)
change_stype!(o, -1)
end
o
end
function Sparse{Tv}(A::SparseMatrixCSC, stype::Integer) where Tv<:VTypes
## Check length of input. This should never fail but see #20024
if length(getcolptr(A)) <= size(A, 2)
throw(ArgumentError("length of colptr must be at least size(A,2) + 1 = $(size(A, 2) + 1) but was $(length(getcolptr(A)))"))
end
if nnz(A) > length(rowvals(A))
throw(ArgumentError("length of rowval is $(length(rowvals(A))) but value of colptr requires length to be at least $(nnz(A))"))
end
if nnz(A) > length(nonzeros(A))
throw(ArgumentError("length of nzval is $(length(nonzeros(A))) but value of colptr requires length to be at least $(nnz(A))"))
end
o = allocate_sparse(size(A, 1), size(A, 2), nnz(A), true, true, stype, Tv)
s = unsafe_load(pointer(o))
for i = 1:(size(A, 2) + 1)
unsafe_store!(s.p, getcolptr(A)[i] - 1, i)
end
for i = 1:nnz(A)
unsafe_store!(s.i, rowvals(A)[i] - 1, i)
end
if Tv <: Complex && stype != 0
# Need to remove any non real elements in the diagonal because, in contrast to
# BLAS/LAPACK these are not ignored by CHOLMOD. If even tiny imaginary parts are
# present CHOLMOD will fail with a non-positive definite/zero pivot error.
for j = 1:size(A, 2)
for ip = getcolptr(A)[j]:getcolptr(A)[j + 1] - 1
v = nonzeros(A)[ip]
unsafe_store!(Ptr{Tv}(s.x), rowvals(A)[ip] == j ? Complex(real(v)) : v, ip)
end
end
elseif Tv == eltype(nonzeros(A))
unsafe_copyto!(Ptr{Tv}(s.x), pointer(nonzeros(A)), nnz(A))
else
for i = 1:nnz(A)
unsafe_store!(Ptr{Tv}(s.x), nonzeros(A)[i], i)
end
end
check_sparse(o)
return o
end
# handle promotion
function Sparse(A::SparseMatrixCSC{Tv,SuiteSparse_long}, stype::Integer) where {Tv}
T = promote_type(Tv, Float64)
return Sparse{T}(A, stype)
end
# convert SparseVectors into CHOLMOD Sparse types through a mx1 CSC matrix
Sparse(A::SparseVector) = Sparse(SparseMatrixCSC(A))
function Sparse(A::SparseMatrixCSC)
o = Sparse(A, 0)
# check if array is symmetric and change stype if it is
if ishermitian(o)
change_stype!(o, -1)
end
o
end
Sparse(A::Symmetric{Tv, SparseMatrixCSC{Tv,Ti}}) where {Tv<:Real, Ti} =
Sparse(A.data, A.uplo == 'L' ? -1 : 1)
Sparse(A::Hermitian{Tv,SparseMatrixCSC{Tv,Ti}}) where {Tv, Ti} =
Sparse(A.data, A.uplo == 'L' ? -1 : 1)
Sparse(A::Dense) = dense_to_sparse(A, SuiteSparse_long)
Sparse(L::Factor) = factor_to_sparse!(copy(L))
function Sparse(filename::String)
open(filename) do f
return read_sparse(f, SuiteSparse_long)
end
end
## conversion back to base Julia types
function Matrix{T}(D::Dense{T}) where T
s = unsafe_load(pointer(D))
a = Matrix{T}(undef, s.nrow, s.ncol)
copyto!(a, D)
end
Base.copyto!(dest::Base.PermutedDimsArrays.PermutedDimsArray, src::Dense) = _copy!(dest, src) # ambig
Base.copyto!(dest::Dense{T}, D::Dense{T}) where {T<:VTypes} = _copy!(dest, D)
Base.copyto!(dest::AbstractArray{T}, D::Dense{T}) where {T<:VTypes} = _copy!(dest, D)
Base.copyto!(dest::AbstractArray{T,2}, D::Dense{T}) where {T<:VTypes} = _copy!(dest, D)
Base.copyto!(dest::AbstractArray, D::Dense) = _copy!(dest, D)
function _copy!(dest::AbstractArray, D::Dense{T}) where {T<:VTypes}
require_one_based_indexing(dest)
s = unsafe_load(pointer(D))
n = s.nrow*s.ncol
n <= length(dest) || throw(BoundsError(dest, n))
if s.d == s.nrow && isa(dest, Array)
unsafe_copyto!(pointer(dest), Ptr{T}(s.x), s.d*s.ncol)
else
k = 0
for j = 1:s.ncol
for i = 1:s.nrow
dest[k+=1] = unsafe_load(Ptr{T}(s.x), i + (j - 1)*s.d)
end
end
end
dest
end
Matrix(D::Dense{T}) where {T} = Matrix{T}(D)
function Vector{T}(D::Dense{T}) where T
if size(D, 2) > 1
throw(DimensionMismatch("input must be a vector but had $(size(D, 2)) columns"))
end
copyto!(Vector{T}(undef, size(D, 1)), D)
end
Vector(D::Dense{T}) where {T} = Vector{T}(D)
function _extract_args(s, ::Type{T}) where {T<:VTypes}
return (s.nrow, s.ncol, increment(unsafe_wrap(Array, s.p, (s.ncol + 1,), own = false)),
increment(unsafe_wrap(Array, s.i, (s.nzmax,), own = false)),
copy(unsafe_wrap(Array, Ptr{T}(s.x), (s.nzmax,), own = false)))
end
# Trim extra elements in rowval and nzval left around sometimes by CHOLMOD rutines
function _trim_nz_builder!(m, n, colptr, rowval, nzval)
l = colptr[end] - 1
resize!(rowval, l)
resize!(nzval, l)
SparseMatrixCSC(m, n, colptr, rowval, nzval)
end
function SparseMatrixCSC{Tv,SuiteSparse_long}(A::Sparse{Tv}) where Tv
s = unsafe_load(pointer(A))
if s.stype != 0
throw(ArgumentError("matrix has stype != 0. Convert to matrix " *
"with stype == 0 before converting to SparseMatrixCSC"))
end
args = _extract_args(s, Tv)
s.sorted == 0 && _sort_buffers!(args...);
return _trim_nz_builder!(args...)
end
function Symmetric{Float64,SparseMatrixCSC{Float64,SuiteSparse_long}}(A::Sparse{Float64})
s = unsafe_load(pointer(A))
issymmetric(A) || throw(ArgumentError("matrix is not symmetric"))
args = _extract_args(s, Float64)
s.sorted == 0 && _sort_buffers!(args...)
Symmetric(_trim_nz_builder!(args...), s.stype > 0 ? :U : :L)
end
convert(T::Type{Symmetric{Float64,SparseMatrixCSC{Float64,SuiteSparse_long}}}, A::Sparse{Float64}) = T(A)
function Hermitian{Tv,SparseMatrixCSC{Tv,SuiteSparse_long}}(A::Sparse{Tv}) where Tv<:VTypes
s = unsafe_load(pointer(A))
ishermitian(A) || throw(ArgumentError("matrix is not Hermitian"))
args = _extract_args(s, Tv)
s.sorted == 0 && _sort_buffers!(args...)
Hermitian(_trim_nz_builder!(args...), s.stype > 0 ? :U : :L)
end
convert(T::Type{Hermitian{Tv,SparseMatrixCSC{Tv,SuiteSparse_long}}}, A::Sparse{Tv}) where {Tv<:VTypes} = T(A)
function sparse(A::Sparse{Float64}) # Notice! Cannot be type stable because of stype
s = unsafe_load(pointer(A))
if s.stype == 0
return SparseMatrixCSC{Float64,SuiteSparse_long}(A)
end
Symmetric{Float64,SparseMatrixCSC{Float64,SuiteSparse_long}}(A)
end
function sparse(A::Sparse{ComplexF64}) # Notice! Cannot be type stable because of stype
s = unsafe_load(pointer(A))
if s.stype == 0
return SparseMatrixCSC{ComplexF64,SuiteSparse_long}(A)
end
Hermitian{ComplexF64,SparseMatrixCSC{ComplexF64,SuiteSparse_long}}(A)
end
function sparse(F::Factor)
s = unsafe_load(pointer(F))
if s.is_ll != 0
L = Sparse(F)
A = sparse(L*L')
else
LD = sparse(F.LD)
L, d = getLd!(LD)
A = (L * Diagonal(d)) * L'
end
# no need to sort buffers here, as A isa SparseMatrixCSC
# and it is taken care in sparse
p = get_perm(F)
if p != [1:s.n;]
pinv = Vector{Int}(undef, length(p))
for k = 1:length(p)
pinv[p[k]] = k
end
A = A[pinv,pinv]
end
A
end
sparse(D::Dense) = sparse(Sparse(D))
function sparse(FC::FactorComponent{Tv,:L}) where Tv
F = Factor(FC)
s = unsafe_load(pointer(F))
if s.is_ll == 0
throw(CHOLMODException("sparse: supported only for :LD on LDLt factorizations"))
end
sparse(Sparse(F))
end
sparse(FC::FactorComponent{Tv,:LD}) where {Tv} = sparse(Sparse(Factor(FC)))
# Calculate the offset into the stype field of the cholmod_sparse_struct and
# change the value
const __SPARSE_STYPE_OFFSET = fieldoffset(cholmod_sparse_struct, findfirst(name -> name === :stype, fieldnames(cholmod_sparse_struct))::Int)
function change_stype!(A::Sparse, i::Integer)
unsafe_store!(Ptr{Cint}(pointer(A) + __SPARSE_STYPE_OFFSET), i)
return A
end
free!(A::Dense) = free!(pointer(A))
free!(A::Sparse) = free!(pointer(A))
free!(F::Factor) = free!(pointer(F))
eltype(::Type{Dense{T}}) where {T<:VTypes} = T
eltype(::Type{Factor{T}}) where {T<:VTypes} = T
eltype(::Type{Sparse{T}}) where {T<:VTypes} = T
nnz(F::Factor) = nnz(Sparse(F))
function show(io::IO, F::Factor)
println(io, typeof(F))
showfactor(io, F)
end
function show(io::IO, FC::FactorComponent)
println(io, typeof(FC))
showfactor(io, Factor(FC))
end
function showfactor(io::IO, F::Factor)
s = unsafe_load(pointer(F))
print(io, """
type: $(s.is_ll!=0 ? "LLt" : "LDLt")
method: $(s.is_super!=0 ? "supernodal" : "simplicial")
maxnnz: $(Int(s.nzmax))
nnz: $(nnz(F))
success: $(s.minor == size(F, 1))
""")
end
# getindex not defined for these, so don't use the normal array printer
show(io::IO, ::MIME"text/plain", FC::FactorComponent) = show(io, FC)
show(io::IO, ::MIME"text/plain", F::Factor) = show(io, F)
isvalid(A::Dense) = check_dense(A)
isvalid(A::Sparse) = check_sparse(A)
isvalid(A::Factor) = check_factor(A)
function size(A::Union{Dense,Sparse})
s = unsafe_load(pointer(A))
return (Int(s.nrow), Int(s.ncol))
end