g2o
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
g2o::BlockSolver< Traits > Class Template Reference

Implementation of a solver operating on the blocks of the Hessian. More...

#include <block_solver.h>

Inheritance diagram for g2o::BlockSolver< Traits >:
Inheritance graph
[legend]
Collaboration diagram for g2o::BlockSolver< Traits >:
Collaboration graph
[legend]

Public Types

typedef Traits::PoseMatrixType PoseMatrixType
 
typedef Traits::LandmarkMatrixType LandmarkMatrixType
 
typedef Traits::PoseLandmarkMatrixType PoseLandmarkMatrixType
 
typedef Traits::PoseVectorType PoseVectorType
 
typedef Traits::LandmarkVectorType LandmarkVectorType
 
typedef Traits::PoseHessianType PoseHessianType
 
typedef Traits::LandmarkHessianType LandmarkHessianType
 
typedef Traits::PoseLandmarkHessianType PoseLandmarkHessianType
 
typedef Traits::LinearSolverType LinearSolverType
 

Public Member Functions

 BlockSolver (std::unique_ptr< LinearSolverType > linearSolver)
 
 ~BlockSolver ()
 
virtual bool init (SparseOptimizer *optmizer, bool online=false)
 
virtual bool buildStructure (bool zeroBlocks=false)
 
virtual bool updateStructure (const std::vector< HyperGraph::Vertex * > &vset, const HyperGraph::EdgeSet &edges)
 
virtual bool buildSystem ()
 
virtual bool solve ()
 
virtual bool computeMarginals (SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices)
 
virtual bool setLambda (double lambda, bool backup=false)
 
virtual void restoreDiagonal ()
 
virtual bool supportsSchur ()
 
virtual bool schur ()
 should the solver perform the schur complement or not
 
virtual void setSchur (bool s)
 
LinearSolver< PoseMatrixType > & linearSolver () const
 
virtual void setWriteDebug (bool writeDebug)
 
virtual bool writeDebug () const
 
virtual bool saveHessian (const std::string &fileName) const
 write the hessian to disk using the specified file name
 
virtual void multiplyHessian (double *dest, const double *src) const
 
- Public Member Functions inherited from g2o::BlockSolverBase
virtual ~BlockSolverBase ()
 
- Public Member Functions inherited from g2o::Solver
 Solver ()
 
virtual ~Solver ()
 
double * x ()
 return x, the solution vector
 
const double * x () const
 
double * b ()
 return b, the right hand side of the system
 
const double * b () const
 
size_t vectorSize () const
 return the size of the solution vector (x) and b
 
SparseOptimizeroptimizer () const
 the optimizer (graph) on which the solver works
 
void setOptimizer (SparseOptimizer *optimizer)
 
bool levenberg () const
 the system is Levenberg-Marquardt
 
void setLevenberg (bool levenberg)
 
size_t additionalVectorSpace () const
 
void setAdditionalVectorSpace (size_t s)
 

Static Public Attributes

static const int PoseDim = Traits::PoseDim
 
static const int LandmarkDim = Traits::LandmarkDim
 

Protected Member Functions

void resize (int *blockPoseIndices, int numPoseBlocks, int *blockLandmarkIndices, int numLandmarkBlocks, int totalDim)
 
void deallocate ()
 
- Protected Member Functions inherited from g2o::Solver
void resizeVector (size_t sx)
 

Protected Attributes

std::unique_ptr< SparseBlockMatrix< PoseMatrixType > > _Hpp
 
std::unique_ptr< SparseBlockMatrix< LandmarkMatrixType > > _Hll
 
std::unique_ptr< SparseBlockMatrix< PoseLandmarkMatrixType > > _Hpl
 
std::unique_ptr< SparseBlockMatrix< PoseMatrixType > > _Hschur
 
std::unique_ptr< SparseBlockMatrixDiagonal< LandmarkMatrixType > > _DInvSchur
 
std::unique_ptr< SparseBlockMatrixCCS< PoseLandmarkMatrixType > > _HplCCS
 
std::unique_ptr< SparseBlockMatrixCCS< PoseMatrixType > > _HschurTransposedCCS
 
std::unique_ptr< LinearSolverType_linearSolver
 
std::vector< PoseVectorType_diagonalBackupPose
 
std::vector< LandmarkVectorType_diagonalBackupLandmark
 
bool _doSchur
 
std::unique_ptr< double[], aligned_deleter< double > > _coefficients
 
std::unique_ptr< double[], aligned_deleter< double > > _bschur
 
int _numPoses
 
int _numLandmarks
 
int _sizePoses
 
int _sizeLandmarks
 
- Protected Attributes inherited from g2o::Solver
SparseOptimizer_optimizer
 
double * _x
 
double * _b
 
size_t _xSize
 
size_t _maxXSize
 
bool _isLevenberg
 
size_t _additionalVectorSpace
 

Detailed Description

template<typename Traits>
class g2o::BlockSolver< Traits >

Implementation of a solver operating on the blocks of the Hessian.

Definition at line 103 of file block_solver.h.

Member Typedef Documentation

◆ LandmarkHessianType

template<typename Traits >
typedef Traits::LandmarkHessianType g2o::BlockSolver< Traits >::LandmarkHessianType

Definition at line 114 of file block_solver.h.

◆ LandmarkMatrixType

template<typename Traits >
typedef Traits::LandmarkMatrixType g2o::BlockSolver< Traits >::LandmarkMatrixType

Definition at line 108 of file block_solver.h.

◆ LandmarkVectorType

template<typename Traits >
typedef Traits::LandmarkVectorType g2o::BlockSolver< Traits >::LandmarkVectorType

Definition at line 111 of file block_solver.h.

◆ LinearSolverType

template<typename Traits >
typedef Traits::LinearSolverType g2o::BlockSolver< Traits >::LinearSolverType

