g2o
Public Types | Public Member Functions | Protected Attributes | List of all members
g2o::SparseBlockMatrix< MatrixType > Class Template Reference

Sparse matrix which uses blocks. More...

#include <sparse_block_matrix.h>

Public Types

typedef MatrixType SparseMatrixBlock
 this is the type of the elementary block, it is an Eigen::Matrix. More...
 
typedef std::map< int, SparseMatrixBlock * > IntBlockMap
 

Public Member Functions

int cols () const
 columns of the matrix More...
 
int rows () const
 rows of the matrix More...
 
 SparseBlockMatrix (const int *rbi, const int *cbi, int rb, int cb, bool hasStorage=true)
 
 SparseBlockMatrix ()
 
 ~SparseBlockMatrix ()
 
void clear (bool dealloc=false)
 this zeroes all the blocks. If dealloc=true the blocks are removed from memory More...
 
SparseMatrixBlockblock (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 More...
 
const SparseMatrixBlockblock (int r, int c) const
 returns the block at location r,c More...
 
int rowsOfBlock (int r) const
 how many rows does the block at block-row r has? More...
 
int colsOfBlock (int c) const
 how many cols does the block at block-col c has? More...
 
int rowBaseOfBlock (int r) const
 where does the row at block-row r starts? More...
 
int colBaseOfBlock (int c) const
 where does the col at block-col r starts? More...
 
size_t nonZeros () const
 number of non-zero elements More...
 
size_t nonZeroBlocks () const
 number of allocated blocks More...
 
SparseBlockMatrixclone () const
 deep copy of a sparse-block-matrix; More...
 
SparseBlockMatrixslice (int rmin, int rmax, int cmin, int cmax, bool alloc=true) const
 
template<class MatrixTransposedType >
bool transpose (SparseBlockMatrix< MatrixTransposedType > *&dest) const
 transposes a block matrix, The transposed type should match the argument false on failure More...
 
bool add (SparseBlockMatrix< MatrixType > *&dest) const
 adds the current matrix to the destination More...
 
template<class MatrixResultType , class MatrixFactorType >
bool multiply (SparseBlockMatrix< MatrixResultType > *&dest, const SparseBlockMatrix< MatrixFactorType > *M) const
 dest = (*this) * M More...
 
void multiply (double *&dest, const double *src) const
 dest = (*this) * src More...
 
void multiplySymmetricUpperTriangle (double *&dest, const double *src) const
 
void rightMultiply (double *&dest, const double *src) const
 dest = M * (*this) More...
 
void scale (double a)
 *this *= a More...
 
bool symmPermutation (SparseBlockMatrix< MatrixType > *&dest, const int *pinv, bool onlyUpper=false) const
 
int fillCCS (int *Cp, int *Ci, double *Cx, bool upperTriangle=false) const
 
int fillCCS (double *Cx, bool upperTriangle=false) const
 
void fillBlockStructure (MatrixStructure &ms) const
 exports the non zero blocks in the structure matrix ms More...
 
const std::vector< IntBlockMap > & blockCols () const
 the block matrices per block-column More...
 
std::vector< IntBlockMap > & blockCols ()
 
const std::vector< int > & rowBlockIndices () const
 indices of the row blocks More...
 
std::vector< int > & rowBlockIndices ()
 
const std::vector< int > & colBlockIndices () const
 indices of the column blocks More...
 
std::vector< int > & colBlockIndices ()
 
bool writeOctave (const char *filename, bool upperTriangle=true) const
 
int fillSparseBlockMatrixCCS (SparseBlockMatrixCCS< MatrixType > &blockCCS) const
 
int fillSparseBlockMatrixCCSTransposed (SparseBlockMatrixCCS< MatrixType > &blockCCS) const
 
void takePatternFromHash (SparseBlockMatrixHashMap< MatrixType > &hashMatrix)
 

Protected Attributes

std::vector< int > _rowBlockIndices
 vector of the indices of the blocks along the rows. More...
 
std::vector< int > _colBlockIndices
 
std::vector< IntBlockMap_blockCols
 
bool _hasStorage
 

Detailed Description

template<class MatrixType = MatrixXD>
class g2o::SparseBlockMatrix< MatrixType >

Sparse matrix which uses blocks.

Template class that specifies a sparse block matrix. A block matrix is a sparse matrix made of dense blocks. These blocks cannot have a random pattern, but follow a (variable) grid structure. This structure is specified by a partition of the rows and the columns of the matrix. The blocks are represented by the Eigen::Matrix structure, thus they can be statically or dynamically allocated. For efficiency reasons it is convenient to allocate them statically, when possible. A static block matrix has all blocks of the same size, and the size of the block is specified by the template argument. If this is not the case, and you have different block sizes than you have to use a dynamic-block matrix (default template argument).

Definition at line 61 of file sparse_block_matrix.h.

Member Typedef Documentation

template<class MatrixType = MatrixXD>
typedef std::map<int, SparseMatrixBlock*> g2o::SparseBlockMatrix< MatrixType >::IntBlockMap

Definition at line 72 of file sparse_block_matrix.h.

template<class MatrixType = MatrixXD>
typedef MatrixType g2o::SparseBlockMatrix< MatrixType >::SparseMatrixBlock

this is the type of the elementary block, it is an Eigen::Matrix.

Definition at line 65 of file sparse_block_matrix.h.

Constructor & Destructor Documentation

template<class MatrixType >
g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix ( const int *  rbi,
const int *  cbi,
int  rb,
int  cb,
bool  hasStorage = true 
)

constructs a sparse block matrix having a specific layout

Parameters
rbiarray of int containing the row layout of the blocks. the component i of the array should contain the index of the first row of the block i+1.
rbiarray of int containing the column layout of the blocks. the component i of the array should contain the index of the first col of the block i+1.
rbnumber of row blocks
cbnumber of col blocks
hasStorageset it to true if the matrix "owns" the blocks, thus it deletes it on destruction. if false the matrix is only a "view" over an existing structure.

Definition at line 53 of file sparse_block_matrix.hpp.

53  :
54  _rowBlockIndices(rbi,rbi+rb),
55  _colBlockIndices(cbi,cbi+cb),
56  _blockCols(cb), _hasStorage(hasStorage)
57  {
58  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix ( )
template<class MatrixType >
g2o::SparseBlockMatrix< MatrixType >::~SparseBlockMatrix ( )

Definition at line 85 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_hasStorage, and g2o::SparseBlockMatrix< MatrixType >::clear().

85  {
86  if (_hasStorage)
87  clear(true);
88  }
void clear(bool dealloc=false)
this zeroes all the blocks. If dealloc=true the blocks are removed from memory

Member Function Documentation

template<class MatrixType>
bool g2o::SparseBlockMatrix< MatrixType >::add ( SparseBlockMatrix< MatrixType > *&  dest) const

adds the current matrix to the destination

Definition at line 168 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::block(), and g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), main(), and g2o::BlockSolver< Traits >::solve().

168  {
169  if (! dest){
171  } else {
172  if (! dest->_hasStorage)
173  return false;
174  if(_rowBlockIndices.size()!=dest->_rowBlockIndices.size())
175  return false;
176  if (_colBlockIndices.size()!=dest->_colBlockIndices.size())
177  return false;
178  for (size_t i=0; i<_rowBlockIndices.size(); ++i){
179  if(_rowBlockIndices[i]!=dest->_rowBlockIndices[i])
180  return false;
181  }
182  for (size_t i=0; i<_colBlockIndices.size(); ++i){
183  if(_colBlockIndices[i]!=dest->_colBlockIndices[i])
184  return false;
185  }
186  }
187  for (size_t i=0; i<_blockCols.size(); ++i){
188  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
190  typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* d=dest->block(it->first,i,true);
191  (*d)+=*s;
192  }
193  }
194  return true;
195  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
SparseBlockMatrix< MatrixType >::SparseMatrixBlock * g2o::SparseBlockMatrix< MatrixType >::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

Definition at line 91 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rowsOfBlock().

Referenced by g2o::SparseBlockMatrix< MatrixType >::add(), g2o::BlockSolver< Traits >::buildStructure(), g2o::MarginalCovarianceCholesky::computeCovariance(), g2o::EdgeLabeler::labelEdge(), main(), g2o::SparseBlockMatrix< MatrixType >::multiply(), g2o::BlockSolver< Traits >::restoreDiagonal(), g2o::BlockSolver< Traits >::setLambda(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), testMarginals(), g2o::SparseBlockMatrix< MatrixType >::transpose(), and g2o::BlockSolver< Traits >::updateStructure().

91  {
92  typename SparseBlockMatrix<MatrixType>::IntBlockMap::iterator it =_blockCols[c].find(r);
94  if (it==_blockCols[c].end()){
95  if (!_hasStorage && ! alloc )
96  return 0;
97  else {
98  int rb=rowsOfBlock(r);
99  int cb=colsOfBlock(c);
100  _block=new typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock(rb,cb);
101  _block->setZero();
102  std::pair < typename SparseBlockMatrix<MatrixType>::IntBlockMap::iterator, bool> result
103  =_blockCols[c].insert(std::make_pair(r,_block)); (void) result;
104  assert (result.second);
105  }
106  } else {
107  _block=it->second;
108  }
109  return _block;
110  }
int rowsOfBlock(int r) const
how many rows does the block at block-row r has?
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
const SparseBlockMatrix< MatrixType >::SparseMatrixBlock * g2o::SparseBlockMatrix< MatrixType >::block ( int  r,
int  c 
) const

returns the block at location r,c

Definition at line 113 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols.

113  {
114  typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it =_blockCols[c].find(r);
115  if (it==_blockCols[c].end())
116  return 0;
117  return it->second;
118  }
std::vector< IntBlockMap > _blockCols
template<class MatrixType = MatrixXD>
const std::vector<IntBlockMap>& g2o::SparseBlockMatrix< MatrixType >::blockCols ( ) const
inline
template<class MatrixType = MatrixXD>
std::vector<IntBlockMap>& g2o::SparseBlockMatrix< MatrixType >::blockCols ( )
inline

Definition at line 178 of file sparse_block_matrix.h.

178 { return _blockCols;}
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::clear ( bool  dealloc = false)

this zeroes all the blocks. If dealloc=true the blocks are removed from memory

Definition at line 67 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, and g2o::SparseBlockMatrix< MatrixType >::_hasStorage.

Referenced by g2o::BlockSolver< Traits >::buildSystem(), g2o::BlockSolver< Traits >::init(), main(), g2o::BlockSolver< Traits >::solve(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), and g2o::SparseBlockMatrix< MatrixType >::~SparseBlockMatrix().

67  {
68 # ifdef G2O_OPENMP
69 # pragma omp parallel for default (shared) if (_blockCols.size() > 100)
70 # endif
71  for (int i=0; i < static_cast<int>(_blockCols.size()); ++i) {
72  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
74  if (_hasStorage && dealloc)
75  delete b;
76  else
77  b->setZero();
78  }
79  if (_hasStorage && dealloc)
80  _blockCols[i].clear();
81  }
82  }
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
SparseBlockMatrix< MatrixType > * g2o::SparseBlockMatrix< MatrixType >::clone ( ) const

