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 "eigenpy/eigenpy.hpp"
10 
11 #include <Eigen/Core>
12 #include <Eigen/Eigenvalues>
13 
15 
16 namespace eigenpy
17 {
18 
19  template<typename _MatrixType>
21  : public boost::python::def_visitor< SelfAdjointEigenSolverVisitor<_MatrixType> >
22  {
23 
24  typedef _MatrixType MatrixType;
25  typedef typename MatrixType::Scalar Scalar;
26  typedef Eigen::SelfAdjointEigenSolver<MatrixType> Solver;
27 
28  template<class PyClass>
29  void visit(PyClass& cl) const
30  {
31  namespace bp = boost::python;
32  cl
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<int> >(bp::args("matrix","options"),
37  "Computes eigendecomposition of given matrix"))
38 
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.",
44  bp::return_internal_reference<>())
45 
46  .def("compute",&SelfAdjointEigenSolverVisitor::compute_proxy<MatrixType>,
47  bp::args("self","matrix"),
48  "Computes the eigendecomposition of given matrix.",
49  bp::return_value_policy<bp::reference_existing_object>())
50  .def("compute",(Solver & (Solver::*)(const Eigen::EigenBase<MatrixType> & matrix, int options))&Solver::compute,
51  bp::args("self","matrix","options"),
52  "Computes the eigendecomposition of given matrix.",
53  bp::return_self<>())
54 
56  bp::args("self","matrix"),
57  "Computes eigendecomposition of given matrix using a closed-form algorithm.",
58  bp::return_self<>())
59  .def("computeDirect",(Solver & (Solver::*)(const MatrixType & matrix, int options))&Solver::computeDirect,
60  bp::args("self","matrix","options"),
61  "Computes eigendecomposition of given matrix using a closed-form algorithm.",
62  bp::return_self<>())
63 
64  .def("operatorInverseSqrt",&Solver::operatorInverseSqrt,bp::arg("self"),
65  "Computes the inverse square root of the matrix.")
66  .def("operatorSqrt",&Solver::operatorSqrt,bp::arg("self"),
67  "Computes the inverse square root of the matrix.")
68 
69  .def("info",&Solver::info,bp::arg("self"),
70  "NumericalIssue if the input contains INF or NaN values or overflow occured. Returns Success otherwise.")
71  ;
72  }
73 
74  static void expose()
75  {
76  static const std::string classname = "SelfAdjointEigenSolver" + scalar_name<Scalar>::shortname();
77  expose(classname);
78  }
79 
80  static void expose(const std::string & name)
81  {
82  namespace bp = boost::python;
83  bp::class_<Solver>(name.c_str(),
84  bp::no_init)
86  }
87 
88  private:
89 
90  template<typename MatrixType>
91  static Solver & compute_proxy(Solver & self, const Eigen::EigenBase<MatrixType> & matrix)
92  {
93  return self.compute(matrix);
94  }
95 
96  static Solver & computeDirect_proxy(Solver & self, const MatrixType & matrix)
97  {
98  return self.computeDirect(matrix);
99  }
100 
101  };
102 
103 } // namespace eigenpy
104 
105 #endif // ifndef __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__
106 
107 
static void expose(const std::string &name)
boost::python::object matrix()
Definition: bnpy.cpp:20
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 Sat Apr 17 2021 02:37:59