g2o
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
g2o::LinearSolverDense< MatrixType > Class Template Reference

linear solver using dense cholesky decomposition More...

#include <linear_solver_dense.h>

Inheritance diagram for g2o::LinearSolverDense< MatrixType >:
Inheritance graph
[legend]
Collaboration diagram for g2o::LinearSolverDense< MatrixType >:
Collaboration graph
[legend]

Public Member Functions

 LinearSolverDense ()
 
virtual ~LinearSolverDense ()
 
virtual bool init ()
 
bool solve (const SparseBlockMatrix< MatrixType > &A, double *x, double *b)
 
- Public Member Functions inherited from g2o::LinearSolver< MatrixType >
 LinearSolver ()
 
virtual ~LinearSolver ()
 
virtual bool solveBlocks (double **&blocks, const SparseBlockMatrix< MatrixType > &A)
 
virtual bool solvePattern (SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices, const SparseBlockMatrix< MatrixType > &A)
 
bool writeDebug () const
 write a debug dump of the system matrix if it is not PSD in solve
 
void setWriteDebug (bool b)
 

Protected Attributes

bool _reset
 
MatrixX _H
 
Eigen::LDLT< MatrixX_cholesky
 
- Protected Attributes inherited from g2o::LinearSolver< MatrixType >
bool _writeDebug
 

Additional Inherited Members

- Static Public Member Functions inherited from g2o::LinearSolver< MatrixType >
static void allocateBlocks (const SparseBlockMatrix< MatrixType > &A, double **&blocks)
 allocate block memory structure
 
static void deallocateBlocks (const SparseBlockMatrix< MatrixType > &A, double **&blocks)
 de-allocate the block structure
 
template<typename BlockDerived , typename ScalarDerived >
static void blockToScalarPermutation (const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
 

Detailed Description

template<typename MatrixType>
class g2o::LinearSolverDense< MatrixType >

linear solver using dense cholesky decomposition

Definition at line 46 of file linear_solver_dense.h.

Constructor & Destructor Documentation

◆ LinearSolverDense()

template<typename MatrixType >
g2o::LinearSolverDense< MatrixType >::LinearSolverDense ( )
inline

Definition at line 48 of file linear_solver_dense.h.

48: LinearSolver<MatrixType>(), _reset(true) {}

◆ ~LinearSolverDense()

template<typename MatrixType >
virtual g2o::LinearSolverDense< MatrixType >::~LinearSolverDense ( )
inlinevirtual

Definition at line 50 of file linear_solver_dense.h.

50{}

Member Function Documentation

◆ init()

template<typename MatrixType >
virtual bool g2o::LinearSolverDense< MatrixType >::init ( )
inlinevirtual

init for operating on matrices with a different non-zero pattern like before

Implements g2o::LinearSolver< MatrixType >.

Definition at line 52 of file linear_solver_dense.h.

52 {
53 _reset = true;
54 return true;
55 }

References g2o::LinearSolverDense< MatrixType >::_reset.

◆ solve()

template<typename MatrixType >
bool g2o::LinearSolverDense< MatrixType >::solve ( const SparseBlockMatrix< MatrixType > &  A,
double *  x,
double *  b 
)
inlinevirtual

Assumes that A is the same matrix for several calls. Among other assumptions, the non-zero pattern does not change! If the matrix changes call init() before. solve system Ax = b, x and b have to allocated beforehand!!

Implements g2o::LinearSolver< MatrixType >.

Definition at line 57 of file linear_solver_dense.h.

57 {
58 int n = A.cols();
59 int m = A.cols();
60
61 MatrixX& H = _H;
62 if (H.cols() != n) {
63 H.resize(n, m);
64 _reset = true;
65 }
66 if (_reset) {
67 _reset = false;
68 H.setZero();
69 }
70
71 // copy the sparse block matrix into a dense matrix
72 int c_idx = 0;
73 for (size_t i = 0; i < A.blockCols().size(); ++i) {
74 int c_size = A.colsOfBlock(i);
75 assert(c_idx == A.colBaseOfBlock(i) && "mismatch in block indices");
76
78 A.blockCols()[i];
79 if (col.size() > 0) {
80 typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it;
81 for (it = col.begin(); it != col.end(); ++it) {
82 int r_idx = A.rowBaseOfBlock(it->first);
83 // only the upper triangular block is processed
84 if (it->first <= (int)i) {
85 int r_size = A.rowsOfBlock(it->first);
86 H.block(r_idx, c_idx, r_size, c_size) = *(it->second);
87 if (r_idx != c_idx) // write the lower triangular block
88 H.block(c_idx, r_idx, c_size, r_size) = it->second->transpose();
89 }
90 }
91 }
92
93 c_idx += c_size;
94 }
95
96 // solving via Cholesky decomposition
97 VectorX::MapType xvec(x, m);
98 VectorX::ConstMapType bvec(b, n);
99 _cholesky.compute(H);
100 if (_cholesky.isPositive()) {
101 xvec = _cholesky.solve(bvec);
102 return true;
103 }
104 return false;
105 }
Eigen::LDLT< MatrixX > _cholesky
std::map< int, SparseMatrixBlock * > IntBlockMap
MatrixN< Eigen::Dynamic > MatrixX
Definition eigen_types.h:74

References g2o::LinearSolverDense< MatrixType >::_cholesky, g2o::LinearSolverDense< MatrixType >::_H, g2o::LinearSolverDense< MatrixType >::_reset, g2o::SparseBlockMatrix< MatrixType >::blockCols(), g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::rowsOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::transpose().

Member Data Documentation

◆ _cholesky

template<typename MatrixType >
Eigen::LDLT<MatrixX> g2o::LinearSolverDense< MatrixType >::_cholesky
protected

Definition at line 110 of file linear_solver_dense.h.

Referenced by g2o::LinearSolverDense< MatrixType >::solve().

◆ _H

template<typename MatrixType >
MatrixX g2o::LinearSolverDense< MatrixType >::_H
protected

Definition at line 109 of file linear_solver_dense.h.

Referenced by g2o::LinearSolverDense< MatrixType >::solve().

◆ _reset

template<typename MatrixType >
bool g2o::LinearSolverDense< MatrixType >::_reset
protected

The documentation for this class was generated from the following file: