g2o
Loading...
Searching...
No Matches
linear_solver_eigen.h
Go to the documentation of this file.
1// g2o - General Graph Optimization
2// Copyright (C) 2011 R. Kuemmerle, G. Grisetti, W. Burgard
3// All rights reserved.
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are
7// met:
8//
9// * Redistributions of source code must retain the above copyright notice,
10// this list of conditions and the following disclaimer.
11// * Redistributions in binary form must reproduce the above copyright
12// notice, this list of conditions and the following disclaimer in the
13// documentation and/or other materials provided with the distribution.
14//
15// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
16// IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
17// TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
18// PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
19// HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
20// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
21// TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
22// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
23// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
24// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
25// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26
27#ifndef G2O_LINEAR_SOLVER_EIGEN_H
28#define G2O_LINEAR_SOLVER_EIGEN_H
29
30#include <Eigen/Sparse>
31#include <Eigen/SparseCholesky>
32#include <cassert>
33#include <iostream>
34#include <vector>
35
39#include "g2o/stuff/logger.h"
40#include "g2o/stuff/timeutil.h"
41
42namespace g2o {
43
50template <typename MatrixType>
51class LinearSolverEigen : public LinearSolverCCS<MatrixType> {
52 public:
53 typedef Eigen::SparseMatrix<double, Eigen::ColMajor> SparseMatrix;
54 typedef Eigen::Triplet<double> Triplet;
55 typedef Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic>
57
59 Eigen::SimplicialLLT<SparseMatrix, Eigen::Upper>;
60
66 public:
68
71 const PermutationMatrix& permutation) {
72 m_Pinv = permutation;
73 m_P = permutation.inverse();
74 int size = a.cols();
75 SparseMatrix ap(size, size);
76 ap.selfadjointView<Eigen::Upper>() =
77 a.selfadjointView<UpLo>().twistedBy(m_P);
78 analyzePattern_preordered(ap, false);
79 }
80
81 protected:
82 using CholeskyDecompositionBase::analyzePattern_preordered;
83 };
84
85 public:
87 LinearSolverEigen() : LinearSolverCCS<MatrixType>(), _init(true) {}
88
89 virtual bool init() {
90 _init = true;
91 return true;
92 }
93
94 bool solve(const SparseBlockMatrix<MatrixType>& A, double* x, double* b) {
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);
104 if (globalStats) {
105 globalStats->timeNumericDecomposition = get_monotonic_time() - t;
106 globalStats->choleskyNNZ =
107 _cholesky.matrixL().nestedExpression().nonZeros();
108 }
109 return true;
110 }
111
112 protected:
113 bool _init;
116
117 // compute the cholesky factor
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 }
140
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 }
174 if (globalStats)
176 }
177
179 bool onlyValues) {
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 }
192
200 std::function<void(MarginalCovarianceCholesky&)> compute) {
201 // compute the cholesky factor
202 double t;
203 bool cholState = computeCholesky(A, t);
204 if (!cholState) return false;
205 // compute the inverse blocks
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
218 if (globalStats) {
219 globalStats->choleskyNNZ =
220 _cholesky.matrixL().nestedExpression().nonZeros();
221 }
222 return true;
223 }
224};
225
226} // namespace g2o
227
228#endif
Solver with faster iterating structure for the linear matrix.
void initMatrixStructure(const SparseBlockMatrix< MatrixType > &A)
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix
bool blockOrdering() const
do the AMD ordering on the blocks or on the scalar matrix
Sub-classing Eigen's SimplicialLLT to perform ordering with a given ordering.
void analyzePatternWithPermutation(SparseMatrix &a, const PermutationMatrix &permutation)
use a given permutation for analyzing the pattern of the sparse matrix
linear solver which uses the sparse Cholesky solver from Eigen
void computeSymbolicDecomposition(const SparseBlockMatrix< MatrixType > &A)
bool solveBlocks_impl(const SparseBlockMatrix< MatrixType > &A, std::function< void(MarginalCovarianceCholesky &)> compute)
Eigen::Triplet< double > Triplet
bool computeCholesky(const SparseBlockMatrix< MatrixType > &A, double &t)
void fillSparseMatrix(const SparseBlockMatrix< MatrixType > &A, bool onlyValues)
Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic > PermutationMatrix
CholeskyDecomposition _cholesky
Eigen::SimplicialLLT< SparseMatrix, Eigen::Upper > CholeskyDecompositionBase
bool solve(const SparseBlockMatrix< MatrixType > &A, double *x, double *b)
Eigen::SparseMatrix< double, Eigen::ColMajor > SparseMatrix
bool writeDebug() const
write a debug dump of the system matrix if it is not PSD in solve
static void blockToScalarPermutation(const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
computing the marginal covariance given a cholesky factor (lower triangle of the factor)
void setCholeskyFactor(int n, int *Lp, int *Li, double *Lx, int *permInv)
Sparse matrix which uses blocks.
size_t nonZeros() const
number of non-zero elements
bool writeOctave(const char *filename, bool upperTriangle=true) const
size_t nonZeroBlocks() const
number of allocated blocks
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
int rows() const
rows of the matrix
int cols() const
columns of the matrix
void fillBlockStructure(MatrixStructure &ms) const
exports the non zero blocks in the structure matrix ms
#define G2O_ERROR(...)
Definition logger.h:89
#define G2O_DEBUG(...)
Definition logger.h:90
double get_monotonic_time()
Definition timeutil.cpp:43
statistics about the optimization
Definition batch_stats.h:40
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77
double timeNumericDecomposition
numeric decomposition (0 if not done)
Definition batch_stats.h:58
double timeSymbolicDecomposition
symbolic decomposition (0 if not done)
Definition batch_stats.h:57
size_t choleskyNNZ
number of non-zeros in the cholesky factor
Definition batch_stats.h:75
utility functions for handling time related stuff