Definition at line 116 of file block_solver.h.

◆ PoseHessianType

template<typename Traits >
typedef Traits::PoseHessianType g2o::BlockSolver< Traits >::PoseHessianType

Definition at line 113 of file block_solver.h.

◆ PoseLandmarkHessianType

template<typename Traits >
typedef Traits::PoseLandmarkHessianType g2o::BlockSolver< Traits >::PoseLandmarkHessianType

Definition at line 115 of file block_solver.h.

◆ PoseLandmarkMatrixType

template<typename Traits >
typedef Traits::PoseLandmarkMatrixType g2o::BlockSolver< Traits >::PoseLandmarkMatrixType

Definition at line 109 of file block_solver.h.

◆ PoseMatrixType

template<typename Traits >
typedef Traits::PoseMatrixType g2o::BlockSolver< Traits >::PoseMatrixType

Definition at line 107 of file block_solver.h.

◆ PoseVectorType

template<typename Traits >
typedef Traits::PoseVectorType g2o::BlockSolver< Traits >::PoseVectorType

Definition at line 110 of file block_solver.h.

Constructor & Destructor Documentation

◆ BlockSolver()

template<typename Traits >
g2o::BlockSolver< Traits >::BlockSolver ( std::unique_ptr< LinearSolverType linearSolver)

allocate a block solver ontop of the underlying linear solver. NOTE: The BlockSolver assumes exclusive access to the linear solver and will therefore free the pointer in its destructor.

Definition at line 41 of file block_solver.hpp.

42 : BlockSolverBase(), _linearSolver(std::move(linearSolver)) {
43 // workspace
44 _xSize = 0;
45 _numPoses = 0;
46 _numLandmarks = 0;
47 _sizePoses = 0;
49 _doSchur = true;
50}
std::unique_ptr< LinearSolverType > _linearSolver
LinearSolver< PoseMatrixType > & linearSolver() const
size_t _xSize
Definition solver.h:142

References g2o::BlockSolver< Traits >::_doSchur, g2o::BlockSolver< Traits >::_numLandmarks, g2o::BlockSolver< Traits >::_numPoses, g2o::BlockSolver< Traits >::_sizeLandmarks, g2o::BlockSolver< Traits >::_sizePoses, and g2o::Solver::_xSize.

◆ ~BlockSolver()

template<typename Traits >
g2o::BlockSolver< Traits >::~BlockSolver ( )

Definition at line 107 of file block_solver.hpp.

107{}

Member Function Documentation

◆ buildStructure()

template<typename Traits >
bool g2o::BlockSolver< Traits >::buildStructure ( bool  zeroBlocks = false)
virtual

build the structure of the system

Implements g2o::Solver.

Definition at line 110 of file block_solver.hpp.

