SymmetricBlockMatrix.cpp
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 
21 #include <gtsam/base/cholesky.h>
22 #include <gtsam/base/timing.h>
24 
25 namespace gtsam {
26 
27 /* ************************************************************************* */
29  variableColOffsets_.push_back(0);
31 }
32 
33 /* ************************************************************************* */
35  const SymmetricBlockMatrix& other) {
37  result.variableColOffsets_.resize(other.nBlocks() + 1);
38  for (size_t i = 0; i < result.variableColOffsets_.size(); ++i)
39  result.variableColOffsets_[i] = other.variableColOffsets_[other.blockStart_
40  + i] - other.variableColOffsets_[other.blockStart_];
41  result.matrix_.resize(other.cols(), other.cols());
42  result.assertInvariants();
43  return result;
44 }
45 
46 /* ************************************************************************* */
48  const VerticalBlockMatrix& other) {
50  result.variableColOffsets_.resize(other.nBlocks() + 1);
51  for (size_t i = 0; i < result.variableColOffsets_.size(); ++i)
52  result.variableColOffsets_[i] = other.variableColOffsets_[other.blockStart_
53  + i] - other.variableColOffsets_[other.blockStart_];
54  result.matrix_.resize(other.cols(), other.cols());
55  result.assertInvariants();
56  return result;
57 }
58 
59 /* ************************************************************************* */
61  if (I == J) {
62  return diagonalBlock(I);
63  } else if (I < J) {
64  return aboveDiagonalBlock(I, J);
65  } else {
66  return aboveDiagonalBlock(J, I).transpose();
67  }
68 }
69 
70 /* ************************************************************************* */
72  full().triangularView<Eigen::Upper>() *= -1.0;
73 }
74 
75 /* ************************************************************************* */
77  const auto identity = Matrix::Identity(rows(), rows());
78  full().triangularView<Eigen::Upper>() =
79  selfadjointView().llt().solve(identity).triangularView<Eigen::Upper>();
80 }
81 
82 /* ************************************************************************* */
84  gttic(VerticalBlockMatrix_choleskyPartial);
86  if (!gtsam::choleskyPartial(matrix_, offset(nFrontals) - topleft, topleft)) {
87  throw CholeskyFailed();
88  }
89 }
90 
91 /* ************************************************************************* */
93  gttic(VerticalBlockMatrix_split);
94 
95  // Construct a VerticalBlockMatrix that contains [R Sd]
96  const size_t n1 = offset(nFrontals);
98 
99  // Copy into it.
100  RSd.full() = matrix_.topRows(n1);
101  RSd.full().triangularView<Eigen::StrictlyLower>().setZero();
102 
103  // Take lower-right block of Ab_ to get the remaining factor
104  blockStart() = nFrontals;
105 
106  return RSd;
107 }
108 
109 /* ************************************************************************* */
110 
111 } //\ namespace gtsam
112 
timing.h
Timing utilities.
gtsam::VerticalBlockMatrix::LikeActiveViewOf
static VerticalBlockMatrix LikeActiveViewOf(const VerticalBlockMatrix &rhs)
Definition: VerticalBlockMatrix.cpp:25
gtsam::SymmetricBlockMatrix::assertInvariants
void assertInvariants() const
Definition: SymmetricBlockMatrix.h:361
benchmark.n1
n1
Definition: benchmark.py: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
Eigen::Upper
@ Upper
Definition: Constants.h:211
cholesky.h
Efficient incomplete Cholesky on rank-deficient matrices, todo: constrained Cholesky.
gtsam::VerticalBlockMatrix::full
Block full()
Definition: VerticalBlockMatrix.h:157
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
result
Values result
Definition: OdometryOptimize.cpp:8
gtsam::SymmetricBlockMatrix::invertInPlace
void invertInPlace()
Invert the entire active matrix in place.
Definition: SymmetricBlockMatrix.cpp:76
gtsam::CholeskyFailed
Indicate Cholesky factorization failure.
Definition: ThreadsafeException.h:115
gtsam::SymmetricBlockMatrix::negate
void negate()
Negate the entire active matrix.
Definition: SymmetricBlockMatrix.cpp:71
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
ThreadsafeException.h
I
#define I
Definition: main.h:112
VerticalBlockMatrix.h
A matrix with column blocks of pre-defined sizes. Used in JacobianFactor and GaussianConditional.
SymmetricBlockMatrix.h
Access to matrices via blocks of pre-defined sizes. Used in GaussianFactor and GaussianConditional.
gtsam::SymmetricBlockMatrix::LikeActiveViewOf
static SymmetricBlockMatrix LikeActiveViewOf(const SymmetricBlockMatrix &other)
Definition: SymmetricBlockMatrix.cpp:34
Eigen::StrictlyLower
@ StrictlyLower
Definition: Constants.h:221
gtsam::SymmetricBlockMatrix::rows
DenseIndex rows() const
Row size.
Definition: SymmetricBlockMatrix.h:116
gtsam::SymmetricBlockMatrix::selfadjointView
Eigen::SelfAdjointView< Block, Eigen::Upper > selfadjointView()
Get self adjoint view.
Definition: SymmetricBlockMatrix.h:244
gtsam
traits
Definition: chartTesting.h:28
gtsam::SymmetricBlockMatrix::setZero
void setZero()
Set the entire active matrix zero.
Definition: SymmetricBlockMatrix.h:260
gtsam::SymmetricBlockMatrix::split
VerticalBlockMatrix split(DenseIndex nFrontals)
Definition: SymmetricBlockMatrix.cpp:92
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::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::SymmetricBlockMatrix::block
Matrix block(DenseIndex I, DenseIndex J) const
Definition: SymmetricBlockMatrix.cpp:60
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::aboveDiagonalBlock
constBlock aboveDiagonalBlock(DenseIndex I, DenseIndex J) const
Get block above the diagonal (I, J).
Definition: SymmetricBlockMatrix.h:152
gtsam::SymmetricBlockMatrix::choleskyPartial
void choleskyPartial(DenseIndex nFrontals)
Definition: SymmetricBlockMatrix.cpp:83
gtsam::SymmetricBlockMatrix
Definition: SymmetricBlockMatrix.h:53
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
pybind_wrapper_test_script.other
other
Definition: pybind_wrapper_test_script.py:42
gtsam::SymmetricBlockMatrix::SymmetricBlockMatrix
SymmetricBlockMatrix()
Construct from an empty matrix (asserts that the matrix is empty)
Definition: SymmetricBlockMatrix.cpp:28
gttic
#define gttic(label)
Definition: timing.h:295
gtsam::SymmetricBlockMatrix::variableColOffsets_
FastVector< DenseIndex > variableColOffsets_
the starting columns of each block (0-based)
Definition: SymmetricBlockMatrix.h:62


gtsam
Author(s):
autogenerated on Tue Jun 25 2024 03:03:49