EigenSolver.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2020 INRIA
3  */
4 
5 #ifndef __eigenpy_decomposition_eigen_solver_hpp__
6 #define __eigenpy_decomposition_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<EigenSolverVisitor<_MatrixType> > {
20  typedef _MatrixType MatrixType;
21  typedef typename MatrixType::Scalar Scalar;
22  typedef Eigen::EigenSolver<MatrixType> Solver;
23 
24  template <class PyClass>
25  void visit(PyClass& cl) const {
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"))
32 
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.")
38 
39  .def("compute", &EigenSolverVisitor::compute_proxy<MatrixType>,
40  bp::args("self", "matrix"),
41  "Computes the eigendecomposition of given matrix.",
42  bp::return_self<>())
43  .def("compute",
44  (Solver &
45  (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix, bool)) &
46  Solver::compute,
47  bp::args("self", "matrix", "compute_eigen_vectors"),
48  "Computes the eigendecomposition of given matrix.",
49  bp::return_self<>())
50 
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.",
56  bp::return_self<>())
57 
58  .def("pseudoEigenvalueMatrix", &Solver::pseudoEigenvalueMatrix,
59  bp::arg("self"),
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<>())
65 
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.");
69  }
70 
71  static void expose() {
72  static const std::string classname =
73  "EigenSolver" + scalar_name<Scalar>::shortname();
74  expose(classname);
75  }
76 
77  static void expose(const std::string& name) {
78  bp::class_<Solver>(name.c_str(), bp::no_init).def(EigenSolverVisitor());
79  }
80 
81  private:
82  template <typename MatrixType>
83  static Solver& compute_proxy(Solver& self,
84  const Eigen::EigenBase<MatrixType>& matrix) {
85  return self.compute(matrix);
86  }
87 };
88 
89 } // namespace eigenpy
90 
91 #endif // ifndef __eigenpy_decomposition_eigen_solver_hpp__
void visit(PyClass &cl) const
Definition: EigenSolver.hpp:25
MatrixType::Scalar Scalar
Definition: EigenSolver.hpp:21
static std::string shortname()
static Solver & compute_proxy(Solver &self, const Eigen::EigenBase< MatrixType > &matrix)
Definition: EigenSolver.hpp:83
Eigen::EigenSolver< MatrixType > Solver
Definition: EigenSolver.hpp:22
static void expose(const std::string &name)
Definition: EigenSolver.hpp:77


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