g2o
sparse_block_matrix_ccs.h
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 #ifndef G2O_SPARSE_BLOCK_MATRIX_CCS_H
28 #define G2O_SPARSE_BLOCK_MATRIX_CCS_H
29 
30 #include <vector>
31 #include <cassert>
32 #include <Eigen/Core>
33 
34 #include "g2o/config.h"
35 #include "matrix_operations.h"
36 
37 #include <unordered_map>
38 
39 namespace g2o {
40 
48  template <class MatrixType>
50  {
51  public:
53  typedef MatrixType SparseMatrixBlock;
54 
56  int cols() const {return _colBlockIndices.size() ? _colBlockIndices.back() : 0;}
58  int rows() const {return _rowBlockIndices.size() ? _rowBlockIndices.back() : 0;}
59 
63  struct RowBlock
64  {
65  int row;
66  MatrixType* block;
67  RowBlock() : row(-1), block(0) {}
68  RowBlock(int r, MatrixType* b) : row(r), block(b) {}
69  bool operator<(const RowBlock& other) const { return row < other.row;}
70  };
71  typedef std::vector<RowBlock> SparseColumn;
72 
73  SparseBlockMatrixCCS(const std::vector<int>& rowIndices, const std::vector<int>& colIndices) :
74  _rowBlockIndices(rowIndices), _colBlockIndices(colIndices)
75  {}
76 
78  int rowsOfBlock(int r) const { return r ? _rowBlockIndices[r] - _rowBlockIndices[r-1] : _rowBlockIndices[0] ; }
79 
81  int colsOfBlock(int c) const { return c ? _colBlockIndices[c] - _colBlockIndices[c-1] : _colBlockIndices[0]; }
82 
84  int rowBaseOfBlock(int r) const { return r ? _rowBlockIndices[r-1] : 0 ; }
85 
87  int colBaseOfBlock(int c) const { return c ? _colBlockIndices[c-1] : 0 ; }
88 
90  const std::vector<SparseColumn>& blockCols() const { return _blockCols;}
91  std::vector<SparseColumn>& blockCols() { return _blockCols;}
92 
94  const std::vector<int>& rowBlockIndices() const { return _rowBlockIndices;}
95 
97  const std::vector<int>& colBlockIndices() const { return _colBlockIndices;}
98 
99  void rightMultiply(double*& dest, const double* src) const
100  {
101  int destSize=cols();
102 
103  if (! dest){
104  dest=new double [ destSize ];
105  memset(dest,0, destSize*sizeof(double));
106  }
107 
108  // map the memory by Eigen
109  Eigen::Map<Eigen::VectorXd> destVec(dest, destSize);
110  Eigen::Map<const Eigen::VectorXd> srcVec(src, rows());
111 
112 # ifdef G2O_OPENMP
113 # pragma omp parallel for default (shared) schedule(dynamic, 10)
114 # endif
115  for (int i=0; i < static_cast<int>(_blockCols.size()); ++i){
116  int destOffset = colBaseOfBlock(i);
117  for (typename SparseColumn::const_iterator it = _blockCols[i].begin(); it!=_blockCols[i].end(); ++it) {
118  const SparseMatrixBlock* a = it->block;
119  int srcOffset = rowBaseOfBlock(it->row);
120  // destVec += *a.transpose() * srcVec (according to the sub-vector parts)
121  internal::template atxpy<SparseMatrixBlock>(*a, srcVec, srcOffset, destVec, destOffset);
122  }
123  }
124  }
125 
129  void sortColumns()
130  {
131  for (int i=0; i < static_cast<int>(_blockCols.size()); ++i){
132  std::sort(_blockCols[i].begin(), _blockCols[i].end());
133  }
134  }
135 
139  int fillCCS(int* Cp, int* Ci, double* Cx, bool upperTriangle = false) const
140  {
141  assert(Cp && Ci && Cx && "Target destination is NULL");
142  int nz=0;
143  for (size_t i=0; i<_blockCols.size(); ++i){
144  int cstart=i ? _colBlockIndices[i-1] : 0;
145  int csize=colsOfBlock(i);
146  for (int c=0; c<csize; ++c) {
147  *Cp=nz;
148  for (typename SparseColumn::const_iterator it = _blockCols[i].begin(); it!=_blockCols[i].end(); ++it) {
149  const SparseMatrixBlock* b=it->block;
150  int rstart=it->row ? _rowBlockIndices[it->row-1] : 0;
151 
152  int elemsToCopy = b->rows();
153  if (upperTriangle && rstart == cstart)
154  elemsToCopy = c + 1;
155  for (int r=0; r<elemsToCopy; ++r){
156  *Cx++ = (*b)(r,c);
157  *Ci++ = rstart++;
158  ++nz;
159  }
160  }
161  ++Cp;
162  }
163  }
164  *Cp=nz;
165  return nz;
166  }
167 
172  int fillCCS(double* Cx, bool upperTriangle = false) const
173  {
174  assert(Cx && "Target destination is NULL");
175  double* CxStart = Cx;
176  int cstart = 0;
177  for (size_t i=0; i<_blockCols.size(); ++i){
178  int csize = _colBlockIndices[i] - cstart;
179  for (int c=0; c<csize; ++c) {
180  for (typename SparseColumn::const_iterator it = _blockCols[i].begin(); it!=_blockCols[i].end(); ++it) {
181  const SparseMatrixBlock* b = it->block;
182  int rstart = it->row ? _rowBlockIndices[it->row-1] : 0;
183 
184  int elemsToCopy = b->rows();
185  if (upperTriangle && rstart == cstart)
186  elemsToCopy = c + 1;
187  memcpy(Cx, b->data() + c*b->rows(), elemsToCopy * sizeof(double));
188  Cx += elemsToCopy;
189 
190  }
191  }
192  cstart = _colBlockIndices[i];
193  }
194  return Cx - CxStart;
195  }
196 
197  protected:
198  const std::vector<int>& _rowBlockIndices;
199  const std::vector<int>& _colBlockIndices;
200  std::vector<SparseColumn> _blockCols;
201  };
202 
203 
204 
210  template <class MatrixType>
212  {
213  public:
215  typedef MatrixType SparseMatrixBlock;
216 
218  int cols() const {return _colBlockIndices.size() ? _colBlockIndices.back() : 0;}
220  int rows() const {return _rowBlockIndices.size() ? _rowBlockIndices.back() : 0;}
221 
222  typedef std::unordered_map<int, MatrixType*> SparseColumn;
223 
224  SparseBlockMatrixHashMap(const std::vector<int>& rowIndices, const std::vector<int>& colIndices) :
225  _rowBlockIndices(rowIndices), _colBlockIndices(colIndices)
226  {}
227 
229  int rowsOfBlock(int r) const { return r ? _rowBlockIndices[r] - _rowBlockIndices[r-1] : _rowBlockIndices[0] ; }
230 
232  int colsOfBlock(int c) const { return c ? _colBlockIndices[c] - _colBlockIndices[c-1] : _colBlockIndices[0]; }
233 
235  int rowBaseOfBlock(int r) const { return r ? _rowBlockIndices[r-1] : 0 ; }
236 
238  int colBaseOfBlock(int c) const { return c ? _colBlockIndices[c-1] : 0 ; }
239 
241  const std::vector<SparseColumn>& blockCols() const { return _blockCols;}
242  std::vector<SparseColumn>& blockCols() { return _blockCols;}
243 
245  const std::vector<int>& rowBlockIndices() const { return _rowBlockIndices;}
246 
248  const std::vector<int>& colBlockIndices() const { return _colBlockIndices;}
249 
253  MatrixType* addBlock(int r, int c, bool zeroBlock = false)
254  {
255  assert(c <(int)_blockCols.size() && "accessing column which is not available");
256  SparseColumn& sparseColumn = _blockCols[c];
257  typename SparseColumn::iterator foundIt = sparseColumn.find(r);
258  if (foundIt == sparseColumn.end()) {
259  int rb = rowsOfBlock(r);
260  int cb = colsOfBlock(c);
261  MatrixType* m = new MatrixType(rb, cb);
262  if (zeroBlock)
263  m->setZero();
264  sparseColumn[r] = m;
265  return m;
266  }
267  return foundIt->second;
268  }
269 
270  protected:
271  const std::vector<int>& _rowBlockIndices;
272  const std::vector<int>& _colBlockIndices;
273  std::vector<SparseColumn> _blockCols;
274  };
275 
276 } //end namespace
277 
278 #endif
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
int fillCCS(double *Cx, bool upperTriangle=false) const
const std::vector< int > & colBlockIndices() const
indices of the column blocks
const std::vector< int > & rowBlockIndices() const
indices of the row blocks
void rightMultiply(double *&dest, const double *src) const
std::unordered_map< int, MatrixType * > SparseColumn
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
int rows() const
rows of the matrix
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
int rowsOfBlock(int r) const
how many rows does the block at block-row r has?
std::vector< SparseColumn > & blockCols()
int colBaseOfBlock(int c) const
where does the col at block-col r start?
int rowBaseOfBlock(int r) const
where does the row at block-row r start?
const std::vector< int > & _colBlockIndices
vector of the indices of the blocks along the cols
std::vector< RowBlock > SparseColumn
SparseBlockMatrixHashMap(const std::vector< int > &rowIndices, const std::vector< int > &colIndices)
std::vector< SparseColumn > & blockCols()
int fillCCS(int *Cp, int *Ci, double *Cx, bool upperTriangle=false) const
int cols() const
columns of the matrix
MatrixType * addBlock(int r, int c, bool zeroBlock=false)
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
Sparse matrix which uses blocks.
int colBaseOfBlock(int c) const
where does the col at block-col r start?
std::vector< SparseColumn > _blockCols
the matrices stored in CCS order
int colsOfBlock(int c) const
how many cols does the block at block-col c has?
int cols() const
columns of the matrix
const std::vector< int > & _rowBlockIndices
vector of the indices of the blocks along the rows.
SparseBlockMatrixCCS(const std::vector< int > &rowIndices, const std::vector< int > &colIndices)
MatrixType * block
matrix pointer for the block
std::vector< SparseColumn > _blockCols
the matrices stored in CCS order
const std::vector< int > & rowBlockIndices() const
indices of the row blocks
const std::vector< int > & _rowBlockIndices
vector of the indices of the blocks along the rows.
const std::vector< int > & _colBlockIndices
vector of the indices of the blocks along the cols
int rowsOfBlock(int r) const
how many rows does the block at block-row r has?
int rows() const
rows of the matrix
Sparse matrix which uses blocks based on hash structures.
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
bool operator<(const RowBlock &other) const
int rowBaseOfBlock(int r) const
where does the row at block-row r start?
const std::vector< int > & colBlockIndices() const
indices of the column blocks