SymmetricBlockMatrix.h
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3 * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4 * Atlanta, Georgia 30332-0415
5 * All Rights Reserved
6 * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8 * See LICENSE for the license information
9 
10 * -------------------------------------------------------------------------- */
11 
18 #pragma once
19 
20 #include <gtsam/base/FastVector.h>
21 #include <gtsam/base/Matrix.h>
22 #include <gtsam/base/types.h>
23 #include <gtsam/dllexport.h>
24 #ifdef GTSAM_ENABLE_BOOST_SERIALIZATION
25 #include <boost/serialization/nvp.hpp>
26 #endif
27 #include <cassert>
28 #include <stdexcept>
29 #include <array>
30 
31 namespace boost {
32 namespace serialization {
33 class access;
34 } /* namespace serialization */
35 } /* namespace boost */
36 
37 namespace gtsam {
38 
39  // Forward declarations
40  class VerticalBlockMatrix;
41 
53  class GTSAM_EXPORT SymmetricBlockMatrix
54  {
55  public:
59 
60  protected:
63 
65 
66  public:
69 
71  template<typename CONTAINER>
72  SymmetricBlockMatrix(const CONTAINER& dimensions, bool appendOneDimension = false) :
73  blockStart_(0)
74  {
75  fillOffsets(dimensions.begin(), dimensions.end(), appendOneDimension);
76  matrix_.resize(variableColOffsets_.back(), variableColOffsets_.back());
77  assertInvariants();
78  }
79 
81  template<typename ITERATOR>
82  SymmetricBlockMatrix(ITERATOR firstBlockDim, ITERATOR lastBlockDim, bool appendOneDimension = false) :
83  blockStart_(0)
84  {
85  fillOffsets(firstBlockDim, lastBlockDim, appendOneDimension);
86  matrix_.resize(variableColOffsets_.back(), variableColOffsets_.back());
87  assertInvariants();
88  }
89 
91  template<typename CONTAINER>
92  SymmetricBlockMatrix(const CONTAINER& dimensions, const Matrix& matrix, bool appendOneDimension = false) :
93  blockStart_(0)
94  {
95  matrix_.resize(matrix.rows(), matrix.cols());
96  matrix_.triangularView<Eigen::Upper>() = matrix.triangularView<Eigen::Upper>();
97  fillOffsets(dimensions.begin(), dimensions.end(), appendOneDimension);
98  if(matrix_.rows() != matrix_.cols())
99  throw std::invalid_argument("Requested to create a SymmetricBlockMatrix from a non-square matrix.");
100  if(variableColOffsets_.back() != matrix_.cols())
101  throw std::invalid_argument("Requested to create a SymmetricBlockMatrix with dimensions that do not sum to the total size of the provided matrix.");
102  assertInvariants();
103  }
104 
108  static SymmetricBlockMatrix LikeActiveViewOf(const SymmetricBlockMatrix& other);
109 
113  static SymmetricBlockMatrix LikeActiveViewOf(const VerticalBlockMatrix& other);
114 
116  DenseIndex rows() const { assertInvariants(); return variableColOffsets_.back() - variableColOffsets_[blockStart_]; }
117 
119  DenseIndex cols() const { return rows(); }
120 
122  DenseIndex nBlocks() const { return nActualBlocks() - blockStart_; }
123 
126  return calcIndices(block, block, 1, 1)[2];
127  }
128 
131 
135 
138  return block_(J, J).selfadjointView<Eigen::Upper>();
139  }
140 
143  return block_(J, J).selfadjointView<Eigen::Upper>();
144  }
145 
148  return block_(J, J).diagonal();
149  }
150 
153  assert(I < J);
154  return block_(I, J);
155  }
156 
159  DenseIndex I, DenseIndex J) const {
160  assert(J > I);
161  return block_(I, I, J - I, J - I).selfadjointView<Eigen::Upper>();
162  }
163 
166  DenseIndex J) const {
167  assert(J > I);
168  return block_(I, I, J - I, J - I).triangularView<Eigen::Upper>();
169  }
170 
173  DenseIndex i_endBlock,
174  DenseIndex j_startBlock,
175  DenseIndex j_endBlock) const {
176  assert(i_startBlock < j_startBlock);
177  assert(i_endBlock <= j_startBlock);
178  return block_(i_startBlock, j_startBlock, i_endBlock - i_startBlock,
179  j_endBlock - j_startBlock);
180  }
181 
183  Block aboveDiagonalRange(DenseIndex i_startBlock, DenseIndex i_endBlock,
184  DenseIndex j_startBlock, DenseIndex j_endBlock) {
185  assert(i_startBlock < j_startBlock);
186  assert(i_endBlock <= j_startBlock);
187  return block_(i_startBlock, j_startBlock, i_endBlock - i_startBlock,
188  j_endBlock - j_startBlock);
189  }
190 
194 
196  template <typename XprType>
197  void setDiagonalBlock(DenseIndex I, const XprType& xpr) {
198  block_(I, I).triangularView<Eigen::Upper>() = xpr.template triangularView<Eigen::Upper>();
199  }
200 
202  template <typename XprType>
204  assert(I != J);
205  if (I < J) {
206  block_(I, J) = xpr;
207  } else {
208  block_(J, I) = xpr.transpose();
209  }
210  }
211 
213  template <typename XprType>
215  // TODO(gareth): Eigen won't let us add triangular or self-adjoint views
216  // here, so we do it manually.
217  auto dest = block_(I, I);
218  assert(dest.rows() == xpr.rows());
219  assert(dest.cols() == xpr.cols());
220  for (DenseIndex col = 0; col < dest.cols(); ++col) {
221  for (DenseIndex row = 0; row <= col; ++row) {
222  dest(row, col) += xpr(row, col);
223  }
224  }
225  }
226 
229  template <typename XprType>
231  assert(I != J);
232  if (I < J) {
233  block_(I, J).noalias() += xpr;
234  } else {
235  block_(J, I).noalias() += xpr.transpose();
236  }
237  }
238 
242 
245  return full().selfadjointView<Eigen::Upper>();
246  }
247 
250  return full().selfadjointView<Eigen::Upper>();
251  }
252 
254  template <typename XprType>
255  void setFullMatrix(const XprType& xpr) {
256  full().triangularView<Eigen::Upper>() = xpr.template triangularView<Eigen::Upper>();
257  }
258 
260  void setZero() {
261  full().triangularView<Eigen::Upper>().setZero();
262  }
263 
265  void negate();
266 
268  void invertInPlace();
269 
271 
275  DenseIndex& blockStart() { return blockStart_; }
276 
279  DenseIndex blockStart() const { return blockStart_; }
280 
291  void choleskyPartial(DenseIndex nFrontals);
292 
299 
300  protected:
301 
304  return variableColOffsets_.size();
305  }
306 
309  return nOffsets() - 1;
310  }
311 
314  assert(block >= 0);
315  const DenseIndex actual_index = block + blockStart();
316  assert(actual_index < nOffsets());
317  return variableColOffsets_[actual_index];
318  }
319 
322  DenseIndex blockRows = 1, DenseIndex blockCols = 1) const {
323  const std::array<DenseIndex, 4> indices =
324  calcIndices(iBlock, jBlock, blockRows, blockCols);
325  return matrix_.block(indices[0], indices[1], indices[2], indices[3]);
326  }
327 
329  Block block_(DenseIndex iBlock, DenseIndex jBlock, DenseIndex blockRows = 1,
330  DenseIndex blockCols = 1) {
331  const std::array<DenseIndex, 4> indices =
332  calcIndices(iBlock, jBlock, blockRows, blockCols);
333  return matrix_.block(indices[0], indices[1], indices[2], indices[3]);
334  }
335 
337  constBlock full() const {
338  return block_(0, 0, nBlocks(), nBlocks());
339  }
340 
343  return block_(0, 0, nBlocks(), nBlocks());
344  }
345 
347  std::array<DenseIndex, 4> calcIndices(DenseIndex iBlock, DenseIndex jBlock,
348  DenseIndex blockRows,
349  DenseIndex blockCols) const {
350  assert(blockRows >= 0);
351  assert(blockCols >= 0);
352 
353  // adjust indices to account for start and size of blocks
354  const DenseIndex denseI = offset(iBlock);
355  const DenseIndex denseJ = offset(jBlock);
356  const DenseIndex denseRows = offset(iBlock + blockRows) - denseI;
357  const DenseIndex denseCols = offset(jBlock + blockCols) - denseJ;
358  return {{denseI, denseJ, denseRows, denseCols}};
359  }
360 
361  void assertInvariants() const
362  {
363  assert(matrix_.rows() == matrix_.cols());
364  assert(matrix_.cols() == variableColOffsets_.back());
365  assert(blockStart_ < (DenseIndex)variableColOffsets_.size());
366  }
367 
368  template<typename ITERATOR>
369  void fillOffsets(ITERATOR firstBlockDim, ITERATOR lastBlockDim, bool appendOneDimension)
370  {
371  variableColOffsets_.resize((lastBlockDim-firstBlockDim) + 1 + (appendOneDimension ? 1 : 0));
372  variableColOffsets_[0] = 0;
373  DenseIndex j=0;
374  for(ITERATOR dim=firstBlockDim; dim!=lastBlockDim; ++dim) {
375  variableColOffsets_[j+1] = variableColOffsets_[j] + *dim;
376  ++ j;
377  }
378  if(appendOneDimension)
379  {
380  variableColOffsets_[j+1] = variableColOffsets_[j] + 1;
381  ++ j;
382  }
383  }
384 
385  friend class VerticalBlockMatrix;
386  template<typename SymmetricBlockMatrixType> friend class SymmetricBlockMatrixBlockExpr;
387 
388  private:
389 #ifdef GTSAM_ENABLE_BOOST_SERIALIZATION
390 
391  friend class boost::serialization::access;
392  template<class ARCHIVE>
393  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
394  // Fill in the lower triangle part of the matrix, so boost::serialization won't
395  // complain about uninitialized data with an input_stream_error exception
396  // http://www.boost.org/doc/libs/1_37_0/libs/serialization/doc/exceptions.html#stream_error
397  matrix_.triangularView<Eigen::Lower>() = matrix_.triangularView<Eigen::Upper>().transpose();
398  ar & BOOST_SERIALIZATION_NVP(matrix_);
399  ar & BOOST_SERIALIZATION_NVP(variableColOffsets_);
400  ar & BOOST_SERIALIZATION_NVP(blockStart_);
401  }
402 #endif
403  };
404 
406  class CholeskyFailed;
407 
408 }
dimensions
const std::vector< size_t > dimensions
Definition: testVerticalBlockMatrix.cpp:27
gtsam::SymmetricBlockMatrix::Block
Eigen::Block< Matrix > Block
Definition: SymmetricBlockMatrix.h:57
col
m col(1)
gtsam::SymmetricBlockMatrix::block_
Block block_(DenseIndex iBlock, DenseIndex jBlock, DenseIndex blockRows=1, DenseIndex blockCols=1)
Get an arbitrary block from the matrix. Indices are in block units.
Definition: SymmetricBlockMatrix.h:329
Eigen::Block
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
gtsam::SymmetricBlockMatrix::assertInvariants
void assertInvariants() const
Definition: SymmetricBlockMatrix.h:361
gtsam::SymmetricBlockMatrix::block_
constBlock block_(DenseIndex iBlock, DenseIndex jBlock, DenseIndex blockRows=1, DenseIndex blockCols=1) const
Get an arbitrary block from the matrix. Indices are in block units.
Definition: SymmetricBlockMatrix.h:321
gtsam::SymmetricBlockMatrix::setOffDiagonalBlock
void setOffDiagonalBlock(DenseIndex I, DenseIndex J, const XprType &xpr)
Set an off-diagonal block. Only the upper triangular portion of xpr is evaluated.
Definition: SymmetricBlockMatrix.h:203
FastVector.h
A thin wrapper around std::vector that uses a custom allocator.
Eigen::CwiseBinaryOp
Generic expression where a coefficient-wise binary operator is applied to two expressions.
Definition: CwiseBinaryOp.h:77
gtsam::SymmetricBlockMatrix::full
constBlock full() const
Get the full matrix as a block.
Definition: SymmetricBlockMatrix.h:337
gtsam::SymmetricBlockMatrix::blockStart_
DenseIndex blockStart_
Changes apparent matrix view, see main class comment.
Definition: SymmetricBlockMatrix.h:64
types.h
Typedefs for easier changing of types.
gtsam::SymmetricBlockMatrix::nActualBlocks
DenseIndex nActualBlocks() const
Number of actual blocks in the full matrix.
Definition: SymmetricBlockMatrix.h:308
gtsam::FastVector
std::vector< T, typename internal::FastDefaultVectorAllocator< T >::type > FastVector
Definition: FastVector.h:34
Matrix.h
typedef and functions to augment Eigen's MatrixXd
Eigen::Upper
@ Upper
Definition: Constants.h:211
gtsam::SymmetricBlockMatrix::selfadjointView
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j).
Definition: SymmetricBlockMatrix.h:158
gtsam::SymmetricBlockMatrix::full
Block full()
Get the full matrix as a block.
Definition: SymmetricBlockMatrix.h:342
gtsam::SymmetricBlockMatrix::blockStart
DenseIndex blockStart() const
Definition: SymmetricBlockMatrix.h:279
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
Eigen::SelfAdjointView
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition: SelfAdjointView.h:49
gtsam::SymmetricBlockMatrix::setDiagonalBlock
void setDiagonalBlock(DenseIndex I, const XprType &xpr)
Set a diagonal block. Only the upper triangular portion of xpr is evaluated.
Definition: SymmetricBlockMatrix.h:197
boost
Definition: boostmultiprec.cpp:109
gtsam::SymmetricBlockMatrix::diagonalBlock
Eigen::SelfAdjointView< constBlock, Eigen::Upper > diagonalBlock(DenseIndex J) const
Return the J'th diagonal block as a self adjoint view.
Definition: SymmetricBlockMatrix.h:142
gtsam::Vector
Eigen::VectorXd Vector
Definition: Vector.h:38
gtsam::SymmetricBlockMatrix::aboveDiagonalRange
constBlock aboveDiagonalRange(DenseIndex i_startBlock, DenseIndex i_endBlock, DenseIndex j_startBlock, DenseIndex j_endBlock) const
Get a range [i,j) from the matrix. Indices are in block units.
Definition: SymmetricBlockMatrix.h:172
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
block
m m block(1, 0, 2, 2)<< 4
test_eigen_tensor.indices
indices
Definition: test_eigen_tensor.py:33
gtsam::SymmetricBlockMatrix::SymmetricBlockMatrix
SymmetricBlockMatrix(const CONTAINER &dimensions, const Matrix &matrix, bool appendOneDimension=false)
Construct from a container of the sizes of each vertical block and a pre-prepared matrix.
Definition: SymmetricBlockMatrix.h:92
gtsam::SymmetricBlockMatrix::aboveDiagonalRange
Block aboveDiagonalRange(DenseIndex i_startBlock, DenseIndex i_endBlock, DenseIndex j_startBlock, DenseIndex j_endBlock)
Get a range [i,j) from the matrix. Indices are in block units.
Definition: SymmetricBlockMatrix.h:183
gtsam::SymmetricBlockMatrix::matrix_
Matrix matrix_
The full matrix.
Definition: SymmetricBlockMatrix.h:61
J
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
gtsam::VerticalBlockMatrix
Definition: VerticalBlockMatrix.h:42
Eigen::internal::negate
T negate(const T &x)
Definition: packetmath_test_shared.h:24
j
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
I
#define I
Definition: main.h:112
gtsam::row
const MATRIX::ConstRowXpr row(const MATRIX &A, size_t j)
Definition: base/Matrix.h:221
gtsam::SymmetricBlockMatrix::setFullMatrix
void setFullMatrix(const XprType &xpr)
Set the entire active matrix. Only reads the upper triangular part of xpr.
Definition: SymmetricBlockMatrix.h:255
gtsam::SymmetricBlockMatrix::SymmetricBlockMatrix
SymmetricBlockMatrix(ITERATOR firstBlockDim, ITERATOR lastBlockDim, bool appendOneDimension=false)
Construct from iterator over the sizes of each vertical block.
Definition: SymmetricBlockMatrix.h:82
gtsam::SymmetricBlockMatrix::nBlocks
DenseIndex nBlocks() const
Block count.
Definition: SymmetricBlockMatrix.h:122
Eigen::Lower
@ Lower
Definition: Constants.h:209
gtsam::SymmetricBlockMatrix::rows
DenseIndex rows() const
Row size.
Definition: SymmetricBlockMatrix.h:116
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
offset
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate offset
Definition: gnuplot_common_settings.hh:64
gtsam::SymmetricBlockMatrix::This
SymmetricBlockMatrix This
Definition: SymmetricBlockMatrix.h:56
gtsam::SymmetricBlockMatrix::nOffsets
DenseIndex nOffsets() const
Number of offsets in the full matrix.
Definition: SymmetricBlockMatrix.h:303
gtsam::SymmetricBlockMatrix::triangularView
Eigen::TriangularView< constBlock, Eigen::Upper > triangularView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j) as a triangular view.
Definition: SymmetricBlockMatrix.h:165
gtsam::SymmetricBlockMatrix::calcIndices
std::array< DenseIndex, 4 > calcIndices(DenseIndex iBlock, DenseIndex jBlock, DenseIndex blockRows, DenseIndex blockCols) const
Compute the indices into the underlying matrix for a given block.
Definition: SymmetricBlockMatrix.h:347
gtsam::SymmetricBlockMatrix::updateDiagonalBlock
void updateDiagonalBlock(DenseIndex I, const XprType &xpr)
Increment the diagonal block by the values in xpr. Only reads the upper triangular part of xpr.
Definition: SymmetricBlockMatrix.h:214
gtsam::SymmetricBlockMatrix::selfadjointView
Eigen::SelfAdjointView< Block, Eigen::Upper > selfadjointView()
Get self adjoint view.
Definition: SymmetricBlockMatrix.h:244
gtsam::SymmetricBlockMatrix::cols
DenseIndex cols() const
Column size.
Definition: SymmetricBlockMatrix.h:119
Eigen::CwiseBinaryOp::cols
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: CwiseBinaryOp.h:125
gtsam
traits
Definition: SFMdata.h:40
gtsam::SymmetricBlockMatrix::setZero
void setZero()
Set the entire active matrix zero.
Definition: SymmetricBlockMatrix.h:260
gtsam::choleskyPartial
bool choleskyPartial(Matrix &ABC, size_t nFrontal, size_t topleft)
Definition: base/cholesky.cpp:107
gtsam::DenseIndex
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:103
gtsam::SymmetricBlockMatrix::updateOffDiagonalBlock
void updateOffDiagonalBlock(DenseIndex I, DenseIndex J, const XprType &xpr)
Definition: SymmetricBlockMatrix.h:230
gtsam::SymmetricBlockMatrix::diagonal
Vector diagonal(DenseIndex J) const
Get the diagonal of the J'th diagonal block.
Definition: SymmetricBlockMatrix.h:147
gtsam::SymmetricBlockMatrix::diagonalBlock
Eigen::SelfAdjointView< Block, Eigen::Upper > diagonalBlock(DenseIndex J)
Return the J'th diagonal block as a self adjoint view.
Definition: SymmetricBlockMatrix.h:137
gtsam::split
void split(const G &g, const PredecessorMap< KEY > &tree, G &Ab1, G &Ab2)
Definition: graph-inl.h:245
gtsam::SymmetricBlockMatrix::getDim
DenseIndex getDim(DenseIndex block) const
Number of dimensions for variable on this diagonal block.
Definition: SymmetricBlockMatrix.h:125
setZero
v setZero(3)
gtsam::SymmetricBlockMatrix::blockStart
DenseIndex & blockStart()
Definition: SymmetricBlockMatrix.h:275
gtsam::SymmetricBlockMatrix::offset
DenseIndex offset(DenseIndex block) const
Get an offset for a block index (in the active view).
Definition: SymmetricBlockMatrix.h:313
gtsam::SymmetricBlockMatrix::SymmetricBlockMatrix
SymmetricBlockMatrix(const CONTAINER &dimensions, bool appendOneDimension=false)
Construct from a container of the sizes of each block.
Definition: SymmetricBlockMatrix.h:72
gtsam::SymmetricBlockMatrix::aboveDiagonalBlock
constBlock aboveDiagonalBlock(DenseIndex I, DenseIndex J) const
Get block above the diagonal (I, J).
Definition: SymmetricBlockMatrix.h:152
gtsam::SymmetricBlockMatrix::constBlock
Eigen::Block< const Matrix > constBlock
Definition: SymmetricBlockMatrix.h:58
gtsam::SymmetricBlockMatrix::selfadjointView
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView() const
Get self adjoint view.
Definition: SymmetricBlockMatrix.h:249
Eigen::TriangularView
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:187
gtsam::SymmetricBlockMatrix
Definition: SymmetricBlockMatrix.h:53
pybind_wrapper_test_script.other
other
Definition: pybind_wrapper_test_script.py:42
gtsam::SymmetricBlockMatrix::fillOffsets
void fillOffsets(ITERATOR firstBlockDim, ITERATOR lastBlockDim, bool appendOneDimension)
Definition: SymmetricBlockMatrix.h:369
gtsam::SymmetricBlockMatrix::variableColOffsets_
FastVector< DenseIndex > variableColOffsets_
the starting columns of each block (0-based)
Definition: SymmetricBlockMatrix.h:62
Eigen::CwiseBinaryOp::rows
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: CwiseBinaryOp.h:120


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:05:22