Medial Code Documentation
Loading...
Searching...
No Matches
CholmodSupport.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_CHOLMODSUPPORT_H
11#define EIGEN_CHOLMODSUPPORT_H
12
13namespace Eigen {
14
15namespace internal {
16
17template<typename Scalar, typename CholmodType>
18void cholmod_configure_matrix(CholmodType& mat)
19{
20 if (internal::is_same<Scalar,float>::value)
21 {
22 mat.xtype = CHOLMOD_REAL;
23 mat.dtype = CHOLMOD_SINGLE;
24 }
25 else if (internal::is_same<Scalar,double>::value)
26 {
27 mat.xtype = CHOLMOD_REAL;
28 mat.dtype = CHOLMOD_DOUBLE;
29 }
30 else if (internal::is_same<Scalar,std::complex<float> >::value)
31 {
32 mat.xtype = CHOLMOD_COMPLEX;
33 mat.dtype = CHOLMOD_SINGLE;
34 }
35 else if (internal::is_same<Scalar,std::complex<double> >::value)
36 {
37 mat.xtype = CHOLMOD_COMPLEX;
38 mat.dtype = CHOLMOD_DOUBLE;
39 }
40 else
41 {
42 eigen_assert(false && "Scalar type not supported by CHOLMOD");
43 }
44}
45
46} // namespace internal
47
51template<typename _Scalar, int _Options, typename _StorageIndex>
52cholmod_sparse viewAsCholmod(SparseMatrix<_Scalar,_Options,_StorageIndex>& mat)
53{
54 cholmod_sparse res;
55 res.nzmax = mat.nonZeros();
56 res.nrow = mat.rows();;
57 res.ncol = mat.cols();
58 res.p = mat.outerIndexPtr();
59 res.i = mat.innerIndexPtr();
60 res.x = mat.valuePtr();
61 res.z = 0;
62 res.sorted = 1;
63 if(mat.isCompressed())
64 {
65 res.packed = 1;
66 res.nz = 0;
67 }
68 else
69 {
70 res.packed = 0;
71 res.nz = mat.innerNonZeroPtr();
72 }
73
74 res.dtype = 0;
75 res.stype = -1;
76
77 if (internal::is_same<_StorageIndex,int>::value)
78 {
79 res.itype = CHOLMOD_INT;
80 }
81 else if (internal::is_same<_StorageIndex,SuiteSparse_long>::value)
82 {
83 res.itype = CHOLMOD_LONG;
84 }
85 else
86 {
87 eigen_assert(false && "Index type not supported yet");
88 }
89
90 // setup res.xtype
91 internal::cholmod_configure_matrix<_Scalar>(res);
92
93 res.stype = 0;
94
95 return res;
96}
97
98template<typename _Scalar, int _Options, typename _Index>
99const cholmod_sparse viewAsCholmod(const SparseMatrix<_Scalar,_Options,_Index>& mat)
100{
101 cholmod_sparse res = viewAsCholmod(mat.const_cast_derived());
102 return res;
103}
104
107template<typename _Scalar, int _Options, typename _Index, unsigned int UpLo>
108cholmod_sparse viewAsCholmod(const SparseSelfAdjointView<const SparseMatrix<_Scalar,_Options,_Index>, UpLo>& mat)
109{
110 cholmod_sparse res = viewAsCholmod(mat.matrix().const_cast_derived());
111
112 if(UpLo==Upper) res.stype = 1;
113 if(UpLo==Lower) res.stype = -1;
114
115 return res;
116}
117
120template<typename Derived>
121cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat)
122{
123 EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
124 typedef typename Derived::Scalar Scalar;
125
126 cholmod_dense res;
127 res.nrow = mat.rows();
128 res.ncol = mat.cols();
129 res.nzmax = res.nrow * res.ncol;
130 res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
131 res.x = (void*)(mat.derived().data());
132 res.z = 0;
133
134 internal::cholmod_configure_matrix<Scalar>(res);
135
136 return res;
137}
138
141template<typename Scalar, int Flags, typename StorageIndex>
142MappedSparseMatrix<Scalar,Flags,StorageIndex> viewAsEigen(cholmod_sparse& cm)
143{
144 return MappedSparseMatrix<Scalar,Flags,StorageIndex>
145 (cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol],
146 static_cast<StorageIndex*>(cm.p), static_cast<StorageIndex*>(cm.i),static_cast<Scalar*>(cm.x) );
147}
148
149enum CholmodMode {
150 CholmodAuto, CholmodSimplicialLLt, CholmodSupernodalLLt, CholmodLDLt
151};
152
153
159template<typename _MatrixType, int _UpLo, typename Derived>
160class CholmodBase : public SparseSolverBase<Derived>
161{
162 protected:
164 using Base::derived;
165 using Base::m_isInitialized;
166 public:
167 typedef _MatrixType MatrixType;
168 enum { UpLo = _UpLo };
169 typedef typename MatrixType::Scalar Scalar;
170 typedef typename MatrixType::RealScalar RealScalar;
171 typedef MatrixType CholMatrixType;
172 typedef typename MatrixType::StorageIndex StorageIndex;
173 enum {
174 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
175 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
176 };
177
178 public:
179
181 : m_cholmodFactor(0), m_info(Success)
182 {
183 m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0);
184 cholmod_start(&m_cholmod);
185 }
186
187 explicit CholmodBase(const MatrixType& matrix)
188 : m_cholmodFactor(0), m_info(Success)
189 {
190 m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0);
191 cholmod_start(&m_cholmod);
192 compute(matrix);
193 }
194
196 {
197 if(m_cholmodFactor)
198 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
199 cholmod_finish(&m_cholmod);
200 }
201
202 inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
203 inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
204
211 {
212 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
213 return m_info;
214 }
215
217 Derived& compute(const MatrixType& matrix)
218 {
219 analyzePattern(matrix);
220 factorize(matrix);
221 return derived();
222 }
223
230 void analyzePattern(const MatrixType& matrix)
231 {
232 if(m_cholmodFactor)
233 {
234 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
235 m_cholmodFactor = 0;
236 }
237 cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
238 m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
239
240 this->m_isInitialized = true;
241 this->m_info = Success;
242 m_analysisIsOk = true;
243 m_factorizationIsOk = false;
244 }
245
252 void factorize(const MatrixType& matrix)
253 {
254 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
255 cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
256 cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod);
257
258 // If the factorization failed, minor is the column at which it did. On success minor == n.
259 this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
260 m_factorizationIsOk = true;
261 }
262
265 cholmod_common& cholmod() { return m_cholmod; }
266
267 #ifndef EIGEN_PARSED_BY_DOXYGEN
269 template<typename Rhs,typename Dest>
270 void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
271 {
272 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
273 const Index size = m_cholmodFactor->n;
274 EIGEN_UNUSED_VARIABLE(size);
275 eigen_assert(size==b.rows());
276
277 // note: cd stands for Cholmod Dense
278 Rhs& b_ref(b.const_cast_derived());
279 cholmod_dense b_cd = viewAsCholmod(b_ref);
280 cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
281 if(!x_cd)
282 {
283 this->m_info = NumericalIssue;
284 return;
285 }
286 // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
287 dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
288 cholmod_free_dense(&x_cd, &m_cholmod);
289 }
290
292 template<typename RhsScalar, int RhsOptions, typename RhsIndex, typename DestScalar, int DestOptions, typename DestIndex>
293 void _solve_impl(const SparseMatrix<RhsScalar,RhsOptions,RhsIndex> &b, SparseMatrix<DestScalar,DestOptions,DestIndex> &dest) const
294 {
295 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
296 const Index size = m_cholmodFactor->n;
297 EIGEN_UNUSED_VARIABLE(size);
298 eigen_assert(size==b.rows());
299
300 // note: cs stands for Cholmod Sparse
301 cholmod_sparse b_cs = viewAsCholmod(b);
302 cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
303 if(!x_cs)
304 {
305 this->m_info = NumericalIssue;
306 return;
307 }
308 // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
309 dest = viewAsEigen<DestScalar,DestOptions,DestIndex>(*x_cs);
310 cholmod_free_sparse(&x_cs, &m_cholmod);
311 }
312 #endif // EIGEN_PARSED_BY_DOXYGEN
313
314
324 Derived& setShift(const RealScalar& offset)
325 {
326 m_shiftOffset[0] = offset;
327 return derived();
328 }
329
330 template<typename Stream>
331 void dumpMemory(Stream& /*s*/)
332 {}
333
334 protected:
335 mutable cholmod_common m_cholmod;
336 cholmod_factor* m_cholmodFactor;
337 RealScalar m_shiftOffset[2];
338 mutable ComputationInfo m_info;
339 int m_factorizationIsOk;
340 int m_analysisIsOk;
341};
342
363template<typename _MatrixType, int _UpLo = Lower>
364class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT<_MatrixType, _UpLo> >
365{
367 using Base::m_cholmod;
368
369 public:
370
371 typedef _MatrixType MatrixType;
372
373 CholmodSimplicialLLT() : Base() { init(); }
374
375 CholmodSimplicialLLT(const MatrixType& matrix) : Base()
376 {
377 init();
378 this->compute(matrix);
379 }
380
382 protected:
383 void init()
384 {
385 m_cholmod.final_asis = 0;
386 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
387 m_cholmod.final_ll = 1;
388 }
389};
390
391
412template<typename _MatrixType, int _UpLo = Lower>
413class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT<_MatrixType, _UpLo> >
414{
416 using Base::m_cholmod;
417
418 public:
419
420 typedef _MatrixType MatrixType;
421
422 CholmodSimplicialLDLT() : Base() { init(); }
423
424 CholmodSimplicialLDLT(const MatrixType& matrix) : Base()
425 {
426 init();
427 this->compute(matrix);
428 }
429
431 protected:
432 void init()
433 {
434 m_cholmod.final_asis = 1;
435 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
436 }
437};
438
459template<typename _MatrixType, int _UpLo = Lower>
460class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT<_MatrixType, _UpLo> >
461{
463 using Base::m_cholmod;
464
465 public:
466
467 typedef _MatrixType MatrixType;
468
469 CholmodSupernodalLLT() : Base() { init(); }
470
471 CholmodSupernodalLLT(const MatrixType& matrix) : Base()
472 {
473 init();
474 this->compute(matrix);
475 }
476
478 protected:
479 void init()
480 {
481 m_cholmod.final_asis = 1;
482 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
483 }
484};
485
508template<typename _MatrixType, int _UpLo = Lower>
509class CholmodDecomposition : public CholmodBase<_MatrixType, _UpLo, CholmodDecomposition<_MatrixType, _UpLo> >
510{
512 using Base::m_cholmod;
513
514 public:
515
516 typedef _MatrixType MatrixType;
517
518 CholmodDecomposition() : Base() { init(); }
519
520 CholmodDecomposition(const MatrixType& matrix) : Base()
521 {
522 init();
523 this->compute(matrix);
524 }
525
527
528 void setMode(CholmodMode mode)
529 {
530 switch(mode)
531 {
532 case CholmodAuto:
533 m_cholmod.final_asis = 1;
534 m_cholmod.supernodal = CHOLMOD_AUTO;
535 break;
536 case CholmodSimplicialLLt:
537 m_cholmod.final_asis = 0;
538 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
539 m_cholmod.final_ll = 1;
540 break;
541 case CholmodSupernodalLLt:
542 m_cholmod.final_asis = 1;
543 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
544 break;
545 case CholmodLDLt:
546 m_cholmod.final_asis = 1;
547 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
548 break;
549 default:
550 break;
551 }
552 }
553 protected:
554 void init()
555 {
556 m_cholmod.final_asis = 1;
557 m_cholmod.supernodal = CHOLMOD_AUTO;
558 }
559};
560
561} // end namespace Eigen
562
563#endif // EIGEN_CHOLMODSUPPORT_H
The base class for the direct Cholesky factorization of Cholmod.
Definition CholmodSupport.h:161
cholmod_common & cholmod()
Returns a reference to the Cholmod's configuration structure to get a full control over the performed...
Definition CholmodSupport.h:265
Derived & compute(const MatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix.
Definition CholmodSupport.h:217
void analyzePattern(const MatrixType &matrix)
Performs a symbolic decomposition on the sparsity pattern of matrix.
Definition CholmodSupport.h:230
void factorize(const MatrixType &matrix)
Performs a numeric decomposition of matrix.
Definition CholmodSupport.h:252
ComputationInfo info() const
Reports whether previous computation was successful.
Definition CholmodSupport.h:210
Derived & setShift(const RealScalar &offset)
Sets the shift parameter that will be used to adjust the diagonal coefficients during the numerical f...
Definition CholmodSupport.h:324
A general Cholesky factorization and solver based on Cholmod.
Definition CholmodSupport.h:510
A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod.
Definition CholmodSupport.h:414
A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.
Definition CholmodSupport.h:365
A supernodal Cholesky (LLT) factorization and solver based on Cholmod.
Definition CholmodSupport.h:461
Pseudo expression representing a solving operation.
Definition Solve.h:63
A base class for sparse solvers.
Definition SparseSolverBase.h:54
ComputationInfo
Enum for reporting the status of a computation.
Definition Constants.h:430
@ NumericalIssue
The provided data did not satisfy the prerequisites.
Definition Constants.h:434
@ Success
Computation was successful.
Definition Constants.h:432
Definition inference.c:32