class mrpt::math::CSparseMatrix::CholeskyDecomp¶
Auxiliary class to hold the results of a Cholesky factorization of a sparse matrix.
This implementation does not allow updating/downdating.
Usage example:
CSparseMatrix SM(100,100); SM.insert_entry(i,j, val); ... SM.compressFromTriplet(); // Do Cholesky decomposition: CSparseMatrix::CholeskyDecomp CD(SM); CD.get_inverse(); ...
Only the upper triangular part of the input matrix is accessed.
This class was initially adapted from “robotvision”, by Hauke Strasdat, Steven Lovegrove and Andrew J. Davison. See http://www.openslam.org/robotvision.html
This class designed to be “uncopiable”.
See also:
The main class: CSparseMatrix
#include <mrpt/math/CSparseMatrix.h> class CholeskyDecomp { public: // construction CholeskyDecomp(const CSparseMatrix& A); CholeskyDecomp(const CholeskyDecomp& A); // methods CholeskyDecomp& operator = (const CholeskyDecomp&); CMatrixDouble get_L() const; void get_L(CMatrixDouble& out_L) const; template <class VECTOR> VECTOR backsub(const VECTOR& b) const; void backsub(const CVectorDouble& b, CVectorDouble& result_x) const; void backsub(const double* b, double* result, const size_t N) const; void update(const CSparseMatrix& new_SM); };
Construction¶
CholeskyDecomp(const CSparseMatrix& A)
Constructor from a square definite-positive sparse matrix A, which can be use to solve Ax=b The actual Cholesky decomposition takes places in this constructor.
Constructor from a square semidefinite-positive sparse matrix.
Only the upper triangular part of the matrix is accessed.
The actual Cholesky decomposition takes places in this constructor.
Parameters:
std::runtime_error |
On non-square input matrix. |
On non-definite-positive matrix as input. |
|
std::runtime_error |
On non-square input matrix. |
On non-semidefinite-positive matrix as input. |
Methods¶
CMatrixDouble get_L() const
Return the L matrix (L*L’ = M), as a dense matrix.
void get_L(CMatrixDouble& out_L) const
Return the L matrix (L*L’ = M), as a dense matrix.
template <class VECTOR> VECTOR backsub(const VECTOR& b) const
Return the vector from a back-substitution step that solves: Ux=b.
void backsub(const CVectorDouble& b, CVectorDouble& result_x) const
Return the vector from a back-substitution step that solves: Ux=b.
Vectors can be Eigen::VectorXd or mrpt::math::CVectorDouble
void backsub(const double* b, double* result, const size_t N) const
overload for double pointers which assume the user has reserved the output memory for result
Return the vector from a back-substitution step that solves: Ux=b.
void update(const CSparseMatrix& new_SM)
Update the Cholesky factorization from an updated vesion of the original input, square definite-positive sparse matrix.
NOTE: This new matrix MUST HAVE exactly the same sparse structure than the original one.