110 {
111 assert(_optimizer);
112
113 size_t sparseDim = 0;
114 _numPoses = 0;
115 _numLandmarks = 0;
116 _sizePoses = 0;
117 _sizeLandmarks = 0;
118 int* blockPoseIndices = new int[_optimizer->indexMapping().size()];
119 int* blockLandmarkIndices = new int[_optimizer->indexMapping().size()];
120
121 for (size_t i = 0; i < _optimizer->indexMapping().size(); ++i) {
122 OptimizableGraph::Vertex* v = _optimizer->indexMapping()[i];
123 int dim = v->dimension();
124 if (!v->marginalized()) {
125 v->setColInHessian(_sizePoses);
126 _sizePoses += dim;
127 blockPoseIndices[_numPoses] = _sizePoses;
128 ++_numPoses;
129 } else {
130 v->setColInHessian(_sizeLandmarks);
131 _sizeLandmarks += dim;
132 blockLandmarkIndices[_numLandmarks] = _sizeLandmarks;
134 }
135 sparseDim += dim;
136 }
137 resize(blockPoseIndices, _numPoses, blockLandmarkIndices, _numLandmarks,
138 sparseDim);
139 delete[] blockLandmarkIndices;
140 delete[] blockPoseIndices;
141
142 // allocate the diagonal on Hpp and Hll
143 int poseIdx = 0;
144 int landmarkIdx = 0;
145 for (size_t i = 0; i < _optimizer->indexMapping().size(); ++i) {
146 OptimizableGraph::Vertex* v = _optimizer->indexMapping()[i];
147 if (!v->marginalized()) {
148 // assert(poseIdx == v->hessianIndex());
149 PoseMatrixType* m = _Hpp->block(poseIdx, poseIdx, true);
150 if (zeroBlocks) m->setZero();
151 v->mapHessianMemory(m->data());
152 ++poseIdx;
153 } else {
154 LandmarkMatrixType* m = _Hll->block(landmarkIdx, landmarkIdx, true);
155 if (zeroBlocks) m->setZero();
156 v->mapHessianMemory(m->data());
157 ++landmarkIdx;
158 }
159 }
160 assert(poseIdx == _numPoses && landmarkIdx == _numLandmarks);
161
162 // temporary structures for building the pattern of the Schur complement
163 SparseBlockMatrixHashMap<PoseMatrixType>* schurMatrixLookup = 0;
164 if (_doSchur) {
165 schurMatrixLookup = new SparseBlockMatrixHashMap<PoseMatrixType>(
166 _Hschur->rowBlockIndices(), _Hschur->colBlockIndices());
167 schurMatrixLookup->blockCols().resize(_Hschur->blockCols().size());
168 }
169
170 // here we assume that the landmark indices start after the pose ones
171 // create the structure in Hpp, Hll and in Hpl
172 for (SparseOptimizer::EdgeContainer::const_iterator it =
173 _optimizer->activeEdges().begin();
174 it != _optimizer->activeEdges().end(); ++it) {
175 OptimizableGraph::Edge* e = *it;
176
177 for (size_t viIdx = 0; viIdx < e->vertices().size(); ++viIdx) {
178 OptimizableGraph::Vertex* v1 =
179 (OptimizableGraph::Vertex*)e->vertex(viIdx);
180 int ind1 = v1->hessianIndex();
181 if (ind1 == -1) continue;
182 int indexV1Bak = ind1;
183 for (size_t vjIdx = viIdx + 1; vjIdx < e->vertices().size(); ++vjIdx) {
184 OptimizableGraph::Vertex* v2 =
185 (OptimizableGraph::Vertex*)e->vertex(vjIdx);
186 int ind2 = v2->hessianIndex();
187 if (ind2 == -1) continue;
188 ind1 = indexV1Bak;
189 bool transposedBlock = ind1 > ind2;
190 if (transposedBlock) { // make sure, we allocate the upper triangle
191 // block
192 std::swap(ind1, ind2);
193 }
194 if (!v1->marginalized() && !v2->marginalized()) {
195 PoseMatrixType* m = _Hpp->block(ind1, ind2, true);
196 if (zeroBlocks) m->setZero();
197 e->mapHessianMemory(m->data(), viIdx, vjIdx, transposedBlock);
198 if (_Hschur) { // assume this is only needed in case we solve with
199 // the schur complement
200 schurMatrixLookup->addBlock(ind1, ind2);
201 }
202 } else if (v1->marginalized() && v2->marginalized()) {
203 // RAINER hmm.... should we ever reach this here????
205 _Hll->block(ind1 - _numPoses, ind2 - _numPoses, true);
206 if (zeroBlocks) m->setZero();
207 e->mapHessianMemory(m->data(), viIdx, vjIdx, false);
208 } else {
209 if (v1->marginalized()) {
210 PoseLandmarkMatrixType* m = _Hpl->block(
211 v2->hessianIndex(), v1->hessianIndex() - _numPoses, true);
212 if (zeroBlocks) m->setZero();
213 e->mapHessianMemory(
214 m->data(), viIdx, vjIdx,
215 true); // transpose the block before writing to it
216 } else {
217 PoseLandmarkMatrixType* m = _Hpl->block(
218 v1->hessianIndex(), v2->hessianIndex() - _numPoses, true);
219 if (zeroBlocks) m->setZero();
220 e->mapHessianMemory(m->data(), viIdx, vjIdx,
221 false); // directly the block
222 }
223 }
224 }
225 }
226 }
227
228 if (!_doSchur) {
229 delete schurMatrixLookup;
230 return true;
231 }
232
233 _DInvSchur->diagonal().resize(landmarkIdx);
234 _Hpl->fillSparseBlockMatrixCCS(*_HplCCS);
235
236 for (OptimizableGraph::Vertex* v : _optimizer->indexMapping()) {
237 if (v->marginalized()) {
238 const HyperGraph::EdgeSet& vedges = v->edges();
239 for (HyperGraph::EdgeSet::const_iterator it1 = vedges.begin();
240 it1 != vedges.end(); ++it1) {
241 for (size_t i = 0; i < (*it1)->vertices().size(); ++i) {
242 OptimizableGraph::Vertex* v1 =
243 (OptimizableGraph::Vertex*)(*it1)->vertex(i);
244 if (v1->hessianIndex() == -1 || v1 == v) continue;
245 for (HyperGraph::EdgeSet::const_iterator it2 = vedges.begin();
246 it2 != vedges.end(); ++it2) {
247 for (size_t j = 0; j < (*it2)->vertices().size(); ++j) {
248 OptimizableGraph::Vertex* v2 =
249 (OptimizableGraph::Vertex*)(*it2)->vertex(j);
250 if (v2->hessianIndex() == -1 || v2 == v) continue;
251 int i1 = v1->hessianIndex();
252 int i2 = v2->hessianIndex();
253 if (i1 <= i2) {
254 schurMatrixLookup->addBlock(i1, i2);
255 }
256 }
257 }
258 }
259 }
260 }
261 }
262
263 _Hschur->takePatternFromHash(*schurMatrixLookup);
264 delete schurMatrixLookup;
265 _Hschur->fillSparseBlockMatrixCCSTransposed(*_HschurTransposedCCS);
266
267 return true;
268}
void resize(int *blockPoseIndices, int numPoseBlocks, int *blockLandmarkIndices, int numLandmarkBlocks, int totalDim)
std::unique_ptr< SparseBlockMatrix< PoseMatrixType > > _Hschur
std::unique_ptr< SparseBlockMatrixCCS< PoseLandmarkMatrixType > > _HplCCS
std::unique_ptr< SparseBlockMatrix< LandmarkMatrixType > > _Hll
std::unique_ptr< SparseBlockMatrix< PoseMatrixType > > _Hpp
Traits::PoseMatrixType PoseMatrixType
Traits::PoseLandmarkMatrixType PoseLandmarkMatrixType
std::unique_ptr< SparseBlockMatrixDiagonal< LandmarkMatrixType > > _DInvSchur
std::unique_ptr< SparseBlockMatrixCCS< PoseMatrixType > > _HschurTransposedCCS
std::unique_ptr< SparseBlockMatrix< PoseLandmarkMatrixType > > _Hpl
Traits::LandmarkMatrixType LandmarkMatrixType
std::set< Edge * > EdgeSet
SparseOptimizer * _optimizer
Definition solver.h:139
const EdgeContainer & activeEdges() const
the edges active in the current optimization
const VertexContainer & indexMapping() const
the index mapping of the vertices

