.. _program_listing_file_include_eigenpy_decompositions_EigenSolver.hpp: Program Listing for File EigenSolver.hpp ======================================== |exhale_lsh| :ref:`Return to documentation for file ` (``include/eigenpy/decompositions/EigenSolver.hpp``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp /* * Copyright 2020 INRIA */ #ifndef __eigenpy_decompositions_eigen_solver_hpp__ #define __eigenpy_decompositions_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 EigenSolverVisitor : public boost::python::def_visitor > { typedef _MatrixType MatrixType; typedef typename MatrixType::Scalar Scalar; typedef Eigen::EigenSolver Solver; template void visit(PyClass& cl) const { cl.def(bp::init<>("Default constructor")) .def(bp::init( bp::arg("size"), "Default constructor with memory preallocation")) .def(bp::init >( bp::args("matrix", "compute_eigen_vectors"), "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.") .def("compute", &EigenSolverVisitor::compute_proxy, bp::args("self", "matrix"), "Computes the eigendecomposition of given matrix.", bp::return_self<>()) .def("compute", (Solver & (Solver::*)(const Eigen::EigenBase& matrix, bool)) & Solver::compute, bp::args("self", "matrix", "compute_eigen_vectors"), "Computes the eigendecomposition of given matrix.", bp::return_self<>()) .def("getMaxIterations", &Solver::getMaxIterations, bp::arg("self"), "Returns the maximum number of iterations.") .def("setMaxIterations", &Solver::setMaxIterations, bp::args("self", "max_iter"), "Sets the maximum number of iterations allowed.", bp::return_self<>()) .def("pseudoEigenvalueMatrix", &Solver::pseudoEigenvalueMatrix, bp::arg("self"), "Returns the block-diagonal matrix in the " "pseudo-eigendecomposition.") .def("pseudoEigenvectors", &Solver::pseudoEigenvectors, bp::arg("self"), "Returns the pseudo-eigenvectors of given matrix.", bp::return_internal_reference<>()) .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 = "EigenSolver" + scalar_name::shortname(); expose(classname); } static void expose(const std::string& name) { bp::class_(name.c_str(), bp::no_init) .def(EigenSolverVisitor()) .def(IdVisitor()); } private: template static Solver& compute_proxy(Solver& self, const Eigen::EigenBase& matrix) { return self.compute(matrix); } }; } // namespace eigenpy #endif // ifndef __eigenpy_decompositions_eigen_solver_hpp__