eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
SelfAdjointEigenSolver.hpp
1/*
2 * Copyright 2020-2024 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__
6#define __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__
7
8#include <Eigen/Core>
9#include <Eigen/Eigenvalues>
10
11#include "eigenpy/eigen-to-python.hpp"
12#include "eigenpy/eigenpy.hpp"
13#include "eigenpy/utils/scalar-name.hpp"
14
15namespace eigenpy {
16
17template <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 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> VectorType;
25
26 template <class PyClass>
27 void visit(PyClass& cl) const {
28 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
29 .def(bp::init<Eigen::DenseIndex>(
30 bp::args("self", "size"),
31 "Default constructor with memory preallocation"))
32 .def(bp::init<MatrixType, bp::optional<int>>(
33 bp::args("self", "matrix", "options"),
34 "Computes eigendecomposition of given matrix"))
35
36 .def(
37 "eigenvalues",
38 +[](const Solver& c) -> const VectorType& {
39 return c.eigenvalues();
40 },
41 "Returns the eigenvalues of given matrix.",
42 bp::return_internal_reference<>())
43 .def(
44 "eigenvectors",
45 +[](const Solver& c) -> const MatrixType& {
46 return c.eigenvectors();
47 },
48 "Returns the eigenvectors of given matrix.",
49 bp::return_internal_reference<>())
50
51 .def("compute",
52 &SelfAdjointEigenSolverVisitor::compute_proxy<MatrixType>,
53 bp::args("self", "matrix"),
54 "Computes the eigendecomposition of given matrix.",
55 bp::return_value_policy<bp::reference_existing_object>())
56 .def("compute",
57 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix,
58 int options)) &
59 Solver::compute,
60 bp::args("self", "matrix", "options"),
61 "Computes the eigendecomposition of given matrix.",
62 bp::return_self<>())
63
64 .def("computeDirect",
65 &SelfAdjointEigenSolverVisitor::computeDirect_proxy,
66 bp::args("self", "matrix"),
67 "Computes eigendecomposition of given matrix using a closed-form "
68 "algorithm.",
69 bp::return_self<>())
70 .def("computeDirect",
71 (Solver & (Solver::*)(const MatrixType& matrix, int options)) &
72 Solver::computeDirect,
73 bp::args("self", "matrix", "options"),
74 "Computes eigendecomposition of given matrix using a closed-form "
75 "algorithm.",
76 bp::return_self<>())
77
78 .def("operatorInverseSqrt", &Solver::operatorInverseSqrt,
79 bp::arg("self"), "Computes the inverse square root of the matrix.")
80 .def("operatorSqrt", &Solver::operatorSqrt, bp::arg("self"),
81 "Computes the inverse square root of the matrix.")
82
83 .def("info", &Solver::info, bp::arg("self"),
84 "NumericalIssue if the input contains INF or NaN values or "
85 "overflow occured. Returns Success otherwise.");
86 }
87
88 static void expose() {
89 static const std::string classname =
90 "SelfAdjointEigenSolver" + scalar_name<Scalar>::shortname();
91 expose(classname);
92 }
93
94 static void expose(const std::string& name) {
95 bp::class_<Solver>(name.c_str(), bp::no_init)
98 }
99
100 private:
101 template <typename MatrixType>
102 static Solver& compute_proxy(Solver& self,
103 const Eigen::EigenBase<MatrixType>& matrix) {
104 return self.compute(matrix);
105 }
106
107 static Solver& computeDirect_proxy(Solver& self, const MatrixType& matrix) {
108 return self.computeDirect(matrix);
109 }
110};
111
112} // namespace eigenpy
113
114#endif // ifndef __eigenpy_decompositions_self_adjoint_eigen_solver_hpp__
void expose()
Call the expose function of a given type T.
Definition expose.hpp:23
Add the Python method id to retrieving a unique id for a given object exposed with Boost....
Definition id.hpp:18