deep copy of a sparse-block-matrix;

Definition at line 122 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, and g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock().

122  {
124  for (size_t i=0; i<_blockCols.size(); ++i){
125  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
127  ret->_blockCols[i].insert(std::make_pair(it->first, b));
128  }
129  }
130  ret->_hasStorage=true;
131  return ret;
132  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock ( int  c) const
inline
template<class MatrixType = MatrixXD>
const std::vector<int>& g2o::SparseBlockMatrix< MatrixType >::colBlockIndices ( ) const
inline
template<class MatrixType = MatrixXD>
std::vector<int>& g2o::SparseBlockMatrix< MatrixType >::colBlockIndices ( )
inline

Definition at line 186 of file sparse_block_matrix.h.

186 { return _colBlockIndices;}
std::vector< int > _colBlockIndices
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::cols ( ) const
inline
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::colsOfBlock ( int  c) const
inline
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure ( MatrixStructure ms) const

exports the non zero blocks in the structure matrix ms

Definition at line 519 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::MatrixStructure::Aii, g2o::MatrixStructure::alloc(), g2o::MatrixStructure::Ap, g2o::MatrixStructure::m, and g2o::SparseBlockMatrix< MatrixType >::nonZeroBlocks().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), g2o::LinearSolverCholmodOnline< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverCSparse< MatrixType >::computeSymbolicDecomposition(), and g2o::LinearSolverCholmod< MatrixType >::computeSymbolicDecomposition().

