5 #ifndef __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__
6 #define __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__
9 #include <Eigen/Eigenvalues>
17 template <
typename _MatrixType>
19 :
public boost::python::def_visitor<
20 SelfAdjointEigenSolverVisitor<_MatrixType> > {
22 typedef typename MatrixType::Scalar
Scalar;
23 typedef Eigen::SelfAdjointEigenSolver<MatrixType>
Solver;
25 template <
class PyClass>
27 cl.def(bp::init<>(bp::arg(
"self"),
"Default constructor"))
28 .def(bp::init<Eigen::DenseIndex>(
29 bp::args(
"self",
"size"),
30 "Default constructor with memory preallocation"))
32 bp::args(
"self",
"matrix",
"options"),
33 "Computes eigendecomposition of given matrix"))
35 .def(
"eigenvalues", &Solver::eigenvalues, bp::arg(
"self"),
36 "Returns the eigenvalues of given matrix.",
37 bp::return_internal_reference<>())
38 .def(
"eigenvectors", &Solver::eigenvectors, bp::arg(
"self"),
39 "Returns the eigenvectors of given matrix.",
40 bp::return_internal_reference<>())
43 &SelfAdjointEigenSolverVisitor::compute_proxy<MatrixType>,
44 bp::args(
"self",
"matrix"),
45 "Computes the eigendecomposition of given matrix.",
46 bp::return_value_policy<bp::reference_existing_object>())
48 (
Solver & (
Solver::*)(
const Eigen::EigenBase<MatrixType>& matrix,
51 bp::args(
"self",
"matrix",
"options"),
52 "Computes the eigendecomposition of given matrix.",
57 bp::args(
"self",
"matrix"),
58 "Computes eigendecomposition of given matrix using a closed-form "
63 Solver::computeDirect,
64 bp::args(
"self",
"matrix",
"options"),
65 "Computes eigendecomposition of given matrix using a closed-form "
69 .def(
"operatorInverseSqrt", &Solver::operatorInverseSqrt,
70 bp::arg(
"self"),
"Computes the inverse square root of the matrix.")
71 .def(
"operatorSqrt", &Solver::operatorSqrt, bp::arg(
"self"),
72 "Computes the inverse square root of the matrix.")
74 .def(
"info", &Solver::info, bp::arg(
"self"),
75 "NumericalIssue if the input contains INF or NaN values or "
76 "overflow occured. Returns Success otherwise.");
80 static const std::string classname =
86 bp::class_<Solver>(
name.c_str(), bp::no_init)
92 template <
typename MatrixType>
94 const Eigen::EigenBase<MatrixType>& matrix) {
95 return self.compute(matrix);
99 return self.computeDirect(matrix);
105 #endif // ifndef __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__