$darkmode
Eigen  5.0.1-dev
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 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
20 template <typename Scalar>
21 struct cholmod_configure_matrix;
22 
23 template <>
24 struct cholmod_configure_matrix<double> {
25  template <typename CholmodType>
26  static void run(CholmodType& mat) {
27  mat.xtype = CHOLMOD_REAL;
28  mat.dtype = CHOLMOD_DOUBLE;
29  }
30 };
31 
32 template <>
33 struct cholmod_configure_matrix<std::complex<double> > {
34  template <typename CholmodType>
35  static void run(CholmodType& mat) {
36  mat.xtype = CHOLMOD_COMPLEX;
37  mat.dtype = CHOLMOD_DOUBLE;
38  }
39 };
40 
41 // Other scalar types are not yet supported by Cholmod
42 // template<> struct cholmod_configure_matrix<float> {
43 // template<typename CholmodType>
44 // static void run(CholmodType& mat) {
45 // mat.xtype = CHOLMOD_REAL;
46 // mat.dtype = CHOLMOD_SINGLE;
47 // }
48 // };
49 //
50 // template<> struct cholmod_configure_matrix<std::complex<float> > {
51 // template<typename CholmodType>
52 // static void run(CholmodType& mat) {
53 // mat.xtype = CHOLMOD_COMPLEX;
54 // mat.dtype = CHOLMOD_SINGLE;
55 // }
56 // };
57 
58 } // namespace internal
59 
63 template <typename Scalar_, int Options_, typename StorageIndex_>
65  cholmod_sparse res;
66  res.nzmax = mat.nonZeros();
67  res.nrow = mat.rows();
68  res.ncol = mat.cols();
69  res.p = mat.outerIndexPtr();
70  res.i = mat.innerIndexPtr();
71  res.x = mat.valuePtr();
72  res.z = 0;
73  res.sorted = 1;
74  if (mat.isCompressed()) {
75  res.packed = 1;
76  res.nz = 0;
77  } else {
78  res.packed = 0;
79  res.nz = mat.innerNonZeroPtr();
80  }
81 
82  res.dtype = 0;
83  res.stype = -1;
84 
85  if (internal::is_same<StorageIndex_, int>::value) {
86  res.itype = CHOLMOD_INT;
87  } else if (internal::is_same<StorageIndex_, SuiteSparse_long>::value) {
88  res.itype = CHOLMOD_LONG;
89  } else {
90  eigen_assert(false && "Index type not supported yet");
91  }
92 
93  // setup res.xtype
94  internal::cholmod_configure_matrix<Scalar_>::run(res);
95 
96  res.stype = 0;
97 
98  return res;
99 }
100 
101 template <typename Scalar_, int Options_, typename Index_>
102 const cholmod_sparse viewAsCholmod(const SparseMatrix<Scalar_, Options_, Index_>& mat) {
103  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<Scalar_, Options_, Index_> >(mat.const_cast_derived()));
104  return res;
105 }
106 
107 template <typename Scalar_, int Options_, typename Index_>
108 const cholmod_sparse viewAsCholmod(const SparseVector<Scalar_, Options_, Index_>& mat) {
109  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<Scalar_, Options_, Index_> >(mat.const_cast_derived()));
110  return res;
111 }
112 
115 template <typename Scalar_, int Options_, typename Index_, unsigned int UpLo>
117  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<Scalar_, Options_, Index_> >(mat.matrix().const_cast_derived()));
118 
119  if (UpLo == Upper) res.stype = 1;
120  if (UpLo == Lower) res.stype = -1;
121  // swap stype for rowmajor matrices (only works for real matrices)
122  EIGEN_STATIC_ASSERT((Options_ & RowMajorBit) == 0 || NumTraits<Scalar_>::IsComplex == 0,
123  THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
124  if (Options_ & RowMajorBit) res.stype *= -1;
125 
126  return res;
127 }
128 
131 template <typename Derived>
132 cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat) {
133  EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags & RowMajorBit) == 0,
134  THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
135  typedef typename Derived::Scalar Scalar;
136 
137  cholmod_dense res;
138  res.nrow = mat.rows();
139  res.ncol = mat.cols();
140  res.nzmax = res.nrow * res.ncol;
141  res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
142  res.x = (void*)(mat.derived().data());
143  res.z = 0;
144 
145  internal::cholmod_configure_matrix<Scalar>::run(res);
146 
147  return res;
148 }
149 
152 template <typename Scalar, typename StorageIndex>
155  cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol], static_cast<StorageIndex*>(cm.p),
156  static_cast<StorageIndex*>(cm.i), static_cast<Scalar*>(cm.x));
157 }
158 
161 template <typename Scalar, typename StorageIndex>
164  cm.n, cm.n, static_cast<StorageIndex*>(cm.p)[cm.n], static_cast<StorageIndex*>(cm.p),
165  static_cast<StorageIndex*>(cm.i), static_cast<Scalar*>(cm.x));
166 }
167 
168 namespace internal {
169 
170 // template specializations for int and long that call the correct cholmod method
171 
172 #define EIGEN_CHOLMOD_SPECIALIZE0(ret, name) \
173  template <typename StorageIndex_> \
174  inline ret cm_##name(cholmod_common& Common) { \
175  return cholmod_##name(&Common); \
176  } \
177  template <> \
178  inline ret cm_##name<SuiteSparse_long>(cholmod_common & Common) { \
179  return cholmod_l_##name(&Common); \
180  }
181 
182 #define EIGEN_CHOLMOD_SPECIALIZE1(ret, name, t1, a1) \
183  template <typename StorageIndex_> \
184  inline ret cm_##name(t1& a1, cholmod_common& Common) { \
185  return cholmod_##name(&a1, &Common); \
186  } \
187  template <> \
188  inline ret cm_##name<SuiteSparse_long>(t1 & a1, cholmod_common & Common) { \
189  return cholmod_l_##name(&a1, &Common); \
190  }
191 
192 EIGEN_CHOLMOD_SPECIALIZE0(int, start)
193 EIGEN_CHOLMOD_SPECIALIZE0(int, finish)
194 
195 EIGEN_CHOLMOD_SPECIALIZE1(int, free_factor, cholmod_factor*, L)
196 EIGEN_CHOLMOD_SPECIALIZE1(int, free_dense, cholmod_dense*, X)
197 EIGEN_CHOLMOD_SPECIALIZE1(int, free_sparse, cholmod_sparse*, A)
198 
199 EIGEN_CHOLMOD_SPECIALIZE1(cholmod_factor*, analyze, cholmod_sparse, A)
200 EIGEN_CHOLMOD_SPECIALIZE1(cholmod_sparse*, factor_to_sparse, cholmod_factor, L)
201 
202 template <typename StorageIndex_>
203 inline cholmod_dense* cm_solve(int sys, cholmod_factor& L, cholmod_dense& B, cholmod_common& Common) {
204  return cholmod_solve(sys, &L, &B, &Common);
205 }
206 template <>
207 inline cholmod_dense* cm_solve<SuiteSparse_long>(int sys, cholmod_factor& L, cholmod_dense& B, cholmod_common& Common) {
208  return cholmod_l_solve(sys, &L, &B, &Common);
209 }
210 
211 template <typename StorageIndex_>
212 inline cholmod_sparse* cm_spsolve(int sys, cholmod_factor& L, cholmod_sparse& B, cholmod_common& Common) {
213  return cholmod_spsolve(sys, &L, &B, &Common);
214 }
215 template <>
216 inline cholmod_sparse* cm_spsolve<SuiteSparse_long>(int sys, cholmod_factor& L, cholmod_sparse& B,
217  cholmod_common& Common) {
218  return cholmod_l_spsolve(sys, &L, &B, &Common);
219 }
220 
221 template <typename StorageIndex_>
222 inline int cm_factorize_p(cholmod_sparse* A, double beta[2], StorageIndex_* fset, std::size_t fsize, cholmod_factor* L,
223  cholmod_common& Common) {
224  return cholmod_factorize_p(A, beta, fset, fsize, L, &Common);
225 }
226 template <>
227 inline int cm_factorize_p<SuiteSparse_long>(cholmod_sparse* A, double beta[2], SuiteSparse_long* fset,
228  std::size_t fsize, cholmod_factor* L, cholmod_common& Common) {
229  return cholmod_l_factorize_p(A, beta, fset, fsize, L, &Common);
230 }
231 
232 #undef EIGEN_CHOLMOD_SPECIALIZE0
233 #undef EIGEN_CHOLMOD_SPECIALIZE1
234 
235 } // namespace internal
236 
237 enum CholmodMode { CholmodAuto, CholmodSimplicialLLt, CholmodSupernodalLLt, CholmodLDLt };
238 
244 template <typename MatrixType_, int UpLo_, typename Derived>
245 class CholmodBase : public SparseSolverBase<Derived> {
246  protected:
248  using Base::derived;
249  using Base::m_isInitialized;
250 
251  public:
252  typedef MatrixType_ MatrixType;
253  enum { UpLo = UpLo_ };
254  typedef typename MatrixType::Scalar Scalar;
255  typedef typename MatrixType::RealScalar RealScalar;
256  typedef MatrixType CholMatrixType;
257  typedef typename MatrixType::StorageIndex StorageIndex;
258  enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
259 
260  public:
261  CholmodBase() : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) {
262  EIGEN_STATIC_ASSERT((internal::is_same<double, RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
263  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
264  internal::cm_start<StorageIndex>(m_cholmod);
265  }
266 
267  explicit CholmodBase(const MatrixType& matrix)
268  : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) {
269  EIGEN_STATIC_ASSERT((internal::is_same<double, RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
270  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
271  internal::cm_start<StorageIndex>(m_cholmod);
272  compute(matrix);
273  }
274 
275  ~CholmodBase() {
276  if (m_cholmodFactor) internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
277  internal::cm_finish<StorageIndex>(m_cholmod);
278  }
279 
280  inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
281  inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
282 
289  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
290  return m_info;
291  }
292 
294  Derived& compute(const MatrixType& matrix) {
295  analyzePattern(matrix);
296  factorize(matrix);
297  return derived();
298  }
299 
306  void analyzePattern(const MatrixType& matrix) {
307  if (m_cholmodFactor) {
308  internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
309  m_cholmodFactor = 0;
310  }
311  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
312  m_cholmodFactor = internal::cm_analyze<StorageIndex>(A, m_cholmod);
313 
314  this->m_isInitialized = true;
315  this->m_info = Success;
316  m_analysisIsOk = true;
317  m_factorizationIsOk = false;
318  }
319 
327  void factorize(const MatrixType& matrix) {
328  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
329  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
330  internal::cm_factorize_p<StorageIndex>(&A, m_shiftOffset, 0, 0, m_cholmodFactor, m_cholmod);
331 
332  // If the factorization failed, either the input matrix was zero (so m_cholmodFactor == nullptr), or minor is the
333  // column at which it failed. On success minor == n.
334  this->m_info =
335  (m_cholmodFactor != nullptr && m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
336  m_factorizationIsOk = true;
337  }
338 
341  cholmod_common& cholmod() { return m_cholmod; }
342 
343 #ifndef EIGEN_PARSED_BY_DOXYGEN
344 
345  template <typename Rhs, typename Dest>
346  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
347  eigen_assert(m_factorizationIsOk &&
348  "The decomposition is not in a valid state for solving, you must first call either compute() or "
349  "symbolic()/numeric()");
350  const Index size = m_cholmodFactor->n;
351  EIGEN_UNUSED_VARIABLE(size);
352  eigen_assert(size == b.rows());
353 
354  // Cholmod needs column-major storage without inner-stride, which corresponds to the default behavior of Ref.
356 
357  cholmod_dense b_cd = viewAsCholmod(b_ref);
358  cholmod_dense* x_cd = internal::cm_solve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cd, m_cholmod);
359  if (!x_cd) {
360  this->m_info = NumericalIssue;
361  return;
362  }
363  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
364  // NOTE Actually, the copy can be avoided by calling cholmod_solve2 instead of cholmod_solve
365  dest = Matrix<Scalar, Dest::RowsAtCompileTime, Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),
366  b.rows(), b.cols());
367  internal::cm_free_dense<StorageIndex>(x_cd, m_cholmod);
368  }
369 
371  template <typename RhsDerived, typename DestDerived>
372  void _solve_impl(const SparseMatrixBase<RhsDerived>& b, SparseMatrixBase<DestDerived>& dest) const {
373  eigen_assert(m_factorizationIsOk &&
374  "The decomposition is not in a valid state for solving, you must first call either compute() or "
375  "symbolic()/numeric()");
376  const Index size = m_cholmodFactor->n;
377  EIGEN_UNUSED_VARIABLE(size);
378  eigen_assert(size == b.rows());
379 
380  // note: cs stands for Cholmod Sparse
381  Ref<SparseMatrix<typename RhsDerived::Scalar, ColMajor, typename RhsDerived::StorageIndex> > b_ref(
382  b.const_cast_derived());
383  cholmod_sparse b_cs = viewAsCholmod(b_ref);
384  cholmod_sparse* x_cs = internal::cm_spsolve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cs, m_cholmod);
385  if (!x_cs) {
386  this->m_info = NumericalIssue;
387  return;
388  }
389  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
390  // NOTE cholmod_spsolve in fact just calls the dense solver for blocks of 4 columns at a time (similar to Eigen's
391  // sparse solver)
392  dest.derived() = viewAsEigen<typename DestDerived::Scalar, typename DestDerived::StorageIndex>(*x_cs);
393  internal::cm_free_sparse<StorageIndex>(x_cs, m_cholmod);
394  }
395 #endif // EIGEN_PARSED_BY_DOXYGEN
396 
406  Derived& setShift(const RealScalar& offset) {
407  m_shiftOffset[0] = double(offset);
408  return derived();
409  }
410 
412  Scalar determinant() const {
413  using std::exp;
414  return exp(logDeterminant());
415  }
416 
418  Scalar logDeterminant() const {
419  using numext::real;
420  using std::log;
421  eigen_assert(m_factorizationIsOk &&
422  "The decomposition is not in a valid state for solving, you must first call either compute() or "
423  "symbolic()/numeric()");
424 
425  RealScalar logDet = 0;
426  Scalar* x = static_cast<Scalar*>(m_cholmodFactor->x);
427  if (m_cholmodFactor->is_super) {
428  // Supernodal factorization stored as a packed list of dense column-major blocks,
429  // as described by the following structure:
430 
431  // super[k] == index of the first column of the j-th super node
432  StorageIndex* super = static_cast<StorageIndex*>(m_cholmodFactor->super);
433  // pi[k] == offset to the description of row indices
434  StorageIndex* pi = static_cast<StorageIndex*>(m_cholmodFactor->pi);
435  // px[k] == offset to the respective dense block
436  StorageIndex* px = static_cast<StorageIndex*>(m_cholmodFactor->px);
437 
438  Index nb_super_nodes = m_cholmodFactor->nsuper;
439  for (Index k = 0; k < nb_super_nodes; ++k) {
440  StorageIndex ncols = super[k + 1] - super[k];
441  StorageIndex nrows = pi[k + 1] - pi[k];
442 
443  Map<const Array<Scalar, 1, Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows + 1));
444  logDet += sk.real().log().sum();
445  }
446  } else {
447  // Simplicial factorization stored as standard CSC matrix.
448  StorageIndex* p = static_cast<StorageIndex*>(m_cholmodFactor->p);
449  Index size = m_cholmodFactor->n;
450  for (Index k = 0; k < size; ++k) logDet += log(real(x[p[k]]));
451  }
452  if (m_cholmodFactor->is_ll) logDet *= 2.0;
453  return logDet;
454  }
455 
456  template <typename Stream>
457  void dumpMemory(Stream& /*s*/) {}
458 
459  protected:
460  mutable cholmod_common m_cholmod;
461  cholmod_factor* m_cholmodFactor;
462  double m_shiftOffset[2];
463  mutable ComputationInfo m_info;
464  int m_factorizationIsOk;
465  int m_analysisIsOk;
466 };
467 
491 template <typename MatrixType_, int UpLo_ = Lower>
492 class CholmodSimplicialLLT : public CholmodBase<MatrixType_, UpLo_, CholmodSimplicialLLT<MatrixType_, UpLo_> > {
494  using Base::m_cholmod;
495 
496  public:
497  typedef MatrixType_ MatrixType;
498  typedef typename MatrixType::Scalar Scalar;
499  typedef typename MatrixType::RealScalar RealScalar;
500  typedef typename MatrixType::StorageIndex StorageIndex;
503 
504  CholmodSimplicialLLT() : Base() { init(); }
505 
506  CholmodSimplicialLLT(const MatrixType& matrix) : Base() {
507  init();
508  this->compute(matrix);
509  }
510 
512 
514  inline MatrixL matrixL() const { return viewAsEigen<Scalar, StorageIndex>(*Base::m_cholmodFactor); }
515 
517  inline MatrixU matrixU() const { return matrixL().adjoint(); }
518 
519  protected:
520  void init() {
521  m_cholmod.final_asis = 0;
522  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
523  m_cholmod.final_ll = 1;
524  }
525 };
526 
550 template <typename MatrixType_, int UpLo_ = Lower>
551 class CholmodSimplicialLDLT : public CholmodBase<MatrixType_, UpLo_, CholmodSimplicialLDLT<MatrixType_, UpLo_> > {
553  using Base::m_cholmod;
554 
555  public:
556  typedef MatrixType_ MatrixType;
557  typedef typename MatrixType::Scalar Scalar;
558  typedef typename MatrixType::RealScalar RealScalar;
559  typedef typename MatrixType::StorageIndex StorageIndex;
563 
564  CholmodSimplicialLDLT() : Base() { init(); }
565 
566  CholmodSimplicialLDLT(const MatrixType& matrix) : Base() {
567  init();
568  this->compute(matrix);
569  }
570 
572 
574  inline VectorType vectorD() const {
575  auto cholmodL = viewAsEigen<Scalar, StorageIndex>(*Base::m_cholmodFactor);
576 
577  VectorType D{cholmodL.rows()};
578 
579  for (Index k = 0; k < cholmodL.outerSize(); ++k) {
580  typename decltype(cholmodL)::InnerIterator it{cholmodL, k};
581  D(k) = it.value();
582  }
583 
584  return D;
585  }
586 
588  inline MatrixL matrixL() const { return viewAsEigen<Scalar, StorageIndex>(*Base::m_cholmodFactor); }
589 
591  inline MatrixU matrixU() const { return matrixL().adjoint(); }
592 
593  protected:
594  void init() {
595  m_cholmod.final_asis = 1;
596  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
597  }
598 };
599 
623 template <typename MatrixType_, int UpLo_ = Lower>
624 class CholmodSupernodalLLT : public CholmodBase<MatrixType_, UpLo_, CholmodSupernodalLLT<MatrixType_, UpLo_> > {
626  using Base::m_cholmod;
627 
628  public:
629  typedef MatrixType_ MatrixType;
630  typedef typename MatrixType::Scalar Scalar;
631  typedef typename MatrixType::RealScalar RealScalar;
632  typedef typename MatrixType::StorageIndex StorageIndex;
633 
634  CholmodSupernodalLLT() : Base() { init(); }
635 
636  CholmodSupernodalLLT(const MatrixType& matrix) : Base() {
637  init();
638  this->compute(matrix);
639  }
640 
642 
644  inline MatrixType matrixL() const {
645  // Convert Cholmod factor's supernodal storage format to Eigen's CSC storage format
646  cholmod_sparse* cholmodL = internal::cm_factor_to_sparse(*Base::m_cholmodFactor, m_cholmod);
647  MatrixType L = viewAsEigen<Scalar, StorageIndex>(*cholmodL);
648  internal::cm_free_sparse<StorageIndex>(cholmodL, m_cholmod);
649 
650  return L;
651  }
652 
654  inline MatrixType matrixU() const { return matrixL().adjoint(); }
655 
656  protected:
657  void init() {
658  m_cholmod.final_asis = 1;
659  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
660  }
661 };
662 
688 template <typename MatrixType_, int UpLo_ = Lower>
689 class CholmodDecomposition : public CholmodBase<MatrixType_, UpLo_, CholmodDecomposition<MatrixType_, UpLo_> > {
691  using Base::m_cholmod;
692 
693  public:
694  typedef MatrixType_ MatrixType;
695 
696  CholmodDecomposition() : Base() { init(); }
697 
698  CholmodDecomposition(const MatrixType& matrix) : Base() {
699  init();
700  this->compute(matrix);
701  }
702 
704 
705  void setMode(CholmodMode mode) {
706  switch (mode) {
707  case CholmodAuto:
708  m_cholmod.final_asis = 1;
709  m_cholmod.supernodal = CHOLMOD_AUTO;
710  break;
711  case CholmodSimplicialLLt:
712  m_cholmod.final_asis = 0;
713  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
714  m_cholmod.final_ll = 1;
715  break;
716  case CholmodSupernodalLLt:
717  m_cholmod.final_asis = 1;
718  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
719  break;
720  case CholmodLDLt:
721  m_cholmod.final_asis = 1;
722  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
723  break;
724  default:
725  break;
726  }
727  }
728 
729  protected:
730  void init() {
731  m_cholmod.final_asis = 1;
732  m_cholmod.supernodal = CHOLMOD_AUTO;
733  }
734 };
735 
736 } // end namespace Eigen
737 
738 #endif // EIGEN_CHOLMODSUPPORT_H
cholmod_common & cholmod()
Definition: CholmodSupport.h:341
constexpr Derived & derived()
Definition: EigenBase.h:49
void analyzePattern(const MatrixType &matrix)
Definition: CholmodSupport.h:306
const AdjointReturnType adjoint() const
Definition: TriangularMatrix.h:224
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
MatrixL matrixL() const
Definition: CholmodSupport.h:588
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Map< const SparseMatrix< Scalar, ColMajor, StorageIndex > > viewAsEigen(cholmod_sparse &cm)
Definition: CholmodSupport.h:153
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
MatrixType matrixL() const
Definition: CholmodSupport.h:644
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:52
Definition: BFloat16.h:231
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
A supernodal Cholesky (LLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:624
const unsigned int RowMajorBit
Definition: Constants.h:70
MatrixU matrixU() const
Definition: CholmodSupport.h:517
Definition: Constants.h:211
Scalar logDeterminant() const
Definition: CholmodSupport.h:418
MatrixL matrixL() const
Definition: CholmodSupport.h:514
Definition: Constants.h:442
Convenience specialization of Stride to specify only an inner stride See class Map for some examples...
Definition: Stride.h:93
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)
Definition: CholmodSupport.h:64
MatrixType matrixU() const
Definition: CholmodSupport.h:654
A general Cholesky factorization and solver based on Cholmod.
Definition: CholmodSupport.h:689
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
MatrixU matrixU() const
Definition: CholmodSupport.h:591
Definition: Constants.h:213
The base class for the direct Cholesky factorization of Cholmod.
Definition: CholmodSupport.h:245
Derived & compute(const MatrixType &matrix)
Definition: CholmodSupport.h:294
Definition: Constants.h:440
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:264
A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:492
VectorType vectorD() const
Definition: CholmodSupport.h:574
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)
void factorize(const MatrixType &matrix)
Definition: CholmodSupport.h:327
A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:551
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:166
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: CholmodSupport.h:288
Scalar value() const
Definition: CoreIterators.h:48
constexpr Index rows() const noexcept
Definition: EigenBase.h:59
constexpr Index cols() const noexcept
Definition: EigenBase.h:61
ComputationInfo
Definition: Constants.h:438
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_exp_op< typename Derived::Scalar >, const Derived > exp(const Eigen::ArrayBase< Derived > &x)
Derived & setShift(const RealScalar &offset)
Definition: CholmodSupport.h:406
An InnerIterator allows to loop over the element of any matrix expression.
Definition: CoreIterators.h:37
Scalar determinant() const
Definition: CholmodSupport.h:412