References g2o::SparseBlockMatrixHashMap< MatrixType >::addBlock(), g2o::SparseBlockMatrixHashMap< MatrixType >::blockCols(), g2o::OptimizableGraph::Vertex::dimension(), g2o::OptimizableGraph::Vertex::hessianIndex(), g2o::OptimizableGraph::Vertex::mapHessianMemory(), g2o::OptimizableGraph::Edge::mapHessianMemory(), g2o::OptimizableGraph::Vertex::marginalized(), g2o::OptimizableGraph::Vertex::setColInHessian(), g2o::HyperGraph::Edge::vertex(), and g2o::HyperGraph::Edge::vertices().

◆ buildSystem()

template<typename Traits >
bool g2o::BlockSolver< Traits >::buildSystem ( )
virtual

build the current system

Implements g2o::Solver.

Definition at line 492 of file block_solver.hpp.

492 {
493 // clear b vector
494#ifdef G2O_OPENMP
495#pragma omp parallel for default( \
496 shared) if (_optimizer->indexMapping().size() > 1000)
497#endif
498 for (int i = 0; i < static_cast<int>(_optimizer->indexMapping().size());
499 ++i) {
500 OptimizableGraph::Vertex* v = _optimizer->indexMapping()[i];
501 assert(v);
502 v->clearQuadraticForm();
503 }
504 _Hpp->clear();
505 if (_doSchur) {
506 _Hll->clear();
507 _Hpl->clear();
508 }
509
510 // resetting the terms for the pairwise constraints
511 // built up the current system by storing the Hessian blocks in the edges and
512 // vertices
513#ifndef G2O_OPENMP
514 // no threading, we do not need to copy the workspace
515 JacobianWorkspace& jacobianWorkspace = _optimizer->jacobianWorkspace();
516#else
517 // if running with threads need to produce copies of the workspace for each
518 // thread
519 JacobianWorkspace jacobianWorkspace = _optimizer->jacobianWorkspace();
520#pragma omp parallel for default(shared) firstprivate( \
521 jacobianWorkspace) if (_optimizer->activeEdges().size() > 100)
522#endif
523 for (int k = 0; k < static_cast<int>(_optimizer->activeEdges().size()); ++k) {
524 OptimizableGraph::Edge* e = _optimizer->activeEdges()[k];
525 e->linearizeOplus(
526 jacobianWorkspace); // jacobian of the nodes' oplus (manifold)
527 e->constructQuadraticForm();
528#ifndef NDEBUG
529 for (size_t i = 0; i < e->vertices().size(); ++i) {
530 const OptimizableGraph::Vertex* v =
531 static_cast<const OptimizableGraph::Vertex*>(e->vertex(i));
532 if (!v->fixed()) {
533 bool hasANan = arrayHasNaN(jacobianWorkspace.workspaceForVertex(i),
534 e->dimension() * v->dimension());
535 if (hasANan) {
536 G2O_WARN(
537 "buildSystem(): NaN within Jacobian for edge {} for vertex {}",
538 static_cast<void*>(e), i);
539 break;
540 }
541 }
542 }
543#endif
544 }
545
546 // flush the current system in a sparse block matrix
547#ifdef G2O_OPENMP
548#pragma omp parallel for default( \
549 shared) if (_optimizer->indexMapping().size() > 1000)
550#endif
551 for (int i = 0; i < static_cast<int>(_optimizer->indexMapping().size());
552 ++i) {
553 OptimizableGraph::Vertex* v = _optimizer->indexMapping()[i];
554 int iBase = v->colInHessian();
555 if (v->marginalized()) iBase += _sizePoses;
556 v->copyB(_b + iBase);
557 }
558
559 return false;
560}
double * _b
Definition solver.h:141
#define G2O_WARN(...)
Definition logger.h:88
bool arrayHasNaN(const double *array, int size, int *nanIndex=0)
Definition misc.h:168
JacobianWorkspace & jacobianWorkspace()
the workspace for storing the Jacobians of the graph

References g2o::arrayHasNaN(), g2o::OptimizableGraph::Vertex::clearQuadraticForm(), g2o::OptimizableGraph::Vertex::colInHessian(), g2o::OptimizableGraph::Edge::constructQuadraticForm(), g2o::OptimizableGraph::Vertex::copyB(), g2o::OptimizableGraph::Vertex::dimension(), g2o::OptimizableGraph::Edge::dimension(), g2o::OptimizableGraph::Vertex::fixed(), G2O_WARN, g2o::OptimizableGraph::Edge::linearizeOplus(), g2o::OptimizableGraph::Vertex::marginalized(), g2o::HyperGraph::Edge::vertex(), g2o::HyperGraph::Edge::vertices(), and g2o::JacobianWorkspace::workspaceForVertex().

◆ computeMarginals()

template<typename Traits >
bool g2o::BlockSolver< Traits >::computeMarginals ( SparseBlockMatrix< MatrixX > &  spinv,
const std::vector< std::pair< int, int > > &  blockIndices 
)
virtual

computes the block diagonal elements of the pattern specified in the input and stores them in given SparseBlockMatrix

Implements g2o::Solver.

Definition at line 479 of file block_solver.hpp.

481 {
482 double t = get_monotonic_time();
483 bool ok = _linearSolver->solvePattern(spinv, blockIndices, *_Hpp);
484 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
485 if (globalStats) {
486 globalStats->timeMarginals = get_monotonic_time() - t;
487 }
488 return ok;
489}
double get_monotonic_time()
Definition timeutil.cpp:43
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77

References g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), and g2o::G2OBatchStatistics::timeMarginals.

◆ deallocate()

template<typename Traits >
void g2o::BlockSolver< Traits >::deallocate ( )
protected

Definition at line 93 of file block_solver.hpp.

93 {
94 _Hpp.reset();
95 _Hll.reset();
96 _Hpl.reset();
97 _Hschur.reset();
98 _DInvSchur.reset();
99 _coefficients.reset();
100 _bschur.reset();
101
102 _HplCCS.reset();
103 _HschurTransposedCCS.reset();
104}
std::unique_ptr< double[], aligned_deleter< double > > _bschur
std::unique_ptr< double[], aligned_deleter< double > > _coefficients

