10 #ifndef EIGEN_CHOLMODSUPPORT_H
11 #define EIGEN_CHOLMODSUPPORT_H
19 template<
typename Scalar>
struct cholmod_configure_matrix;
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;
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;
59 template<
typename Scalar_,
int Options_,
typename StorageIndex_>
63 res.nzmax =
mat.nonZeros();
66 res.p =
mat.outerIndexPtr();
67 res.i =
mat.innerIndexPtr();
71 if(
mat.isCompressed())
79 res.nz =
mat.innerNonZeroPtr();
85 if (internal::is_same<StorageIndex_,int>::value)
87 res.itype = CHOLMOD_INT;
89 else if (internal::is_same<StorageIndex_,SuiteSparse_long>::value)
91 res.itype = CHOLMOD_LONG;
99 internal::cholmod_configure_matrix<Scalar_>::run(
res);
106 template<
typename Scalar_,
int Options_,
typename Index_>
113 template<
typename Scalar_,
int Options_,
typename Index_>
122 template<
typename Scalar_,
int Options_,
typename Index_,
unsigned int UpLo>
138 template<
typename Derived>
142 typedef typename Derived::Scalar Scalar;
152 internal::cholmod_configure_matrix<Scalar>::run(
res);
159 template<
typename Scalar,
int Flags,
typename StorageIndex>
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) );
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); }
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); }
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); }
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); }
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); }
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); }
199 #undef EIGEN_CHOLMOD_SPECIALIZE0
200 #undef EIGEN_CHOLMOD_SPECIALIZE1
215 template<
typename MatrixType_,
int UpLo_,
typename Derived>
225 typedef typename MatrixType::Scalar
Scalar;
239 EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
241 internal::cm_start<StorageIndex>(
m_cholmod);
247 EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
249 internal::cm_start<StorageIndex>(
m_cholmod);
257 internal::cm_finish<StorageIndex>(
m_cholmod);
295 cholmod_sparse
A =
viewAsCholmod(matrix.template selfadjointView<UpLo>());
313 cholmod_sparse
A =
viewAsCholmod(matrix.template selfadjointView<UpLo>());
325 #ifndef EIGEN_PARSED_BY_DOXYGEN
327 template<
typename Rhs,
typename Dest>
330 eigen_assert(
m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
348 internal::cm_free_dense<StorageIndex>(x_cd,
m_cholmod);
352 template<
typename RhsDerived,
typename DestDerived>
353 void _solve_impl(
const SparseMatrixBase<RhsDerived> &
b, SparseMatrixBase<DestDerived> &dest)
const
355 eigen_assert(
m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
361 Ref<SparseMatrix<typename RhsDerived::Scalar,ColMajor,typename RhsDerived::StorageIndex> > b_ref(
b.const_cast_derived());
371 dest.derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
372 internal::cm_free_sparse<StorageIndex>(x_cs,
m_cholmod);
404 eigen_assert(
m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
421 for (
Index k=0; k < nb_super_nodes; ++k)
427 logDet += sk.real().log().sum();
443 template<
typename Stream>
478 template<
typename MatrixType_,
int UpLo_ = Lower>
501 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
529 template<
typename MatrixType_,
int UpLo_ = Lower>
552 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
578 template<
typename MatrixType_,
int UpLo_ = Lower>
601 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
629 template<
typename MatrixType_,
int UpLo_ = Lower>
659 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
664 m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
668 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
const ExpReturnType exp() const
const LogReturnType log() const
#define EIGEN_CHOLMOD_SPECIALIZE1(ret, name, t1, a1)
#define EIGEN_CHOLMOD_SPECIALIZE0(ret, name)
RealReturnType real() const
#define EIGEN_UNUSED_VARIABLE(var)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT(X, MSG)
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)
Scalar logDeterminant() const
StorageIndex rows() const
void dumpMemory(Stream &)
cholmod_factor * m_cholmodFactor
Derived & setShift(const RealScalar &offset)
MatrixType::StorageIndex StorageIndex
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.
A matrix or vector expression mapping an existing array of data.
Base class for all dense matrices, vectors, and expressions.
MatrixBase< Derived > & matrix()
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
static ConstMapType Map(const Scalar *data)
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
A matrix or vector expression mapping an existing expression.
A versatible sparse matrix representation.
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
A base class for sparse solvers.
const unsigned int RowMajorBit
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)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)
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)
Derived & const_cast_derived() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.