SelfAdjointEigenSolver.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2020 INRIA
3  */
4 
5 #ifndef __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__
6 #define __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__
7 
8 #include <Eigen/Core>
9 #include <Eigen/Eigenvalues>
10 
12 #include "eigenpy/eigenpy.hpp"
14 
15 namespace eigenpy {
16 
17 template <typename _MatrixType>
19  : public boost::python::def_visitor<
20  SelfAdjointEigenSolverVisitor<_MatrixType> > {
21  typedef _MatrixType MatrixType;
22  typedef typename MatrixType::Scalar Scalar;
23  typedef Eigen::SelfAdjointEigenSolver<MatrixType> Solver;
24 
25  template <class PyClass>
26  void visit(PyClass& cl) const {
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"))
31  .def(bp::init<MatrixType, bp::optional<int> >(
32  bp::args("self", "matrix", "options"),
33  "Computes eigendecomposition of given matrix"))
34 
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<>())
41 
42  .def("compute",
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>())
47  .def("compute",
48  (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix,
49  int options)) &
50  Solver::compute,
51  bp::args("self", "matrix", "options"),
52  "Computes the eigendecomposition of given matrix.",
53  bp::return_self<>())
54 
55  .def("computeDirect",
57  bp::args("self", "matrix"),
58  "Computes eigendecomposition of given matrix using a closed-form "
59  "algorithm.",
60  bp::return_self<>())
61  .def("computeDirect",
62  (Solver & (Solver::*)(const MatrixType& matrix, int options)) &
63  Solver::computeDirect,
64  bp::args("self", "matrix", "options"),
65  "Computes eigendecomposition of given matrix using a closed-form "
66  "algorithm.",
67  bp::return_self<>())
68 
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.")
73 
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.");
77  }
78 
79  static void expose() {
80  static const std::string classname =
81  "SelfAdjointEigenSolver" + scalar_name<Scalar>::shortname();
82  expose(classname);
83  }
84 
85  static void expose(const std::string& name) {
86  bp::class_<Solver>(name.c_str(), bp::no_init)
88  }
89 
90  private:
91  template <typename MatrixType>
92  static Solver& compute_proxy(Solver& self,
93  const Eigen::EigenBase<MatrixType>& matrix) {
94  return self.compute(matrix);
95  }
96 
97  static Solver& computeDirect_proxy(Solver& self, const MatrixType& matrix) {
98  return self.computeDirect(matrix);
99  }
100 };
101 
102 } // namespace eigenpy
103 
104 #endif // ifndef __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__
static void expose(const std::string &name)
static Solver & compute_proxy(Solver &self, const Eigen::EigenBase< MatrixType > &matrix)
static std::string shortname()
Eigen::SelfAdjointEigenSolver< MatrixType > Solver
static Solver & computeDirect_proxy(Solver &self, const MatrixType &matrix)


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Fri Jun 2 2023 02:10:26