-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathlinear_solver_eigen.h
237 lines (209 loc) · 8.89 KB
/
linear_solver_eigen.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
// g2o - General Graph Optimization
// Copyright (C) 2011 R. Kuemmerle, G. Grisetti, W. Burgard
// All rights reserved.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// * Redistributions of source code must retain the above copyright notice,
// this list of conditions and the following disclaimer.
// * Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
// IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
// TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
// PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
// HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
// TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#ifndef G2O_LINEAR_SOLVER_EIGEN_H
#define G2O_LINEAR_SOLVER_EIGEN_H
#include <Eigen/Sparse>
#include <Eigen/SparseCholesky>
#include "../core/linear_solver.h"
#include "../core/batch_stats.h"
#include "../stuff/timeutil.h"
#include "../core/eigen_types.h"
#include <iostream>
#include <vector>
namespace g2o {
/**
* \brief 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, I guess.
*/
template <typename MatrixType>
class LinearSolverEigen: public LinearSolver<MatrixType>
{
public:
typedef Eigen::SparseMatrix<double, Eigen::ColMajor> SparseMatrix;
typedef Eigen::Triplet<double> Triplet;
typedef Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic, SparseMatrix::StorageIndex> PermutationMatrix;
/**
* \brief Sub-classing Eigen's SimplicialLDLT to perform ordering with a given ordering
*/
class CholeskyDecomposition : public Eigen::SimplicialLDLT<SparseMatrix, Eigen::Upper>
{
public:
CholeskyDecomposition() : Eigen::SimplicialLDLT<SparseMatrix, Eigen::Upper>() {}
using Eigen::SimplicialLDLT< SparseMatrix, Eigen::Upper>::analyzePattern_preordered;
void analyzePatternWithPermutation(SparseMatrix& a, const PermutationMatrix& permutation)
{
m_Pinv = permutation;
m_P = permutation.inverse();
int size = a.cols();
SparseMatrix ap(size, size);
ap.selfadjointView<Eigen::Upper>() = a.selfadjointView<UpLo>().twistedBy(m_P);
analyzePattern_preordered(ap, true);
}
};
public:
LinearSolverEigen() :
LinearSolver<MatrixType>(),
_init(true), _blockOrdering(false), _writeDebug(false)
{
}
virtual ~LinearSolverEigen()
{
}
virtual bool init()
{
_init = true;
return true;
}
bool solve(const SparseBlockMatrix<MatrixType>& A, double* x, double* b)
{
if (_init)
_sparseMatrix.resize(A.rows(), A.cols());
fillSparseMatrix(A, !_init);
if (_init) // compute the symbolic composition once
computeSymbolicDecomposition(A);
_init = false;
double t=get_monotonic_time();
_cholesky.factorize(_sparseMatrix);
if (_cholesky.info() != Eigen::Success) { // the matrix is not positive definite
if (_writeDebug) {
std::cerr << "Cholesky failure, writing debug.txt (Hessian loadable by Octave)" << std::endl;
A.writeOctave("debug.txt");
}
return false;
}
// Solving the system
VectorXD::MapType xx(x, _sparseMatrix.cols());
VectorXD::ConstMapType bb(b, _sparseMatrix.cols());
xx = _cholesky.solve(bb);
G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
if (globalStats) {
globalStats->timeNumericDecomposition = get_monotonic_time() - t;
globalStats->choleskyNNZ = _cholesky.matrixL().nestedExpression().nonZeros() + _sparseMatrix.cols(); // the elements of D
}
return true;
}
//! do the AMD ordering on the blocks or on the scalar matrix
bool blockOrdering() const { return _blockOrdering;}
void setBlockOrdering(bool blockOrdering) { _blockOrdering = blockOrdering;}
//! write a debug dump of the system matrix if it is not SPD in solve
virtual bool writeDebug() const { return _writeDebug;}
virtual void setWriteDebug(bool b) { _writeDebug = b;}
protected:
bool _init;
bool _blockOrdering;
bool _writeDebug;
SparseMatrix _sparseMatrix;
CholeskyDecomposition _cholesky;
/**
* compute the symbolic decompostion 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.
*/
void computeSymbolicDecomposition(const SparseBlockMatrix<MatrixType>& A)
{
double t=get_monotonic_time();
if (! _blockOrdering) {
_cholesky.analyzePattern(_sparseMatrix);
} else {
// block ordering with the Eigen Interface
// This is really ugly currently, as it calls internal functions from Eigen
// and modifies the SparseMatrix class
Eigen::PermutationMatrix<Eigen::Dynamic,Eigen::Dynamic> blockP;
{
// prepare a block structure matrix for calling AMD
std::vector<Triplet> triplets;
for (size_t c = 0; c < A.blockCols().size(); ++c){
const typename SparseBlockMatrix<MatrixType>::IntBlockMap& column = A.blockCols()[c];
for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it = column.begin(); it != column.end(); ++it) {
const int& r = it->first;
if (r > static_cast<int>(c)) // only upper triangle
break;
triplets.push_back(Triplet(r, c, 0.));
}
}
// call the AMD ordering on the block matrix.
// Relies on Eigen's internal stuff, probably bad idea
SparseMatrix auxBlockMatrix(A.blockCols().size(), A.blockCols().size());
auxBlockMatrix.setFromTriplets(triplets.begin(), triplets.end());
typename CholeskyDecomposition::CholMatrixType C;
C = auxBlockMatrix.selfadjointView<Eigen::Upper>();
Eigen::internal::minimum_degree_ordering(C, blockP);
}
int rows = A.rows();
assert(rows == A.cols() && "Matrix A is not square");
// Adapt the block permutation to the scalar matrix
PermutationMatrix scalarP;
scalarP.resize(rows);
int scalarIdx = 0;
for (int i = 0; i < blockP.size(); ++i) {
const int& p = blockP.indices()(i);
int base = A.colBaseOfBlock(p);
int nCols = A.colsOfBlock(p);
for (int j = 0; j < nCols; ++j)
scalarP.indices()(scalarIdx++) = base++;
}
assert(scalarIdx == rows && "did not completely fill the permutation matrix");
// analyze with the scalar permutation
_cholesky.analyzePatternWithPermutation(_sparseMatrix, scalarP);
}
G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
if (globalStats)
globalStats->timeSymbolicDecomposition = get_monotonic_time() - t;
}
void fillSparseMatrix(const SparseBlockMatrix<MatrixType>& A, bool onlyValues)
{
if (onlyValues) {
A.fillCCS(_sparseMatrix.valuePtr(), true);
} else {
// create from triplet structure
std::vector<Triplet> triplets;
triplets.reserve(A.nonZeros());
for (size_t c = 0; c < A.blockCols().size(); ++c) {
int colBaseOfBlock = A.colBaseOfBlock(c);
const typename SparseBlockMatrix<MatrixType>::IntBlockMap& column = A.blockCols()[c];
for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it = column.begin(); it != column.end(); ++it) {
int rowBaseOfBlock = A.rowBaseOfBlock(it->first);
const MatrixType& m = *(it->second);
for (int cc = 0; cc < m.cols(); ++cc) {
int aux_c = colBaseOfBlock + cc;
for (int rr = 0; rr < m.rows(); ++rr) {
int aux_r = rowBaseOfBlock + rr;
if (aux_r > aux_c)
break;
triplets.push_back(Triplet(aux_r, aux_c, m(rr, cc)));
}
}
}
}
_sparseMatrix.setFromTriplets(triplets.begin(), triplets.end());
}
}
};
} // end namespace
#endif