520  {
521  int n = _colBlockIndices.size();
522  int nzMax = (int)nonZeroBlocks();
523 
524  ms.alloc(n, nzMax);
525  ms.m = _rowBlockIndices.size();
526 
527  int nz = 0;
528  int* Cp = ms.Ap;
529  int* Ci = ms.Aii;
530  for (int i = 0; i < static_cast<int>(_blockCols.size()); ++i){
531  *Cp = nz;
532  const int& c = i;
533  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it) {
534  const int& r = it->first;
535  if (r <= c) {
536  *Ci++ = r;
537  ++nz;
538  }
539  }
540  Cp++;
541  }
542  *Cp=nz;
543  assert(nz <= nzMax);
544  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
size_t nonZeroBlocks() const
number of allocated blocks
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
int g2o::SparseBlockMatrix< MatrixType >::fillCCS ( int *  Cp,
int *  Ci,
double *  Cx,
bool  upperTriangle = false 
) const

fill the CCS arrays of a matrix, arrays have to be allocated beforehand

Definition at line 489 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt(), and g2o::LinearSolverEigen< MatrixType >::fillSparseMatrix().

490  {
491  assert(Cp && Ci && Cx && "Target destination is NULL");
492  int nz=0;
493  for (size_t i=0; i<_blockCols.size(); ++i){
494  int cstart=i ? _colBlockIndices[i-1] : 0;
495  int csize=colsOfBlock(i);
496  for (int c=0; c<csize; ++c) {
497  *Cp=nz;
498  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
499  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* b=it->second;
500  int rstart=it->first ? _rowBlockIndices[it->first-1] : 0;
501 
502  int elemsToCopy = b->rows();
503  if (upperTriangle && rstart == cstart)
504  elemsToCopy = c + 1;
505  for (int r=0; r<elemsToCopy; ++r){
506  *Cx++ = (*b)(r,c);
507  *Ci++ = rstart++;
508  ++nz;
509  }
510  }
511  ++Cp;
512  }
513  }
514  *Cp=nz;
515  return nz;
516  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
int g2o::SparseBlockMatrix< MatrixType >::fillCCS ( double *  Cx,
bool  upperTriangle = false 
) const

