Eigen::CholmodBase< MatrixType_, UpLo_, Derived > Class Template Reference

The base class for the direct Cholesky factorization of Cholmod. More...

+ Inheritance diagram for Eigen::CholmodBase< MatrixType_, UpLo_, Derived >:

Public Types

enum  { UpLo }
 
enum  {
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef MatrixType CholMatrixType
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

void analyzePattern (const MatrixType &matrix)
 
cholmod_common & cholmod ()
 
 CholmodBase ()
 
 CholmodBase (const MatrixType &matrix)
 
StorageIndex cols () const
 
Derived & compute (const MatrixType &matrix)
 
Scalar determinant () const
 
template<typename Stream >
void dumpMemory (Stream &)
 
void factorize (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
Scalar logDeterminant () const
 
StorageIndex rows () const
 
Derived & setShift (const RealScalar &offset)
 
 ~CholmodBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Protected Types

typedef SparseSolverBase< Derived > Base
 

Protected Member Functions

Derived & derived ()
 
const Derived & derived () const
 

Protected Attributes

int m_analysisIsOk
 
cholmod_common m_cholmod
 
cholmod_factor * m_cholmodFactor
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
double m_shiftOffset [2]
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Detailed Description

template<typename MatrixType_, int UpLo_, typename Derived>
class Eigen::CholmodBase< MatrixType_, UpLo_, Derived >

The base class for the direct Cholesky factorization of Cholmod.

See also
class CholmodSupernodalLLT, class CholmodSimplicialLDLT, class CholmodSimplicialLLT

Definition at line 216 of file CholmodSupport.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef SparseSolverBase<Derived> Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::Base
protected

Definition at line 219 of file CholmodSupport.h.

◆ CholMatrixType

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholMatrixType

Definition at line 227 of file CholmodSupport.h.

◆ MatrixType

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType_ Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::MatrixType

Definition at line 223 of file CholmodSupport.h.

◆ RealScalar

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::RealScalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::RealScalar

Definition at line 226 of file CholmodSupport.h.

◆ Scalar

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::Scalar

Definition at line 225 of file CholmodSupport.h.

◆ StorageIndex

template<typename MatrixType_ , int UpLo_, typename Derived >
typedef MatrixType::StorageIndex Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::StorageIndex

Definition at line 228 of file CholmodSupport.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Derived >
anonymous enum
Enumerator
UpLo 

Definition at line 224 of file CholmodSupport.h.

224 { UpLo = UpLo_ };

◆ anonymous enum

template<typename MatrixType_ , int UpLo_, typename Derived >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 229 of file CholmodSupport.h.

229  {
230  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
231  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
232  };

Constructor & Destructor Documentation

◆ CholmodBase() [1/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholmodBase ( )
inline

Definition at line 236 of file CholmodSupport.h.

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  }
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
ComputationInfo m_info
cholmod_factor * m_cholmodFactor
cholmod_common m_cholmod
@ Success
Definition: Constants.h:446

◆ CholmodBase() [2/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::CholmodBase ( const MatrixType matrix)
inlineexplicit

Definition at line 244 of file CholmodSupport.h.

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  }
Derived & compute(const MatrixType &matrix)

◆ ~CholmodBase()

template<typename MatrixType_ , int UpLo_, typename Derived >
Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::~CholmodBase ( )
inline

Definition at line 253 of file CholmodSupport.h.

254  {
255  if(m_cholmodFactor)
256  internal::cm_free_factor<StorageIndex>(m_cholmodFactor, m_cholmod);
257  internal::cm_finish<StorageIndex>(m_cholmod);
258  }

Member Function Documentation

◆ analyzePattern()

template<typename MatrixType_ , int UpLo_, typename Derived >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::analyzePattern ( const MatrixType matrix)
inline

Performs a symbolic decomposition on the sparsity pattern of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()

Definition at line 288 of file CholmodSupport.h.

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  }
MatrixXcf A
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)

◆ cholmod()

template<typename MatrixType_ , int UpLo_, typename Derived >
cholmod_common& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::cholmod ( )
inline

Returns a reference to the Cholmod's configuration structure to get a full control over the performed operations. See the Cholmod user guide for details.

Definition at line 323 of file CholmodSupport.h.

323 { return m_cholmod; }

◆ cols()

template<typename MatrixType_ , int UpLo_, typename Derived >
StorageIndex Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::cols ( ) const
inline

Definition at line 260 of file CholmodSupport.h.

260 { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }

◆ compute()

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::compute ( const MatrixType matrix)
inline

Computes the sparse Cholesky decomposition of matrix

Definition at line 275 of file CholmodSupport.h.

