10 #ifndef EIGEN_CHOLMODSUPPORT_H 11 #define EIGEN_CHOLMODSUPPORT_H 20 template<
typename CholmodType>
22 mat.xtype = CHOLMOD_REAL;
23 mat.dtype = CHOLMOD_DOUBLE;
28 template<
typename CholmodType>
30 mat.xtype = CHOLMOD_COMPLEX;
31 mat.dtype = CHOLMOD_DOUBLE;
57 template<
typename _Scalar,
int _Options,
typename _StorageIndex>
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();
69 if(
mat.isCompressed())
77 res.nz =
mat.innerNonZeroPtr();
85 res.itype = CHOLMOD_INT;
89 res.itype = CHOLMOD_LONG;
104 template<
typename _Scalar,
int _Options,
typename _Index>
111 template<
typename _Scalar,
int _Options,
typename _Index>
120 template<
typename _Scalar,
int _Options,
typename _Index,
unsigned int UpLo>
125 if(UpLo==
Upper) res.stype = 1;
126 if(UpLo==
Lower) res.stype = -1;
133 template<
typename Derived>
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());
154 template<
typename Scalar,
int Flags,
typename StorageIndex>
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) );
172 template<
typename _MatrixType,
int _UpLo,
typename Derived>
178 using Base::m_isInitialized;
181 enum { UpLo = _UpLo };
187 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
188 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
194 : m_cholmodFactor(0), m_info(
Success), m_factorizationIsOk(false), m_analysisIsOk(false)
197 m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
198 cholmod_start(&m_cholmod);
202 : m_cholmodFactor(0), m_info(
Success), m_factorizationIsOk(false), m_analysisIsOk(false)
205 m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
206 cholmod_start(&m_cholmod);
213 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
214 cholmod_finish(&m_cholmod);
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); }
227 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
234 analyzePattern(matrix);
249 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
252 cholmod_sparse
A =
viewAsCholmod(matrix.template selfadjointView<UpLo>());
253 m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
255 this->m_isInitialized =
true;
257 m_analysisIsOk =
true;
258 m_factorizationIsOk =
false;
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);
275 m_factorizationIsOk =
true;
280 cholmod_common&
cholmod() {
return m_cholmod; }
282 #ifndef EIGEN_PARSED_BY_DOXYGEN 284 template<
typename Rhs,
typename Dest>
287 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
296 cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
304 cholmod_free_dense(&x_cd, &m_cholmod);
308 template<
typename RhsDerived,
typename DestDerived>
311 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
319 cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
326 dest.
derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
327 cholmod_free_sparse(&x_cs, &m_cholmod);
329 #endif // EIGEN_PARSED_BY_DOXYGEN 343 m_shiftOffset[0] = double(offset);
359 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
361 RealScalar logDet = 0;
362 Scalar *
x =
static_cast<Scalar*
>(m_cholmodFactor->x);
363 if (m_cholmodFactor->is_super)
369 StorageIndex *super =
static_cast<StorageIndex*
>(m_cholmodFactor->super);
371 StorageIndex *pi =
static_cast<StorageIndex*
>(m_cholmodFactor->pi);
373 StorageIndex *
px =
static_cast<StorageIndex*
>(m_cholmodFactor->px);
375 Index nb_super_nodes = m_cholmodFactor->nsuper;
376 for (
Index k=0; k < nb_super_nodes; ++k)
378 StorageIndex ncols = super[k + 1] - super[k];
379 StorageIndex nrows = pi[k + 1] - pi[k];
382 logDet += sk.real().log().sum();
388 StorageIndex *
p =
static_cast<StorageIndex*
>(m_cholmodFactor->p);
391 logDet +=
log(
real( x[p[k]] ));
393 if (m_cholmodFactor->is_ll)
398 template<
typename Stream>
405 double m_shiftOffset[2];
433 template<
typename _MatrixType,
int _UpLo = Lower>
437 using Base::m_cholmod;
455 m_cholmod.final_asis = 0;
456 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
457 m_cholmod.final_ll = 1;
484 template<
typename _MatrixType,
int _UpLo = Lower>
488 using Base::m_cholmod;
506 m_cholmod.final_asis = 1;
507 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
533 template<
typename _MatrixType,
int _UpLo = Lower>
537 using Base::m_cholmod;
555 m_cholmod.final_asis = 1;
556 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
584 template<
typename _MatrixType,
int _UpLo = Lower>
588 using Base::m_cholmod;
609 m_cholmod.final_asis = 1;
610 m_cholmod.supernodal = CHOLMOD_AUTO;
613 m_cholmod.final_asis = 0;
614 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
615 m_cholmod.final_ll = 1;
618 m_cholmod.final_asis = 1;
619 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
622 m_cholmod.final_asis = 1;
623 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
632 m_cholmod.final_asis = 1;
633 m_cholmod.supernodal = CHOLMOD_AUTO;
639 #endif // EIGEN_CHOLMODSUPPORT_H
void setMode(CholmodMode mode)
StorageIndex rows() const
CholmodBase< _MatrixType, _UpLo, CholmodSupernodalLLT > Base
Scalar determinant() const
CholmodBase< _MatrixType, _UpLo, CholmodDecomposition > Base
A versatible sparse matrix representation.
EIGEN_DEVICE_FUNC const ExpReturnType exp() const
SparseSolverBase< Derived > Base
A matrix or vector expression mapping an existing array of data.
cholmod_factor * m_cholmodFactor
void factorize(const MatrixType &matrix)
double logDeterminant(const typename BAYESTREE::sharedClique &clique)
A base class for sparse solvers.
MappedSparseMatrix< Scalar, Flags, StorageIndex > viewAsEigen(cholmod_sparse &cm)
EIGEN_DEVICE_FUNC const LogReturnType log() const
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate offset
Namespace containing all symbols from the Eigen library.
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
CholmodBase(const MatrixType &matrix)
A supernodal Cholesky (LLT) factorization and solver based on Cholmod.
MatrixType CholMatrixType
CholmodDecomposition(const MatrixType &matrix)
const unsigned int RowMajorBit
CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLDLT > Base
MatrixType::StorageIndex StorageIndex
cholmod_common & cholmod()
StorageIndex cols() const
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Base class of any sparse matrices or sparse expressions.
void dumpMemory(Stream &)
Convenience specialization of Stride to specify only an inner stride See class Map for some examples...
void analyzePattern(const MatrixType &matrix)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
RealScalar RealScalar * px
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLLT > Base
A general Cholesky factorization and solver based on Cholmod.
NumTraits< Scalar >::Real RealScalar
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< _Scalar, _Options, _StorageIndex > > mat)
The base class for the direct Cholesky factorization of Cholmod.
ComputationInfo info() const
Reports whether previous computation was successful.
A matrix or vector expression mapping an existing expression.
A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod.
MatrixType::RealScalar RealScalar
A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod.
Derived & compute(const MatrixType &matrix)
const Derived & derived() const
Derived & setShift(const RealScalar &offset)
detail::initimpl::constructor< Args... > init()
Binds an existing constructor taking arguments Args...
static ConstMapType Map(const Scalar *data)
MatrixType::Scalar Scalar
CholmodSimplicialLLT(const MatrixType &matrix)
CholmodSupernodalLLT(const MatrixType &matrix)
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
void run(Expr &expr, Dev &dev)
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Base class for all dense matrices, vectors, and expressions.
Scalar logDeterminant() const
#define EIGEN_UNUSED_VARIABLE(var)
void _solve_impl(const SparseMatrixBase< RhsDerived > &b, SparseMatrixBase< DestDerived > &dest) const
CholmodSimplicialLDLT(const MatrixType &matrix)
SparseMatrix< _Scalar, _Options, _StorageIndex > & const_cast_derived() const