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

linear solver which uses the sparse Cholesky solver from Eigen More...

#include <linear_solver_eigen.h>

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

Classes

class  CholeskyDecomposition
 Sub-classing Eigen's SimplicialLLT to perform ordering with a given ordering. More...
 

Public Types

typedef Eigen::SparseMatrix< double, Eigen::ColMajor > SparseMatrix
 
typedef Eigen::Triplet< double > Triplet
 
typedef Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic > PermutationMatrix
 
using CholeskyDecompositionBase = Eigen::SimplicialLLT< SparseMatrix, Eigen::Upper >
 

Public Member Functions

 LinearSolverEigen ()
 
virtual bool init ()
 
bool solve (const SparseBlockMatrix< MatrixType > &A, double *x, double *b)
 
- Public Member Functions inherited from g2o::LinearSolverCCS< MatrixType >
 LinearSolverCCS ()
 
 ~LinearSolverCCS ()
 
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 blockOrdering () const
 do the AMD ordering on the blocks or on the scalar matrix
 
void setBlockOrdering (bool blockOrdering)
 
- Public Member Functions inherited from g2o::LinearSolver< MatrixType >
 LinearSolver ()
 
virtual ~LinearSolver ()
 
bool writeDebug () const
 write a debug dump of the system matrix if it is not PSD in solve
 
void setWriteDebug (bool b)
 

Public Attributes

 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
 

Protected Member Functions

bool computeCholesky (const SparseBlockMatrix< MatrixType > &A, double &t)
 
void computeSymbolicDecomposition (const SparseBlockMatrix< MatrixType > &A)
 
void fillSparseMatrix (const SparseBlockMatrix< MatrixType > &A, bool onlyValues)
 
bool solveBlocks_impl (const SparseBlockMatrix< MatrixType > &A, std::function< void(MarginalCovarianceCholesky &)> compute)
 
- Protected Member Functions inherited from g2o::LinearSolverCCS< MatrixType >
void initMatrixStructure (const SparseBlockMatrix< MatrixType > &A)
 

Protected Attributes

bool _init
 
SparseMatrix _sparseMatrix
 
CholeskyDecomposition _cholesky
 
- Protected Attributes inherited from g2o::LinearSolverCCS< MatrixType >
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix
 
bool _blockOrdering
 
- 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::LinearSolverEigen< MatrixType >

linear solver which uses the sparse Cholesky solver from Eigen

Has no dependencies except Eigen. Hence, should compile almost everywhere without to much issues. Performance should be similar to CSparse.

Definition at line 51 of file linear_solver_eigen.h.

Member Typedef Documentation

◆ CholeskyDecompositionBase

template<typename MatrixType >
using g2o::LinearSolverEigen< MatrixType >::CholeskyDecompositionBase = Eigen::SimplicialLLT<SparseMatrix, Eigen::Upper>

Definition at line 58 of file linear_solver_eigen.h.

◆ PermutationMatrix

template<typename MatrixType >
typedef Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic> g2o::LinearSolverEigen< MatrixType >::PermutationMatrix

Definition at line 56 of file linear_solver_eigen.h.

◆ SparseMatrix

template<typename MatrixType >
typedef Eigen::SparseMatrix<double, Eigen::ColMajor> g2o::LinearSolverEigen< MatrixType >::SparseMatrix

Definition at line 53 of file linear_solver_eigen.h.

◆ Triplet

template<typename MatrixType >
typedef Eigen::Triplet<double> g2o::LinearSolverEigen< MatrixType >::Triplet

Definition at line 54 of file linear_solver_eigen.h.

Constructor & Destructor Documentation

◆ LinearSolverEigen()

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

Definition at line 87 of file linear_solver_eigen.h.

87: LinearSolverCCS<MatrixType>(), _init(true) {}

Member Function Documentation

◆ computeCholesky()

template<typename MatrixType >
bool g2o::LinearSolverEigen< MatrixType >::computeCholesky ( const SparseBlockMatrix< MatrixType > &  A,
double &  t 
)
inlineprotected

Definition at line 118 of file linear_solver_eigen.h.

118 {
119 // perform some operations only once upon init
120 if (_init) _sparseMatrix.resize(A.rows(), A.cols());
123 _init = false;
124
125 t = get_monotonic_time();
126 _cholesky.factorize(_sparseMatrix);
127 if (_cholesky.info() !=
128 Eigen::Success) { // the matrix is not positive definite
129 if (this->writeDebug()) {
130 G2O_ERROR(
131 "Cholesky failure, writing debug.txt (Hessian loadable by Octave)");
132 A.writeOctave("debug.txt");
133 } else {
134 G2O_DEBUG("Cholesky failure");
135 }
136 return false;
137 }
138 return true;
139 }
void computeSymbolicDecomposition(const SparseBlockMatrix< MatrixType > &A)
void fillSparseMatrix(const SparseBlockMatrix< MatrixType > &A, bool onlyValues)
CholeskyDecomposition _cholesky
bool writeDebug() const
write a debug dump of the system matrix if it is not PSD in solve
#define G2O_ERROR(...)
Definition logger.h:89
#define G2O_DEBUG(...)
Definition logger.h:90
double get_monotonic_time()
Definition timeutil.cpp:43

