g2o
Loading...
Searching...
No Matches
linear_solver.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_H
28#define G2O_LINEAR_SOLVER_H
29
30#include <cassert>
31#include <functional>
32
34#include "sparse_block_matrix.h"
36
37namespace g2o {
38
46template <typename MatrixType>
48 public:
50 virtual ~LinearSolver() {}
51
56 virtual bool init() = 0;
57
64 virtual bool solve(const SparseBlockMatrix<MatrixType>& A, double* x,
65 double* b) = 0;
66
71 virtual bool solveBlocks(double**& blocks,
73 (void)blocks;
74 (void)A;
75 return false;
76 }
77
82 virtual bool solvePattern(
84 const std::vector<std::pair<int, int> >& blockIndices,
86 (void)spinv;
87 (void)blockIndices;
88 (void)A;
89 return false;
90 }
91
93 bool writeDebug() const { return _writeDebug; }
94 void setWriteDebug(bool b) { _writeDebug = b; }
95
98 double**& blocks) {
99 blocks = new double*[A.rows()];
100 double** block = blocks;
101 for (size_t i = 0; i < A.rowBlockIndices().size(); ++i) {
102 int dim = A.rowsOfBlock(i) * A.colsOfBlock(i);
103 *block = new double[dim];
104 block++;
105 }
106 }
107
110 double**& blocks) {
111 for (size_t i = 0; i < A.rowBlockIndices().size(); ++i) {
112 delete[] blocks[i];
113 }
114 delete[] blocks;
115 blocks = nullptr;
116 }
117
121 template <typename BlockDerived, typename ScalarDerived>
124 const Eigen::MatrixBase<BlockDerived>& p,
125 const Eigen::MatrixBase<ScalarDerived>& scalar /* output */) {
126 int n = A.cols();
127 Eigen::MatrixBase<ScalarDerived>& scalarPermutation =
128 const_cast<Eigen::MatrixBase<ScalarDerived>&>(scalar);
129 if (scalarPermutation.size() == 0) scalarPermutation.derived().resize(n);
130 if (scalarPermutation.size() < n) scalarPermutation.derived().resize(2 * n);
131 size_t scalarIdx = 0;
132 for (size_t i = 0; i < A.colBlockIndices().size(); ++i) {
133 int base = A.colBaseOfBlock(p(i));
134 int nCols = A.colsOfBlock(p(i));
135 for (int j = 0; j < nCols; ++j) {
136 scalarPermutation(scalarIdx++) = base++;
137 }
138 }
139 assert((int)scalarIdx == n);
140 }
141
142 protected:
144};
145
149template <typename MatrixType>
150class LinearSolverCCS : public LinearSolver<MatrixType> {
151 public:
153 : LinearSolver<MatrixType>(), _ccsMatrix(0), _blockOrdering(true) {}
155
156 virtual bool solveBlocks(double**& blocks,
158 auto compute = [&](MarginalCovarianceCholesky& mcc) {
159 if (!blocks) LinearSolverCCS<MatrixType>::allocateBlocks(A, blocks);
160 mcc.computeCovariance(blocks, A.rowBlockIndices());
161 };
162 return solveBlocks_impl(A, compute);
163 }
164
165 virtual bool solvePattern(
167 const std::vector<std::pair<int, int> >& blockIndices,
169 auto compute = [&](MarginalCovarianceCholesky& mcc) {
170 mcc.computeCovariance(spinv, A.rowBlockIndices(), blockIndices);
171 };
172 return solveBlocks_impl(A, compute);
173 }
174
176 bool blockOrdering() const { return _blockOrdering; }
178
179 protected:
182
189
190 virtual bool solveBlocks_impl(
192 std::function<void(MarginalCovarianceCholesky&)> compute) = 0;
193};
194
195} // namespace g2o
196
197#endif
Solver with faster iterating structure for the linear matrix.
void initMatrixStructure(const SparseBlockMatrix< MatrixType > &A)
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix
virtual bool solvePattern(SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices, const SparseBlockMatrix< MatrixType > &A)
virtual bool solveBlocks(double **&blocks, const SparseBlockMatrix< MatrixType > &A)
void setBlockOrdering(bool blockOrdering)
virtual bool solveBlocks_impl(const SparseBlockMatrix< MatrixType > &A, std::function< void(MarginalCovarianceCholesky &)> compute)=0
bool blockOrdering() const
do the AMD ordering on the blocks or on the scalar matrix
basic solver for Ax = b
static void allocateBlocks(const SparseBlockMatrix< MatrixType > &A, double **&blocks)
allocate block memory structure
virtual bool solveBlocks(double **&blocks, const SparseBlockMatrix< MatrixType > &A)
virtual ~LinearSolver()
static void deallocateBlocks(const SparseBlockMatrix< MatrixType > &A, double **&blocks)
de-allocate the block structure
void setWriteDebug(bool b)
bool writeDebug() const
write a debug dump of the system matrix if it is not PSD in solve
virtual bool solve(const SparseBlockMatrix< MatrixType > &A, double *x, double *b)=0
virtual bool init()=0
static void blockToScalarPermutation(const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
virtual bool solvePattern(SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices, const SparseBlockMatrix< MatrixType > &A)
computing the marginal covariance given a cholesky factor (lower triangle of the factor)
Sparse matrix which uses blocks.
Sparse matrix which uses blocks.
int colBaseOfBlock(int c) const
where does the col at block-col r starts?
int rowsOfBlock(int r) const
how many rows does the block at block-row r has?
int fillSparseBlockMatrixCCS(SparseBlockMatrixCCS< MatrixType > &blockCCS) const
const std::vector< int > & colBlockIndices() const
indices of the column blocks
int rows() const
rows of the matrix
int cols() const
columns of the matrix
const std::vector< int > & rowBlockIndices() const
indices of the row blocks
int colsOfBlock(int c) const
how many cols does the block at block-col c has?