fill the CCS arrays of a matrix, arrays have to be allocated beforehand. This function only writes the values and assumes that column and row structures have already been written.

Definition at line 464 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rows().

465  {
466  assert(Cx && "Target destination is NULL");
467  double* CxStart = Cx;
468  for (size_t i=0; i<_blockCols.size(); ++i){
469  int cstart=i ? _colBlockIndices[i-1] : 0;
470  int csize=colsOfBlock(i);
471  for (int c=0; c<csize; ++c) {
472  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
473  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* b=it->second;
474  int rstart=it->first ? _rowBlockIndices[it->first-1] : 0;
475 
476  int elemsToCopy = b->rows();
477  if (upperTriangle && rstart == cstart)
478  elemsToCopy = c + 1;
479  memcpy(Cx, b->data() + c*b->rows(), elemsToCopy * sizeof(double));
480  Cx += elemsToCopy;
481 
482  }
483  }
484  }
485  return Cx - CxStart;
486  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType>
int g2o::SparseBlockMatrix< MatrixType >::fillSparseBlockMatrixCCS ( SparseBlockMatrixCCS< MatrixType > &  blockCCS) const

copy into CCS structure

Returns
number of processed blocks, -1 on error

Definition at line 591 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrixCCS< MatrixType >::blockCols(), and g2o::SparseBlockMatrix< MatrixType >::blockCols().

Referenced by g2o::BlockSolver< Traits >::buildStructure(), g2o::SparseBlockMatrix< PoseMatrixType >::colBlockIndices(), and g2o::LinearSolverCCS< MatrixType >::initMatrixStructure().

592  {
593  blockCCS.blockCols().resize(blockCols().size());
594  int numblocks = 0;
595  for (size_t i = 0; i < blockCols().size(); ++i) {
596  const IntBlockMap& row = blockCols()[i];
597  typename SparseBlockMatrixCCS<MatrixType>::SparseColumn& dest = blockCCS.blockCols()[i];
598  dest.clear();
599  dest.reserve(row.size());
600  for (typename IntBlockMap::const_iterator it = row.begin(); it != row.end(); ++it) {
601  dest.push_back(typename SparseBlockMatrixCCS<MatrixType>::RowBlock(it->first, it->second));
602  ++numblocks;
603  }
604  }
605  return numblocks;
606  }
std::vector< RowBlock > SparseColumn
std::map< int, SparseMatrixBlock * > IntBlockMap
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
template<class MatrixType>
int g2o::SparseBlockMatrix< MatrixType >::fillSparseBlockMatrixCCSTransposed ( SparseBlockMatrixCCS< MatrixType > &  blockCCS) const

copy as transposed into a CCS structure

Returns
number of processed blocks, -1 on error

Definition at line 609 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrixCCS< MatrixType >::blockCols(), and g2o::SparseBlockMatrix< MatrixType >::blockCols().

Referenced by g2o::BlockSolver< Traits >::buildStructure(), and g2o::SparseBlockMatrix< PoseMatrixType >::colBlockIndices().

610  {
611  blockCCS.blockCols().clear();
612  blockCCS.blockCols().resize(_rowBlockIndices.size());
613  int numblocks = 0;
614  for (size_t i = 0; i < blockCols().size(); ++i) {
615  const IntBlockMap& row = blockCols()[i];
616  for (typename IntBlockMap::const_iterator it = row.begin(); it != row.end(); ++it) {
617  typename SparseBlockMatrixCCS<MatrixType>::SparseColumn& dest = blockCCS.blockCols()[it->first];
618  dest.push_back(typename SparseBlockMatrixCCS<MatrixType>::RowBlock(i, it->second));
619  ++numblocks;
620  }
621  }
622  return numblocks;
623  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< RowBlock > SparseColumn
std::map< int, SparseMatrixBlock * > IntBlockMap
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
template<class MatrixType >
template<class MatrixResultType , class MatrixFactorType >
bool g2o::SparseBlockMatrix< MatrixType >::multiply ( SparseBlockMatrix< MatrixResultType > *&  dest,
const SparseBlockMatrix< MatrixFactorType > *  M 
) const

dest = (*this) * M

Definition at line 199 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::block(), g2o::SparseBlockMatrix< MatrixType >::cols(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), and main().

199  {
200  // sanity check
201  if (_colBlockIndices.size()!=M->_rowBlockIndices.size())
202  return false;
203  for (size_t i=0; i<_colBlockIndices.size(); ++i){
204  if (_colBlockIndices[i]!=M->_rowBlockIndices[i])
205  return false;
206  }
207  if (! dest) {
208  dest=new SparseBlockMatrix<MatrixResultType>(&_rowBlockIndices[0],&(M->_colBlockIndices[0]), _rowBlockIndices.size(), M->_colBlockIndices.size() );
209  }
210  if (! dest->_hasStorage)
211  return false;
212  for (size_t i=0; i<M->_blockCols.size(); ++i){
213  for (typename SparseBlockMatrix<MatrixFactorType>::IntBlockMap::const_iterator it=M->_blockCols[i].begin(); it!=M->_blockCols[i].end(); ++it){
214  // look for a non-zero block in a row of column it
215  int colM=i;
216  const typename SparseBlockMatrix<MatrixFactorType>::SparseMatrixBlock *b=it->second;
217  typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator rbt=_blockCols[it->first].begin();
218  while(rbt!=_blockCols[it->first].end()){
219  //int colA=it->first;
220  int rowA=rbt->first;
221  typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock *a=rbt->second;
222  typename SparseBlockMatrix<MatrixResultType>::SparseMatrixBlock *c=dest->block(rowA,colM,true);
223  assert (c->rows()==a->rows());
224  assert (c->cols()==b->cols());
225  ++rbt;
226  (*c)+=(*a)*(*b);
227  }
228  }
229  }
230  return false;
231  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::multiply ( double *&  dest,
const double *  src 
) const

dest = (*this) * src

Definition at line 234 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::cols(), and g2o::SparseBlockMatrix< MatrixType >::rows().

234  {
235  if (! dest){
236  dest=new double [_rowBlockIndices[_rowBlockIndices.size()-1] ];
237  memset(dest,0, _rowBlockIndices[_rowBlockIndices.size()-1]*sizeof(double));
238  }
239 
240  // map the memory by Eigen
241  Eigen::Map<VectorXD> destVec(dest, rows());
242  const Eigen::Map<const VectorXD> srcVec(src, cols());
243 
244  for (size_t i=0; i<_blockCols.size(); ++i){
245  int srcOffset = i ? _colBlockIndices[i-1] : 0;
246 
247  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
248  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* a=it->second;
249  int destOffset = it->first ? _rowBlockIndices[it->first - 1] : 0;
250  // destVec += *a * srcVec (according to the sub-vector parts)
251  internal::template axpy<typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock>(*a, srcVec, srcOffset, destVec, destOffset);
252  }
253  }
254  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
int cols() const
columns of the matrix
int rows() const
rows of the matrix
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::multiplySymmetricUpperTriangle ( double *&  dest,
const double *  src 
) const

compute dest = (*this) * src However, assuming that this is a symmetric matrix where only the upper triangle is stored

Definition at line 257 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock().

258  {
259  if (! dest){
260  dest=new double [_rowBlockIndices[_rowBlockIndices.size()-1] ];
261  memset(dest,0, _rowBlockIndices[_rowBlockIndices.size()-1]*sizeof(double));
262  }
263 
264  // map the memory by Eigen
265  Eigen::Map<VectorXD> destVec(dest, rows());
266  const Eigen::Map<const VectorXD> srcVec(src, cols());
267 
268  for (size_t i=0; i<_blockCols.size(); ++i){
269  int srcOffset = colBaseOfBlock(i);
270  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
271  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* a=it->second;
272  int destOffset = rowBaseOfBlock(it->first);
273  if (destOffset > srcOffset) // only upper triangle
274  break;
275  // destVec += *a * srcVec (according to the sub-vector parts)
276  internal::template axpy<typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock>(*a, srcVec, srcOffset, destVec, destOffset);
277  if (destOffset < srcOffset)
278  internal::template atxpy<typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock>(*a, srcVec, destOffset, destVec, srcOffset);
279  }
280  }
281  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
int cols() const
columns of the matrix
int rows() const
rows of the matrix
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
int colBaseOfBlock(int c) const
where does the col at block-col r starts?
int rowBaseOfBlock(int r) const
where does the row at block-row r starts?
template<class MatrixType >
size_t g2o::SparseBlockMatrix< MatrixType >::nonZeroBlocks ( ) const

number of allocated blocks

Definition at line 353 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols.

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure(), and g2o::SparseBlockMatrix< MatrixType >::nonZeros().

353  {
354  size_t count=0;
355  for (size_t i=0; i<_blockCols.size(); ++i)
356  count+=_blockCols[i].size();
357  return count;
358  }
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
size_t g2o::SparseBlockMatrix< MatrixType >::nonZeros ( ) const

number of non-zero elements

Definition at line 361 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::nonZeroBlocks(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt(), g2o::LinearSolverCholmod< MatrixType >::fillCholmodExt(), g2o::LinearSolverCSparse< MatrixType >::fillCSparse(), and g2o::LinearSolverEigen< MatrixType >::fillSparseMatrix().

361  {
362  if (MatrixType::SizeAtCompileTime != Eigen::Dynamic) {
363  size_t nnz = nonZeroBlocks() * MatrixType::SizeAtCompileTime;
364  return nnz;
365  } else {
366  size_t count=0;
367  for (size_t i=0; i<_blockCols.size(); ++i){
368  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
369  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* a=it->second;
370  count += a->cols()*a->rows();
371  }
372  }
373  return count;
374  }
375  }
size_t nonZeroBlocks() const
number of allocated blocks
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::rightMultiply ( double *&  dest,
const double *  src 
) const

dest = M * (*this)

Definition at line 284 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock().

284  {
285  int destSize=cols();
286 
287  if (! dest){
288  dest=new double [ destSize ];
289  memset(dest,0, destSize*sizeof(double));
290  }
291 
292  // map the memory by Eigen
293  Eigen::Map<VectorXD> destVec(dest, destSize);
294  Eigen::Map<const VectorXD> srcVec(src, rows());
295 
296 # ifdef G2O_OPENMP
297 # pragma omp parallel for default (shared) schedule(dynamic, 10)
298 # endif
299  for (int i=0; i < static_cast<int>(_blockCols.size()); ++i){
300  int destOffset = colBaseOfBlock(i);
301  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin();
302  it!=_blockCols[i].end();
303  ++it){
304  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* a=it->second;
305  int srcOffset = rowBaseOfBlock(it->first);
306  // destVec += *a.transpose() * srcVec (according to the sub-vector parts)
307  internal::template atxpy<typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock>(*a, srcVec, srcOffset, destVec, destOffset);
308  }
309  }
310 
311  }
int cols() const
columns of the matrix
int rows() const
rows of the matrix
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
int colBaseOfBlock(int c) const
where does the col at block-col r starts?
int rowBaseOfBlock(int r) const
where does the row at block-row r starts?
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock ( int  r) const
inline
template<class MatrixType = MatrixXD>
const std::vector<int>& g2o::SparseBlockMatrix< MatrixType >::rowBlockIndices ( ) const
inline
template<class MatrixType = MatrixXD>
std::vector<int>& g2o::SparseBlockMatrix< MatrixType >::rowBlockIndices ( )
inline

Definition at line 182 of file sparse_block_matrix.h.

182 { return _rowBlockIndices;}
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::rows ( ) const
inline
template<class MatrixType = MatrixXD>
int g2o::SparseBlockMatrix< MatrixType >::rowsOfBlock ( int  r) const
inline
template<class MatrixType >
void g2o::SparseBlockMatrix< MatrixType >::scale ( double  a)

*this *= a

Definition at line 314 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols.

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock().

314  {
315  for (size_t i=0; i<_blockCols.size(); ++i){
316  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
318  *a *= a_;
319  }
320  }
321  }
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
SparseBlockMatrix< MatrixType > * g2o::SparseBlockMatrix< MatrixType >::slice ( int  rmin,
int  rmax,
int  cmin,
int  cmax,
bool  alloc = true 
) const