◆ init()

template<typename Traits >
bool g2o::BlockSolver< Traits >::init ( SparseOptimizer optimizer,
bool  online = false 
)
virtual

initialize the solver, called once before the first iteration

Implements g2o::Solver.

Definition at line 604 of file block_solver.hpp.

604 {
606 if (!online) {
607 if (_Hpp) _Hpp->clear();
608 if (_Hpl) _Hpl->clear();
609 if (_Hll) _Hll->clear();
610 }
611 _linearSolver->init();
612 return true;
613}
SparseOptimizer * optimizer() const
the optimizer (graph) on which the solver works
Definition solver.h:108

References g2o::SparseOptimizer::clear().

◆ linearSolver()

template<typename Traits >
LinearSolver< PoseMatrixType > & g2o::BlockSolver< Traits >::linearSolver ( ) const
inline

◆ multiplyHessian()

template<typename Traits >
virtual void g2o::BlockSolver< Traits >::multiplyHessian ( double *  dest,
const double *  src 
) const
inlinevirtual

compute dest = H * src

Implements g2o::BlockSolverBase.

Definition at line 149 of file block_solver.h.

149 {
150 _Hpp->multiplySymmetricUpperTriangle(dest, src);
151 }

References g2o::BlockSolver< Traits >::_Hpp.

◆ resize()

template<typename Traits >
void g2o::BlockSolver< Traits >::resize ( int *  blockPoseIndices,
int  numPoseBlocks,
int *  blockLandmarkIndices,
int  numLandmarkBlocks,
int  totalDim 
)
protected

Definition at line 53 of file block_solver.hpp.

55 {
56 deallocate();
57
58 resizeVector(s);
59
60 if (_doSchur) {
61 // the following two are only used in schur
62 assert(_sizePoses > 0 && "allocating with wrong size");
63 _coefficients.reset(allocate_aligned<double>(s));
64 _bschur.reset(allocate_aligned<double>(_sizePoses));
65 }
66
67 _Hpp = std::make_unique<PoseHessianType>(blockPoseIndices, blockPoseIndices,
68 numPoseBlocks, numPoseBlocks);
69 if (_doSchur) {
70 _Hschur = std::make_unique<PoseHessianType>(
71 blockPoseIndices, blockPoseIndices, numPoseBlocks, numPoseBlocks);
72 _Hll = std::make_unique<LandmarkHessianType>(
73 blockLandmarkIndices, blockLandmarkIndices, numLandmarkBlocks,
74 numLandmarkBlocks);
76 std::make_unique<SparseBlockMatrixDiagonal<LandmarkMatrixType>>(
77 _Hll->colBlockIndices());
78 _Hpl = std::make_unique<PoseLandmarkHessianType>(
79 blockPoseIndices, blockLandmarkIndices, numPoseBlocks,
80 numLandmarkBlocks);
81 _HplCCS = std::make_unique<SparseBlockMatrixCCS<PoseLandmarkMatrixType>>(
82 _Hpl->rowBlockIndices(), _Hpl->colBlockIndices());
84 std::make_unique<SparseBlockMatrixCCS<PoseMatrixType>>(
85 _Hschur->colBlockIndices(), _Hschur->rowBlockIndices());
86#ifdef G2O_OPENMP
87 _coefficientsMutex.resize(numPoseBlocks);
88#endif
89 }
90}
void resizeVector(size_t sx)
Definition solver.cpp:50

◆ restoreDiagonal()

template<typename Traits >
void g2o::BlockSolver< Traits >::restoreDiagonal ( )
virtual

restore a previously made backup of the diagonal

Implements g2o::Solver.

Definition at line 588 of file block_solver.hpp.

588 {
589 assert((int)_diagonalBackupPose.size() == _numPoses &&
590 "Mismatch in dimensions");
591 assert((int)_diagonalBackupLandmark.size() == _numLandmarks &&
592 "Mismatch in dimensions");
593 for (int i = 0; i < _numPoses; ++i) {
594 PoseMatrixType* b = _Hpp->block(i, i);
595 b->diagonal() = _diagonalBackupPose[i];
596 }
597 for (int i = 0; i < _numLandmarks; ++i) {
598 LandmarkMatrixType* b = _Hll->block(i, i);
599 b->diagonal() = _diagonalBackupLandmark[i];
600 }
601}
std::vector< LandmarkVectorType > _diagonalBackupLandmark
std::vector< PoseVectorType > _diagonalBackupPose
double * b()
return b, the right hand side of the system
Definition solver.h:101

◆ saveHessian()

template<typename Traits >
bool g2o::BlockSolver< Traits >::saveHessian ( const std::string &  ) const
virtual

write the hessian to disk using the specified file name

Implements g2o::Solver.

Definition at line 621 of file block_solver.hpp.

621 {
622 return _Hpp->writeOctave(fileName.c_str(), true);
623}

◆ schur()

template<typename Traits >
virtual bool g2o::BlockSolver< Traits >::schur ( )
inlinevirtual

should the solver perform the schur complement or not

Implements g2o::Solver.

Definition at line 139 of file block_solver.h.

139{ return _doSchur; }

References g2o::BlockSolver< Traits >::_doSchur.

◆ setLambda()

template<typename Traits >
bool g2o::BlockSolver< Traits >::setLambda ( double  lambda,
bool  backup = false 
)
virtual

update the system while performing Levenberg, i.e., modifying the diagonal components of A by doing += lambda along the main diagonal of the Matrix. Note that this function may be called with a positive and a negative lambda. The latter is used to undo a former modification. If backup is true, then the solver should store a backup of the diagonal, which can be restored by restoreDiagonal()