276  {
277  analyzePattern(matrix);
278  factorize(matrix);
279  return derived();
280  }
void factorize(const MatrixType &matrix)
void analyzePattern(const MatrixType &matrix)

◆ derived() [1/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

Definition at line 83 of file SparseSolverBase.h.

83 { return *static_cast<Derived*>(this); }

◆ derived() [2/2]

template<typename MatrixType_ , int UpLo_, typename Derived >
const Derived& Eigen::SparseSolverBase< Derived >::derived
inlineprotected

Definition at line 84 of file SparseSolverBase.h.

84 { return *static_cast<const Derived*>(this); }

◆ determinant()

template<typename MatrixType_ , int UpLo_, typename Derived >
Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::determinant ( ) const
inline
Returns
the determinant of the underlying matrix from the current factorization

Definition at line 393 of file CholmodSupport.h.

394  {
395  using std::exp;
396  return exp(logDeterminant());
397  }
const ExpReturnType exp() const
Scalar logDeterminant() const
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_exp_op< typename Derived::Scalar >, const Derived > exp(const Eigen::ArrayBase< Derived > &x)

◆ dumpMemory()

template<typename MatrixType_ , int UpLo_, typename Derived >
template<typename Stream >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::dumpMemory ( Stream &  )
inline

Definition at line 444 of file CholmodSupport.h.

445  {}

◆ factorize()

template<typename MatrixType_ , int UpLo_, typename Derived >
void Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::factorize ( const MatrixType matrix)
inline

Performs a numeric decomposition of matrix

The given matrix must have the same sparsity pattern as the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()

Definition at line 310 of file CholmodSupport.h.

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  }
#define eigen_assert(x)
Definition: Macros.h:902
@ NumericalIssue
Definition: Constants.h:448

◆ info()

template<typename MatrixType_ , int UpLo_, typename Derived >
ComputationInfo Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the matrix.appears to be negative.

Definition at line 268 of file CholmodSupport.h.

269  {
270  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
271  return m_info;
272  }

◆ logDeterminant()

template<typename MatrixType_ , int UpLo_, typename Derived >
Scalar Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::logDeterminant ( ) const
inline
Returns
the log determinant of the underlying matrix from the current factorization

Definition at line 400 of file CholmodSupport.h.

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  }
const LogReturnType log() const
RealReturnType real() const
float * p
MatrixType::Scalar Scalar
MatrixType::RealScalar RealScalar
MatrixType::StorageIndex StorageIndex
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_log_op< typename Derived::Scalar >, const Derived > log(const Eigen::ArrayBase< Derived > &x)

◆ rows()

template<typename MatrixType_ , int UpLo_, typename Derived >
StorageIndex Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::rows ( ) const
inline

Definition at line 261 of file CholmodSupport.h.

261 { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }

◆ setShift()

template<typename MatrixType_ , int UpLo_, typename Derived >
Derived& Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::setShift ( const RealScalar offset)
inline

Sets the shift parameter that will be used to adjust the diagonal coefficients during the numerical factorization.

During the numerical factorization, an offset term is added to the diagonal coefficients:
d_ii = offset + d_ii

The default is offset=0.

Returns
a reference to *this.

Definition at line 386 of file CholmodSupport.h.

387  {
388  m_shiftOffset[0] = double(offset);
389  return derived();
390  }

Member Data Documentation

◆ m_analysisIsOk

template<typename MatrixType_ , int UpLo_, typename Derived >
int Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_analysisIsOk
protected

Definition at line 453 of file CholmodSupport.h.

◆ m_cholmod

template<typename MatrixType_ , int UpLo_, typename Derived >
cholmod_common Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmod
mutableprotected

Definition at line 448 of file CholmodSupport.h.

◆ m_cholmodFactor

template<typename MatrixType_ , int UpLo_, typename Derived >
cholmod_factor* Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_cholmodFactor
protected

Definition at line 449 of file CholmodSupport.h.

◆ m_factorizationIsOk

template<typename MatrixType_ , int UpLo_, typename Derived >
int Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_factorizationIsOk
protected

Definition at line 452 of file CholmodSupport.h.

◆ m_info

template<typename MatrixType_ , int UpLo_, typename Derived >
ComputationInfo Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_info
mutableprotected

Definition at line 451 of file CholmodSupport.h.

◆ m_isInitialized

template<typename MatrixType_ , int UpLo_, typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

Definition at line 123 of file SparseSolverBase.h.

◆ m_shiftOffset

template<typename MatrixType_ , int UpLo_, typename Derived >
double Eigen::CholmodBase< MatrixType_, UpLo_, Derived >::m_shiftOffset[2]
protected

Definition at line 450 of file CholmodSupport.h.


The documentation for this class was generated from the following file: