31template <
class MatrixType>
35 : _rowBlockIndices(rbi, rbi + rb),
36 _colBlockIndices(cbi, cbi + cb),
38 _hasStorage(hasStorage) {}
40template <
class MatrixType>
42 : _blockCols(0), _hasStorage(true) {}
44template <
class MatrixType>
47#pragma omp parallel for default(shared) if (_blockCols.size() > 100)
49 for (
int i = 0; i < static_cast<int>(_blockCols.size()); ++i) {
51 it = _blockCols[i].begin();
52 it != _blockCols[i].end(); ++it) {
54 if (_hasStorage && dealloc)
59 if (_hasStorage && dealloc) _blockCols[i].
clear();
63template <
class MatrixType>
65 if (_hasStorage) clear(
true);
68template <
class MatrixType>
72 _blockCols[c].find(r);
74 if (it == _blockCols[c].end()) {
75 if (!_hasStorage && !alloc)
78 int rb = rowsOfBlock(r);
79 int cb = colsOfBlock(c);
83 std::pair<typename SparseBlockMatrix<MatrixType>::IntBlockMap::iterator,
85 result = _blockCols[c].insert(std::make_pair(r, _block));
87 assert(result.second);
95template <
class MatrixType>
99 _blockCols[c].find(r);
100 if (it == _blockCols[c].end())
return nullptr;
104template <
class MatrixType>
108 _rowBlockIndices.size(), _colBlockIndices.size());
109 for (
size_t i = 0; i < _blockCols.size(); ++i) {
111 it = _blockCols[i].begin();
112 it != _blockCols[i].end(); ++it) {
116 ret->
_blockCols[i].insert(std::make_pair(it->first, b));
123template <
class MatrixType>
124template <
class MatrixTransposedType>
127 for (
size_t i = 0; i < _blockCols.size(); ++i) {
129 it = _blockCols[i].begin();
130 it != _blockCols[i].end(); ++it) {
133 dest.
block(i, it->first,
true);
139template <
class MatrixType>
140template <
class MatrixTransposedType>
146 for (
size_t i = 0; i < _rowBlockIndices.size(); ++i) {
149 for (
size_t i = 0; i < _colBlockIndices.size(); ++i) {
153 transpose_internal(dest);
157template <
class MatrixType>
158template <
class MatrixTransposedType>
159std::unique_ptr<SparseBlockMatrix<MatrixTransposedType>>
161 auto dest = std::make_unique<SparseBlockMatrix<MatrixTransposedType>>(
162 &_colBlockIndices[0], &_rowBlockIndices[0], _colBlockIndices.size(),
163 _rowBlockIndices.size());
164 transpose_internal(*dest);
168template <
class MatrixType>
171 for (
size_t i = 0; i < _blockCols.size(); ++i) {
173 it = _blockCols[i].begin();
174 it != _blockCols[i].end(); ++it) {
183template <
class MatrixType>
189 for (
size_t i = 0; i < _rowBlockIndices.size(); ++i) {
192 for (
size_t i = 0; i < _colBlockIndices.size(); ++i) {
200template <
class MatrixType>
201std::unique_ptr<SparseBlockMatrix<MatrixType>>
203 auto a = std::make_unique<SparseBlockMatrix>(
204 &_rowBlockIndices[0], &_colBlockIndices[0], _rowBlockIndices.size(),
205 _colBlockIndices.size());
210template <
class MatrixType>
211template <
class MatrixResultType,
class MatrixFactorType>
217 for (
size_t i = 0; i < _colBlockIndices.size(); ++i) {
226 for (
size_t i = 0; i < M->
_blockCols.size(); ++i) {
228 MatrixFactorType>::IntBlockMap::const_iterator it =
236 _blockCols[it->first].begin();
237 while (rbt != _blockCols[it->first].end()) {
239 int rowA = rbt->first;
243 dest->
block(rowA, colM,
true);
254template <
class MatrixType>
256 const double* src)
const {
258 dest =
new double[_rowBlockIndices[_rowBlockIndices.size() - 1]];
260 _rowBlockIndices[_rowBlockIndices.size() - 1] *
sizeof(
double));
264 Eigen::Map<VectorX> destVec(dest, rows());
265 const Eigen::Map<const VectorX> srcVec(src, cols());
267 for (
size_t i = 0; i < _blockCols.size(); ++i) {
268 int srcOffset = i ? _colBlockIndices[i - 1] : 0;
271 it = _blockCols[i].begin();
272 it != _blockCols[i].end(); ++it) {
275 int destOffset = it->first ? _rowBlockIndices[it->first - 1] : 0;
277 internal::template axpy<
279 *a, srcVec, srcOffset, destVec, destOffset);
284template <
class MatrixType>
286 double*& dest,
const double* src)
const {
288 dest =
new double[_rowBlockIndices[_rowBlockIndices.size() - 1]];
290 _rowBlockIndices[_rowBlockIndices.size() - 1] *
sizeof(
double));
294 Eigen::Map<VectorX> destVec(dest, rows());
295 const Eigen::Map<const VectorX> srcVec(src, cols());
297 for (
size_t i = 0; i < _blockCols.size(); ++i) {
298 int srcOffset = colBaseOfBlock(i);
300 it = _blockCols[i].begin();
301 it != _blockCols[i].end(); ++it) {
304 int destOffset = rowBaseOfBlock(it->first);
305 if (destOffset > srcOffset)
308 internal::template axpy<
310 *a, srcVec, srcOffset, destVec, destOffset);
311 if (destOffset < srcOffset)
312 internal::template atxpy<
314 *a, srcVec, destOffset, destVec, srcOffset);
319template <
class MatrixType>
321 const double* src)
const {
322 int destSize = cols();
325 dest =
new double[destSize];
326 memset(dest, 0, destSize *
sizeof(
double));
330 Eigen::Map<VectorX> destVec(dest, destSize);
331 Eigen::Map<const VectorX> srcVec(src, rows());
334#pragma omp parallel for default(shared) schedule(dynamic, 10)
336 for (
int i = 0; i < static_cast<int>(_blockCols.size()); ++i) {
337 int destOffset = colBaseOfBlock(i);
339 it = _blockCols[i].begin();
340 it != _blockCols[i].end(); ++it) {
343 int srcOffset = rowBaseOfBlock(it->first);
345 internal::template atxpy<
347 *a, srcVec, srcOffset, destVec, destOffset);
352template <
class MatrixType>
354 for (
size_t i = 0; i < _blockCols.size(); ++i) {
356 it = _blockCols[i].begin();
357 it != _blockCols[i].end(); ++it) {
364template <
class MatrixType>
366 int rmin,
int rmax,
int cmin,
int cmax,
bool alloc)
const {
370 rowIdx[0] = rowsOfBlock(rmin);
371 for (
int i = 1; i < m; ++i) {
372 rowIdx[i] = rowIdx[i - 1] + rowsOfBlock(rmin + i);
376 colIdx[0] = colsOfBlock(cmin);
377 for (
int i = 1; i < n; ++i) {
378 colIdx[i] = colIdx[i - 1] + colsOfBlock(cmin + i);
381 for (
int i = 0; i < n; ++i) {
384 it = _blockCols[mc].begin();
385 it != _blockCols[mc].end(); ++it) {
386 if (it->first >= rmin && it->first < rmax) {
392 s->
_blockCols[i].insert(std::make_pair(it->first - rmin, b));
400template <
class MatrixType>
403 for (
size_t i = 0; i < _blockCols.size(); ++i) count += _blockCols[i].size();
407template <
class MatrixType>
409 if (MatrixType::SizeAtCompileTime != Eigen::Dynamic) {
410 size_t nnz = nonZeroBlocks() * MatrixType::SizeAtCompileTime;
414 for (
size_t i = 0; i < _blockCols.size(); ++i) {
416 it = _blockCols[i].begin();
417 it != _blockCols[i].end(); ++it) {
427template <
class MatrixType>
439 for (
size_t i = 0; i < m.
blockCols().size(); ++i) {
445 os <<
"BLOCK: " << it->first <<
" " << i << std::endl;
446 os << *b << std::endl;
452template <
class MatrixType>
455 bool upperTriangle)
const {
457 size_t n = _rowBlockIndices.size();
459 std::vector<int> blockSizes(_rowBlockIndices.size());
460 blockSizes[0] = _rowBlockIndices[0];
461 for (
size_t i = 1; i < n; ++i) {
462 blockSizes[i] = _rowBlockIndices[i] - _rowBlockIndices[i - 1];
465 std::vector<int> pBlockIndices(_rowBlockIndices.size());
466 for (
size_t i = 0; i < n; ++i) {
467 pBlockIndices[pinv[i]] = blockSizes[i];
469 for (
size_t i = 1; i < n; ++i) {
470 pBlockIndices[i] += pBlockIndices[i - 1];
479 for (
size_t i = 0; i < n; ++i) {
486 for (
size_t i = 0; i < n; ++i) {
489 it = _blockCols[i].begin();
490 it != _blockCols[i].end(); ++it) {
491 int pj = pinv[it->first];
497 if (!upperTriangle || pj <= pi) {
498 b = dest->
block(pj, pi,
true);
503 b = dest->
block(pi, pj,
true);
515template <
class MatrixType>
517 bool upperTriangle)
const {
518 assert(Cx &&
"Target destination is NULL");
519 double* CxStart = Cx;
520 for (
size_t i = 0; i < _blockCols.size(); ++i) {
521 int cstart = i ? _colBlockIndices[i - 1] : 0;
522 int csize = colsOfBlock(i);
523 for (
int c = 0; c < csize; ++c) {
525 it = _blockCols[i].begin();
526 it != _blockCols[i].end(); ++it) {
529 int rstart = it->first ? _rowBlockIndices[it->first - 1] : 0;
531 int elemsToCopy = b->
rows();
532 if (upperTriangle && rstart == cstart) elemsToCopy = c + 1;
533 memcpy(Cx, b->data() + c * b->
rows(), elemsToCopy *
sizeof(
double));
541template <
class MatrixType>
543 bool upperTriangle)
const {
544 assert(Cp && Ci && Cx &&
"Target destination is NULL");
546 for (
size_t i = 0; i < _blockCols.size(); ++i) {
547 int cstart = i ? _colBlockIndices[i - 1] : 0;
548 int csize = colsOfBlock(i);
549 for (
int c = 0; c < csize; ++c) {
552 it = _blockCols[i].begin();
553 it != _blockCols[i].end(); ++it) {
556 int rstart = it->first ? _rowBlockIndices[it->first - 1] : 0;
558 int elemsToCopy = b->
rows();
559 if (upperTriangle && rstart == cstart) elemsToCopy = c + 1;
560 for (
int r = 0; r < elemsToCopy; ++r) {
573template <
class MatrixType>
576 ms.
alloc(_colBlockIndices.size(), nonZeroBlocks());
577 ms.
m = _rowBlockIndices.size();
578 fillBlockStructure(ms.
Ap, ms.
Aii);
581template <
class MatrixType>
584 for (
int c = 0; c < static_cast<int>(_blockCols.size()); ++c) {
586 for (
auto it = _blockCols[c].begin(); it != _blockCols[c].end(); ++it) {
587 const int& r = it->first;
596 assert(nz <=
static_cast<int>(nonZeroBlocks()));
599template <
class MatrixType>
601 bool upperTriangle)
const {
602 std::string name = filename;
603 std::string::size_type lastDot = name.find_last_of(
'.');
604 if (lastDot != std::string::npos) name.resize(lastDot);
606 std::vector<TripletEntry> entries;
607 for (
size_t i = 0; i < _blockCols.size(); ++i) {
610 it = _blockCols[i].begin();
611 it != _blockCols[i].end(); ++it) {
612 const int& r = it->first;
613 const MatrixType& m = *(it->second);
614 for (
int cc = 0; cc < m.cols(); ++cc)
615 for (
int rr = 0; rr < m.rows(); ++rr) {
616 int aux_r = rowBaseOfBlock(r) + rr;
617 int aux_c = colBaseOfBlock(c) + cc;
618 entries.push_back(
TripletEntry(aux_r, aux_c, m(rr, cc)));
619 if (upperTriangle && r != c) {
620 entries.push_back(
TripletEntry(aux_c, aux_r, m(rr, cc)));
626 int nz = entries.size();
629 std::ofstream fout(filename);
630 fout <<
"# name: " << name << std::endl;
631 fout <<
"# type: sparse matrix" << std::endl;
632 fout <<
"# nnz: " << nz << std::endl;
633 fout <<
"# rows: " << rows() << std::endl;
634 fout <<
"# columns: " << cols() << std::endl;
635 fout << std::setprecision(9) << std::fixed << std::endl;
637 for (std::vector<TripletEntry>::const_iterator it = entries.begin();
638 it != entries.end(); ++it) {
640 fout << entry.
r + 1 <<
" " << entry.
c + 1 <<
" " << entry.
x << std::endl;
645template <
class MatrixType>
648 blockCCS.
blockCols().resize(blockCols().size());
650 for (
size_t i = 0; i < blockCols().size(); ++i) {
655 dest.reserve(row.size());
656 for (
typename IntBlockMap::const_iterator it = row.begin(); it != row.end();
659 it->first, it->second));
666template <
class MatrixType>
670 blockCCS.
blockCols().resize(_rowBlockIndices.size());
672 for (
size_t i = 0; i < blockCols().size(); ++i) {
674 for (
typename IntBlockMap::const_iterator it = row.begin(); it != row.end();
686template <
class MatrixType>
691 typedef std::pair<int, MatrixType*> SparseColumnPair;
694 for (
size_t i = 0; i < hashMatrix.
blockCols().size(); ++i) {
696 HashSparseColumn& column = hashMatrix.
blockCols()[i];
697 if (column.size() == 0)
continue;
698 std::vector<SparseColumnPair> sparseRowSorted;
699 sparseRowSorted.reserve(column.size());
700 for (
typename HashSparseColumn::const_iterator it = column.begin();
701 it != column.end(); ++it)
702 sparseRowSorted.push_back(*it);
703 std::sort(sparseRowSorted.begin(), sparseRowSorted.end(),
706 HashSparseColumn aux;
707 std::swap(aux, column);
710 destColumnMap.insert(sparseRowSorted[0]);
711 for (
size_t j = 1; j < sparseRowSorted.size(); ++j) {
716 destColumnMap.insert(hint, sparseRowSorted[j]);
representing the structure of a matrix in column compressed structure (only the upper triangular part...
int * Aii
row indices of A, of size nz = Ap [n]
int m
A is m-by-n. m must be >= 0.
void alloc(int n_, int nz)
int * Ap
column pointers for A, of size n+1
Sparse matrix which uses blocks.
std::vector< RowBlock > SparseColumn
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
Sparse matrix which uses blocks based on hash structures.
std::unordered_map< int, MatrixType * > SparseColumn
const std::vector< SparseColumn > & blockCols() const
the block matrices per block-column
Sparse matrix which uses blocks.
size_t nonZeros() const
number of non-zero elements
void takePatternFromHash(SparseBlockMatrixHashMap< MatrixType > &hashMatrix)
bool symmPermutation(SparseBlockMatrix< MatrixType > *&dest, const int *pinv, bool onlyUpper=false) const
bool add(SparseBlockMatrix< MatrixType > &dest) const
adds the current matrix to the destination
bool transpose(SparseBlockMatrix< MatrixTransposedType > &dest) const
bool writeOctave(const char *filename, bool upperTriangle=true) const
int fillSparseBlockMatrixCCSTransposed(SparseBlockMatrixCCS< MatrixType > &blockCCS) const
size_t nonZeroBlocks() const
number of allocated blocks
int fillCCS(int *Cp, int *Ci, double *Cx, bool upperTriangle=false) const
const std::vector< IntBlockMap > & blockCols() const
the block matrices per block-column
int fillSparseBlockMatrixCCS(SparseBlockMatrixCCS< MatrixType > &blockCCS) const
std::unique_ptr< SparseBlockMatrix< MatrixType > > added() const
SparseBlockMatrix * slice(int rmin, int rmax, int cmin, int cmax, bool alloc=true) const
void rightMultiply(double *&dest, const double *src) const
dest = M * (*this)
void multiplySymmetricUpperTriangle(double *&dest, const double *src) const
const std::vector< int > & colBlockIndices() const
indices of the column blocks
std::unique_ptr< SparseBlockMatrix< MatrixTransposedType > > transposed() const
std::map< int, SparseMatrixBlock * > IntBlockMap
int rows() const
rows of the matrix
SparseMatrixBlock * block(int r, int c, bool alloc=false)
SparseBlockMatrix * clone() const
deep copy of a sparse-block-matrix;
std::vector< int > _rowBlockIndices
MatrixType SparseMatrixBlock
this is the type of the elementary block, it is an Eigen::Matrix.
bool multiply(SparseBlockMatrix< MatrixResultType > *&dest, const SparseBlockMatrix< MatrixFactorType > *M) const
dest = (*this) * M
int cols() const
columns of the matrix
std::vector< int > _colBlockIndices
void add_internal(SparseBlockMatrix< MatrixType > &dest) const
void fillBlockStructure(MatrixStructure &ms) const
exports the non zero blocks in the structure matrix ms
std::vector< IntBlockMap > _blockCols
const std::vector< int > & rowBlockIndices() const
indices of the row blocks
void clear(bool dealloc=false)
void transpose_internal(SparseBlockMatrix< MatrixTransposedType > &dest) const
void scale(double a)
*this *= a
std::ostream & operator<<(std::ostream &os, const G2OBatchStatistics &st)