accelerate.hpp
Go to the documentation of this file.
1 /*
2  * Copyright 2024 INRIA
3  */
4 
5 #ifndef __eigenpy_decomposition_sparse_accelerate_accelerate_hpp__
6 #define __eigenpy_decomposition_sparse_accelerate_accelerate_hpp__
7 
8 #include "eigenpy/eigenpy.hpp"
11 
12 #include <Eigen/AccelerateSupport>
13 
14 namespace eigenpy {
15 
16 template <typename AccelerateDerived>
17 struct AccelerateImplVisitor : public boost::python::def_visitor<
18  AccelerateImplVisitor<AccelerateDerived> > {
19  typedef AccelerateDerived Solver;
20 
21  typedef typename AccelerateDerived::MatrixType MatrixType;
22  typedef typename MatrixType::Scalar Scalar;
23  typedef typename MatrixType::RealScalar RealScalar;
25  typedef typename MatrixType::StorageIndex StorageIndex;
26 
27  template <class PyClass>
28  void visit(PyClass &cl) const {
29  cl
30 
31  .def("analyzePattern", &Solver::analyzePattern,
32  bp::args("self", "matrix"),
33  "Performs a symbolic decomposition on the sparcity of matrix.\n"
34  "This function is particularly useful when solving for several "
35  "problems having the same structure.")
36 
39 
40  .def("compute",
41  (Solver & (Solver::*)(const MatrixType &matrix)) & Solver::compute,
42  bp::args("self", "matrix"),
43  "Computes the sparse Cholesky decomposition of a given matrix.",
44  bp::return_self<>())
45 
46  .def("factorize", &Solver::factorize, bp::args("self", "matrix"),
47  "Performs a numeric decomposition of a given matrix.\n"
48  "The given matrix must has the same sparcity than the matrix on "
49  "which the symbolic decomposition has been performed.\n"
50  "See also analyzePattern().")
51 
52  .def("info", &Solver::info, bp::arg("self"),
53  "NumericalIssue if the input contains INF or NaN values or "
54  "overflow occured. Returns Success otherwise.")
55 
56  .def("setOrder", &Solver::setOrder, bp::arg("self"), "Set order");
57  }
58 
59  static void expose(const std::string &name, const std::string &doc = "") {
60  bp::class_<Solver, boost::noncopyable>(name.c_str(), doc.c_str(),
61  bp::no_init)
62  .def(AccelerateImplVisitor())
63 
64  .def(bp::init<>(bp::arg("self"), "Default constructor"))
65  .def(bp::init<MatrixType>(bp::args("self", "matrix"),
66  "Constructs and performs the "
67  "factorization from a given matrix."));
68  }
69 };
70 
71 } // namespace eigenpy
72 
73 #endif // ifndef __eigenpy_decomposition_sparse_accelerate_accelerate_hpp__
eigenpy::SparseSolverBaseVisitor
Definition: decompositions/sparse/SparseSolverBase.hpp:16
eigenpy::AccelerateImplVisitor::Scalar
MatrixType::Scalar Scalar
Definition: accelerate.hpp:22
eigenpy::AccelerateImplVisitor::StorageIndex
MatrixType::StorageIndex StorageIndex
Definition: accelerate.hpp:25
eigenpy::AccelerateImplVisitor::CholMatrixType
MatrixType CholMatrixType
Definition: accelerate.hpp:24
eigenpy
Definition: alignment.hpp:14
eigenpy::EigenBaseVisitor
Definition: EigenBase.hpp:13
eigenpy::AccelerateImplVisitor::expose
static void expose(const std::string &name, const std::string &doc="")
Definition: accelerate.hpp:59
eigenpy::AccelerateImplVisitor::RealScalar
MatrixType::RealScalar RealScalar
Definition: accelerate.hpp:23
setup.name
name
Definition: setup.in.py:179
eigenpy.hpp
eigenpy::AccelerateImplVisitor::Solver
AccelerateDerived Solver
Definition: accelerate.hpp:19
eigenpy::AccelerateImplVisitor
Definition: accelerate.hpp:17
eigenpy::AccelerateImplVisitor::MatrixType
AccelerateDerived::MatrixType MatrixType
Definition: accelerate.hpp:21
eigenpy::AccelerateImplVisitor::visit
void visit(PyClass &cl) const
Definition: accelerate.hpp:28
SparseSolverBase.hpp
EigenBase.hpp


eigenpy
Author(s): Justin Carpentier, Nicolas Mansard
autogenerated on Fri Jun 14 2024 02:15:58