5 #ifndef __eigenpy_decomposition_sparse_cholmod_cholmod_supernodal_llt_hpp__
6 #define __eigenpy_decomposition_sparse_cholmod_cholmod_supernodal_llt_hpp__
14 template <
typename MatrixType_,
int UpLo_ = Eigen::Lower>
16 :
public boost::python::def_visitor<
17 CholmodSupernodalLLTVisitor<MatrixType_, UpLo_> > {
19 typedef typename MatrixType::Scalar
Scalar;
22 typedef Eigen::CholmodSupernodalLLT<MatrixType_, UpLo_>
Solver;
24 template <
class PyClass>
29 .def(bp::init<>(bp::arg(
"self"),
"Default constructor"))
30 .def(bp::init<MatrixType>(bp::args(
"self",
"matrix"),
31 "Constructs and performs the LLT "
32 "factorization from a given matrix."))
38 static const std::string classname =
44 bp::class_<Solver, boost::noncopyable>(
46 "A supernodal direct Cholesky (LLT) factorization and solver based on "
48 "This class allows to solve for A.X = B sparse linear problems via a "
49 "supernodal LL^T Cholesky factorization using the Cholmod library."
50 "This supernodal variant performs best on dense enough problems, e.g., "
51 "3D FEM, or very high order 2D FEM."
52 "The sparse matrix A must be selfadjoint and positive definite. The "
53 "vectors or matrices X and B can be either dense or sparse.",