eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
IterativeSolverBase.hpp
1/*
2 * Copyright 2017 CNRS
3 * Copyright 2025 INRIA
4 */
5
6#ifndef __eigenpy_solvers_iterative_solver_base_hpp__
7#define __eigenpy_solvers_iterative_solver_base_hpp__
8
9#include "eigenpy/fwd.hpp"
10#include "eigenpy/solvers/SparseSolverBase.hpp"
11
12namespace eigenpy {
13
14template <typename IterativeSolver>
15struct IterativeSolverVisitor : public boost::python::def_visitor<
16 IterativeSolverVisitor<IterativeSolver>> {
17 typedef typename IterativeSolver::MatrixType MatrixType;
18 typedef typename IterativeSolver::Preconditioner Preconditioner;
19 typedef Eigen::VectorXd VectorType;
20
21 template <class PyClass>
22 void visit(PyClass& cl) const {
23 typedef IterativeSolver IS;
24
26 .def("error", &IS::error,
27 "Returns the tolerance error reached during the last solve.\n"
28 "It is a close approximation of the true relative residual error "
29 "|Ax-b|/|b|.")
30 .def("info", &IS::info,
31 "Returns success if the iterations converged, and NoConvergence "
32 "otherwise.")
33 .def(
34 "iterations", &IS::iterations,
35 "Returns the number of iterations performed during the last solve.")
36 .def("maxIterations", &IS::maxIterations,
37 "Returns the max number of iterations.\n"
38 "It is either the value setted by setMaxIterations or, by "
39 "default, twice the number of columns of the matrix.")
40 .def("setMaxIterations", &IS::setMaxIterations,
41 "Sets the max number of iterations.\n"
42 "Default is twice the number of columns of the matrix.",
43 bp::return_value_policy<bp::reference_existing_object>())
44 .def("tolerance", &IS::tolerance,
45 "Returns he tolerance threshold used by the stopping criteria.")
46 .def("setTolerance", &IS::setTolerance,
47 "Sets the tolerance threshold used by the stopping criteria.\n"
48 "This value is used as an upper bound to the relative residual "
49 "error: |Ax-b|/|b|. The default value is the machine precision.",
50 bp::return_value_policy<bp::reference_existing_object>())
51 .def("analyzePattern", &analyzePattern, bp::arg("A"),
52 "Initializes the iterative solver for the sparsity pattern of the "
53 "matrix A for further solving Ax=b problems.\n"
54 "Currently, this function mostly calls analyzePattern on the "
55 "preconditioner.\n"
56 "In the future we might, for instance, implement column "
57 "reordering for faster matrix vector products.",
58 bp::return_value_policy<bp::reference_existing_object>())
59 .def("factorize", &factorize, bp::arg("A"),
60 "Initializes the iterative solver with the numerical values of "
61 "the matrix A for further solving Ax=b problems.\n"
62 "Currently, this function mostly calls factorize on the "
63 "preconditioner.",
64 bp::return_value_policy<bp::reference_existing_object>())
65 .def("compute", &compute, bp::arg("A"),
66 "Initializes the iterative solver with the numerical values of "
67 "the matrix A for further solving Ax=b problems.\n"
68 "Currently, this function mostly calls factorize on the "
69 "preconditioner.\n"
70 "In the future we might, for instance, implement column "
71 "reordering for faster matrix vector products.",
72 bp::return_value_policy<bp::reference_existing_object>())
73 .def("solveWithGuess", &solveWithGuess, bp::args("b", "x0"),
74 "Returns the solution x of Ax = b using the current decomposition "
75 "of A and x0 as an initial solution.")
76 .def("preconditioner",
77 (Preconditioner & (IS::*)(void)) & IS::preconditioner,
78 "Returns a read-write reference to the preconditioner for custom "
79 "configuration.",
80 bp::return_internal_reference<>());
81 }
82
83 private:
84 static IterativeSolver& factorize(IterativeSolver& self,
85 const MatrixType& m) {
86 return self.factorize(m);
87 }
88
89 static IterativeSolver& compute(IterativeSolver& self, const MatrixType& m) {
90 return self.compute(m);
91 }
92
93 static IterativeSolver& analyzePattern(IterativeSolver& self,
94 const MatrixType& m) {
95 return self.analyzePattern(m);
96 }
97
98 static VectorType solveWithGuess(IterativeSolver& self,
99 const Eigen::VectorXd& b,
100 const Eigen::VectorXd& x0) {
101 return self.solveWithGuess(b, x0);
102 }
103};
104
105} // namespace eigenpy
106
107#endif // ifndef __eigenpy_solvers_iterative_solver_base_hpp__