g2o
Loading...
Searching...
No Matches
linear_solver_dense.h
Go to the documentation of this file.
1// g2o - General Graph Optimization
2// Copyright (C) 2011 H. Strasdat
3// Copyright (C) 2012 R. Kümmerle
4// All rights reserved.
5//
6// Redistribution and use in source and binary forms, with or without
7// modification, are permitted provided that the following conditions are
8// met:
9//
10// * Redistributions of source code must retain the above copyright notice,
11// this list of conditions and the following disclaimer.
12// * Redistributions in binary form must reproduce the above copyright
13// notice, this list of conditions and the following disclaimer in the
14// documentation and/or other materials provided with the distribution.
15//
16// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
17// IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
18// TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
19// PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
20// HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
21// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
22// TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
23// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
24// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
25// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
26// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
27
28#ifndef G2O_LINEAR_SOLVER_DENSE_H
29#define G2O_LINEAR_SOLVER_DENSE_H
30
31#include <Eigen/Cholesky>
32#include <Eigen/Core>
33#include <cassert>
34#include <utility>
35#include <vector>
36
39
40namespace g2o {
41
45template <typename MatrixType>
46class LinearSolverDense : public LinearSolver<MatrixType> {
47 public:
48 LinearSolverDense() : LinearSolver<MatrixType>(), _reset(true) {}
49
50 virtual ~LinearSolverDense() {}
51
52 virtual bool init() {
53 _reset = true;
54 return true;
55 }
56
57 bool solve(const SparseBlockMatrix<MatrixType>& A, double* x, double* b) {
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) {
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 }
106
107 protected:
108 bool _reset;
110 Eigen::LDLT<MatrixX> _cholesky;
111};
112
113} // namespace g2o
114
115#endif
linear solver using dense cholesky decomposition
Eigen::LDLT< MatrixX > _cholesky
bool solve(const SparseBlockMatrix< MatrixType > &A, double *x, double *b)
basic solver for Ax = b
Sparse matrix which uses blocks.
bool transpose(SparseBlockMatrix< MatrixTransposedType > &dest) const
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?
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
std::map< int, SparseMatrixBlock * > IntBlockMap
int cols() const
columns of the matrix
int rowBaseOfBlock(int r) const
where does the row at block-row r starts?
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
MatrixN< Eigen::Dynamic > MatrixX
Definition eigen_types.h:74