returns a view or a copy of the block matrix

Parameters
rminstarting block row
rmaxending block row
cminstarting block col
cmaxending block col
allocif true it makes a deep copy, if false it creates a view.

Definition at line 324 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), g2o::SparseBlockMatrix< MatrixType >::rowsOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock().

324  {
325  int m=rmax-rmin;
326  int n=cmax-cmin;
327  int rowIdx [m];
328  rowIdx[0] = rowsOfBlock(rmin);
329  for (int i=1; i<m; ++i){
330  rowIdx[i]=rowIdx[i-1]+rowsOfBlock(rmin+i);
331  }
332 
333  int colIdx [n];
334  colIdx[0] = colsOfBlock(cmin);
335  for (int i=1; i<n; ++i){
336  colIdx[i]=colIdx[i-1]+colsOfBlock(cmin+i);
337  }
338  typename SparseBlockMatrix<MatrixType>::SparseBlockMatrix* s=new SparseBlockMatrix(rowIdx, colIdx, m, n, true);
339  for (int i=0; i<n; ++i){
340  int mc=cmin+i;
341  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[mc].begin(); it!=_blockCols[mc].end(); ++it){
342  if (it->first >= rmin && it->first < rmax){
343  typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* b = alloc ? new typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock (* (it->second) ) : it->second;
344  s->_blockCols[i].insert(std::make_pair(it->first-rmin, b));
345  }
346  }
347  }
348  s->_hasStorage=alloc;
349  return s;
350  }
int rowsOfBlock(int r) const
how many rows does the block at block-row r has?
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType>
bool g2o::SparseBlockMatrix< MatrixType >::symmPermutation ( SparseBlockMatrix< MatrixType > *&  dest,
const int *  pinv,
bool  onlyUpper = false 
) const

writes in dest a block permutaton specified by pinv.

Parameters
pinvarray such that new_block[i] = old_block[pinv[i]]

