eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
IncompleteLUT.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_solvers_incomplete_lut_hpp__
6#define __eigenpy_solvers_incomplete_lut_hpp__
7
8#include "eigenpy/eigenpy.hpp"
9#include "eigenpy/utils/scalar-name.hpp"
10
11namespace eigenpy {
12
13template <typename _MatrixType>
15 : public boost::python::def_visitor<IncompleteLUTVisitor<_MatrixType>> {
16 typedef _MatrixType MatrixType;
17 typedef typename MatrixType::Scalar Scalar;
18 typedef typename MatrixType::RealScalar RealScalar;
19 typedef Eigen::IncompleteLUT<Scalar> Solver;
20 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
21 DenseVectorXs;
22 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
23 MatrixType::Options>
24 DenseMatrixXs;
25
26 template <class PyClass>
27 void visit(PyClass& cl) const {
28 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
29 .def(bp::init<MatrixType>(bp::args("self", "matrix"),
30 "Constructs and performs the LDLT "
31 "factorization from a given matrix."))
32 .def(bp::init<const MatrixType&, RealScalar, int>(
33 (bp::arg("matrix"),
34 bp::arg("droptol") = Eigen::NumTraits<Scalar>::dummy_precision(),
35 bp::arg("fillfactor") = 10),
36 "Constructs an incomplete LU factorization from a given matrix."))
37
38 .def("rows", &Solver::rows, bp::arg("self"),
39 "Returns the number of rows of the matrix.")
40 .def("cols", &Solver::cols, bp::arg("self"),
41 "Returns the number of cols of the matrix.")
42
43 .def("info", &Solver::info, bp::arg("self"),
44 "Reports whether previous computation was successful.")
45
46 .def(
47 "analyzePattern",
48 +[](Solver& self, const MatrixType& amat) {
49 self.analyzePattern(amat);
50 },
51 bp::arg("matrix"))
52 .def(
53 "factorize",
54 +[](Solver& self, const MatrixType& amat) { self.factorize(amat); },
55 bp::arg("matrix"))
56 .def(
57 "compute",
58 +[](Solver& self, const MatrixType& amat) { self.compute(amat); },
59 bp::arg("matrix"))
60
61 .def("setDroptol", &Solver::setDroptol, bp::arg("self"))
62 .def("setFillfactor", &Solver::setFillfactor, bp::arg("self"))
63
64 .def(
65 "solve",
66 +[](const Solver& self, const Eigen::Ref<DenseVectorXs const>& b)
67 -> DenseVectorXs { return self.solve(b); },
68 bp::arg("b"),
69 "Returns the solution x of A x = b using the current decomposition "
70 "of A, where b is a right hand side vector.")
71 .def(
72 "solve",
73 +[](const Solver& self, const Eigen::Ref<DenseMatrixXs const>& B)
74 -> DenseMatrixXs { return self.solve(B); },
75 bp::arg("b"),
76 "Returns the solution X of A X = B using the current decomposition "
77 "of A where B is a right hand side matrix.")
78 .def(
79 "solve",
80 +[](const Solver& self, const MatrixType& B) -> MatrixType {
81 DenseMatrixXs B_dense = DenseMatrixXs(B);
82 DenseMatrixXs X_dense = self.solve(B_dense);
83 return MatrixType(X_dense.sparseView());
84 },
85 bp::arg("b"),
86 "Returns the solution X of A X = B using the current decomposition "
87 "of A where B is a right hand side matrix.");
88 }
89
90 static void expose() {
91 static const std::string classname =
92 "IncompleteLUT_" + scalar_name<Scalar>::shortname();
93 expose(classname);
94 }
95
96 static void expose(const std::string& name) {
97 bp::class_<Solver, boost::noncopyable>(
98 name.c_str(),
99 "Incomplete LU factorization with dual-threshold strategy.",
100 bp::no_init)
102 .def(IdVisitor<Solver>());
103 }
104};
105
106} // namespace eigenpy
107
108#endif // ifndef __eigenpy_solvers_incomplete_lut_hpp__
void expose()
Call the expose function of a given type T.
Definition expose.hpp:23
Add the Python method id to retrieving a unique id for a given object exposed with Boost....
Definition id.hpp:18