1+ module GBCHOLMOD
2+ import Base: (\ ), getproperty, show, size
3+ using LinearAlgebra
4+ import LinearAlgebra: (\ ),
5+ cholesky, cholesky!, det, diag, ishermitian, isposdef,
6+ issuccess, issymmetric, ldlt, ldlt!, logdet
7+ using .. SuiteSparseGraphBLAS: AbstractGBMatrix, AbstractGBVector, unsafepack!, unsafeunpack!, GBMatrix,
8+ GBVector, AbstractGBArray, LibGraphBLAS, ColMajor, sparsitystatus,
9+ _sizedjlmalloc, increment!, isshallow, nnz, tempunpack!, storedeltype
10+ using SuiteSparseGraphBLAS
11+
12+ using StorageOrders
13+
14+ import .. increment, .. increment!, .. decrement, .. decrement!
15+
16+ using SuiteSparse. LibSuiteSparse
17+ using SuiteSparse. CHOLMOD
18+ using SuiteSparse. CHOLMOD: VTypes, cholesky, cholesky!, ldlt!, Factor, FactorComponent,
19+ spsolve, CHOLMOD_A, change_stype!
20+
21+ function CHOLMOD. Sparse (A:: AbstractGBMatrix{T} , stype:: Integer ) where T
22+ colptr, rowval, nzval, repack! = tempunpack! (A, SuiteSparseGraphBLAS. Sparse (); order = ColMajor ())
23+ nzval = ! (T <: VTypes ) ? promote_type (T, Float64).(nzval) : nzval
24+ C = try
25+ CHOLMOD. Sparse (size (A)... , colptr, rowval, nzval, stype)
26+ catch
27+ rethrow ()
28+ finally
29+ repack! ()
30+ end
31+ if ishermitian (C)
32+ change_stype! (C, - 1 )
33+ end
34+ return C
35+ end
36+ CHOLMOD. Sparse (A:: AbstractGBMatrix ) = CHOLMOD. Sparse (A, 0 )
37+ CHOLMOD. Sparse (A:: Symmetric{Tv, AbstractGBMatrix{Tv}} ) where {Tv<: Real } =
38+ CHOLMOD. Sparse (A. data, A. uplo == ' L' ? - 1 : 1 )
39+ CHOLMOD. Sparse (A:: Hermitian{Tv,AbstractGBMatrix{Tv}} ) where {Tv} =
40+ CHOLMOD. Sparse (A. data, A. uplo == ' L' ? - 1 : 1 )
41+
42+ function CHOLMOD. Dense (A:: AbstractGBMatrix )
43+ x, repack! = tempunpack! (A, SuiteSparseGraphBLAS. Dense (); order = ColMajor ())
44+ return try
45+ CHOLMOD. Dense (x)
46+ catch
47+ rethrow ()
48+ finally
49+ repack! ()
50+ end
51+ end
52+
53+ # TODO : Dense -> GBMatrix
54+
55+ LinearAlgebra. cholesky! (F:: Factor , A:: Union {AbstractGBMatrix{T},
56+ AbstractGBMatrix{Complex{T}},
57+ Symmetric{T,AbstractGBMatrix{T}},
58+ Hermitian{Complex{T},AbstractGBMatrix{Complex{T}}},
59+ Hermitian{T,AbstractGBMatrix{T}}};
60+ shift = 0.0 , check:: Bool = true ) where {T<: Real } =
61+ CHOLMOD. cholesky! (F, CHOLMOD. Sparse (A); shift = shift, check = check)
62+
63+ LinearAlgebra. cholesky (A:: Union {AbstractGBMatrix{T}, AbstractGBMatrix{Complex{T}},
64+ Symmetric{T,AbstractGBMatrix{T}},
65+ Hermitian{Complex{T},AbstractGBMatrix{Complex{T}}},
66+ Hermitian{T,AbstractGBMatrix{T}}};
67+ kws... ) where {T<: Real } = CHOLMOD. cholesky (CHOLMOD. Sparse (A); kws... )
68+
69+ LinearAlgebra. ldlt! (F:: Factor , A:: Union {AbstractGBMatrix{T},
70+ AbstractGBMatrix{Complex{T}},
71+ Symmetric{T,AbstractGBMatrix{T}},
72+ Hermitian{Complex{T},AbstractGBMatrix{Complex{T}}},
73+ Hermitian{T,AbstractGBMatrix{T}}};
74+ shift = 0.0 , check:: Bool = true ) where {T<: Real } =
75+ CHOLMOD. ldlt! (F, CHOLMOD. Sparse (A), shift = shift, check = check)
76+
77+ LinearAlgebra. ldlt (A:: Union {AbstractGBMatrix{T},AbstractGBMatrix{Complex{T}},
78+ Symmetric{T,AbstractGBMatrix{T}},
79+ Hermitian{Complex{T},AbstractGBMatrix{Complex{T}}},
80+ Hermitian{T,AbstractGBMatrix{T}}};
81+ kws... ) where {T<: Real } = CHOLMOD. ldlt (CHOLMOD. Sparse (A); kws... )
82+
83+ function (\ )(L:: FactorComponent , B:: AbstractGBArray )
84+ sparse (L\ CHOLMOD. Sparse (B,0 ))
85+ end
86+
87+ (\ )(L:: FactorComponent , B:: Adjoint{<:Any,<:AbstractGBMatrix} ) = L \ copy (B)
88+ (\ )(L:: FactorComponent , B:: Transpose{<:Any,<:AbstractGBMatrix} ) = L \ copy (B)
89+ \ (adjL:: Adjoint{<:Any,<:FactorComponent} , B:: Union{VecOrMat,AbstractGBArray} ) = (L = adjL. parent; adjoint (L)\ B)
90+ (\ )(L:: Factor , B:: S ) where {S<: AbstractGBMatrix } = S (spsolve (CHOLMOD_A, L, CHOLMOD. Sparse (B, 0 )))
91+ (\ )(L:: Factor , B:: Adjoint{<:Any,<:AbstractGBMatrix} ) = L \ copy (B)
92+ (\ )(L:: Factor , B:: Transpose{<:Any,<:AbstractGBMatrix} ) = L \ copy (B)
93+ (\ )(L:: Factor , B:: V ) where {V<: AbstractGBVector } = V (spsolve (CHOLMOD_A, L, CHOLMOD. Sparse (B)))
94+ \ (adjL:: Adjoint{<:Any,<:Factor} , B:: AbstractGBArray ) = (L = adjL. parent; \ (adjoint (L), CHOLMOD. Sparse (B)))
95+
96+ # TODO : 1681 -> 1693
97+ # TODO : Correct Adjoint copy for GBMatrices
98+ # Dense -> GBMatrix
99+ # Sparse -> GBMatrix
100+ end
0 commit comments