References g2o::LinearSolverEigen< MatrixType >::_cholesky, g2o::LinearSolverEigen< MatrixType >::_init, g2o::LinearSolverEigen< MatrixType >::_sparseMatrix, g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverEigen< MatrixType >::fillSparseMatrix(), G2O_DEBUG, G2O_ERROR, g2o::get_monotonic_time(), g2o::SparseBlockMatrix< MatrixType >::rows(), g2o::LinearSolver< MatrixType >::writeDebug(), and g2o::SparseBlockMatrix< MatrixType >::writeOctave().

Referenced by g2o::LinearSolverEigen< MatrixType >::solve(), and g2o::LinearSolverEigen< MatrixType >::solveBlocks_impl().

◆ computeSymbolicDecomposition()

template<typename MatrixType >
void g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition ( const SparseBlockMatrix< MatrixType > &  A)
inlineprotected

compute the symbolic decomposition of the matrix only once. Since A has the same pattern in all the iterations, we only compute the fill-in reducing ordering once and re-use for all the following iterations.

Definition at line 147 of file linear_solver_eigen.h.

147 {
148 double t = get_monotonic_time();
149 if (!this->blockOrdering()) {
150 _cholesky.analyzePattern(_sparseMatrix);
151 } else {
152 assert(A.rows() == A.cols() && "Matrix A is not square");
153 // block ordering with the Eigen Interface
154 Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic> blockP;
155 {
156 SparseMatrix auxBlockMatrix(A.blockCols().size(), A.blockCols().size());
157 auxBlockMatrix.resizeNonZeros(A.nonZeroBlocks());
158 // fill the CCS structure of the Eigen SparseMatrix
159 A.fillBlockStructure(auxBlockMatrix.outerIndexPtr(),
160 auxBlockMatrix.innerIndexPtr());
161 // determine ordering by AMD
162 using Ordering = Eigen::AMDOrdering<SparseMatrix::StorageIndex>;
163 Ordering ordering;
164 ordering(auxBlockMatrix, blockP);
165 }
166
167 // Adapt the block permutation to the scalar matrix
168 PermutationMatrix scalarP(A.rows());
169 this->blockToScalarPermutation(A, blockP.indices(), scalarP.indices());
170 // analyze with the scalar permutation
172 }
173 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
174 if (globalStats)
175 globalStats->timeSymbolicDecomposition = get_monotonic_time() - t;
176 }
bool blockOrdering() const
do the AMD ordering on the blocks or on the scalar matrix
void analyzePatternWithPermutation(SparseMatrix &a, const PermutationMatrix &permutation)
use a given permutation for analyzing the pattern of the sparse matrix
Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic > PermutationMatrix
Eigen::SparseMatrix< double, Eigen::ColMajor > SparseMatrix
static void blockToScalarPermutation(const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77

References g2o::LinearSolverEigen< MatrixType >::_cholesky, g2o::LinearSolverEigen< MatrixType >::_sparseMatrix, g2o::LinearSolverEigen< MatrixType >::CholeskyDecomposition::analyzePatternWithPermutation(), g2o::SparseBlockMatrix< MatrixType >::blockCols(), g2o::LinearSolverCCS< MatrixType >::blockOrdering(), g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure(), g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), g2o::SparseBlockMatrix< MatrixType >::nonZeroBlocks(), g2o::SparseBlockMatrix< MatrixType >::rows(), and g2o::G2OBatchStatistics::timeSymbolicDecomposition.

Referenced by g2o::LinearSolverEigen< MatrixType >::computeCholesky().

◆ fillSparseMatrix()

template<typename MatrixType >
void g2o::LinearSolverEigen< MatrixType >::fillSparseMatrix ( const SparseBlockMatrix< MatrixType > &  A,
bool  onlyValues 
)
inlineprotected

Definition at line 178 of file linear_solver_eigen.h.

179 {
180 if (onlyValues) {
181 this->_ccsMatrix->fillCCS(_sparseMatrix.valuePtr(), true);
182 return;
183 }
184 this->initMatrixStructure(A);
185 _sparseMatrix.resizeNonZeros(A.nonZeros());
186 int nz = this->_ccsMatrix->fillCCS(_sparseMatrix.outerIndexPtr(),
187 _sparseMatrix.innerIndexPtr(),
188 _sparseMatrix.valuePtr(), true);
189 (void)nz;
190 assert(nz <= static_cast<int>(_sparseMatrix.data().size()));
191 }
void initMatrixStructure(const SparseBlockMatrix< MatrixType > &A)
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix

