|
g2o
|
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. | |
| typedef std::map< int, SparseMatrixBlock * > | IntBlockMap |
Public Member Functions | |
| int | cols () const |
| columns of the matrix | |
| int | rows () const |
| rows of the matrix | |
| SparseBlockMatrix (const int *rbi, const int *cbi, int rb, int cb, bool hasStorage=true) | |
| SparseBlockMatrix () | |
| ~SparseBlockMatrix () | |
| void | clear (bool dealloc=false) |
| SparseMatrixBlock * | block (int r, int c, bool alloc=false) |
| const SparseMatrixBlock * | block (int r, int c) const |
| returns the block at location r,c | |
| 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? | |
| int | rowBaseOfBlock (int r) const |
| where does the row at block-row r starts? | |
| int | colBaseOfBlock (int c) const |
| where does the col at block-col r starts? | |
| size_t | nonZeros () const |
| number of non-zero elements | |
| size_t | nonZeroBlocks () const |
| number of allocated blocks | |
| SparseBlockMatrix * | clone () const |
| deep copy of a sparse-block-matrix; | |
| SparseBlockMatrix * | slice (int rmin, int rmax, int cmin, int cmax, bool alloc=true) const |
| template<class MatrixTransposedType > | |
| bool | transpose (SparseBlockMatrix< MatrixTransposedType > &dest) const |
| template<class MatrixTransposedType > | |
| std::unique_ptr< SparseBlockMatrix< MatrixTransposedType > > | transposed () const |
| bool | add (SparseBlockMatrix< MatrixType > &dest) const |
| adds the current matrix to the destination | |
| std::unique_ptr< SparseBlockMatrix< MatrixType > > | added () const |
| template<class MatrixResultType , class MatrixFactorType > | |
| bool | multiply (SparseBlockMatrix< MatrixResultType > *&dest, const SparseBlockMatrix< MatrixFactorType > *M) const |
| dest = (*this) * M | |
| void | multiply (double *&dest, const double *src) const |
| dest = (*this) * src | |
| void | multiplySymmetricUpperTriangle (double *&dest, const double *src) const |
| void | rightMultiply (double *&dest, const double *src) const |
| dest = M * (*this) | |
| void | scale (double a) |
| *this *= a | |
| 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 | |
| void | fillBlockStructure (int *Cp, int *Ci) const |
| exports the non zero blocks into a column compressed structure | |
| const std::vector< IntBlockMap > & | blockCols () const |
| the block matrices per block-column | |
| std::vector< IntBlockMap > & | blockCols () |
| const std::vector< int > & | rowBlockIndices () const |
| indices of the row blocks | |
| std::vector< int > & | rowBlockIndices () |
| const std::vector< int > & | colBlockIndices () const |
| indices of the column blocks | |
| 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 |
| std::vector< int > | _colBlockIndices |
| std::vector< IntBlockMap > | _blockCols |
| bool | _hasStorage |
Private Member Functions | |
| template<class MatrixTransposedType > | |
| void | transpose_internal (SparseBlockMatrix< MatrixTransposedType > &dest) const |
| void | add_internal (SparseBlockMatrix< MatrixType > &dest) const |
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 64 of file sparse_block_matrix.h.
| typedef std::map<int, SparseMatrixBlock*> g2o::SparseBlockMatrix< MatrixType >::IntBlockMap |
Definition at line 78 of file sparse_block_matrix.h.
| typedef MatrixType g2o::SparseBlockMatrix< MatrixType >::SparseMatrixBlock |
this is the type of the elementary block, it is an Eigen::Matrix.
Definition at line 67 of file sparse_block_matrix.h.
| 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
| rbi | array 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. |
| cbi | array 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. |
| rb | number of row blocks |
| cb | number of col blocks |
| hasStorage | set 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 32 of file sparse_block_matrix.hpp.
| g2o::SparseBlockMatrix< MatrixType >::SparseBlockMatrix | ( | ) |
Definition at line 41 of file sparse_block_matrix.hpp.
| g2o::SparseBlockMatrix< MatrixType >::~SparseBlockMatrix | ( | ) |
Definition at line 64 of file sparse_block_matrix.hpp.
| bool g2o::SparseBlockMatrix< MatrixType >::add | ( | SparseBlockMatrix< MatrixType > & | dest | ) | const |
adds the current matrix to the destination
Definition at line 184 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, and g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
|
private |
Definition at line 169 of file sparse_block_matrix.hpp.
| std::unique_ptr< SparseBlockMatrix< MatrixType > > g2o::SparseBlockMatrix< MatrixType >::added | ( | ) | const |
Definition at line 202 of file sparse_block_matrix.hpp.
| const SparseBlockMatrix< MatrixType >::SparseMatrixBlock * g2o::SparseBlockMatrix< MatrixType >::block | ( | int | r, |
| int | c | ||
| ) | const |
returns the block at location r,c
Definition at line 97 of file sparse_block_matrix.hpp.
| 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 70 of file sparse_block_matrix.hpp.
Referenced by g2o::MarginalCovarianceCholesky::computeCovariance(), g2o::Slam2DViewer::draw(), g2o::EdgeLabeler::labelEdge(), main(), main(), g2o::SolverSLAM2DLinear::solveOrientation(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), testMarginals(), g2o::SparseBlockMatrix< MatrixType >::transpose_internal(), and g2o::SparseOptimizerIncremental::updateInitialization().
|
inline |
Definition at line 212 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_blockCols.
|
inline |
the block matrices per block-column
Definition at line 211 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_blockCols.
Referenced by g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition(), g2o::operator<<(), g2o::LinearSolverDense< MatrixType >::solve(), g2o::LinearSolverPCG< MatrixType >::solve(), and g2o::SparseOptimizerIncremental::updateInitialization().
| 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 45 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::clear().
Referenced by g2o::SparseBlockMatrix< MatrixType >::clear(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), g2o::SparseOptimizerIncremental::updateInitialization(), and g2o::SparseOptimizerIncremental::~SparseOptimizerIncremental().
| SparseBlockMatrix< MatrixType > * g2o::SparseBlockMatrix< MatrixType >::clone | ( | ) | const |
deep copy of a sparse-block-matrix;
Definition at line 105 of file sparse_block_matrix.hpp.
|
inline |
where does the col at block-col r starts?
Definition at line 129 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::MarginalCovarianceCholesky::computeCovariance(), g2o::LinearSolverCholmodOnline< MatrixType >::computeSymbolicDecomposition(), and g2o::LinearSolverDense< MatrixType >::solve().
|
inline |
Definition at line 220 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices.
|
inline |
indices of the column blocks
Definition at line 219 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::LinearSolverCCS< MatrixType >::initMatrixStructure(), g2o::operator<<(), g2o::LinearSolverPCG< MatrixType >::solve(), and g2o::SparseOptimizerIncremental::updateInitialization().
|
inline |
columns of the matrix
Definition at line 70 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::LinearSolverEigen< MatrixType >::computeCholesky(), g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt(), g2o::LinearSolverCholmod< MatrixType >::fillCholmodExt(), g2o::LinearSolverCSparse< MatrixType >::fillCSparse(), g2o::SparseBlockMatrix< MatrixType >::nonZeros(), g2o::LinearSolverDense< MatrixType >::solve(), g2o::LinearSolverPCG< MatrixType >::solve(), and g2o::SparseBlockMatrix< MatrixType >::symmPermutation().
|
inline |
how many cols does the block at block-col c has?
Definition at line 118 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::allocateBlocks(), g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::LinearSolverCholmodOnline< MatrixType >::computeSymbolicDecomposition(), and g2o::LinearSolverDense< MatrixType >::solve().
| void g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure | ( | int * | Cp, |
| int * | Ci | ||
| ) | const |
exports the non zero blocks into a column compressed structure
Definition at line 582 of file sparse_block_matrix.hpp.
| void g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure | ( | MatrixStructure & | ms | ) | const |
exports the non zero blocks in the structure matrix ms
Definition at line 574 of file sparse_block_matrix.hpp.
References g2o::MatrixStructure::Aii, g2o::MatrixStructure::alloc(), g2o::MatrixStructure::Ap, and g2o::MatrixStructure::m.
Referenced by g2o::LinearSolverCholmodOnline< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverCholmod< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverCSparse< MatrixType >::computeSymbolicDecomposition(), and g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition().
| 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 516 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::rows().
| 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 542 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::rows().
Referenced by g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), and g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt().
| int g2o::SparseBlockMatrix< MatrixType >::fillSparseBlockMatrixCCS | ( | SparseBlockMatrixCCS< MatrixType > & | blockCCS | ) | const |
copy into CCS structure
Definition at line 646 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrixCCS< MatrixType >::blockCols().
Referenced by g2o::LinearSolverCCS< MatrixType >::initMatrixStructure().
| int g2o::SparseBlockMatrix< MatrixType >::fillSparseBlockMatrixCCSTransposed | ( | SparseBlockMatrixCCS< MatrixType > & | blockCCS | ) | const |
copy as transposed into a CCS structure
Definition at line 667 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrixCCS< MatrixType >::blockCols().
| void g2o::SparseBlockMatrix< MatrixType >::multiply | ( | double *& | dest, |
| const double * | src | ||
| ) | const |
dest = (*this) * src
Definition at line 255 of file sparse_block_matrix.hpp.
| bool g2o::SparseBlockMatrix< MatrixType >::multiply | ( | SparseBlockMatrix< MatrixResultType > *& | dest, |
| const SparseBlockMatrix< MatrixFactorType > * | M | ||
| ) | const |
dest = (*this) * M
Definition at line 212 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::_blockCols, g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, and g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
| 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 285 of file sparse_block_matrix.hpp.
| size_t g2o::SparseBlockMatrix< MatrixType >::nonZeroBlocks | ( | ) | const |
number of allocated blocks
Definition at line 401 of file sparse_block_matrix.hpp.
Referenced by g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition().
| size_t g2o::SparseBlockMatrix< MatrixType >::nonZeros | ( | ) | const |
number of non-zero elements
Definition at line 408 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::cols(), and g2o::SparseBlockMatrix< MatrixType >::rows().
Referenced by g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt(), g2o::LinearSolverCholmod< MatrixType >::fillCholmodExt(), g2o::LinearSolverCSparse< MatrixType >::fillCSparse(), and g2o::LinearSolverEigen< MatrixType >::fillSparseMatrix().
| void g2o::SparseBlockMatrix< MatrixType >::rightMultiply | ( | double *& | dest, |
| const double * | src | ||
| ) | const |
dest = M * (*this)
Definition at line 320 of file sparse_block_matrix.hpp.
|
inline |
where does the row at block-row r starts?
Definition at line 124 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
Referenced by g2o::MarginalCovarianceCholesky::computeCovariance(), and g2o::LinearSolverDense< MatrixType >::solve().
|
inline |
Definition at line 216 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
|
inline |
indices of the row blocks
Definition at line 215 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::allocateBlocks(), g2o::LinearSolver< MatrixType >::deallocateBlocks(), g2o::LinearSolverCCS< MatrixType >::initMatrixStructure(), g2o::operator<<(), g2o::LinearSolverPCG< MatrixType >::solve(), g2o::LinearSolverCCS< MatrixType >::solveBlocks(), g2o::LinearSolverCCS< MatrixType >::solvePattern(), and g2o::SparseOptimizerIncremental::updateInitialization().
|
inline |
rows of the matrix
Definition at line 74 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::allocateBlocks(), g2o::LinearSolverEigen< MatrixType >::computeCholesky(), g2o::SparseOptimizerIncremental::computeCholeskyUpdate(), g2o::LinearSolverEigen< MatrixType >::computeSymbolicDecomposition(), g2o::SparseBlockMatrix< MatrixType >::fillCCS(), g2o::SparseBlockMatrix< MatrixType >::fillCCS(), g2o::LinearSolverCholmodOnline< MatrixType >::fillCholmodExt(), g2o::LinearSolverCholmod< MatrixType >::fillCholmodExt(), g2o::LinearSolverCSparse< MatrixType >::fillCSparse(), g2o::SparseBlockMatrix< MatrixType >::nonZeros(), g2o::LinearSolverPCG< MatrixType >::solve(), and g2o::SparseBlockMatrix< MatrixType >::symmPermutation().
|
inline |
how many rows does the block at block-row r has?
Definition at line 112 of file sparse_block_matrix.h.
References g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
Referenced by g2o::LinearSolver< MatrixType >::allocateBlocks(), and g2o::LinearSolverDense< MatrixType >::solve().
| void g2o::SparseBlockMatrix< MatrixType >::scale | ( | double | a | ) |
*this *= a
Definition at line 353 of file sparse_block_matrix.hpp.
| 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
| rmin | starting block row |
| rmax | ending block row |
| cmin | starting block col |
| cmax | ending block col |
| alloc | if true it makes a deep copy, if false it creates a view. |
Definition at line 365 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::_blockCols, and g2o::SparseBlockMatrix< MatrixType >::_hasStorage.
| bool g2o::SparseBlockMatrix< MatrixType >::symmPermutation | ( | SparseBlockMatrix< MatrixType > *& | dest, |
| const int * | pinv, | ||
| bool | onlyUpper = false |
||
| ) | const |
writes in dest a block permutaton specified by pinv.
| pinv | array such that new_block[i] = old_block[pinv[i]] |
Definition at line 453 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices, g2o::SparseBlockMatrix< MatrixType >::block(), g2o::SparseBlockMatrix< MatrixType >::clear(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::rows(), and g2o::SparseBlockMatrix< MatrixType >::transpose().
| 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 687 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrixHashMap< MatrixType >::blockCols().
| 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 141 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::_colBlockIndices, g2o::SparseBlockMatrix< MatrixType >::_hasStorage, and g2o::SparseBlockMatrix< MatrixType >::_rowBlockIndices.
Referenced by g2o::LinearSolverDense< MatrixType >::solve(), and g2o::SparseBlockMatrix< MatrixType >::symmPermutation().
|
private |
Definition at line 125 of file sparse_block_matrix.hpp.
References g2o::SparseBlockMatrix< MatrixType >::block().
| std::unique_ptr< SparseBlockMatrix< MatrixTransposedType > > g2o::SparseBlockMatrix< MatrixType >::transposed | ( | ) | const |
Definition at line 160 of file sparse_block_matrix.hpp.
| bool g2o::SparseBlockMatrix< MatrixType >::writeOctave | ( | const char * | filename, |
| bool | upperTriangle = true |
||
| ) | const |
write the content of this matrix to a stream loadable by Octave
| upperTriangle | does this matrix store only the upper triangular blocks |
Definition at line 600 of file sparse_block_matrix.hpp.
References g2o::TripletEntry::c, g2o::TripletEntry::r, and g2o::TripletEntry::x.
Referenced by g2o::LinearSolverEigen< MatrixType >::computeCholesky().
|
protected |
array of maps of blocks. The index of the array represent a block column of the matrix and the block column is stored as a map row_block -> matrix_block_ptr.
Definition at line 257 of file sparse_block_matrix.h.
Referenced by g2o::SparseBlockMatrix< g2o::MatrixX >::block(), g2o::SparseBlockMatrix< MatrixType >::blockCols(), g2o::SparseBlockMatrix< MatrixType >::blockCols(), g2o::SparseBlockMatrix< MatrixType >::multiply(), and g2o::SparseBlockMatrix< MatrixType >::slice().
|
protected |
vector of the indices of the blocks along the cols
Definition at line 253 of file sparse_block_matrix.h.
Referenced by g2o::SparseBlockMatrix< MatrixType >::add(), g2o::SparseBlockMatrix< MatrixType >::colBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::colBlockIndices(), g2o::SparseBlockMatrix< MatrixType >::colBlockIndices(), g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::SparseBlockMatrix< MatrixType >::colsOfBlock(), g2o::SparseBlockMatrix< MatrixType >::multiply(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), and g2o::SparseBlockMatrix< MatrixType >::transpose().
|
protected |
Definition at line 258 of file sparse_block_matrix.h.
Referenced by g2o::SparseBlockMatrix< MatrixType >::add(), g2o::SparseBlockMatrix< MatrixType >::multiply(), g2o::SparseBlockMatrix< MatrixType >::slice(), and g2o::SparseBlockMatrix< MatrixType >::transpose().
|
protected |
vector of the indices of the blocks along the rows.
Definition at line 250 of file sparse_block_matrix.h.
Referenced by g2o::SparseBlockMatrix< MatrixType >::add(), g2o::SparseBlockMatrix< MatrixType >::multiply(), g2o::SparseBlockMatrix< MatrixType >::rowBaseOfBlock(), g2o::SparseBlockMatrix< MatrixType >::rowBlockIndices(), g2o::SparseBlockMatrix< MatrixType >::rowBlockIndices(), g2o::SparseBlockMatrix< MatrixType >::rows(), g2o::SparseBlockMatrix< MatrixType >::rowsOfBlock(), g2o::SparseBlockMatrix< g2o::MatrixX >::slice(), g2o::SparseBlockMatrix< g2o::MatrixX >::symmPermutation(), g2o::SparseBlockMatrix< MatrixType >::symmPermutation(), and g2o::SparseBlockMatrix< MatrixType >::transpose().