Implements g2o::Solver.

Definition at line 563 of file block_solver.hpp.

563 {
564 if (backup) {
567 }
568#ifdef G2O_OPENMP
569#pragma omp parallel for default(shared) if (_numPoses > 100)
570#endif
571 for (int i = 0; i < _numPoses; ++i) {
572 PoseMatrixType* b = _Hpp->block(i, i);
573 if (backup) _diagonalBackupPose[i] = b->diagonal();
574 b->diagonal().array() += lambda;
575 }
576#ifdef G2O_OPENMP
577#pragma omp parallel for default(shared) if (_numLandmarks > 100)
578#endif
579 for (int i = 0; i < _numLandmarks; ++i) {
580 LandmarkMatrixType* b = _Hll->block(i, i);
581 if (backup) _diagonalBackupLandmark[i] = b->diagonal();
582 b->diagonal().array() += lambda;
583 }
584 return true;
585}

◆ setSchur()

template<typename Traits >
virtual void g2o::BlockSolver< Traits >::setSchur ( bool  s)
inlinevirtual

Implements g2o::Solver.

Definition at line 140 of file block_solver.h.

140{ _doSchur = s; }

References g2o::BlockSolver< Traits >::_doSchur.

Referenced by g2o::SparseOptimizerIncremental::initSolver().

◆ setWriteDebug()

template<typename Traits >
void g2o::BlockSolver< Traits >::setWriteDebug ( bool  )
virtual

write debug output of the Hessian if system is not positive definite

Implements g2o::Solver.

Definition at line 616 of file block_solver.hpp.

616 {
617 _linearSolver->setWriteDebug(writeDebug);
618}
virtual bool writeDebug() const

◆ solve()

template<typename Traits >
bool g2o::BlockSolver< Traits >::solve ( )
virtual

solve Ax = b

Implements g2o::Solver.

Definition at line 331 of file block_solver.hpp.

331 {
332 if (!_doSchur) {
333 double t = get_monotonic_time();
334 bool ok = _linearSolver->solve(*_Hpp, _x, _b);
335 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
336 if (globalStats) {
337 globalStats->timeLinearSolver = get_monotonic_time() - t;
338 globalStats->hessianDimension = globalStats->hessianPoseDimension =
339 _Hpp->cols();
340 }
341 return ok;
342 }
343
344 // schur thing
345
346 // backup the coefficient matrix
347 double t = get_monotonic_time();
348
349 // _Hschur = _Hpp, but keeping the pattern of _Hschur
350 _Hschur->clear();
351 _Hpp->add(*_Hschur);
352
353 //_DInvSchur->clear();
354 memset(_coefficients.get(), 0, _sizePoses * sizeof(double));
355#ifdef G2O_OPENMP
356#pragma omp parallel for default(shared) schedule(dynamic, 10)
357#endif
358 for (int landmarkIndex = 0;
359 landmarkIndex < static_cast<int>(_Hll->blockCols().size());
360 ++landmarkIndex) {
362 marginalizeColumn = _Hll->blockCols()[landmarkIndex];
363 assert(marginalizeColumn.size() == 1 &&
364 "more than one block in _Hll column");
365
366 // calculate inverse block for the landmark
367 const LandmarkMatrixType* D = marginalizeColumn.begin()->second;
368 assert(D && D->rows() == D->cols() && "Error in landmark matrix");
369 LandmarkMatrixType& Dinv = _DInvSchur->diagonal()[landmarkIndex];
370 Dinv = D->inverse();
371
372 LandmarkVectorType db(D->rows());
373 for (int j = 0; j < D->rows(); ++j) {
374 db[j] = _b[_Hll->rowBaseOfBlock(landmarkIndex) + _sizePoses + j];
375 }
376 db = Dinv * db;
377
378 assert((size_t)landmarkIndex < _HplCCS->blockCols().size() &&
379 "Index out of bounds");
381 landmarkColumn = _HplCCS->blockCols()[landmarkIndex];
382
383 for (typename SparseBlockMatrixCCS<
384 PoseLandmarkMatrixType>::SparseColumn::const_iterator it_outer =
385 landmarkColumn.begin();
386 it_outer != landmarkColumn.end(); ++it_outer) {
387 int i1 = it_outer->row;
388
389 const PoseLandmarkMatrixType* Bi = it_outer->block;
390 assert(Bi);
391
392 PoseLandmarkMatrixType BDinv = (*Bi) * (Dinv);
393 assert(_HplCCS->rowBaseOfBlock(i1) < _sizePoses && "Index out of bounds");
394 typename PoseVectorType::MapType Bb(
395 &_coefficients[_HplCCS->rowBaseOfBlock(i1)], Bi->rows());
396#ifdef G2O_OPENMP
397 ScopedOpenMPMutex mutexLock(&_coefficientsMutex[i1]);
398#endif
399 Bb.noalias() += (*Bi) * db;
400
401 assert(i1 >= 0 &&
402 i1 < static_cast<int>(_HschurTransposedCCS->blockCols().size()) &&
403 "Index out of bounds");
404 typename SparseBlockMatrixCCS<PoseMatrixType>::SparseColumn::iterator
405 targetColumnIt = _HschurTransposedCCS->blockCols()[i1].begin();
406
407 typename SparseBlockMatrixCCS<PoseLandmarkMatrixType>::RowBlock aux(i1,
408 0);
409 typename SparseBlockMatrixCCS<
410 PoseLandmarkMatrixType>::SparseColumn::const_iterator it_inner =
411 lower_bound(landmarkColumn.begin(), landmarkColumn.end(), aux);
412 for (; it_inner != landmarkColumn.end(); ++it_inner) {
413 int i2 = it_inner->row;
414 const PoseLandmarkMatrixType* Bj = it_inner->block;
415 assert(Bj);
416 while (targetColumnIt->row < i2 /*&& targetColumnIt != _HschurTransposedCCS->blockCols()[i1].end()*/)
417 ++targetColumnIt;
418 assert(targetColumnIt != _HschurTransposedCCS->blockCols()[i1].end() &&
419 targetColumnIt->row == i2 &&
420 "invalid iterator, something wrong with the matrix structure");
421 PoseMatrixType* Hi1i2 = targetColumnIt->block; //_Hschur->block(i1,i2);
422 assert(Hi1i2);
423 (*Hi1i2).noalias() -= BDinv * Bj->transpose();
424 }
425 }
426 }
427
428 // _bschur = _b for calling solver, and not touching _b
429 memcpy(_bschur.get(), _b, _sizePoses * sizeof(double));
430 for (int i = 0; i < _sizePoses; ++i) {
431 _bschur[i] -= _coefficients[i];
432 }
433
434 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
435 if (globalStats) {
436 globalStats->timeSchurComplement = get_monotonic_time() - t;
437 }
438
439 t = get_monotonic_time();
440 bool solvedPoses = _linearSolver->solve(*_Hschur, _x, _bschur.get());
441 if (globalStats) {
442 globalStats->timeLinearSolver = get_monotonic_time() - t;
443 globalStats->hessianPoseDimension = _Hpp->cols();
444 globalStats->hessianLandmarkDimension = _Hll->cols();
445 globalStats->hessianDimension = globalStats->hessianPoseDimension +
446 globalStats->hessianLandmarkDimension;
447 }
448
449 if (!solvedPoses) return false;
450
451 // _x contains the solution for the poses, now applying it to the landmarks to
452 // get the new part of the solution;
453 double* xp = _x;
454 double* cp = _coefficients.get();
455
456 double* xl = _x + _sizePoses;
457 double* cl = _coefficients.get() + _sizePoses;
458 double* bl = _b + _sizePoses;
459
460 // cp = -xp
461 for (int i = 0; i < _sizePoses; ++i) cp[i] = -xp[i];
462
463 // cl = bl
464 memcpy(cl, bl, _sizeLandmarks * sizeof(double));
465
466 // cl = bl - Bt * xp
467 // Bt->multiply(cl, cp);
468 _HplCCS->rightMultiply(cl, cp);
469
470 // xl = Dinv * cl
471 memset(xl, 0, _sizeLandmarks * sizeof(double));
472 _DInvSchur->multiply(xl, cl);
473 //_DInvSchur->rightMultiply(xl,cl);
474
475 return true;
476}
Traits::LandmarkVectorType LandmarkVectorType
double * _x
Definition solver.h:140
std::vector< RowBlock > SparseColumn
std::map< int, SparseMatrixBlock * > IntBlockMap

