g2o
block_solver.hpp
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 #include "sparse_optimizer.h"
28 #include <Eigen/LU>
29 #include <fstream>
30 #include <iomanip>
31 
32 #include "g2o/stuff/timeutil.h"
33 #include "g2o/stuff/macros.h"
34 #include "g2o/stuff/misc.h"
35 
36 namespace g2o {
37 
38 template <typename Traits>
41  _linearSolver(linearSolver)
42 {
43  // workspace
44  _Hpp=0;
45  _Hll=0;
46  _Hpl=0;
47  _HplCCS = 0;
49  _Hschur=0;
50  _DInvSchur=0;
51  _coefficients=0;
52  _bschur = 0;
53  _xSize=0;
54  _numPoses=0;
55  _numLandmarks=0;
56  _sizePoses=0;
58  _doSchur=true;
59 }
60 
61 template <typename Traits>
62 void BlockSolver<Traits>::resize(int* blockPoseIndices, int numPoseBlocks,
63  int* blockLandmarkIndices, int numLandmarkBlocks,
64  int s)
65 {
66  deallocate();
67 
68  resizeVector(s);
69 
70  if (_doSchur) {
71  // the following two are only used in schur
72  assert(_sizePoses > 0 && "allocating with wrong size");
73  _coefficients = new double [s];
74  _bschur = new double[_sizePoses];
75  }
76 
77  _Hpp=new PoseHessianType(blockPoseIndices, blockPoseIndices, numPoseBlocks, numPoseBlocks);
78  if (_doSchur) {
79  _Hschur=new PoseHessianType(blockPoseIndices, blockPoseIndices, numPoseBlocks, numPoseBlocks);
80  _Hll=new LandmarkHessianType(blockLandmarkIndices, blockLandmarkIndices, numLandmarkBlocks, numLandmarkBlocks);
82  _Hpl=new PoseLandmarkHessianType(blockPoseIndices, blockLandmarkIndices, numPoseBlocks, numLandmarkBlocks);
83  _HplCCS = new SparseBlockMatrixCCS<PoseLandmarkMatrixType>(_Hpl->rowBlockIndices(), _Hpl->colBlockIndices());
85 #ifdef G2O_OPENMP
86  _coefficientsMutex.resize(numPoseBlocks);
87 #endif
88  }
89 }
90 
91 template <typename Traits>
93 {
94  delete _Hpp;
95  _Hpp=0;
96  delete _Hll;
97  _Hll=0;
98  delete _Hpl;
99  _Hpl = 0;
100  delete _Hschur;
101  _Hschur=0;
102  delete _DInvSchur;
103  _DInvSchur=0;
104  delete[] _coefficients;
105  _coefficients = 0;
106  delete[] _bschur;
107  _bschur = 0;
108  delete _HplCCS;
109  _HplCCS = 0;
110  delete _HschurTransposedCCS;
112 }
113 
114 template <typename Traits>
116 {
117  delete _linearSolver;
118  deallocate();
119 }
120 
121 template <typename Traits>
123 {
124  assert(_optimizer);
125 
126  size_t sparseDim = 0;
127  _numPoses=0;
128  _numLandmarks=0;
129  _sizePoses=0;
130  _sizeLandmarks=0;
131  int* blockPoseIndices = new int[_optimizer->indexMapping().size()];
132  int* blockLandmarkIndices = new int[_optimizer->indexMapping().size()];
133 
134  for (size_t i = 0; i < _optimizer->indexMapping().size(); ++i) {
136  int dim = v->dimension();
137  if (! v->marginalized()){
139  _sizePoses+=dim;
140  blockPoseIndices[_numPoses]=_sizePoses;
141  ++_numPoses;
142  } else {
144  _sizeLandmarks+=dim;
145  blockLandmarkIndices[_numLandmarks]=_sizeLandmarks;
146  ++_numLandmarks;
147  }
148  sparseDim += dim;
149  }
150  resize(blockPoseIndices, _numPoses, blockLandmarkIndices, _numLandmarks, sparseDim);
151  delete[] blockLandmarkIndices;
152  delete[] blockPoseIndices;
153 
154  // allocate the diagonal on Hpp and Hll
155  int poseIdx = 0;
156  int landmarkIdx = 0;
157  for (size_t i = 0; i < _optimizer->indexMapping().size(); ++i) {
159  if (! v->marginalized()){
160  //assert(poseIdx == v->hessianIndex());
161  PoseMatrixType* m = _Hpp->block(poseIdx, poseIdx, true);
162  if (zeroBlocks)
163  m->setZero();
164  v->mapHessianMemory(m->data());
165  ++poseIdx;
166  } else {
167  LandmarkMatrixType* m = _Hll->block(landmarkIdx, landmarkIdx, true);
168  if (zeroBlocks)
169  m->setZero();
170  v->mapHessianMemory(m->data());
171  ++landmarkIdx;
172  }
173  }
174  assert(poseIdx == _numPoses && landmarkIdx == _numLandmarks);
175 
176  // temporary structures for building the pattern of the Schur complement
177  SparseBlockMatrixHashMap<PoseMatrixType>* schurMatrixLookup = 0;
178  if (_doSchur) {
180  schurMatrixLookup->blockCols().resize(_Hschur->blockCols().size());
181  }
182 
183  // here we assume that the landmark indices start after the pose ones
184  // create the structure in Hpp, Hll and in Hpl
185  for (SparseOptimizer::EdgeContainer::const_iterator it=_optimizer->activeEdges().begin(); it!=_optimizer->activeEdges().end(); ++it){
186  OptimizableGraph::Edge* e = *it;
187 
188  for (size_t viIdx = 0; viIdx < e->vertices().size(); ++viIdx) {
190  int ind1 = v1->hessianIndex();
191  if (ind1 == -1)
192  continue;
193  int indexV1Bak = ind1;
194  for (size_t vjIdx = viIdx + 1; vjIdx < e->vertices().size(); ++vjIdx) {
196  int ind2 = v2->hessianIndex();
197  if (ind2 == -1)
198  continue;
199  ind1 = indexV1Bak;
200  bool transposedBlock = ind1 > ind2;
201  if (transposedBlock){ // make sure, we allocate the upper triangle block
202  std::swap(ind1, ind2);
203  }
204  if (! v1->marginalized() && !v2->marginalized()){
205  PoseMatrixType* m = _Hpp->block(ind1, ind2, true);
206  if (zeroBlocks)
207  m->setZero();
208  e->mapHessianMemory(m->data(), viIdx, vjIdx, transposedBlock);
209  if (_Hschur) {// assume this is only needed in case we solve with the schur complement
210  schurMatrixLookup->addBlock(ind1, ind2);
211  }
212  } else if (v1->marginalized() && v2->marginalized()){
213  // RAINER hmm.... should we ever reach this here????
214  LandmarkMatrixType* m = _Hll->block(ind1-_numPoses, ind2-_numPoses, true);
215  if (zeroBlocks)
216  m->setZero();
217  e->mapHessianMemory(m->data(), viIdx, vjIdx, false);
218  } else {
219  if (v1->marginalized()){
221  if (zeroBlocks)
222  m->setZero();
223  e->mapHessianMemory(m->data(), viIdx, vjIdx, true); // transpose the block before writing to it
224  } else {
226  if (zeroBlocks)
227  m->setZero();
228  e->mapHessianMemory(m->data(), viIdx, vjIdx, false); // directly the block
229  }
230  }
231  }
232  }
233  }
234 
235  if (! _doSchur)
236  return true;
237 
238  _DInvSchur->diagonal().resize(landmarkIdx);
240 
241  for (size_t i = 0; i < _optimizer->indexMapping().size(); ++i) {
243  if (v->marginalized()){
244  const HyperGraph::EdgeSet& vedges=v->edges();
245  for (HyperGraph::EdgeSet::const_iterator it1=vedges.begin(); it1!=vedges.end(); ++it1){
246  for (size_t i=0; i<(*it1)->vertices().size(); ++i)
247  {
248  OptimizableGraph::Vertex* v1= (OptimizableGraph::Vertex*) (*it1)->vertex(i);
249  if (v1->hessianIndex()==-1 || v1==v)
250  continue;
251  for (HyperGraph::EdgeSet::const_iterator it2=vedges.begin(); it2!=vedges.end(); ++it2){
252  for (size_t j=0; j<(*it2)->vertices().size(); ++j)
253  {
254  OptimizableGraph::Vertex* v2= (OptimizableGraph::Vertex*) (*it2)->vertex(j);
255  if (v2->hessianIndex()==-1 || v2==v)
256  continue;
257  int i1=v1->hessianIndex();
258  int i2=v2->hessianIndex();
259  if (i1<=i2) {
260  schurMatrixLookup->addBlock(i1, i2);
261  }
262  }
263  }
264  }
265  }
266  }
267  }
268 
269  _Hschur->takePatternFromHash(*schurMatrixLookup);
270  delete schurMatrixLookup;
272 
273  return true;
274 }
275 
276 template <typename Traits>
277 bool BlockSolver<Traits>::updateStructure(const std::vector<HyperGraph::Vertex*>& vset, const HyperGraph::EdgeSet& edges)
278 {
279  for (std::vector<HyperGraph::Vertex*>::const_iterator vit = vset.begin(); vit != vset.end(); ++vit) {
280  OptimizableGraph::Vertex* v = static_cast<OptimizableGraph::Vertex*>(*vit);
281  int dim = v->dimension();
282  if (! v->marginalized()){
284  _sizePoses+=dim;
285  _Hpp->rowBlockIndices().push_back(_sizePoses);
286  _Hpp->colBlockIndices().push_back(_sizePoses);
288  ++_numPoses;
289  int ind = v->hessianIndex();
290  PoseMatrixType* m = _Hpp->block(ind, ind, true);
291  v->mapHessianMemory(m->data());
292  } else {
293  std::cerr << "updateStructure(): Schur not supported" << std::endl;
294  abort();
295  }
296  }
298 
299  for (HyperGraph::EdgeSet::const_iterator it = edges.begin(); it != edges.end(); ++it) {
300  OptimizableGraph::Edge* e = static_cast<OptimizableGraph::Edge*>(*it);
301 
302  for (size_t viIdx = 0; viIdx < e->vertices().size(); ++viIdx) {
304  int ind1 = v1->hessianIndex();
305  int indexV1Bak = ind1;
306  if (ind1 == -1)
307  continue;
308  for (size_t vjIdx = viIdx + 1; vjIdx < e->vertices().size(); ++vjIdx) {
310  int ind2 = v2->hessianIndex();
311  if (ind2 == -1)
312  continue;
313  ind1 = indexV1Bak;
314  bool transposedBlock = ind1 > ind2;
315  if (transposedBlock) // make sure, we allocate the upper triangular block
316  std::swap(ind1, ind2);
317 
318  if (! v1->marginalized() && !v2->marginalized()) {
319  PoseMatrixType* m = _Hpp->block(ind1, ind2, true);
320  e->mapHessianMemory(m->data(), viIdx, vjIdx, transposedBlock);
321  } else {
322  std::cerr << __PRETTY_FUNCTION__ << ": not supported" << std::endl;
323  }
324  }
325  }
326 
327  }
328 
329  return true;
330 }
331 
332 template <typename Traits>
334  //cerr << __PRETTY_FUNCTION__ << endl;
335  if (! _doSchur){
336  double t=get_monotonic_time();
337  bool ok = _linearSolver->solve(*_Hpp, _x, _b);
339  if (globalStats) {
340  globalStats->timeLinearSolver = get_monotonic_time() - t;
341  globalStats->hessianDimension = globalStats->hessianPoseDimension = _Hpp->cols();
342  }
343  return ok;
344  }
345 
346  // schur thing
347 
348  // backup the coefficient matrix
349  double t=get_monotonic_time();
350 
351  // _Hschur = _Hpp, but keeping the pattern of _Hschur
352  _Hschur->clear();
353  _Hpp->add(_Hschur);
354 
355  //_DInvSchur->clear();
356  memset (_coefficients, 0, _sizePoses*sizeof(double));
357 # ifdef G2O_OPENMP
358 # pragma omp parallel for default (shared) schedule(dynamic, 10)
359 # endif
360  for (int landmarkIndex = 0; landmarkIndex < static_cast<int>(_Hll->blockCols().size()); ++landmarkIndex) {
361  const typename SparseBlockMatrix<LandmarkMatrixType>::IntBlockMap& marginalizeColumn = _Hll->blockCols()[landmarkIndex];
362  assert(marginalizeColumn.size() == 1 && "more than one block in _Hll column");
363 
364  // calculate inverse block for the landmark
365  const LandmarkMatrixType * D = marginalizeColumn.begin()->second;
366  assert (D && D->rows()==D->cols() && "Error in landmark matrix");
367  LandmarkMatrixType& Dinv = _DInvSchur->diagonal()[landmarkIndex];
368  Dinv = D->inverse();
369 
370  LandmarkVectorType db(D->rows());
371  for (int j=0; j<D->rows(); ++j) {
372  db[j]=_b[_Hll->rowBaseOfBlock(landmarkIndex) + _sizePoses + j];
373  }
374  db=Dinv*db;
375 
376  assert((size_t)landmarkIndex < _HplCCS->blockCols().size() && "Index out of bounds");
377  const typename SparseBlockMatrixCCS<PoseLandmarkMatrixType>::SparseColumn& landmarkColumn = _HplCCS->blockCols()[landmarkIndex];
378 
379  for (typename SparseBlockMatrixCCS<PoseLandmarkMatrixType>::SparseColumn::const_iterator it_outer = landmarkColumn.begin();
380  it_outer != landmarkColumn.end(); ++it_outer) {
381  int i1 = it_outer->row;
382 
383  const PoseLandmarkMatrixType* Bi = it_outer->block;
384  assert(Bi);
385 
386  PoseLandmarkMatrixType BDinv = (*Bi)*(Dinv);
387  assert(_HplCCS->rowBaseOfBlock(i1) < _sizePoses && "Index out of bounds");
388  typename PoseVectorType::MapType Bb(&_coefficients[_HplCCS->rowBaseOfBlock(i1)], Bi->rows());
389 # ifdef G2O_OPENMP
390  ScopedOpenMPMutex mutexLock(&_coefficientsMutex[i1]);
391 # endif
392  Bb.noalias() += (*Bi)*db;
393 
394  assert(i1 >= 0 && i1 < static_cast<int>(_HschurTransposedCCS->blockCols().size()) && "Index out of bounds");
396 
398  typename SparseBlockMatrixCCS<PoseLandmarkMatrixType>::SparseColumn::const_iterator it_inner = lower_bound(landmarkColumn.begin(), landmarkColumn.end(), aux);
399  for (; it_inner != landmarkColumn.end(); ++it_inner) {
400  int i2 = it_inner->row;
401  const PoseLandmarkMatrixType* Bj = it_inner->block;
402  assert(Bj);
403  while (targetColumnIt->row < i2 /*&& targetColumnIt != _HschurTransposedCCS->blockCols()[i1].end()*/)
404  ++targetColumnIt;
405  assert(targetColumnIt != _HschurTransposedCCS->blockCols()[i1].end() && targetColumnIt->row == i2 && "invalid iterator, something wrong with the matrix structure");
406  PoseMatrixType* Hi1i2 = targetColumnIt->block;//_Hschur->block(i1,i2);
407  assert(Hi1i2);
408  (*Hi1i2).noalias() -= BDinv*Bj->transpose();
409  }
410  }
411  }
412  //cerr << "Solve [marginalize] = " << get_monotonic_time()-t << endl;
413 
414  // _bschur = _b for calling solver, and not touching _b
415  memcpy(_bschur, _b, _sizePoses * sizeof(double));
416  for (int i=0; i<_sizePoses; ++i){
417  _bschur[i]-=_coefficients[i];
418  }
419 
421  if (globalStats){
422  globalStats->timeSchurComplement = get_monotonic_time() - t;
423  }
424 
425  t=get_monotonic_time();
426  bool solvedPoses = _linearSolver->solve(*_Hschur, _x, _bschur);
427  if (globalStats) {
428  globalStats->timeLinearSolver = get_monotonic_time() - t;
429  globalStats->hessianPoseDimension = _Hpp->cols();
430  globalStats->hessianLandmarkDimension = _Hll->cols();
431  globalStats->hessianDimension = globalStats->hessianPoseDimension + globalStats->hessianLandmarkDimension;
432  }
433  //cerr << "Solve [decompose and solve] = " << get_monotonic_time()-t << endl;
434 
435  if (! solvedPoses)
436  return false;
437 
438  // _x contains the solution for the poses, now applying it to the landmarks to get the new part of the
439  // solution;
440  double* xp = _x;
441  double* cp = _coefficients;
442 
443  double* xl=_x+_sizePoses;
444  double* cl=_coefficients + _sizePoses;
445  double* bl=_b+_sizePoses;
446 
447  // cp = -xp
448  for (int i=0; i<_sizePoses; ++i)
449  cp[i]=-xp[i];
450 
451  // cl = bl
452  memcpy(cl,bl,_sizeLandmarks*sizeof(double));
453 
454  // cl = bl - Bt * xp
455  //Bt->multiply(cl, cp);
456  _HplCCS->rightMultiply(cl, cp);
457 
458  // xl = Dinv * cl
459  memset(xl,0, _sizeLandmarks*sizeof(double));
460  _DInvSchur->multiply(xl,cl);
461  //_DInvSchur->rightMultiply(xl,cl);
462  //cerr << "Solve [landmark delta] = " << get_monotonic_time()-t << endl;
463 
464  return true;
465 }
466 
467 
468 template <typename Traits>
469 bool BlockSolver<Traits>::computeMarginals(SparseBlockMatrix<MatrixXD>& spinv, const std::vector<std::pair<int, int> >& blockIndices)
470 {
471  double t = get_monotonic_time();
472  bool ok = _linearSolver->solvePattern(spinv, blockIndices, *_Hpp);
474  if (globalStats) {
475  globalStats->timeMarginals = get_monotonic_time() - t;
476  }
477  return ok;
478 }
479 
480 template <typename Traits>
482 {
483  // clear b vector
484 # ifdef G2O_OPENMP
485 # pragma omp parallel for default (shared) if (_optimizer->indexMapping().size() > 1000)
486 # endif
487  for (int i = 0; i < static_cast<int>(_optimizer->indexMapping().size()); ++i) {
489  assert(v);
490  v->clearQuadraticForm();
491  }
492  _Hpp->clear();
493  if (_doSchur) {
494  _Hll->clear();
495  _Hpl->clear();
496  }
497 
498  // resetting the terms for the pairwise constraints
499  // built up the current system by storing the Hessian blocks in the edges and vertices
500 # ifndef G2O_OPENMP
501  // no threading, we do not need to copy the workspace
502  JacobianWorkspace& jacobianWorkspace = _optimizer->jacobianWorkspace();
503 # else
504  // if running with threads need to produce copies of the workspace for each thread
505  JacobianWorkspace jacobianWorkspace = _optimizer->jacobianWorkspace();
506 # pragma omp parallel for default (shared) firstprivate(jacobianWorkspace) if (_optimizer->activeEdges().size() > 100)
507 # endif
508  for (int k = 0; k < static_cast<int>(_optimizer->activeEdges().size()); ++k) {
510  e->linearizeOplus(jacobianWorkspace); // jacobian of the nodes' oplus (manifold)
512 # ifndef NDEBUG
513  for (size_t i = 0; i < e->vertices().size(); ++i) {
514  const OptimizableGraph::Vertex* v = static_cast<const OptimizableGraph::Vertex*>(e->vertex(i));
515  if (! v->fixed()) {
516  bool hasANan = arrayHasNaN(jacobianWorkspace.workspaceForVertex(i), e->dimension() * v->dimension());
517  if (hasANan) {
518  std::cerr << "buildSystem(): NaN within Jacobian for edge " << e << " for vertex " << i << std::endl;
519  break;
520  }
521  }
522  }
523 # endif
524  }
525 
526  // flush the current system in a sparse block matrix
527 # ifdef G2O_OPENMP
528 # pragma omp parallel for default (shared) if (_optimizer->indexMapping().size() > 1000)
529 # endif
530  for (int i = 0; i < static_cast<int>(_optimizer->indexMapping().size()); ++i) {
532  int iBase = v->colInHessian();
533  if (v->marginalized())
534  iBase+=_sizePoses;
535  v->copyB(_b+iBase);
536  }
537 
538  return 0;
539 }
540 
541 
542 template <typename Traits>
543 bool BlockSolver<Traits>::setLambda(double lambda, bool backup)
544 {
545  if (backup) {
548  }
549 # ifdef G2O_OPENMP
550 # pragma omp parallel for default (shared) if (_numPoses > 100)
551 # endif
552  for (int i = 0; i < _numPoses; ++i) {
553  PoseMatrixType *b=_Hpp->block(i,i);
554  if (backup)
555  _diagonalBackupPose[i] = b->diagonal();
556  b->diagonal().array() += lambda;
557  }
558 # ifdef G2O_OPENMP
559 # pragma omp parallel for default (shared) if (_numLandmarks > 100)
560 # endif
561  for (int i = 0; i < _numLandmarks; ++i) {
563  if (backup)
564  _diagonalBackupLandmark[i] = b->diagonal();
565  b->diagonal().array() += lambda;
566  }
567  return true;
568 }
569 
570 template <typename Traits>
572 {
573  assert((int) _diagonalBackupPose.size() == _numPoses && "Mismatch in dimensions");
574  assert((int) _diagonalBackupLandmark.size() == _numLandmarks && "Mismatch in dimensions");
575  for (int i = 0; i < _numPoses; ++i) {
576  PoseMatrixType *b=_Hpp->block(i,i);
577  b->diagonal() = _diagonalBackupPose[i];
578  }
579  for (int i = 0; i < _numLandmarks; ++i) {
581  b->diagonal() = _diagonalBackupLandmark[i];
582  }
583 }
584 
585 template <typename Traits>
587 {
589  if (! online) {
590  if (_Hpp)
591  _Hpp->clear();
592  if (_Hpl)
593  _Hpl->clear();
594  if (_Hll)
595  _Hll->clear();
596  }
597  _linearSolver->init();
598  return true;
599 }
600 
601 template <typename Traits>
603 {
604  _linearSolver->setWriteDebug(writeDebug);
605 }
606 
607 template <typename Traits>
608 bool BlockSolver<Traits>::saveHessian(const std::string& fileName) const
609 {
610  return _Hpp->writeOctave(fileName.c_str(), true);
611 }
612 
613 } // end namespace
double * _b
Definition: solver.h:138
const std::vector< int > & rowBlockIndices() const
indices of the row blocks
SparseBlockMatrix< PoseMatrixType > * _Hschur
Definition: block_solver.h:153
Traits::PoseHessianType PoseHessianType
Definition: block_solver.h:107
double get_monotonic_time()
Definition: timeutil.cpp:113
int fillSparseBlockMatrixCCS(SparseBlockMatrixCCS< MatrixType > &blockCCS) const
std::vector< LandmarkVectorType, Eigen::aligned_allocator< LandmarkVectorType > > _diagonalBackupLandmark
Definition: block_solver.h:162
statistics about the optimization
Definition: batch_stats.h:40
Traits::LandmarkMatrixType LandmarkMatrixType
Definition: block_solver.h:102
void setColInHessian(int c)
set the row of this vertex in the Hessian
const Vertex * vertex(size_t i) const
Definition: hyper_graph.h:186
virtual bool solve(const SparseBlockMatrix< MatrixType > &A, double *x, double *b)=0
#define __PRETTY_FUNCTION__
Definition: macros.h:89
some general case utility functions
virtual bool solvePattern(SparseBlockMatrix< MatrixXD > &spinv, const std::vector< std::pair< int, int > > &blockIndices, const SparseBlockMatrix< MatrixType > &A)
Definition: linear_solver.h:71
virtual bool saveHessian(const std::string &fileName) const
write the hessian to disk using the specified file name
int dimension() const
returns the dimensions of the error function
int cols() const
columns of the matrix
virtual void setWriteDebug(bool writeDebug)
SparseBlockMatrix< PoseMatrixType > * _Hpp
Definition: block_solver.h:149
Traits::LandmarkHessianType LandmarkHessianType
Definition: block_solver.h:108
virtual int copyB(double *b_) const =0
void rightMultiply(double *&dest, const double *src) const
virtual bool init()=0
double * workspaceForVertex(int vertexIndex)
bool writeOctave(const char *filename, bool upperTriangle=true) const
SparseOptimizer * _optimizer
Definition: solver.h:136
SparseBlockMatrixDiagonal< LandmarkMatrixType > * _DInvSchur
Definition: block_solver.h:154
SparseMatrixBlock * block(int r, int c, bool alloc=false)
returns the block at location r,c. if alloc=true he block is created if it does not exist ...
virtual bool init(SparseOptimizer *optmizer, bool online=false)
void resize(int *blockPoseIndices, int numPoseBlocks, int *blockLandmarkIndices, int numLandmarkBlocks, int totalDim)
lock a mutex within a scope
Definition: openmp_mutex.h:84
size_t _xSize
Definition: solver.h:139
const EdgeContainer & activeEdges() const
the edges active in the current optimization
int hessianIndex() const
temporary index of this node in the parameter vector obtained from linearization
SparseBlockMatrixCCS< PoseLandmarkMatrixType > * _HplCCS
Definition: block_solver.h:156
SparseBlockMatrix< LandmarkMatrixType > * _Hll
Definition: block_solver.h:150
void clear(bool dealloc=false)
this zeroes all the blocks. If dealloc=true the blocks are removed from memory
virtual void linearizeOplus(JacobianWorkspace &jacobianWorkspace)=0
Traits::PoseMatrixType PoseMatrixType
Definition: block_solver.h:101
double * _coefficients
Definition: block_solver.h:170
double * _x
Definition: solver.h:137
bool arrayHasNaN(const double *array, int size, int *nanIndex=0)
Definition: misc.h:177
size_t hessianPoseDimension
dimension of the pose matrix in Schur
Definition: batch_stats.h:69
utility functions for handling time related stuff
JacobianWorkspace & jacobianWorkspace()
the workspace for storing the Jacobians of the graph
const DiagonalVector & diagonal() const
the block matrices per block-column
Traits::LinearSolverType LinearSolverType
Definition: block_solver.h:110
const std::vector< int > & colBlockIndices() const
indices of the column blocks
std::set< Edge * > EdgeSet
Definition: hyper_graph.h:135
int dimension() const
dimension of the estimated state belonging to this node
virtual void setWriteDebug(bool)
Definition: linear_solver.h:80
virtual bool buildStructure(bool zeroBlocks=false)
void takePatternFromHash(SparseBlockMatrixHashMap< MatrixType > &hashMatrix)
MatrixType * addBlock(int r, int c, bool zeroBlock=false)
const VertexContainer & vertices() const
Definition: hyper_graph.h:178
virtual bool setLambda(double lambda, bool backup=false)
const EdgeSet & edges() const
returns the set of hyper-edges that are leaving/entering in this vertex
Definition: hyper_graph.h:151
SparseBlockMatrixCCS< PoseMatrixType > * _HschurTransposedCCS
Definition: block_solver.h:157
Traits::PoseLandmarkHessianType PoseLandmarkHessianType
Definition: block_solver.h:109
virtual void constructQuadraticForm()=0
std::vector< PoseVectorType, Eigen::aligned_allocator< PoseVectorType > > _diagonalBackupPose
Definition: block_solver.h:161
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
virtual bool updateStructure(const std::vector< HyperGraph::Vertex * > &vset, const HyperGraph::EdgeSet &edges)
int colInHessian() const
get the row of this vertex in the Hessian
size_t hessianLandmarkDimension
dimension of the landmark matrix in Schur
Definition: batch_stats.h:70
Traits::PoseLandmarkMatrixType PoseLandmarkMatrixType
Definition: block_solver.h:103
int fillSparseBlockMatrixCCSTransposed(SparseBlockMatrixCCS< MatrixType > &blockCCS) const
Traits::LandmarkVectorType LandmarkVectorType
Definition: block_solver.h:105
static G2OBatchStatistics * globalStats()
Definition: batch_stats.h:73
virtual bool solve()
const VertexContainer & indexMapping() const
the index mapping of the vertices
virtual bool computeMarginals(SparseBlockMatrix< MatrixXD > &spinv, const std::vector< std::pair< int, int > > &blockIndices)
void multiply(double *&dest, const double *src) const
void resizeVector(size_t sx)
Definition: solver.cpp:46
SparseOptimizer * optimizer() const
the optimizer (graph) on which the solver works
Definition: solver.h:106
virtual bool buildSystem()
virtual void mapHessianMemory(double *d)=0
virtual void clearQuadraticForm()=0
bool marginalized() const
true => this node is marginalized out during the optimization
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
A general case Vertex for optimization.
SparseBlockMatrix< PoseLandmarkMatrixType > * _Hpl
Definition: block_solver.h:151
size_t hessianDimension
rows / cols of the Hessian
Definition: batch_stats.h:68
double * b()
return b, the right hand side of the system
Definition: solver.h:99
bool fixed() const
true => this node is fixed during the optimization
LinearSolver< PoseMatrixType > * _linearSolver
Definition: block_solver.h:159
double timeSchurComplement
compute schur complement (0 if not done)
Definition: batch_stats.h:54
provide memory workspace for computing the Jacobians
Sparse matrix which uses blocks based on hash structures.
virtual void mapHessianMemory(double *d, int i, int j, bool rowMajor)=0
bool add(SparseBlockMatrix< MatrixType > *&dest) const
adds the current matrix to the destination
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
virtual bool writeDebug() const
Definition: block_solver.h:137
BlockSolver(LinearSolverType *linearSolver)
int rowBaseOfBlock(int r) const
where does the row at block-row r start?
int rowBaseOfBlock(int r) const
where does the row at block-row r starts?
base for the block solvers with some basic function interfaces
Definition: block_solver.h:82
Sparse matrix which uses blocks.
double timeLinearSolver
time for solving, excluding Schur setup
Definition: batch_stats.h:60
double timeMarginals
computing the inverse elements (solve blocks) and thus the marginal covariances
Definition: batch_stats.h:65
virtual void restoreDiagonal()