MOAB
4.9.3pre
|
The base class for the direct Cholesky factorization of Cholmod. More...
#include <CholmodSupport.h>
Public Types | |
enum | { UpLo = _UpLo } |
enum | { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime } |
typedef _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType | CholMatrixType |
typedef MatrixType::StorageIndex | StorageIndex |
Public Member Functions | |
CholmodBase () | |
CholmodBase (const MatrixType &matrix) | |
~CholmodBase () | |
StorageIndex | cols () const |
StorageIndex | rows () const |
ComputationInfo | info () const |
Reports whether previous computation was successful. | |
Derived & | compute (const MatrixType &matrix) |
void | analyzePattern (const MatrixType &matrix) |
void | factorize (const MatrixType &matrix) |
cholmod_common & | cholmod () |
template<typename Rhs , typename Dest > | |
void | _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const |
template<typename RhsScalar , int RhsOptions, typename RhsIndex , typename DestScalar , int DestOptions, typename DestIndex > | |
void | _solve_impl (const SparseMatrix< RhsScalar, RhsOptions, RhsIndex > &b, SparseMatrix< DestScalar, DestOptions, DestIndex > &dest) const |
Derived & | setShift (const RealScalar &offset) |
Scalar | determinant () const |
Scalar | logDeterminant () const |
template<typename Stream > | |
void | dumpMemory (Stream &) |
Protected Types | |
typedef SparseSolverBase< Derived > | Base |
Protected Attributes | |
cholmod_common | m_cholmod |
cholmod_factor * | m_cholmodFactor |
RealScalar | m_shiftOffset [2] |
ComputationInfo | m_info |
int | m_factorizationIsOk |
int | m_analysisIsOk |
The base class for the direct Cholesky factorization of Cholmod.
Definition at line 160 of file CholmodSupport.h.
typedef SparseSolverBase<Derived> Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::Base [protected] |
Reimplemented in Eigen::CholmodDecomposition< _MatrixType, _UpLo >, Eigen::CholmodSupernodalLLT< _MatrixType, _UpLo >, Eigen::CholmodSimplicialLDLT< _MatrixType, _UpLo >, and Eigen::CholmodSimplicialLLT< _MatrixType, _UpLo >.
Definition at line 163 of file CholmodSupport.h.
typedef MatrixType Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::CholMatrixType |
Definition at line 171 of file CholmodSupport.h.
typedef _MatrixType Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::MatrixType |
Reimplemented in Eigen::CholmodDecomposition< _MatrixType, _UpLo >, Eigen::CholmodSupernodalLLT< _MatrixType, _UpLo >, Eigen::CholmodSimplicialLDLT< _MatrixType, _UpLo >, and Eigen::CholmodSimplicialLLT< _MatrixType, _UpLo >.
Definition at line 167 of file CholmodSupport.h.
typedef MatrixType::RealScalar Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::RealScalar |
Definition at line 170 of file CholmodSupport.h.
typedef MatrixType::Scalar Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::Scalar |
Definition at line 169 of file CholmodSupport.h.
typedef MatrixType::StorageIndex Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::StorageIndex |
Definition at line 172 of file CholmodSupport.h.
anonymous enum |
anonymous enum |
Definition at line 173 of file CholmodSupport.h.
{ ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::CholmodBase | ( | ) | [inline] |
Definition at line 180 of file CholmodSupport.h.
: m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) { m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); cholmod_start(&m_cholmod); }
Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::CholmodBase | ( | const MatrixType & | matrix | ) | [inline, explicit] |
Definition at line 187 of file CholmodSupport.h.
: m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) { m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); cholmod_start(&m_cholmod); compute(matrix); }
Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::~CholmodBase | ( | ) | [inline] |
Definition at line 195 of file CholmodSupport.h.
{ if(m_cholmodFactor) cholmod_free_factor(&m_cholmodFactor, &m_cholmod); cholmod_finish(&m_cholmod); }
void Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::_solve_impl | ( | const MatrixBase< Rhs > & | b, |
MatrixBase< Dest > & | dest | ||
) | const [inline] |
Definition at line 270 of file CholmodSupport.h.
{ eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); const Index size = m_cholmodFactor->n; EIGEN_UNUSED_VARIABLE(size); eigen_assert(size==b.rows()); // Cholmod needs column-major stoarge without inner-stride, which corresponds to the default behavior of Ref. Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b.derived()); cholmod_dense b_cd = viewAsCholmod(b_ref); cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod); if(!x_cd) { this->m_info = NumericalIssue; return; } // TODO optimize this copy by swapping when possible (be careful with alignment, etc.) dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols()); cholmod_free_dense(&x_cd, &m_cholmod); }
void Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::_solve_impl | ( | const SparseMatrix< RhsScalar, RhsOptions, RhsIndex > & | b, |
SparseMatrix< DestScalar, DestOptions, DestIndex > & | dest | ||
) | const [inline] |
Definition at line 294 of file CholmodSupport.h.
{ eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); const Index size = m_cholmodFactor->n; EIGEN_UNUSED_VARIABLE(size); eigen_assert(size==b.rows()); // note: cs stands for Cholmod Sparse cholmod_sparse b_cs = viewAsCholmod(b); cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod); if(!x_cs) { this->m_info = NumericalIssue; return; } // TODO optimize this copy by swapping when possible (be careful with alignment, etc.) dest = viewAsEigen<DestScalar,DestOptions,DestIndex>(*x_cs); cholmod_free_sparse(&x_cs, &m_cholmod); }
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.
Definition at line 230 of file CholmodSupport.h.
{ if(m_cholmodFactor) { cholmod_free_factor(&m_cholmodFactor, &m_cholmod); m_cholmodFactor = 0; } cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>()); m_cholmodFactor = cholmod_analyze(&A, &m_cholmod); this->m_isInitialized = true; this->m_info = Success; m_analysisIsOk = true; m_factorizationIsOk = false; }
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 265 of file CholmodSupport.h.
{ return m_cholmod; }
StorageIndex Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::cols | ( | ) | const [inline] |
Definition at line 202 of file CholmodSupport.h.
{ return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
Derived& Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::compute | ( | const MatrixType & | matrix | ) | [inline] |
Computes the sparse Cholesky decomposition of matrix
Definition at line 217 of file CholmodSupport.h.
{ analyzePattern(matrix); factorize(matrix); return derived(); }
Scalar Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::determinant | ( | ) | const [inline] |
Definition at line 332 of file CholmodSupport.h.
{ using std::exp; return exp(logDeterminant()); }
void Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::dumpMemory | ( | Stream & | ) | [inline] |
Definition at line 383 of file CholmodSupport.h.
{}
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.
Definition at line 252 of file CholmodSupport.h.
{ eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>()); cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod); // If the factorization failed, minor is the column at which it did. On success minor == n. this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue); m_factorizationIsOk = true; }
ComputationInfo Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::info | ( | ) | const [inline] |
Reports whether previous computation was successful.
Success
if computation was succesful, NumericalIssue
if the matrix.appears to be negative. Definition at line 210 of file CholmodSupport.h.
{ eigen_assert(m_isInitialized && "Decomposition is not initialized."); return m_info; }
Scalar Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::logDeterminant | ( | ) | const [inline] |
Definition at line 339 of file CholmodSupport.h.
{ using std::log; using numext::real; eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); RealScalar logDet = 0; Scalar *x = static_cast<Scalar*>(m_cholmodFactor->x); if (m_cholmodFactor->is_super) { // Supernodal factorization stored as a packed list of dense column-major blocs, // as described by the following structure: // super[k] == index of the first column of the j-th super node StorageIndex *super = static_cast<StorageIndex*>(m_cholmodFactor->super); // pi[k] == offset to the description of row indices StorageIndex *pi = static_cast<StorageIndex*>(m_cholmodFactor->pi); // px[k] == offset to the respective dense block StorageIndex *px = static_cast<StorageIndex*>(m_cholmodFactor->px); Index nb_super_nodes = m_cholmodFactor->nsuper; for (Index k=0; k < nb_super_nodes; ++k) { StorageIndex ncols = super[k + 1] - super[k]; StorageIndex nrows = pi[k + 1] - pi[k]; Map<const Array<Scalar,1,Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows+1)); logDet += sk.real().log().sum(); } } else { // Simplicial factorization stored as standard CSC matrix. StorageIndex *p = static_cast<StorageIndex*>(m_cholmodFactor->p); Index size = m_cholmodFactor->n; for (Index k=0; k<size; ++k) logDet += log(real( x[p[k]] )); } if (m_cholmodFactor->is_ll) logDet *= 2.0; return logDet; };
StorageIndex Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::rows | ( | ) | const [inline] |
Definition at line 203 of file CholmodSupport.h.
{ return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
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.
*this
. Definition at line 325 of file CholmodSupport.h.
{ m_shiftOffset[0] = offset; return derived(); }
int Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_analysisIsOk [protected] |
Definition at line 392 of file CholmodSupport.h.
cholmod_common Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_cholmod [mutable, protected] |
Definition at line 387 of file CholmodSupport.h.
cholmod_factor* Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_cholmodFactor [protected] |
Definition at line 388 of file CholmodSupport.h.
int Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_factorizationIsOk [protected] |
Definition at line 391 of file CholmodSupport.h.
ComputationInfo Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_info [mutable, protected] |
Definition at line 390 of file CholmodSupport.h.
RealScalar Eigen::CholmodBase< _MatrixType, _UpLo, Derived >::m_shiftOffset[2] [protected] |
Definition at line 389 of file CholmodSupport.h.