Definition at line 399 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::block(), g2o::SparseBlockMatrix< MatrixType >::clear(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::rows(), g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix(), and g2o::SparseBlockMatrix< MatrixType >::transpose().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), and main().

399  {
400  // compute the permuted version of the new row/column layout
401  size_t n=_rowBlockIndices.size();
402  // computed the block sizes
403  int blockSizes[_rowBlockIndices.size()];
404  blockSizes[0]=_rowBlockIndices[0];
405  for (size_t i=1; i<n; ++i){
406  blockSizes[i]=_rowBlockIndices[i]-_rowBlockIndices[i-1];
407  }
408  // permute them
409  int pBlockIndices[_rowBlockIndices.size()];
410  for (size_t i=0; i<n; ++i){
411  pBlockIndices[pinv[i]]=blockSizes[i];
412  }
413  for (size_t i=1; i<n; ++i){
414  pBlockIndices[i]+=pBlockIndices[i-1];
415  }
416  // allocate C, or check the structure;
417  if (! dest){
418  dest=new SparseBlockMatrix(pBlockIndices, pBlockIndices, n, n);
419  } else {
420  if (dest->_rowBlockIndices.size()!=n)
421  return false;
422  if (dest->_colBlockIndices.size()!=n)
423  return false;
424  for (size_t i=0; i<n; ++i){
425  if (dest->_rowBlockIndices[i]!=pBlockIndices[i])
426  return false;
427  if (dest->_colBlockIndices[i]!=pBlockIndices[i])
428  return false;
429  }
430  dest->clear();
431  }
432  // now ready to permute the columns
433  for (size_t i=0; i<n; ++i){
434  //cerr << PVAR(i) << " ";
435  int pi=pinv[i];
436  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin();
437  it!=_blockCols[i].end(); ++it){
438  int pj=pinv[it->first];
439 
440  const typename SparseBlockMatrix<MatrixType>::SparseMatrixBlock* s=it->second;
441 
443  if (! upperTriangle || pj<=pi) {
444  b=dest->block(pj,pi,true);
445  assert(b->cols()==s->cols());
446  assert(b->rows()==s->rows());
447  *b=*s;
448  } else {
449  b=dest->block(pi,pj,true);
450  assert(b);
451  assert(b->rows()==s->cols());
452  assert(b->cols()==s->rows());
453  *b=s->transpose();
454  }
455  }
456  //cerr << endl;
457  // within each row,
458  }
459  return true;
460 
461  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType>
void g2o::SparseBlockMatrix< MatrixType >::takePatternFromHash ( SparseBlockMatrixHashMap< MatrixType > &  hashMatrix)

take over the memory and matrix pattern from a hash matrix. The structure of the hash matrix will be cleared.

Definition at line 626 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::blockCols(), and g2o::SparseBlockMatrixHashMap< MatrixType >::blockCols().

Referenced by g2o::BlockSolver< Traits >::buildStructure(), and g2o::SparseBlockMatrix< PoseMatrixType >::colBlockIndices().

627  {
628  // sort the sparse columns and add them to the map structures by
629  // exploiting that we are inserting a sorted structure
630  typedef std::pair<int, MatrixType*> SparseColumnPair;
631  typedef typename SparseBlockMatrixHashMap<MatrixType>::SparseColumn HashSparseColumn;
632  for (size_t i = 0; i < hashMatrix.blockCols().size(); ++i) {
633  // prepare a temporary vector for sorting
634  HashSparseColumn& column = hashMatrix.blockCols()[i];
635  if (column.size() == 0)
636  continue;
637  std::vector<SparseColumnPair> sparseRowSorted; // temporary structure
638  sparseRowSorted.reserve(column.size());
639  for (typename HashSparseColumn::const_iterator it = column.begin(); it != column.end(); ++it)
640  sparseRowSorted.push_back(*it);
641  std::sort(sparseRowSorted.begin(), sparseRowSorted.end(), CmpPairFirst<int, MatrixType*>());
642  // try to free some memory early
643  HashSparseColumn aux;
644  std::swap(aux, column);
645  // now insert sorted vector to the std::map structure
646  IntBlockMap& destColumnMap = blockCols()[i];
647  destColumnMap.insert(sparseRowSorted[0]);
648  for (size_t j = 1; j < sparseRowSorted.size(); ++j) {
649  typename SparseBlockMatrix<MatrixType>::IntBlockMap::iterator hint = destColumnMap.end();
650  --hint; // cppreference says the element goes after the hint (until C++11)
651  destColumnMap.insert(hint, sparseRowSorted[j]);
652  }
653  }
654  }
std::unordered_map< int, MatrixType * > SparseColumn
std::map< int, SparseMatrixBlock * > IntBlockMap
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
template<class MatrixType >
template<class MatrixTransposedType >
bool g2o::SparseBlockMatrix< MatrixType >::transpose ( SparseBlockMatrix< MatrixTransposedType > *&  dest) const

transposes a block matrix, The transposed type should match the argument false on failure

