.. _program_listing_file_include_eigenpy_decompositions_minres.hpp: Program Listing for File minres.hpp =================================== |exhale_lsh| :ref:`Return to documentation for file ` (``include/eigenpy/decompositions/minres.hpp``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp /* * Copyright 2021 INRIA */ #ifndef __eigenpy_decompositions_minres_hpp__ #define __eigenpy_decompositions_minres_hpp__ #include #include #include #include "eigenpy/eigenpy.hpp" #include "eigenpy/utils/scalar-name.hpp" namespace eigenpy { template struct IterativeSolverBaseVisitor : public boost::python::def_visitor> { typedef _Solver Solver; typedef typename Solver::MatrixType MatrixType; typedef typename Solver::Preconditioner Preconditioner; typedef typename Solver::Scalar Scalar; typedef typename Solver::RealScalar RealScalar; typedef Eigen::Matrix MatrixXs; template void visit(PyClass& cl) const { cl.def("analyzePattern", (Solver & (Solver::*)(const Eigen::EigenBase& matrix)) & Solver::analyzePattern, bp::args("self", "A"), "Initializes the iterative solver for the sparsity pattern of the " "matrix A for further solving Ax=b problems.", bp::return_self<>()) .def( "factorize", (Solver & (Solver::*)(const Eigen::EigenBase& matrix)) & Solver::factorize, bp::args("self", "A"), "Initializes the iterative solver with the numerical values of the " "matrix A for further solving Ax=b problems.", bp::return_self<>()) .def( "compute", (Solver & (Solver::*)(const Eigen::EigenBase& matrix)) & Solver::compute, bp::args("self", "A"), "Initializes the iterative solver with the matrix A for further " "solving Ax=b problems.", bp::return_self<>()) .def("rows", &Solver::rows, bp::arg("self"), "Returns the number of rows.") .def("cols", &Solver::cols, bp::arg("self"), "Returns the number of columns.") .def("tolerance", &Solver::tolerance, bp::arg("self"), "Returns the tolerance threshold used by the stopping criteria.") .def("setTolerance", &Solver::setTolerance, bp::args("self", "tolerance"), "Sets the tolerance threshold used by the stopping criteria.\n" "This value is used as an upper bound to the relative residual " "error: |Ax-b|/|b|.\n" "The default value is the machine precision given by " "NumTraits::epsilon().", bp::return_self<>()) .def("preconditioner", (Preconditioner & (Solver::*)()) & Solver::preconditioner, bp::arg("self"), "Returns a read-write reference to the preconditioner for custom " "configuration.", bp::return_internal_reference<>()) .def("maxIterations", &Solver::maxIterations, bp::arg("self"), "Returns the max number of iterations.\n" "It is either the value setted by setMaxIterations or, by " "default, twice the number of columns of the matrix.") .def("setMaxIterations", &Solver::setMaxIterations, bp::args("self", "max_iterations"), "Sets the max number of iterations.\n" "Default is twice the number of columns of the matrix.", bp::return_self<>()) .def( "iterations", &Solver::iterations, bp::arg("self"), "Returns the number of iterations performed during the last solve.") .def("error", &Solver::error, bp::arg("self"), "Returns the tolerance error reached during the last solve.\n" "It is a close approximation of the true relative residual error " "|Ax-b|/|b|.") .def("info", &Solver::error, bp::arg("info"), "Returns Success if the iterations converged, and NoConvergence " "otherwise.") .def("solveWithGuess", &solveWithGuess, bp::args("self", "b", "x0"), "Returns the solution x of A x = b using the current " "decomposition of A and x0 as an initial solution.") .def( "solve", &solve, bp::args("self", "b"), "Returns the solution x of A x = b using the current decomposition " "of A where b is a right hand side matrix or vector."); } private: template static MatrixOrVector1 solveWithGuess(const Solver& self, const MatrixOrVector1& b, const MatrixOrVector2& guess) { return self.solveWithGuess(b, guess); } template static MatrixOrVector solve(const Solver& self, const MatrixOrVector& mat_or_vec) { MatrixOrVector res = self.solve(mat_or_vec); return res; } }; template struct MINRESSolverVisitor : public boost::python::def_visitor> { typedef _MatrixType MatrixType; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; typedef Eigen::Matrix VectorXs; typedef Eigen::Matrix MatrixXs; typedef Eigen::MINRES Solver; template void visit(PyClass& cl) const { cl.def(bp::init<>(bp::arg("self"), "Default constructor")) .def(bp::init( bp::args("self", "matrix"), "Initialize the solver with matrix A for further Ax=b solving.\n" "This constructor is a shortcut for the default constructor " "followed by a call to compute().")) .def(IterativeSolverBaseVisitor()); } static void expose() { static const std::string classname = "MINRES" + scalar_name::shortname(); expose(classname); } static void expose(const std::string& name) { bp::class_( name.c_str(), "A minimal residual solver for sparse symmetric problems.\n" "This class allows to solve for A.x = b sparse linear problems using " "the MINRES algorithm of Paige and Saunders (1975). The sparse matrix " "A must be symmetric (possibly indefinite). The vectors x and b can be " "either dense or sparse.\n" "The maximal number of iterations and tolerance value can be " "controlled via the setMaxIterations() and setTolerance() methods. The " "defaults are the size of the problem for the maximal number of " "iterations and NumTraits::epsilon() for the tolerance.\n", bp::no_init) .def(MINRESSolverVisitor()) .def(IdVisitor()); } private: template static MatrixOrVector solve(const Solver& self, const MatrixOrVector& vec) { return self.solve(vec); } }; } // namespace eigenpy #endif // ifndef __eigenpy_decompositions_minres_hpp__