CholmodSupport.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename Scalar> struct cholmod_configure_matrix;
20 
21 template<> struct cholmod_configure_matrix<double> {
22  template<typename CholmodType>
23  static void run(CholmodType& mat) {
24  mat.xtype = CHOLMOD_REAL;
25  mat.dtype = CHOLMOD_DOUBLE;
26  }
27 };
28 
29 template<> struct cholmod_configure_matrix<std::complex<double> > {
30  template<typename CholmodType>
31  static void run(CholmodType& mat) {
32  mat.xtype = CHOLMOD_COMPLEX;
33  mat.dtype = CHOLMOD_DOUBLE;
34  }
35 };
36 
37 // Other scalar types are not yet supported by Cholmod
38 // template<> struct cholmod_configure_matrix<float> {
39 // template<typename CholmodType>
40 // static void run(CholmodType& mat) {
41 // mat.xtype = CHOLMOD_REAL;
42 // mat.dtype = CHOLMOD_SINGLE;
43 // }
44 // };
45 //
46 // template<> struct cholmod_configure_matrix<std::complex<float> > {
47 // template<typename CholmodType>
48 // static void run(CholmodType& mat) {
49 // mat.xtype = CHOLMOD_COMPLEX;
50 // mat.dtype = CHOLMOD_SINGLE;
51 // }
52 // };
53 
54 } // namespace internal
55 
59 template<typename Scalar_, int Options_, typename StorageIndex_>
61 {
62  cholmod_sparse res;
63  res.nzmax = mat.nonZeros();
64  res.nrow = mat.rows();
65  res.ncol = mat.cols();
66  res.p = mat.outerIndexPtr();
67  res.i = mat.innerIndexPtr();
68  res.x = mat.valuePtr();
69  res.z = 0;
70  res.sorted = 1;
71  if(mat.isCompressed())
72  {
73  res.packed = 1;
74  res.nz = 0;
75  }
76  else
77  {
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  {
87  res.itype = CHOLMOD_INT;
88  }
89  else if (internal::is_same<StorageIndex_,SuiteSparse_long>::value)
90  {
91  res.itype = CHOLMOD_LONG;
92  }
93  else
94  {
95  eigen_assert(false && "Index type not supported yet");
96  }
97 
98  // setup res.xtype
99  internal::cholmod_configure_matrix<Scalar_>::run(res);
100 
101  res.stype = 0;
102 
103  return res;
104 }
105 
106 template<typename Scalar_, int Options_, typename Index_>
108 {
110  return res;
111 }
112 
113 template<typename Scalar_, int Options_, typename Index_>
115 {
117  return res;
118 }
119 
122 template<typename Scalar_, int Options_, typename Index_, unsigned int UpLo>
124 {
125  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<Scalar_,Options_,Index_> >(mat.matrix().const_cast_derived()));
126 
127  if(UpLo==Upper) res.stype = 1;
128  if(UpLo==Lower) res.stype = -1;
129  // swap stype for rowmajor matrices (only works for real matrices)
130  EIGEN_STATIC_ASSERT((Options_ & RowMajorBit) == 0 || NumTraits<Scalar_>::IsComplex == 0, THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
131  if(Options_ & RowMajorBit) res.stype *=-1;
132 
133  return res;
134 }
135 
138 template<typename Derived>
140 {
141  EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
142  typedef typename Derived::Scalar Scalar;
143 
144  cholmod_dense res;
145  res.nrow = mat.rows();
146  res.ncol = mat.cols();
147  res.nzmax = res.nrow * res.ncol;
148  res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
149  res.x = (void*)(mat.derived().data());
150  res.z = 0;
151 
152  internal::cholmod_configure_matrix<Scalar>::run(res);
153 
154  return res;
155 }
156 
159 template<typename Scalar, int Flags, typename StorageIndex>
161 {
163  (cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol],
164  static_cast<StorageIndex*>(cm.p), static_cast<StorageIndex*>(cm.i),static_cast<Scalar*>(cm.x) );
165 }
166 
167 namespace internal {
168 
169 // template specializations for int and long that call the correct cholmod method
170 
171 #define EIGEN_CHOLMOD_SPECIALIZE0(ret, name) \
172  template<typename StorageIndex_> inline ret cm_ ## name (cholmod_common &Common) { return cholmod_ ## name (&Common); } \
173  template<> inline ret cm_ ## name<SuiteSparse_long> (cholmod_common &Common) { return cholmod_l_ ## name (&Common); }
174 
175 #define EIGEN_CHOLMOD_SPECIALIZE1(ret, name, t1, a1) \
176  template<typename StorageIndex_> inline ret cm_ ## name (t1& a1, cholmod_common &Common) { return cholmod_ ## name (&a1, &Common); } \
177  template<> inline ret cm_ ## name<SuiteSparse_long> (t1& a1, cholmod_common &Common) { return cholmod_l_ ## name (&a1, &Common); }
178 
179 EIGEN_CHOLMOD_SPECIALIZE0(int, start)
180 EIGEN_CHOLMOD_SPECIALIZE0(int, finish)
181 
182 EIGEN_CHOLMOD_SPECIALIZE1(int, free_factor, cholmod_factor*, L)
183 EIGEN_CHOLMOD_SPECIALIZE1(int, free_dense, cholmod_dense*, X)
184 EIGEN_CHOLMOD_SPECIALIZE1(int, free_sparse, cholmod_sparse*, A)
185 
186 EIGEN_CHOLMOD_SPECIALIZE1(cholmod_factor*, analyze, cholmod_sparse, A)
187 
188 template<typename StorageIndex_> inline cholmod_dense* cm_solve (int sys, cholmod_factor& L, cholmod_dense& B, cholmod_common &Common) { return cholmod_solve (sys, &L, &B, &Common); }
189 template<> inline cholmod_dense* cm_solve<SuiteSparse_long> (int sys, cholmod_factor& L, cholmod_dense& B, cholmod_common &Common) { return cholmod_l_solve (sys, &L, &B, &Common); }
190 
191 template<typename StorageIndex_> inline cholmod_sparse* cm_spsolve (int sys, cholmod_factor& L, cholmod_sparse& B, cholmod_common &Common) { return cholmod_spsolve (sys, &L, &B, &Common); }
192 template<> inline cholmod_sparse* cm_spsolve<SuiteSparse_long> (int sys, cholmod_factor& L, cholmod_sparse& B, cholmod_common &Common) { return cholmod_l_spsolve (sys, &L, &B, &Common); }
193 
194 template<typename StorageIndex_>
195 inline int cm_factorize_p (cholmod_sparse* A, double beta[2], StorageIndex_* fset, std::size_t fsize, cholmod_factor* L, cholmod_common &Common) { return cholmod_factorize_p (A, beta, fset, fsize, L, &Common); }
196 template<>
197 inline int cm_factorize_p<SuiteSparse_long> (cholmod_sparse* A, double beta[2], SuiteSparse_long* fset, std::size_t fsize, cholmod_factor* L, cholmod_common &Common) { return cholmod_l_factorize_p (A, beta, fset, fsize, L, &Common); }
198 
199 #undef EIGEN_CHOLMOD_SPECIALIZE0
200 #undef EIGEN_CHOLMOD_SPECIALIZE1
201 
202 } // namespace internal
203 
204 
207 };
208 
209 
215 template<typename MatrixType_, int UpLo_, typename Derived>
216 class CholmodBase : public SparseSolverBase<Derived>
217 {
218  protected:
220  using Base::derived;
221  using Base::m_isInitialized;
222  public:
223  typedef MatrixType_ MatrixType;
224  enum { UpLo = UpLo_ };
225  typedef typename MatrixType::Scalar Scalar;
226  typedef typename MatrixType::RealScalar RealScalar;
228  typedef typename MatrixType::StorageIndex StorageIndex;
229  enum {
230  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
231  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
232  };
233 
234  public:
235 
238  {
239  EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
240  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
241  internal::cm_start<StorageIndex>(m_cholmod);
242  }
243 
244  explicit CholmodBase(const MatrixType& matrix)
246  {
247  EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
248  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
249  internal::cm_start<StorageIndex>(m_cholmod);
250  compute(matrix);
251  }
252 
254  {
255  if(m_cholmodFactor)
256  internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
257  internal::cm_finish<StorageIndex>(m_cholmod);
258  }
259 
260  inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
261  inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
262 
269  {
270  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
271  return m_info;
272  }
273 
275  Derived& compute(const MatrixType& matrix)
276  {
277  analyzePattern(matrix);
278  factorize(matrix);
279  return derived();
280  }
281 
288  void analyzePattern(const MatrixType& matrix)
289  {
290  if(m_cholmodFactor)
291  {
292  internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
293  m_cholmodFactor = 0;
294  }
295  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
296  m_cholmodFactor = internal::cm_analyze<StorageIndex>(A, m_cholmod);
297 
298  this->m_isInitialized = true;
299  this->m_info = Success;
300  m_analysisIsOk = true;
301  m_factorizationIsOk = false;
302  }
303 
310  void factorize(const MatrixType& matrix)
311  {
312  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
313  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
314  internal::cm_factorize_p<StorageIndex>(&A, m_shiftOffset, 0, 0, m_cholmodFactor, m_cholmod);
315 
316  // If the factorization failed, minor is the column at which it did. On success minor == n.
317  this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
318  m_factorizationIsOk = true;
319  }
320 
323  cholmod_common& cholmod() { return m_cholmod; }
324 
325  #ifndef EIGEN_PARSED_BY_DOXYGEN
327  template<typename Rhs,typename Dest>
328  void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
329  {
330  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
331  const Index size = m_cholmodFactor->n;
333  eigen_assert(size==b.rows());
334 
335  // Cholmod needs column-major storage without inner-stride, which corresponds to the default behavior of Ref.
337 
338  cholmod_dense b_cd = viewAsCholmod(b_ref);
339  cholmod_dense* x_cd = internal::cm_solve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cd, m_cholmod);
340  if(!x_cd)
341  {
342  this->m_info = NumericalIssue;
343  return;
344  }
345  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
346  // NOTE Actually, the copy can be avoided by calling cholmod_solve2 instead of cholmod_solve
347  dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
348  internal::cm_free_dense<StorageIndex>(x_cd, m_cholmod);
349  }
350 
352  template<typename RhsDerived, typename DestDerived>
353  void _solve_impl(const SparseMatrixBase<RhsDerived> &b, SparseMatrixBase<DestDerived> &dest) const
354  {
355  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
356  const Index size = m_cholmodFactor->n;
358  eigen_assert(size==b.rows());
359 
360  // note: cs stands for Cholmod Sparse
361  Ref<SparseMatrix<typename RhsDerived::Scalar,ColMajor,typename RhsDerived::StorageIndex> > b_ref(b.const_cast_derived());
362  cholmod_sparse b_cs = viewAsCholmod(b_ref);
363  cholmod_sparse* x_cs = internal::cm_spsolve<StorageIndex>(CHOLMOD_A, *m_cholmodFactor, b_cs, m_cholmod);
364  if(!x_cs)
365  {
366  this->m_info = NumericalIssue;
367  return;
368  }
369  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
370  // NOTE cholmod_spsolve in fact just calls the dense solver for blocks of 4 columns at a time (similar to Eigen's sparse solver)
371  dest.derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
372  internal::cm_free_sparse<StorageIndex>(x_cs, m_cholmod);
373  }
374  #endif // EIGEN_PARSED_BY_DOXYGEN
375 
376 
386  Derived& setShift(const RealScalar& offset)
387  {
388  m_shiftOffset[0] = double(offset);
389  return derived();
390  }
391 
394  {
395  using std::exp;
396  return exp(logDeterminant());
397  }
398 
401  {
402  using std::log;
403  using numext::real;
404  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
405 
406  RealScalar logDet = 0;
407  Scalar *x = static_cast<Scalar*>(m_cholmodFactor->x);
408  if (m_cholmodFactor->is_super)
409  {
410  // Supernodal factorization stored as a packed list of dense column-major blocs,
411  // as described by the following structure:
412 
413  // super[k] == index of the first column of the j-th super node
414  StorageIndex *super = static_cast<StorageIndex*>(m_cholmodFactor->super);
415  // pi[k] == offset to the description of row indices
416  StorageIndex *pi = static_cast<StorageIndex*>(m_cholmodFactor->pi);
417  // px[k] == offset to the respective dense block
418  StorageIndex *px = static_cast<StorageIndex*>(m_cholmodFactor->px);
419 
420  Index nb_super_nodes = m_cholmodFactor->nsuper;
421  for (Index k=0; k < nb_super_nodes; ++k)
422  {
423  StorageIndex ncols = super[k + 1] - super[k];
424  StorageIndex nrows = pi[k + 1] - pi[k];
425 
426  Map<const Array<Scalar,1,Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows+1));
427  logDet += sk.real().log().sum();
428  }
429  }
430  else
431  {
432  // Simplicial factorization stored as standard CSC matrix.
433  StorageIndex *p = static_cast<StorageIndex*>(m_cholmodFactor->p);
435  for (Index k=0; k<size; ++k)
436  logDet += log(real( x[p[k]] ));
437  }
438  if (m_cholmodFactor->is_ll)
439  logDet *= 2.0;
440  return logDet;
441  }
442 
443  template<typename Stream>
444  void dumpMemory(Stream& /*s*/)
445  {}
446 
447  protected:
448  mutable cholmod_common m_cholmod;
449  cholmod_factor* m_cholmodFactor;
450  double m_shiftOffset[2];
454 };
455 
478 template<typename MatrixType_, int UpLo_ = Lower>
479 class CholmodSimplicialLLT : public CholmodBase<MatrixType_, UpLo_, CholmodSimplicialLLT<MatrixType_, UpLo_> >
480 {
482  using Base::m_cholmod;
483 
484  public:
485 
486  typedef MatrixType_ MatrixType;
487 
489 
491  {
492  init();
493  this->compute(matrix);
494  }
495 
497  protected:
498  void init()
499  {
500  m_cholmod.final_asis = 0;
501  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
502  m_cholmod.final_ll = 1;
503  }
504 };
505 
506 
529 template<typename MatrixType_, int UpLo_ = Lower>
530 class CholmodSimplicialLDLT : public CholmodBase<MatrixType_, UpLo_, CholmodSimplicialLDLT<MatrixType_, UpLo_> >
531 {
533  using Base::m_cholmod;
534 
535  public:
536 
537  typedef MatrixType_ MatrixType;
538 
540 
542  {
543  init();
544  this->compute(matrix);
545  }
546 
548  protected:
549  void init()
550  {
551  m_cholmod.final_asis = 1;
552  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
553  }
554 };
555 
578 template<typename MatrixType_, int UpLo_ = Lower>
579 class CholmodSupernodalLLT : public CholmodBase<MatrixType_, UpLo_, CholmodSupernodalLLT<MatrixType_, UpLo_> >
580 {
582  using Base::m_cholmod;
583 
584  public:
585 
586  typedef MatrixType_ MatrixType;
587 
589 
591  {
592  init();
593  this->compute(matrix);
594  }
595 
597  protected:
598  void init()
599  {
600  m_cholmod.final_asis = 1;
601  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
602  }
603 };
604 
629 template<typename MatrixType_, int UpLo_ = Lower>
630 class CholmodDecomposition : public CholmodBase<MatrixType_, UpLo_, CholmodDecomposition<MatrixType_, UpLo_> >
631 {
633  using Base::m_cholmod;
634 
635  public:
636 
637  typedef MatrixType_ MatrixType;
638 
640 
642  {
643  init();
644  this->compute(matrix);
645  }
646 
648 
649  void setMode(CholmodMode mode)
650  {
651  switch(mode)
652  {
653  case CholmodAuto:
654  m_cholmod.final_asis = 1;
655  m_cholmod.supernodal = CHOLMOD_AUTO;
656  break;
658  m_cholmod.final_asis = 0;
659  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
660  m_cholmod.final_ll = 1;
661  break;
663  m_cholmod.final_asis = 1;
664  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
665  break;
666  case CholmodLDLt:
667  m_cholmod.final_asis = 1;
668  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
669  break;
670  default:
671  break;
672  }
673  }
674  protected:
675  void init()
676  {
677  m_cholmod.final_asis = 1;
678  m_cholmod.supernodal = CHOLMOD_AUTO;
679  }
680 };
681 
682 } // end namespace Eigen
683 
684 #endif // EIGEN_CHOLMODSUPPORT_H
const ExpReturnType exp() const
const LogReturnType log() const
Array< int, 3, 1 > b
#define EIGEN_CHOLMOD_SPECIALIZE1(ret, name, t1, a1)
#define EIGEN_CHOLMOD_SPECIALIZE0(ret, name)
RealReturnType real() const
MatrixXcf A
MatrixXf B
MatrixXd L
Definition: LLT_example.cpp:6
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:957
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
float * p
The base class for the direct Cholesky factorization of Cholmod.
MatrixType CholMatrixType
void factorize(const MatrixType &matrix)
ComputationInfo info() const
Reports whether previous computation was successful.
cholmod_common & cholmod()
Scalar determinant() const
SparseSolverBase< Derived > Base
CholmodBase(const MatrixType &matrix)
MatrixType::Scalar Scalar
MatrixType::RealScalar RealScalar
StorageIndex cols() const
Derived & compute(const MatrixType &matrix)
ComputationInfo m_info
Scalar logDeterminant() const
StorageIndex rows() const
void dumpMemory(Stream &)
cholmod_factor * m_cholmodFactor
Derived & setShift(const RealScalar &offset)
cholmod_common m_cholmod
MatrixType::StorageIndex StorageIndex
MatrixType_ MatrixType
void analyzePattern(const MatrixType &matrix)
A general Cholesky factorization and solver based on Cholmod.
CholmodDecomposition(const MatrixType &matrix)
CholmodBase< MatrixType_, UpLo_, CholmodDecomposition > Base
void setMode(CholmodMode mode)
A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod.
CholmodBase< MatrixType_, UpLo_, CholmodSimplicialLDLT > Base
CholmodSimplicialLDLT(const MatrixType &matrix)
A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.
CholmodSimplicialLLT(const MatrixType &matrix)
CholmodBase< MatrixType_, UpLo_, CholmodSimplicialLLT > Base
A supernodal Cholesky (LLT) factorization and solver based on Cholmod.
CholmodBase< MatrixType_, UpLo_, CholmodSupernodalLLT > Base
CholmodSupernodalLLT(const MatrixType &matrix)
Convenience specialization of Stride to specify only an inner stride See class Map for some examples.
Definition: Stride.h:102
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:98
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
MatrixBase< Derived > & matrix()
Definition: MatrixBase.h:319
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:285
A versatible sparse matrix representation.
Definition: SparseMatrix.h:125
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
A base class for sparse solvers.
a sparse vector class
Definition: SparseVector.h:68
ComputationInfo
Definition: Constants.h:444
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
const unsigned int RowMajorBit
Definition: Constants.h:68
cholmod_sparse * cm_spsolve(int sys, cholmod_factor &L, cholmod_sparse &B, cholmod_common &Common)
cholmod_sparse * cm_spsolve< SuiteSparse_long >(int sys, cholmod_factor &L, cholmod_sparse &B, cholmod_common &Common)
cholmod_dense * cm_solve< SuiteSparse_long >(int sys, cholmod_factor &L, cholmod_dense &B, cholmod_common &Common)
int cm_factorize_p< SuiteSparse_long >(cholmod_sparse *A, double beta[2], SuiteSparse_long *fset, std::size_t fsize, cholmod_factor *L, cholmod_common &Common)
cholmod_dense * cm_solve(int sys, cholmod_factor &L, cholmod_dense &B, cholmod_common &Common)
int cm_factorize_p(cholmod_sparse *A, double beta[2], StorageIndex_ *fset, std::size_t fsize, cholmod_factor *L, cholmod_common &Common)
: InteropHeaders
Definition: Core:139
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)
@ CholmodSimplicialLLt
@ CholmodSupernodalLLt
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
Map< SparseMatrix< Scalar, Flags, StorageIndex > > viewAsEigen(cholmod_sparse &cm)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_exp_op< typename Derived::Scalar >, const Derived > exp(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)
Definition: BFloat16.h:222
Derived & const_cast_derived() const
Definition: EigenBase.h:54
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231