.. _program_listing_file__tmp_ws_src_eigenpy_include_eigenpy_decompositions_SelfAdjointEigenSolver.hpp: Program Listing for File SelfAdjointEigenSolver.hpp =================================================== |exhale_lsh| :ref:`Return to documentation for file ` (``/tmp/ws/src/eigenpy/include/eigenpy/decompositions/SelfAdjointEigenSolver.hpp``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp /* * Copyright 2020 INRIA */ #ifndef __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__ #define __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__ #include #include #include "eigenpy/eigen-to-python.hpp" #include "eigenpy/eigenpy.hpp" #include "eigenpy/utils/scalar-name.hpp" namespace eigenpy { template struct SelfAdjointEigenSolverVisitor : public boost::python::def_visitor< SelfAdjointEigenSolverVisitor<_MatrixType> > { typedef _MatrixType MatrixType; typedef typename MatrixType::Scalar Scalar; typedef Eigen::SelfAdjointEigenSolver Solver; template void visit(PyClass& cl) const { cl.def(bp::init<>(bp::arg("self"), "Default constructor")) .def(bp::init( bp::args("self", "size"), "Default constructor with memory preallocation")) .def(bp::init >( bp::args("self", "matrix", "options"), "Computes eigendecomposition of given matrix")) .def("eigenvalues", &Solver::eigenvalues, bp::arg("self"), "Returns the eigenvalues of given matrix.", bp::return_internal_reference<>()) .def("eigenvectors", &Solver::eigenvectors, bp::arg("self"), "Returns the eigenvectors of given matrix.", bp::return_internal_reference<>()) .def("compute", &SelfAdjointEigenSolverVisitor::compute_proxy, bp::args("self", "matrix"), "Computes the eigendecomposition of given matrix.", bp::return_value_policy()) .def("compute", (Solver & (Solver::*)(const Eigen::EigenBase& matrix, int options)) & Solver::compute, bp::args("self", "matrix", "options"), "Computes the eigendecomposition of given matrix.", bp::return_self<>()) .def("computeDirect", &SelfAdjointEigenSolverVisitor::computeDirect_proxy, bp::args("self", "matrix"), "Computes eigendecomposition of given matrix using a closed-form " "algorithm.", bp::return_self<>()) .def("computeDirect", (Solver & (Solver::*)(const MatrixType& matrix, int options)) & Solver::computeDirect, bp::args("self", "matrix", "options"), "Computes eigendecomposition of given matrix using a closed-form " "algorithm.", bp::return_self<>()) .def("operatorInverseSqrt", &Solver::operatorInverseSqrt, bp::arg("self"), "Computes the inverse square root of the matrix.") .def("operatorSqrt", &Solver::operatorSqrt, bp::arg("self"), "Computes the inverse square root of the matrix.") .def("info", &Solver::info, bp::arg("self"), "NumericalIssue if the input contains INF or NaN values or " "overflow occured. Returns Success otherwise."); } static void expose() { static const std::string classname = "SelfAdjointEigenSolver" + scalar_name::shortname(); expose(classname); } static void expose(const std::string& name) { bp::class_(name.c_str(), bp::no_init) .def(SelfAdjointEigenSolverVisitor()); } private: template static Solver& compute_proxy(Solver& self, const Eigen::EigenBase& matrix) { return self.compute(matrix); } static Solver& computeDirect_proxy(Solver& self, const MatrixType& matrix) { return self.computeDirect(matrix); } }; } // namespace eigenpy #endif // ifndef __eigenpy_decomposition_self_adjoint_eigen_solver_hpp__