eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
IncompleteCholesky.hpp
1/*
2 * Copyright 2025 INRIA
3 */
4
5#ifndef __eigenpy_solvers_incomplete_cholesky_hpp__
6#define __eigenpy_solvers_incomplete_cholesky_hpp__
7
8#include "eigenpy/eigenpy.hpp"
9#include "eigenpy/utils/scalar-name.hpp"
10
11namespace eigenpy {
12
13template <typename _MatrixType>
14struct IncompleteCholeskyVisitor : public boost::python::def_visitor<
15 IncompleteCholeskyVisitor<_MatrixType>> {
16 typedef _MatrixType MatrixType;
17 typedef typename MatrixType::Scalar Scalar;
18 typedef typename MatrixType::RealScalar RealScalar;
19 typedef Eigen::IncompleteCholesky<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 typedef Eigen::SparseMatrix<Scalar, Eigen::ColMajor> FactorType;
27 typedef Eigen::Matrix<RealScalar, Eigen::Dynamic, 1> VectorRx;
28 typedef Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic>
29 PermutationType;
30
31 template <class PyClass>
32 void visit(PyClass& cl) const {
33 cl.def(bp::init<>(bp::arg("self"), "Default constructor"))
34 .def(bp::init<MatrixType>(bp::args("self", "matrix"),
35 "Constructs and performs the LDLT "
36 "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("setInitialShift", &Solver::setInitialShift,
47 bp::args("self", "shift"), "Set the initial shift parameter.")
48
49 .def(
50 "analyzePattern",
51 +[](Solver& self, const MatrixType& amat) {
52 self.analyzePattern(amat);
53 },
54 bp::arg("matrix"))
55 .def(
56 "factorize",
57 +[](Solver& self, const MatrixType& amat) { self.factorize(amat); },
58 bp::arg("matrix"))
59 .def(
60 "compute",
61 +[](Solver& self, const MatrixType& amat) { self.compute(amat); },
62 bp::arg("matrix"))
63
64 .def(
65 "matrixL",
66 +[](const Solver& self) -> const FactorType& {
67 return self.matrixL();
68 },
69 bp::return_value_policy<bp::copy_const_reference>())
70 .def(
71 "scalingS",
72 +[](const Solver& self) -> const VectorRx& {
73 return self.scalingS();
74 },
75 bp::return_value_policy<bp::copy_const_reference>())
76 .def(
77 "permutationP",
78 +[](const Solver& self) -> const PermutationType& {
79 return self.permutationP();
80 },
81 bp::return_value_policy<bp::copy_const_reference>())
82
83 .def(
84 "solve",
85 +[](const Solver& self, const Eigen::Ref<DenseVectorXs const>& b)
86 -> DenseVectorXs { return self.solve(b); },
87 bp::arg("b"),
88 "Returns the solution x of A x = b using the current decomposition "
89 "of A, where b is a right hand side vector.")
90 .def(
91 "solve",
92 +[](const Solver& self, const Eigen::Ref<DenseMatrixXs const>& B)
93 -> DenseMatrixXs { return self.solve(B); },
94 bp::arg("b"),
95 "Returns the solution X of A X = B using the current decomposition "
96 "of A where B is a right hand side matrix.")
97 .def(
98 "solve",
99 +[](const Solver& self, const MatrixType& B) -> MatrixType {
100 DenseMatrixXs B_dense = DenseMatrixXs(B);
101 DenseMatrixXs X_dense = self.solve(B_dense);
102 return MatrixType(X_dense.sparseView());
103 },
104 bp::arg("b"),
105 "Returns the solution X of A X = B using the current decomposition "
106 "of A where B is a right hand side matrix.");
107 }
108
109 static void expose() {
110 static const std::string classname =
111 "IncompleteCholesky_" + scalar_name<Scalar>::shortname();
112 expose(classname);
113 }
114
115 static void expose(const std::string& name) {
116 bp::class_<Solver, boost::noncopyable>(name.c_str(), "Incomplete Cholesky.",
117 bp::no_init)
119 .def(IdVisitor<Solver>());
120 }
121};
122
123} // namespace eigenpy
124
125#endif // ifndef __eigenpy_solvers_incomplete_cholesky_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