5 #ifndef __eigenpy_decompositions_eigen_solver_hpp__
6 #define __eigenpy_decompositions_eigen_solver_hpp__
9 #include <Eigen/Eigenvalues>
17 template <
typename _MatrixType>
19 :
public boost::python::def_visitor<EigenSolverVisitor<_MatrixType> > {
21 typedef typename MatrixType::Scalar
Scalar;
22 typedef Eigen::EigenSolver<MatrixType>
Solver;
24 template <
class PyClass>
26 cl.def(bp::init<>(
"Default constructor"))
27 .def(bp::init<Eigen::DenseIndex>(
28 bp::arg(
"size"),
"Default constructor with memory preallocation"))
29 .def(bp::init<
MatrixType, bp::optional<bool> >(
30 bp::args(
"matrix",
"compute_eigen_vectors"),
31 "Computes eigendecomposition of given matrix"))
33 .def(
"eigenvalues", &Solver::eigenvalues, bp::arg(
"self"),
34 "Returns the eigenvalues of given matrix.",
35 bp::return_internal_reference<>())
36 .def(
"eigenvectors", &Solver::eigenvectors, bp::arg(
"self"),
37 "Returns the eigenvectors of given matrix.")
39 .def(
"compute", &EigenSolverVisitor::compute_proxy<MatrixType>,
40 bp::args(
"self",
"matrix"),
41 "Computes the eigendecomposition of given matrix.",
45 (
Solver::*)(
const Eigen::EigenBase<MatrixType>& matrix,
bool)) &
47 bp::args(
"self",
"matrix",
"compute_eigen_vectors"),
48 "Computes the eigendecomposition of given matrix.",
51 .def(
"getMaxIterations", &Solver::getMaxIterations, bp::arg(
"self"),
52 "Returns the maximum number of iterations.")
53 .def(
"setMaxIterations", &Solver::setMaxIterations,
54 bp::args(
"self",
"max_iter"),
55 "Sets the maximum number of iterations allowed.",
58 .def(
"pseudoEigenvalueMatrix", &Solver::pseudoEigenvalueMatrix,
60 "Returns the block-diagonal matrix in the "
61 "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 "
68 "overflow occured. Returns Success otherwise.");
72 static const std::string classname =
78 bp::class_<Solver>(
name.c_str(), bp::no_init)
84 template <
typename MatrixType>
86 const Eigen::EigenBase<MatrixType>& matrix) {
87 return self.compute(matrix);
93 #endif // ifndef __eigenpy_decompositions_eigen_solver_hpp__