References g2o::LinearSolverCCS< MatrixType >::_ccsMatrix, g2o::LinearSolverEigen< MatrixType >::_sparseMatrix, g2o::LinearSolverCCS< MatrixType >::initMatrixStructure(), and g2o::SparseBlockMatrix< MatrixType >::nonZeros().

Referenced by g2o::LinearSolverEigen< MatrixType >::computeCholesky().

◆ init()

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

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

Implements g2o::LinearSolver< MatrixType >.

Definition at line 89 of file linear_solver_eigen.h.

89 {
90 _init = true;
91 return true;
92 }

References g2o::LinearSolverEigen< MatrixType >::_init.

Referenced by g2o::SolverSLAM2DLinear::solveOrientation().

◆ solve()

template<typename MatrixType >
bool g2o::LinearSolverEigen< 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 94 of file linear_solver_eigen.h.

94 {
95 double t;
96 bool cholState = computeCholesky(A, t);
97 if (!cholState) return false;
98
99 // Solving the system
100 VectorX::MapType xx(x, _sparseMatrix.cols());
101 VectorX::ConstMapType bb(b, _sparseMatrix.cols());
102 xx = _cholesky.solve(bb);
103 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
104 if (globalStats) {
105 globalStats->timeNumericDecomposition = get_monotonic_time() - t;
106 globalStats->choleskyNNZ =
107 _cholesky.matrixL().nestedExpression().nonZeros();
108 }
109 return true;
110 }
bool computeCholesky(const SparseBlockMatrix< MatrixType > &A, double &t)

References g2o::LinearSolverEigen< MatrixType >::_cholesky, g2o::LinearSolverEigen< MatrixType >::_sparseMatrix, g2o::G2OBatchStatistics::choleskyNNZ, g2o::LinearSolverEigen< MatrixType >::computeCholesky(), g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), and g2o::G2OBatchStatistics::timeNumericDecomposition.

◆ solveBlocks_impl()

template<typename MatrixType >
bool g2o::LinearSolverEigen< MatrixType >::solveBlocks_impl ( const SparseBlockMatrix< MatrixType > &  A,
std::function< void(MarginalCovarianceCholesky &)>  compute 
)
inlineprotectedvirtual

Implementation of the general parts for computing the inverse blocks of the linear system matrix. Here we call a function to do the underlying computation.

Implements g2o::LinearSolverCCS< MatrixType >.

Definition at line 198 of file linear_solver_eigen.h.

200 {
201 // compute the cholesky factor
202 double t;
203 bool cholState = computeCholesky(A, t);
204 if (!cholState) return false;
205 // compute the inverse blocks
206 MarginalCovarianceCholesky mcc;
207 mcc.setCholeskyFactor(
208 _cholesky.matrixL().rows(),
209 const_cast<int*>(
210 _cholesky.matrixL().nestedExpression().outerIndexPtr()),
211 const_cast<int*>(
212 _cholesky.matrixL().nestedExpression().innerIndexPtr()),
213 const_cast<double*>(_cholesky.matrixL().nestedExpression().valuePtr()),
214 const_cast<int*>(_cholesky.permutationP().indices().data()));
215 compute(mcc); // call the desired computation on the mcc object
216 // book keeping statistics
217 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
218 if (globalStats) {
219 globalStats->choleskyNNZ =
220 _cholesky.matrixL().nestedExpression().nonZeros();
221 }
222 return true;
223 }

References g2o::LinearSolverEigen< MatrixType >::_cholesky, g2o::G2OBatchStatistics::choleskyNNZ, g2o::LinearSolverEigen< MatrixType >::computeCholesky(), g2o::G2OBatchStatistics::globalStats(), and g2o::MarginalCovarianceCholesky::setCholeskyFactor().

Member Data Documentation

◆ _cholesky

template<typename MatrixType >
CholeskyDecomposition g2o::LinearSolverEigen< MatrixType >::_cholesky
protected

◆ _init

template<typename MatrixType >
bool g2o::LinearSolverEigen< MatrixType >::_init
protected

◆ _sparseMatrix

template<typename MatrixType >
SparseMatrix g2o::LinearSolverEigen< MatrixType >::_sparseMatrix
protected

◆ EIGEN_MAKE_ALIGNED_OPERATOR_NEW

template<typename MatrixType >
g2o::LinearSolverEigen< MatrixType >::EIGEN_MAKE_ALIGNED_OPERATOR_NEW

Definition at line 86 of file linear_solver_eigen.h.


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