GCC Code Coverage Report


Directory: ./
File: include/eigenpy/decompositions/minres.hpp
Date: 2025-04-25 14:39:22
Exec Total Coverage
Lines: 38 40 95.0%
Branches: 39 80 48.8%

Line Branch Exec Source
1 /*
2 * Copyright 2021 INRIA
3 */
4
5 #ifndef __eigenpy_decompositions_minres_hpp__
6 #define __eigenpy_decompositions_minres_hpp__
7
8 #include <Eigen/Core>
9 #include <iostream>
10 #include <unsupported/Eigen/IterativeSolvers>
11
12 #include "eigenpy/eigenpy.hpp"
13 #include "eigenpy/utils/scalar-name.hpp"
14
15 namespace eigenpy {
16 template <typename _Solver>
17 struct IterativeSolverBaseVisitor
18 : public boost::python::def_visitor<IterativeSolverBaseVisitor<_Solver>> {
19 typedef _Solver Solver;
20 typedef typename Solver::MatrixType MatrixType;
21 typedef typename Solver::Preconditioner Preconditioner;
22 typedef typename Solver::Scalar Scalar;
23 typedef typename Solver::RealScalar RealScalar;
24
25 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
26 MatrixType::Options>
27 MatrixXs;
28
29 template <class PyClass>
30 12 void visit(PyClass& cl) const {
31
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 cl.def("analyzePattern",
32 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
33 Solver::analyzePattern,
34 bp::args("self", "A"),
35 "Initializes the iterative solver for the sparsity pattern of the "
36 "matrix A for further solving Ax=b problems.",
37 12 bp::return_self<>())
38
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
12 .def(
39 "factorize",
40 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
41 Solver::factorize,
42 bp::args("self", "A"),
43 "Initializes the iterative solver with the numerical values of the "
44 "matrix A for further solving Ax=b problems.",
45 12 bp::return_self<>())
46
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
12 .def(
47 "compute",
48 (Solver & (Solver::*)(const Eigen::EigenBase<MatrixType>& matrix)) &
49 Solver::compute,
50 bp::args("self", "A"),
51 "Initializes the iterative solver with the matrix A for further "
52 "solving Ax=b problems.",
53 12 bp::return_self<>())
54
55
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("rows", &Solver::rows, bp::arg("self"),
56 "Returns the number of rows.")
57
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("cols", &Solver::cols, bp::arg("self"),
58 "Returns the number of columns.")
59
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("tolerance", &Solver::tolerance, bp::arg("self"),
60 "Returns the tolerance threshold used by the stopping criteria.")
61
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
12 .def("setTolerance", &Solver::setTolerance,
62 bp::args("self", "tolerance"),
63 "Sets the tolerance threshold used by the stopping criteria.\n"
64 "This value is used as an upper bound to the relative residual "
65 "error: |Ax-b|/|b|.\n"
66 "The default value is the machine precision given by "
67 "NumTraits<Scalar>::epsilon().",
68 12 bp::return_self<>())
69
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 .def("preconditioner",
70 (Preconditioner & (Solver::*)()) & Solver::preconditioner,
71
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 bp::arg("self"),
72 "Returns a read-write reference to the preconditioner for custom "
73 "configuration.",
74 12 bp::return_internal_reference<>())
75
76
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("maxIterations", &Solver::maxIterations, bp::arg("self"),
77 "Returns the max number of iterations.\n"
78 "It is either the value setted by setMaxIterations or, by "
79 "default, twice the number of columns of the matrix.")
80
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
12 .def("setMaxIterations", &Solver::setMaxIterations,
81 bp::args("self", "max_iterations"),
82 "Sets the max number of iterations.\n"
83 "Default is twice the number of columns of the matrix.",
84 12 bp::return_self<>())
85
86
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 .def(
87
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
24 "iterations", &Solver::iterations, bp::arg("self"),
88 "Returns the number of iterations performed during the last solve.")
89
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("error", &Solver::error, bp::arg("self"),
90 "Returns the tolerance error reached during the last solve.\n"
91 "It is a close approximation of the true relative residual error "
92 "|Ax-b|/|b|.")
93
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
24 .def("info", &Solver::error, bp::arg("info"),
94 "Returns Success if the iterations converged, and NoConvergence "
95 "otherwise.")
96
97
3/6
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 12 times.
✗ Branch 8 not taken.
24 .def("solveWithGuess", &solveWithGuess<MatrixXs, MatrixXs>,
98 bp::args("self", "b", "x0"),
99 "Returns the solution x of A x = b using the current "
100 "decomposition of A and x0 as an initial solution.")
101
102
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
12 .def(
103 "solve", &solve<MatrixXs>, bp::args("self", "b"),
104 "Returns the solution x of A x = b using the current decomposition "
105 "of A where b is a right hand side matrix or vector.");
106 12 }
107
108 private:
109 template <typename MatrixOrVector1, typename MatrixOrVector2>
110 static MatrixOrVector1 solveWithGuess(const Solver& self,
111 const MatrixOrVector1& b,
112 const MatrixOrVector2& guess) {
113 return self.solveWithGuess(b, guess);
114 }
115
116 template <typename MatrixOrVector>
117 1 static MatrixOrVector solve(const Solver& self,
118 const MatrixOrVector& mat_or_vec) {
119
1/2
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
1 MatrixOrVector res = self.solve(mat_or_vec);
120 1 return res;
121 }
122 };
123
124 template <typename _MatrixType>
125 struct MINRESSolverVisitor
126 : public boost::python::def_visitor<MINRESSolverVisitor<_MatrixType>> {
127 typedef _MatrixType MatrixType;
128 typedef typename MatrixType::Scalar Scalar;
129 typedef typename MatrixType::RealScalar RealScalar;
130 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
131 VectorXs;
132 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
133 MatrixType::Options>
134 MatrixXs;
135 typedef Eigen::MINRES<MatrixType> Solver;
136
137 template <class PyClass>
138 12 void visit(PyClass& cl) const {
139
2/4
✓ Branch 2 taken 12 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 12 times.
✗ Branch 6 not taken.
12 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
140
3/6
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 12 times.
✗ Branch 8 not taken.
24 .def(bp::init<MatrixType>(
141 bp::args("self", "matrix"),
142 "Initialize the solver with matrix A for further Ax=b solving.\n"
143 "This constructor is a shortcut for the default constructor "
144 "followed by a call to compute()."))
145
146
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 .def(IterativeSolverBaseVisitor<Solver>());
147 12 }
148
149 static void expose() {
150 static const std::string classname =
151 "MINRES" + scalar_name<Scalar>::shortname();
152 expose(classname);
153 }
154
155 12 static void expose(const std::string& name) {
156 12 bp::class_<Solver, boost::noncopyable>(
157 name.c_str(),
158 "A minimal residual solver for sparse symmetric problems.\n"
159 "This class allows to solve for A.x = b sparse linear problems using "
160 "the MINRES algorithm of Paige and Saunders (1975). The sparse matrix "
161 "A must be symmetric (possibly indefinite). The vectors x and b can be "
162 "either dense or sparse.\n"
163 "The maximal number of iterations and tolerance value can be "
164 "controlled via the setMaxIterations() and setTolerance() methods. The "
165 "defaults are the size of the problem for the maximal number of "
166 "iterations and NumTraits<Scalar>::epsilon() for the tolerance.\n",
167 bp::no_init)
168
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 .def(MINRESSolverVisitor())
169
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 .def(IdVisitor<Solver>());
170 12 }
171
172 private:
173 template <typename MatrixOrVector>
174 static MatrixOrVector solve(const Solver& self, const MatrixOrVector& vec) {
175 return self.solve(vec);
176 }
177 };
178
179 } // namespace eigenpy
180
181 #endif // ifndef __eigenpy_decompositions_minres_hpp__
182