References g2o::SparseBlockMatrixCCS< MatrixType >::blockCols(), g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), g2o::G2OBatchStatistics::hessianDimension, g2o::G2OBatchStatistics::hessianLandmarkDimension, g2o::G2OBatchStatistics::hessianPoseDimension, g2o::G2OBatchStatistics::timeLinearSolver, and g2o::G2OBatchStatistics::timeSchurComplement.

◆ supportsSchur()

template<typename Traits >
virtual bool g2o::BlockSolver< Traits >::supportsSchur ( )
inlinevirtual

does this solver support the Schur complement for solving a system consisting of poses and landmarks. Re-implement in a derived solver, if your solver supports it.

Reimplemented from g2o::Solver.

Definition at line 138 of file block_solver.h.

138{ return true; }

◆ updateStructure()

template<typename Traits >
bool g2o::BlockSolver< Traits >::updateStructure ( const std::vector< HyperGraph::Vertex * > &  vset,
const HyperGraph::EdgeSet edges 
)
virtual

update the structures for online processing

Implements g2o::Solver.

Definition at line 271 of file block_solver.hpp.

273 {
274 for (std::vector<HyperGraph::Vertex*>::const_iterator vit = vset.begin();
275 vit != vset.end(); ++vit) {
276 OptimizableGraph::Vertex* v = static_cast<OptimizableGraph::Vertex*>(*vit);
277 int dim = v->dimension();
278 if (!v->marginalized()) {
279 v->setColInHessian(_sizePoses);
280 _sizePoses += dim;
281 _Hpp->rowBlockIndices().push_back(_sizePoses);
282 _Hpp->colBlockIndices().push_back(_sizePoses);
283 _Hpp->blockCols().push_back(
285 ++_numPoses;
286 int ind = v->hessianIndex();
287 PoseMatrixType* m = _Hpp->block(ind, ind, true);
288 v->mapHessianMemory(m->data());
289 } else {
290 G2O_ERROR("updateStructure(): Schur not supported");
291 abort();
292 }
293 }
295
296 for (HyperGraph::EdgeSet::const_iterator it = edges.begin();
297 it != edges.end(); ++it) {
298 OptimizableGraph::Edge* e = static_cast<OptimizableGraph::Edge*>(*it);
299
300 for (size_t viIdx = 0; viIdx < e->vertices().size(); ++viIdx) {
301 OptimizableGraph::Vertex* v1 =
302 (OptimizableGraph::Vertex*)e->vertex(viIdx);
303 int ind1 = v1->hessianIndex();
304 int indexV1Bak = ind1;
305 if (ind1 == -1) continue;
306 for (size_t vjIdx = viIdx + 1; vjIdx < e->vertices().size(); ++vjIdx) {
307 OptimizableGraph::Vertex* v2 =
308 (OptimizableGraph::Vertex*)e->vertex(vjIdx);
309 int ind2 = v2->hessianIndex();
310 if (ind2 == -1) continue;
311 ind1 = indexV1Bak;
312 bool transposedBlock = ind1 > ind2;
313 if (transposedBlock) // make sure, we allocate the upper triangular
314 // block
315 std::swap(ind1, ind2);
316
317 if (!v1->marginalized() && !v2->marginalized()) {
318 PoseMatrixType* m = _Hpp->block(ind1, ind2, true);
319 e->mapHessianMemory(m->data(), viIdx, vjIdx, transposedBlock);
320 } else {
321 G2O_ERROR("{}: not supported", __PRETTY_FUNCTION__);
322 }
323 }
324 }
325 }
326
327 return true;
328}
#define G2O_ERROR(...)
Definition logger.h:89
#define __PRETTY_FUNCTION__
Definition macros.h:90

References __PRETTY_FUNCTION__, g2o::OptimizableGraph::Vertex::dimension(), G2O_ERROR, g2o::OptimizableGraph::Vertex::hessianIndex(), g2o::OptimizableGraph::Vertex::mapHessianMemory(), g2o::OptimizableGraph::Edge::mapHessianMemory(), g2o::OptimizableGraph::Vertex::marginalized(), g2o::OptimizableGraph::Vertex::setColInHessian(), g2o::HyperGraph::Edge::vertex(), and g2o::HyperGraph::Edge::vertices().

◆ writeDebug()

template<typename Traits >
virtual bool g2o::BlockSolver< Traits >::writeDebug ( ) const
inlinevirtual

Implements g2o::Solver.

Definition at line 145 of file block_solver.h.

145{ return _linearSolver->writeDebug(); }

References g2o::BlockSolver< Traits >::_linearSolver.

Member Data Documentation

◆ _bschur

template<typename Traits >
std::unique_ptr<double[], aligned_deleter<double> > g2o::BlockSolver< Traits >::_bschur
protected

Definition at line 181 of file block_solver.h.

◆ _coefficients

template<typename Traits >
std::unique_ptr<double[], aligned_deleter<double> > g2o::BlockSolver< Traits >::_coefficients
protected

Definition at line 180 of file block_solver.h.

◆ _diagonalBackupLandmark

template<typename Traits >
std::vector<LandmarkVectorType> g2o::BlockSolver< Traits >::_diagonalBackupLandmark
protected

Definition at line 172 of file block_solver.h.

◆ _diagonalBackupPose

template<typename Traits >
std::vector<PoseVectorType> g2o::BlockSolver< Traits >::_diagonalBackupPose
protected

Definition at line 171 of file block_solver.h.

◆ _DInvSchur

template<typename Traits >
std::unique_ptr<SparseBlockMatrixDiagonal<LandmarkMatrixType> > g2o::BlockSolver< Traits >::_DInvSchur
protected

Definition at line 164 of file block_solver.h.

◆ _doSchur

template<typename Traits >
bool g2o::BlockSolver< Traits >::_doSchur
protected

◆ _Hll

template<typename Traits >
std::unique_ptr<SparseBlockMatrix<LandmarkMatrixType> > g2o::BlockSolver< Traits >::_Hll
protected

Definition at line 160 of file block_solver.h.

◆ _Hpl

template<typename Traits >
std::unique_ptr<SparseBlockMatrix<PoseLandmarkMatrixType> > g2o::BlockSolver< Traits >::_Hpl
protected

Definition at line 161 of file block_solver.h.

◆ _HplCCS

template<typename Traits >
std::unique_ptr<SparseBlockMatrixCCS<PoseLandmarkMatrixType> > g2o::BlockSolver< Traits >::_HplCCS
protected

Definition at line 166 of file block_solver.h.

◆ _Hpp

template<typename Traits >
std::unique_ptr<SparseBlockMatrix<PoseMatrixType> > g2o::BlockSolver< Traits >::_Hpp
protected

Definition at line 159 of file block_solver.h.

Referenced by g2o::BlockSolver< Traits >::multiplyHessian().

◆ _Hschur

template<typename Traits >
std::unique_ptr<SparseBlockMatrix<PoseMatrixType> > g2o::BlockSolver< Traits >::_Hschur
protected

Definition at line 163 of file block_solver.h.

◆ _HschurTransposedCCS

template<typename Traits >
std::unique_ptr<SparseBlockMatrixCCS<PoseMatrixType> > g2o::BlockSolver< Traits >::_HschurTransposedCCS
protected

Definition at line 167 of file block_solver.h.

◆ _linearSolver

template<typename Traits >
std::unique_ptr<LinearSolverType> g2o::BlockSolver< Traits >::_linearSolver
protected

◆ _numLandmarks

template<typename Traits >
int g2o::BlockSolver< Traits >::_numLandmarks
protected

Definition at line 183 of file block_solver.h.

Referenced by g2o::BlockSolver< Traits >::BlockSolver().

◆ _numPoses

template<typename Traits >
int g2o::BlockSolver< Traits >::_numPoses
protected

Definition at line 183 of file block_solver.h.

Referenced by g2o::BlockSolver< Traits >::BlockSolver().

◆ _sizeLandmarks

template<typename Traits >
int g2o::BlockSolver< Traits >::_sizeLandmarks
protected

Definition at line 184 of file block_solver.h.

Referenced by g2o::BlockSolver< Traits >::BlockSolver().

◆ _sizePoses

template<typename Traits >
int g2o::BlockSolver< Traits >::_sizePoses
protected

Definition at line 184 of file block_solver.h.

Referenced by g2o::BlockSolver< Traits >::BlockSolver().

◆ LandmarkDim

template<typename Traits >
const int g2o::BlockSolver< Traits >::LandmarkDim = Traits::LandmarkDim
static

Definition at line 106 of file block_solver.h.

◆ PoseDim

template<typename Traits >
const int g2o::BlockSolver< Traits >::PoseDim = Traits::PoseDim
static

Definition at line 105 of file block_solver.h.


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