5 #ifndef __eigenpy_decomposition_eigen_solver_hpp__ 6 #define __eigenpy_decomposition_eigen_solver_hpp__ 12 #include <Eigen/Eigenvalues> 19 template<
typename _MatrixType>
21 :
public boost::python::def_visitor< EigenSolverVisitor<_MatrixType> >
25 typedef typename MatrixType::Scalar
Scalar;
26 typedef Eigen::EigenSolver<MatrixType>
Solver;
28 template<
class PyClass>
33 .def(bp::init<>(
"Default constructor"))
34 .def(bp::init<Eigen::DenseIndex>(bp::arg(
"size"),
35 "Default constructor with memory preallocation"))
36 .def(bp::init<MatrixType,bp::optional<bool> >(bp::args(
"matrix",
"compute_eigen_vectors"),
37 "Computes eigendecomposition of given matrix"))
39 .def(
"eigenvalues",&Solver::eigenvalues,bp::arg(
"self"),
40 "Returns the eigenvalues of given matrix.",
41 bp::return_internal_reference<>())
42 .def(
"eigenvectors",&Solver::eigenvectors,bp::arg(
"self"),
43 "Returns the eigenvectors of given matrix.")
45 .def(
"compute",&EigenSolverVisitor::compute_proxy<MatrixType>,
46 bp::args(
"self",
"matrix"),
47 "Computes the eigendecomposition of given matrix.",
49 .def(
"compute",(Solver & (Solver::*)(
const Eigen::EigenBase<MatrixType> &
matrix,
bool))&Solver::compute,
50 bp::args(
"self",
"matrix",
"compute_eigen_vectors"),
51 "Computes the eigendecomposition of given matrix.",
54 .def(
"getMaxIterations",&Solver::getMaxIterations,bp::arg(
"self"),
55 "Returns the maximum number of iterations.")
56 .def(
"setMaxIterations",&Solver::setMaxIterations,bp::args(
"self",
"max_iter"),
57 "Sets the maximum number of iterations allowed.",
60 .def(
"pseudoEigenvalueMatrix",&Solver::pseudoEigenvalueMatrix,bp::arg(
"self"),
61 "Returns the block-diagonal matrix in the pseudo-eigendecomposition.")
62 .def(
"pseudoEigenvectors",&Solver::pseudoEigenvectors ,bp::arg(
"self"),
63 "Returns the pseudo-eigenvectors of given matrix.",
64 bp::return_internal_reference<>())
66 .def(
"info",&Solver::info,bp::arg(
"self"),
67 "NumericalIssue if the input contains INF or NaN values or overflow occured. Returns Success otherwise.")
80 bp::class_<Solver>(name.c_str(),
87 template<
typename MatrixType>
90 return self.compute(matrix);
97 #endif // ifndef __eigenpy_decomposition_eigen_solver_hpp__ boost::python::object matrix()
MatrixType::Scalar Scalar
void visit(PyClass &cl) const
static std::string shortname()
static Solver & compute_proxy(Solver &self, const Eigen::EigenBase< MatrixType > &matrix)
Eigen::EigenSolver< MatrixType > Solver
static void expose(const std::string &name)