eigenpy 3.12.0
Bindings between Numpy and Eigen using Boost.Python
Loading...
Searching...
No Matches
SimplicialCholesky.hpp
1/*
2 * Copyright 2024 INRIA
3 */
4
5#ifndef __eigenpy_decompositions_sparse_simplicial_cholesky_hpp__
6#define __eigenpy_decompositions_sparse_simplicial_cholesky_hpp__
7
8#include "eigenpy/eigenpy.hpp"
9#include "eigenpy/eigen/EigenBase.hpp"
10#include "eigenpy/decompositions/sparse/SparseSolverBase.hpp"
11
12#include <Eigen/SparseCholesky>
13
14namespace eigenpy {
15
16template <typename SimplicialDerived>
18 : public boost::python::def_visitor<
19 SimplicialCholeskyVisitor<SimplicialDerived>> {
20 typedef SimplicialDerived Solver;
21
22 typedef typename SimplicialDerived::MatrixType MatrixType;
23 typedef typename MatrixType::Scalar Scalar;
24 typedef typename MatrixType::RealScalar RealScalar;
25
26 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1, MatrixType::Options>
27 DenseVectorXs;
28 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic,
29 MatrixType::Options>
30 DenseMatrixXs;
31
32 template <class PyClass>
33 void visit(PyClass &cl) const {
34 cl.def("analyzePattern", &Solver::analyzePattern,
35 bp::args("self", "matrix"),
36 "Performs a symbolic decomposition on the sparcity of matrix.\n"
37 "This function is particularly useful when solving for several "
38 "problems having the same structure.")
39
42
43 .def("matrixL", &matrixL, bp::arg("self"),
44 "Returns the lower triangular matrix L.")
45 .def("matrixU", &matrixU, bp::arg("self"),
46 "Returns the upper triangular matrix U.")
47
48 .def("compute",
49 (Solver & (Solver::*)(const MatrixType &matrix)) & Solver::compute,
50 bp::args("self", "matrix"),
51 "Computes the sparse Cholesky decomposition of a given matrix.",
52 bp::return_self<>())
53
54 .def("determinant", &Solver::determinant, bp::arg("self"),
55 "Returns the determinant of the underlying matrix from the "
56 "current factorization.")
57
58 .def("factorize", &Solver::factorize, bp::args("self", "matrix"),
59 "Performs a numeric decomposition of a given matrix.\n"
60 "The given matrix must has the same sparcity than the matrix on "
61 "which the symbolic decomposition has been performed.\n"
62 "See also analyzePattern().")
63
64 .def("info", &Solver::info, bp::arg("self"),
65 "NumericalIssue if the input contains INF or NaN values or "
66 "overflow occured. Returns Success otherwise.")
67
68 .def("setShift", &Solver::setShift,
69 (bp::args("self", "offset"), bp::arg("scale") = RealScalar(1)),
70 "Sets the shift parameters that will be used to adjust the "
71 "diagonal coefficients during the numerical factorization.\n"
72 "During the numerical factorization, the diagonal coefficients "
73 "are transformed by the following linear model: d_ii = offset + "
74 "scale * d_ii.\n"
75 "The default is the identity transformation with offset=0, and "
76 "scale=1.",
77 bp::return_self<>())
78
79 .def("permutationP", &Solver::permutationP, bp::arg("self"),
80 "Returns the permutation P.",
81 bp::return_value_policy<bp::copy_const_reference>())
82 .def("permutationPinv", &Solver::permutationPinv, bp::arg("self"),
83 "Returns the inverse P^-1 of the permutation P.",
84 bp::return_value_policy<bp::copy_const_reference>());
85 }
86
87 private:
88 static MatrixType matrixL(const Solver &self) { return self.matrixL(); }
89 static MatrixType matrixU(const Solver &self) { return self.matrixU(); }
90};
91
92} // namespace eigenpy
93
94#endif // ifndef __eigenpy_decompositions_sparse_simplicial_cholesky_hpp__