Definition at line 137 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::block(), and g2o::SparseBlockMatrix< MatrixType >::transpose().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBaseOfBlock(), main(), g2o::LinearSolverDense< MatrixType >::solve(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), and g2o::SparseBlockMatrix< MatrixType >::transpose().

137  {
138  if (! dest){
139  dest=new SparseBlockMatrix<MatrixTransposedType>(&_colBlockIndices[0], &_rowBlockIndices[0], _colBlockIndices.size(), _rowBlockIndices.size());
140  } else {
141  if (! dest->_hasStorage)
142  return false;
143  if(_rowBlockIndices.size()!=dest->_colBlockIndices.size())
144  return false;
145  if (_colBlockIndices.size()!=dest->_rowBlockIndices.size())
146  return false;
147  for (size_t i=0; i<_rowBlockIndices.size(); ++i){
148  if(_rowBlockIndices[i]!=dest->_colBlockIndices[i])
149  return false;
150  }
151  for (size_t i=0; i<_colBlockIndices.size(); ++i){
152  if(_colBlockIndices[i]!=dest->_rowBlockIndices[i])
153  return false;
154  }
155  }
156 
157  for (size_t i=0; i<_blockCols.size(); ++i){
158  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it){
160  typename SparseBlockMatrix<MatrixTransposedType>::SparseMatrixBlock* d=dest->block(i,it->first,true);
161  *d = s->transpose();
162  }
163  }
164  return true;
165  }
std::vector< int > _rowBlockIndices
vector of the indices of the blocks along the rows.
std::vector< int > _colBlockIndices
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
std::vector< IntBlockMap > _blockCols
template<class MatrixType >
bool g2o::SparseBlockMatrix< MatrixType >::writeOctave ( const char *  filename,
bool  upperTriangle = true 
) const

write the content of this matrix to a stream loadable by Octave

Parameters
upperTriangledoes this matrix store only the upper triangular blocks

Definition at line 547 of file sparse_block_matrix.hpp.

References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock(), and g2o::SparseBlockMatrix< MatrixType >::rows().

Referenced by g2o::SparseBlockMatrix< PoseMatrixType >::colBlockIndices(), g2o::BlockSolver< Traits >::saveHessian(), and g2o::LinearSolverEigen< MatrixType >::solve().

548  {
549  std::string name = filename;
550  std::string::size_type lastDot = name.find_last_of('.');
551  if (lastDot != std::string::npos)
552  name = name.substr(0, lastDot);
553 
554  std::vector<TripletEntry> entries;
555  for (size_t i = 0; i<_blockCols.size(); ++i){
556  const int& c = i;
557  for (typename SparseBlockMatrix<MatrixType>::IntBlockMap::const_iterator it=_blockCols[i].begin(); it!=_blockCols[i].end(); ++it) {
558  const int& r = it->first;
559  const MatrixType& m = *(it->second);
560  for (int cc = 0; cc < m.cols(); ++cc)
561  for (int rr = 0; rr < m.rows(); ++rr) {
562  int aux_r = rowBaseOfBlock(r) + rr;
563  int aux_c = colBaseOfBlock(c) + cc;
564  entries.push_back(TripletEntry(aux_r, aux_c, m(rr, cc)));
565  if (upperTriangle && r != c) {
566  entries.push_back(TripletEntry(aux_c, aux_r, m(rr, cc)));
567  }
568  }
569  }
570  }
571 
572  int nz = entries.size();
573  std::sort(entries.begin(), entries.end(), TripletColSort());
574 
575  std::ofstream fout(filename);
576  fout << "# name: " << name << std::endl;
577  fout << "# type: sparse matrix" << std::endl;
578  fout << "# nnz: " << nz << std::endl;
579  fout << "# rows: " << rows() << std::endl;
580  fout << "# columns: " << cols() << std::endl;
581  fout << std::setprecision(9) << std::fixed << std::endl;
582 
583  for (std::vector<TripletEntry>::const_iterator it = entries.begin(); it != entries.end(); ++it) {
584  const TripletEntry& entry = *it;
585  fout << entry.r+1 << " " << entry.c+1 << " " << entry.x << std::endl;
586  }
587  return fout.good();
588  }
int cols() const
columns of the matrix
int rows() const
rows of the matrix
std::vector< IntBlockMap > _blockCols
int colBaseOfBlock(int c) const
where does the col at block-col r starts?
int rowBaseOfBlock(int r) const
where does the row at block-row r starts?

Member Data Documentation

template<class MatrixType = MatrixXD>
std::vector<IntBlockMap> g2o::SparseBlockMatrix< MatrixType >::_blockCols
protected
template<class MatrixType = MatrixXD>
std::vector<int> g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices
protected
template<class MatrixType = MatrixXD>
bool g2o::SparseBlockMatrix< MatrixType >::_hasStorage
protected
template<class MatrixType = MatrixXD>
std::vector<int> g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices
protected

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