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


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Sat Nov 2 2024 02:14:45