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 namespace Eigen {
14 
15 namespace internal {
16 
17 template<typename Scalar> struct cholmod_configure_matrix;
18 
19 template<> struct cholmod_configure_matrix<double> {
20  template<typename CholmodType>
21  static void run(CholmodType& mat) {
22  mat.xtype = CHOLMOD_REAL;
23  mat.dtype = CHOLMOD_DOUBLE;
24  }
25 };
26 
27 template<> struct cholmod_configure_matrix<std::complex<double> > {
28  template<typename CholmodType>
29  static void run(CholmodType& mat) {
30  mat.xtype = CHOLMOD_COMPLEX;
31  mat.dtype = CHOLMOD_DOUBLE;
32  }
33 };
34 
35 // Other scalar types are not yet suppotred by Cholmod
36 // template<> struct cholmod_configure_matrix<float> {
37 // template<typename CholmodType>
38 // static void run(CholmodType& mat) {
39 // mat.xtype = CHOLMOD_REAL;
40 // mat.dtype = CHOLMOD_SINGLE;
41 // }
42 // };
43 //
44 // template<> struct cholmod_configure_matrix<std::complex<float> > {
45 // template<typename CholmodType>
46 // static void run(CholmodType& mat) {
47 // mat.xtype = CHOLMOD_COMPLEX;
48 // mat.dtype = CHOLMOD_SINGLE;
49 // }
50 // };
51 
52 } // namespace internal
53 
57 template<typename _Scalar, int _Options, typename _StorageIndex>
59 {
60  cholmod_sparse res;
61  res.nzmax = mat.nonZeros();
62  res.nrow = mat.rows();
63  res.ncol = mat.cols();
64  res.p = mat.outerIndexPtr();
65  res.i = mat.innerIndexPtr();
66  res.x = mat.valuePtr();
67  res.z = 0;
68  res.sorted = 1;
69  if(mat.isCompressed())
70  {
71  res.packed = 1;
72  res.nz = 0;
73  }
74  else
75  {
76  res.packed = 0;
77  res.nz = mat.innerNonZeroPtr();
78  }
79 
80  res.dtype = 0;
81  res.stype = -1;
82 
84  {
85  res.itype = CHOLMOD_INT;
86  }
88  {
89  res.itype = CHOLMOD_LONG;
90  }
91  else
92  {
93  eigen_assert(false && "Index type not supported yet");
94  }
95 
96  // setup res.xtype
98 
99  res.stype = 0;
100 
101  return res;
102 }
103 
104 template<typename _Scalar, int _Options, typename _Index>
106 {
107  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.const_cast_derived()));
108  return res;
109 }
110 
111 template<typename _Scalar, int _Options, typename _Index>
113 {
114  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.const_cast_derived()));
115  return res;
116 }
117 
120 template<typename _Scalar, int _Options, typename _Index, unsigned int UpLo>
122 {
123  cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.matrix().const_cast_derived()));
124 
125  if(UpLo==Upper) res.stype = 1;
126  if(UpLo==Lower) res.stype = -1;
127 
128  return res;
129 }
130 
133 template<typename Derived>
135 {
136  EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
137  typedef typename Derived::Scalar Scalar;
138 
139  cholmod_dense res;
140  res.nrow = mat.rows();
141  res.ncol = mat.cols();
142  res.nzmax = res.nrow * res.ncol;
143  res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
144  res.x = (void*)(mat.derived().data());
145  res.z = 0;
146 
148 
149  return res;
150 }
151 
154 template<typename Scalar, int Flags, typename StorageIndex>
156 {
158  (cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol],
159  static_cast<StorageIndex*>(cm.p), static_cast<StorageIndex*>(cm.i),static_cast<Scalar*>(cm.x) );
160 }
161 
164 };
165 
166 
172 template<typename _MatrixType, int _UpLo, typename Derived>
173 class CholmodBase : public SparseSolverBase<Derived>
174 {
175  protected:
177  using Base::derived;
178  using Base::m_isInitialized;
179  public:
180  typedef _MatrixType MatrixType;
181  enum { UpLo = _UpLo };
182  typedef typename MatrixType::Scalar Scalar;
185  typedef typename MatrixType::StorageIndex StorageIndex;
186  enum {
187  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
188  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
189  };
190 
191  public:
192 
195  {
196  EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
197  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
198  cholmod_start(&m_cholmod);
199  }
200 
201  explicit CholmodBase(const MatrixType& matrix)
203  {
204  EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
205  m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
206  cholmod_start(&m_cholmod);
207  compute(matrix);
208  }
209 
211  {
212  if(m_cholmodFactor)
213  cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
214  cholmod_finish(&m_cholmod);
215  }
216 
217  inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
218  inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
219 
226  {
227  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
228  return m_info;
229  }
230 
232  Derived& compute(const MatrixType& matrix)
233  {
235  factorize(matrix);
236  return derived();
237  }
238 
246  {
247  if(m_cholmodFactor)
248  {
249  cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
250  m_cholmodFactor = 0;
251  }
252  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
253  m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
254 
255  this->m_isInitialized = true;
256  this->m_info = Success;
257  m_analysisIsOk = true;
258  m_factorizationIsOk = false;
259  }
260 
267  void factorize(const MatrixType& matrix)
268  {
269  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
270  cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
271  cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod);
272 
273  // If the factorization failed, minor is the column at which it did. On success minor == n.
274  this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
275  m_factorizationIsOk = true;
276  }
277 
280  cholmod_common& cholmod() { return m_cholmod; }
281 
282  #ifndef EIGEN_PARSED_BY_DOXYGEN
283 
284  template<typename Rhs,typename Dest>
285  void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
286  {
287  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
288  const Index size = m_cholmodFactor->n;
290  eigen_assert(size==b.rows());
291 
292  // Cholmod needs column-major stoarge without inner-stride, which corresponds to the default behavior of Ref.
294 
295  cholmod_dense b_cd = viewAsCholmod(b_ref);
296  cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
297  if(!x_cd)
298  {
299  this->m_info = NumericalIssue;
300  return;
301  }
302  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
303  dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
304  cholmod_free_dense(&x_cd, &m_cholmod);
305  }
306 
308  template<typename RhsDerived, typename DestDerived>
310  {
311  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
312  const Index size = m_cholmodFactor->n;
314  eigen_assert(size==b.rows());
315 
316  // note: cs stands for Cholmod Sparse
318  cholmod_sparse b_cs = viewAsCholmod(b_ref);
319  cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
320  if(!x_cs)
321  {
322  this->m_info = NumericalIssue;
323  return;
324  }
325  // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
326  dest.derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
327  cholmod_free_sparse(&x_cs, &m_cholmod);
328  }
329  #endif // EIGEN_PARSED_BY_DOXYGEN
330 
331 
341  Derived& setShift(const RealScalar& offset)
342  {
343  m_shiftOffset[0] = double(offset);
344  return derived();
345  }
346 
349  {
350  using std::exp;
351  return exp(logDeterminant());
352  }
353 
356  {
357  using std::log;
358  using numext::real;
359  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
360 
361  RealScalar logDet = 0;
362  Scalar *x = static_cast<Scalar*>(m_cholmodFactor->x);
363  if (m_cholmodFactor->is_super)
364  {
365  // Supernodal factorization stored as a packed list of dense column-major blocs,
366  // as described by the following structure:
367 
368  // super[k] == index of the first column of the j-th super node
369  StorageIndex *super = static_cast<StorageIndex*>(m_cholmodFactor->super);
370  // pi[k] == offset to the description of row indices
371  StorageIndex *pi = static_cast<StorageIndex*>(m_cholmodFactor->pi);
372  // px[k] == offset to the respective dense block
373  StorageIndex *px = static_cast<StorageIndex*>(m_cholmodFactor->px);
374 
375  Index nb_super_nodes = m_cholmodFactor->nsuper;
376  for (Index k=0; k < nb_super_nodes; ++k)
377  {
378  StorageIndex ncols = super[k + 1] - super[k];
379  StorageIndex nrows = pi[k + 1] - pi[k];
380 
381  Map<const Array<Scalar,1,Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows+1));
382  logDet += sk.real().log().sum();
383  }
384  }
385  else
386  {
387  // Simplicial factorization stored as standard CSC matrix.
388  StorageIndex *p = static_cast<StorageIndex*>(m_cholmodFactor->p);
390  for (Index k=0; k<size; ++k)
391  logDet += log(real( x[p[k]] ));
392  }
393  if (m_cholmodFactor->is_ll)
394  logDet *= 2.0;
395  return logDet;
396  };
397 
398  template<typename Stream>
399  void dumpMemory(Stream& /*s*/)
400  {}
401 
402  protected:
403  mutable cholmod_common m_cholmod;
404  cholmod_factor* m_cholmodFactor;
405  double m_shiftOffset[2];
409 };
410 
433 template<typename _MatrixType, int _UpLo = Lower>
434 class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT<_MatrixType, _UpLo> >
435 {
437  using Base::m_cholmod;
438 
439  public:
440 
441  typedef _MatrixType MatrixType;
442 
444 
446  {
447  init();
448  this->compute(matrix);
449  }
450 
452  protected:
453  void init()
454  {
455  m_cholmod.final_asis = 0;
456  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
457  m_cholmod.final_ll = 1;
458  }
459 };
460 
461 
484 template<typename _MatrixType, int _UpLo = Lower>
485 class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT<_MatrixType, _UpLo> >
486 {
488  using Base::m_cholmod;
489 
490  public:
491 
492  typedef _MatrixType MatrixType;
493 
495 
497  {
498  init();
499  this->compute(matrix);
500  }
501 
503  protected:
504  void init()
505  {
506  m_cholmod.final_asis = 1;
507  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
508  }
509 };
510 
533 template<typename _MatrixType, int _UpLo = Lower>
534 class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT<_MatrixType, _UpLo> >
535 {
537  using Base::m_cholmod;
538 
539  public:
540 
541  typedef _MatrixType MatrixType;
542 
544 
546  {
547  init();
548  this->compute(matrix);
549  }
550 
552  protected:
553  void init()
554  {
555  m_cholmod.final_asis = 1;
556  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
557  }
558 };
559 
584 template<typename _MatrixType, int _UpLo = Lower>
585 class CholmodDecomposition : public CholmodBase<_MatrixType, _UpLo, CholmodDecomposition<_MatrixType, _UpLo> >
586 {
588  using Base::m_cholmod;
589 
590  public:
591 
592  typedef _MatrixType MatrixType;
593 
595 
597  {
598  init();
599  this->compute(matrix);
600  }
601 
603 
604  void setMode(CholmodMode mode)
605  {
606  switch(mode)
607  {
608  case CholmodAuto:
609  m_cholmod.final_asis = 1;
610  m_cholmod.supernodal = CHOLMOD_AUTO;
611  break;
613  m_cholmod.final_asis = 0;
614  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
615  m_cholmod.final_ll = 1;
616  break;
618  m_cholmod.final_asis = 1;
619  m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
620  break;
621  case CholmodLDLt:
622  m_cholmod.final_asis = 1;
623  m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
624  break;
625  default:
626  break;
627  }
628  }
629  protected:
630  void init()
631  {
632  m_cholmod.final_asis = 1;
633  m_cholmod.supernodal = CHOLMOD_AUTO;
634  }
635 };
636 
637 } // end namespace Eigen
638 
639 #endif // EIGEN_CHOLMODSUPPORT_H
Eigen::CholmodBase
The base class for the direct Cholesky factorization of Cholmod.
Definition: CholmodSupport.h:173
Eigen::CholmodBase::m_factorizationIsOk
int m_factorizationIsOk
Definition: CholmodSupport.h:407
Eigen::viewAsEigen
MappedSparseMatrix< Scalar, Flags, StorageIndex > viewAsEigen(cholmod_sparse &cm)
Definition: CholmodSupport.h:155
Eigen::CholmodBase::factorize
void factorize(const MatrixType &matrix)
Definition: CholmodSupport.h:267
Eigen::NumericalIssue
@ NumericalIssue
Definition: Constants.h:434
Eigen::CholmodSimplicialLLT::init
void init()
Definition: CholmodSupport.h:453
Eigen::CholmodDecomposition::Base
CholmodBase< _MatrixType, _UpLo, CholmodDecomposition > Base
Definition: CholmodSupport.h:587
Eigen::CholmodSupernodalLLT::~CholmodSupernodalLLT
~CholmodSupernodalLLT()
Definition: CholmodSupport.h:551
Eigen::internal::cholmod_configure_matrix< std::complex< double > >::run
static void run(CholmodType &mat)
Definition: CholmodSupport.h:29
Eigen
Definition: common.h:73
Eigen::CholmodBase::UpLo
@ UpLo
Definition: CholmodSupport.h:181
Eigen::SparseMatrix
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
Eigen::CholmodSimplicialLLT::~CholmodSimplicialLLT
~CholmodSimplicialLLT()
Definition: CholmodSupport.h:451
b
Scalar * b
Definition: cholesky.cpp:56
Eigen::CholmodDecomposition::setMode
void setMode(CholmodMode mode)
Definition: CholmodSupport.h:604
Eigen::CholmodBase::CholmodBase
CholmodBase()
Definition: CholmodSupport.h:193
Eigen::CholmodBase::CholmodBase
CholmodBase(const MatrixType &matrix)
Definition: CholmodSupport.h:201
Eigen::CholmodDecomposition::init
void init()
Definition: CholmodSupport.h:630
Eigen::CholmodSimplicialLLT::CholmodSimplicialLLT
CholmodSimplicialLLT()
Definition: CholmodSupport.h:443
Eigen::CholmodBase::_solve_impl
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition: CholmodSupport.h:285
Eigen::InnerStride
Convenience specialization of Stride to specify only an inner stride See class Map for some examples.
Definition: Stride.h:90
RealScalar
NumTraits< Scalar >::Real RealScalar
Definition: common.h:85
Eigen::CholmodBase::~CholmodBase
~CholmodBase()
Definition: CholmodSupport.h:210
eigen_assert
#define eigen_assert(x)
Definition: Macros.h:579
Eigen::CholmodSupernodalLLT::init
void init()
Definition: CholmodSupport.h:553
Eigen::RowMajorBit
const unsigned int RowMajorBit
Definition: Constants.h:61
Eigen::Upper
@ Upper
Definition: Constants.h:206
Eigen::Success
@ Success
Definition: Constants.h:432
real
float real
Definition: datatypes.h:10
Eigen::CholmodSimplicialLDLT::Base
CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLDLT > Base
Definition: CholmodSupport.h:487
Eigen::CholmodSimplicialLDLT::~CholmodSimplicialLDLT
~CholmodSimplicialLDLT()
Definition: CholmodSupport.h:502
log
const EIGEN_DEVICE_FUNC LogReturnType log() const
Definition: ArrayCwiseUnaryOps.h:105
Eigen::CholmodDecomposition::CholmodDecomposition
CholmodDecomposition(const MatrixType &matrix)
Definition: CholmodSupport.h:596
Scalar
SCALAR Scalar
Definition: common.h:84
exp
const EIGEN_DEVICE_FUNC ExpReturnType exp() const
Definition: ArrayCwiseUnaryOps.h:88
Eigen::CholmodSimplicialLDLT
A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:485
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:102
Eigen::CholmodBase::Base
SparseSolverBase< Derived > Base
Definition: CholmodSupport.h:176
Eigen::CholmodBase::logDeterminant
Scalar logDeterminant() const
Definition: CholmodSupport.h:355
Eigen::CholmodBase::MatrixType
_MatrixType MatrixType
Definition: CholmodSupport.h:180
Eigen::viewAsCholmod
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< _Scalar, _Options, _StorageIndex > > mat)
Definition: CholmodSupport.h:58
Eigen::CholmodDecomposition::MatrixType
_MatrixType MatrixType
Definition: CholmodSupport.h:592
Eigen::CholmodBase::analyzePattern
void analyzePattern(const MatrixType &matrix)
Definition: CholmodSupport.h:245
Eigen::CholmodBase::cholmod
cholmod_common & cholmod()
Definition: CholmodSupport.h:280
Eigen::CholmodBase::m_analysisIsOk
int m_analysisIsOk
Definition: CholmodSupport.h:408
EIGEN_UNUSED_VARIABLE
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:618
Eigen::CholmodSupernodalLLT::MatrixType
_MatrixType MatrixType
Definition: CholmodSupport.h:541
Eigen::CholmodBase::CholMatrixType
MatrixType CholMatrixType
Definition: CholmodSupport.h:184
Eigen::CholmodAuto
@ CholmodAuto
Definition: CholmodSupport.h:163
Eigen::CholmodMode
CholmodMode
Definition: CholmodSupport.h:162
Eigen::internal::cholmod_configure_matrix
Definition: CholmodSupport.h:17
Eigen::CholmodSupernodalLLT::CholmodSupernodalLLT
CholmodSupernodalLLT()
Definition: CholmodSupport.h:543
x
Scalar * x
Definition: level1_cplx_impl.h:89
Eigen::CholmodDecomposition::~CholmodDecomposition
~CholmodDecomposition()
Definition: CholmodSupport.h:602
Eigen::CholmodBase::determinant
Scalar determinant() const
Definition: CholmodSupport.h:348
Eigen::CholmodBase::MaxColsAtCompileTime
@ MaxColsAtCompileTime
Definition: CholmodSupport.h:188
Eigen::CholmodBase::info
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: CholmodSupport.h:225
Eigen::CholmodSimplicialLDLT::CholmodSimplicialLDLT
CholmodSimplicialLDLT(const MatrixType &matrix)
Definition: CholmodSupport.h:496
Eigen::internal::cholmod_configure_matrix< double >::run
static void run(CholmodType &mat)
Definition: CholmodSupport.h:21
Eigen::CholmodBase::m_info
ComputationInfo m_info
Definition: CholmodSupport.h:406
Eigen::CholmodSimplicialLDLT::init
void init()
Definition: CholmodSupport.h:504
Eigen::CholmodSimplicialLDLT::MatrixType
_MatrixType MatrixType
Definition: CholmodSupport.h:492
Eigen::CholmodBase::dumpMemory
void dumpMemory(Stream &)
Definition: CholmodSupport.h:399
Eigen::CholmodSimplicialLDLT::CholmodSimplicialLDLT
CholmodSimplicialLDLT()
Definition: CholmodSupport.h:494
Eigen::Lower
@ Lower
Definition: Constants.h:204
Eigen::CholmodSupernodalLLT
A supernodal Cholesky (LLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:534
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
Eigen::CholmodBase::StorageIndex
MatrixType::StorageIndex StorageIndex
Definition: CholmodSupport.h:185
Eigen::PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::Map
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:587
Eigen::CholmodSimplicialLLT::Base
CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLLT > Base
Definition: CholmodSupport.h:436
Eigen::CholmodBase::m_cholmod
cholmod_common m_cholmod
Definition: CholmodSupport.h:403
Eigen::CholmodSimplicialLLT::CholmodSimplicialLLT
CholmodSimplicialLLT(const MatrixType &matrix)
Definition: CholmodSupport.h:445
Eigen::real
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:545
Eigen::CholmodBase::Scalar
MatrixType::Scalar Scalar
Definition: CholmodSupport.h:182
Eigen::SparseSolverBase
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
Eigen::Ref
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:192
Eigen::CholmodBase::RealScalar
MatrixType::RealScalar RealScalar
Definition: CholmodSupport.h:183
Eigen::internal::traits
Definition: ForwardDeclarations.h:17
EIGEN_STATIC_ASSERT
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:124
Eigen::SparseSolverBase::derived
Derived & derived()
Definition: SparseSolverBase.h:79
Eigen::SparseVector
a sparse vector class
Definition: SparseUtil.h:54
std
Definition: Half.h:150
Eigen::CholmodBase::rows
StorageIndex rows() const
Definition: CholmodSupport.h:218
Eigen::TensorSycl::run
void run(Expr &expr, Dev &dev)
Definition: TensorSyclRun.h:47
Eigen::SparseMatrixBase
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:281
Eigen::CholmodSimplicialLLT
A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.
Definition: CholmodSupport.h:434
Eigen::SparseMatrixBase::derived
const Derived & derived() const
Definition: SparseMatrixBase.h:138
utility::tuple::size
static constexpr size_t size(Tuple< Args... > &)
Provides access to the number of elements in a tuple as a compile-time constant expression.
Definition: TensorSyclTuple.h:143
Eigen::internal::is_same
Definition: Meta.h:63
Eigen::CholmodSupernodalLLt
@ CholmodSupernodalLLt
Definition: CholmodSupport.h:163
Eigen::SparseSelfAdjointView
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:43
internal
Definition: BandTriangularSolver.h:13
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Eigen::CholmodSimplicialLLt
@ CholmodSimplicialLLt
Definition: CholmodSupport.h:163
A
MatrixType A(a, *n, *n, *lda)
Eigen::CholmodBase::_solve_impl
void _solve_impl(const SparseMatrixBase< RhsDerived > &b, SparseMatrixBase< DestDerived > &dest) const
Definition: CholmodSupport.h:309
Eigen::ComputationInfo
ComputationInfo
Definition: Constants.h:430
Eigen::CholmodBase::ColsAtCompileTime
@ ColsAtCompileTime
Definition: CholmodSupport.h:187
Eigen::CholmodBase::compute
Derived & compute(const MatrixType &matrix)
Definition: CholmodSupport.h:232
Eigen::CholmodLDLt
@ CholmodLDLt
Definition: CholmodSupport.h:163
Eigen::CholmodSupernodalLLT::Base
CholmodBase< _MatrixType, _UpLo, CholmodSupernodalLLT > Base
Definition: CholmodSupport.h:536
Eigen::CholmodBase::cols
StorageIndex cols() const
Definition: CholmodSupport.h:217
Eigen::CholmodBase::m_cholmodFactor
cholmod_factor * m_cholmodFactor
Definition: CholmodSupport.h:404
Eigen::CholmodDecomposition
A general Cholesky factorization and solver based on Cholmod.
Definition: CholmodSupport.h:585
Eigen::CholmodSimplicialLLT::MatrixType
_MatrixType MatrixType
Definition: CholmodSupport.h:441
Eigen::MappedSparseMatrix
Sparse matrix.
Definition: MappedSparseMatrix.h:32
Eigen::CholmodBase::m_shiftOffset
double m_shiftOffset[2]
Definition: CholmodSupport.h:405
complex
Definition: datatypes.h:12
mat
else mat
Definition: eigenvalues.cpp:43
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Eigen::CholmodDecomposition::CholmodDecomposition
CholmodDecomposition()
Definition: CholmodSupport.h:594
Eigen::CholmodSupernodalLLT::CholmodSupernodalLLT
CholmodSupernodalLLT(const MatrixType &matrix)
Definition: CholmodSupport.h:545
Eigen::SparseSolverBase::m_isInitialized
bool m_isInitialized
Definition: SparseSolverBase.h:119
px
RealScalar RealScalar * px
Definition: level1_cplx_impl.h:28
Eigen::CholmodBase::setShift
Derived & setShift(const RealScalar &offset)
Definition: CholmodSupport.h:341


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:05:38