JacobianFactor-inl.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 
19 #pragma once
20 
22 
23 #if defined(__GNUC__) && !defined(__clang__) && __GNUC__ >= 13
24 #pragma GCC diagnostic warning "-Wstringop-overread"
25 #endif
26 
27 namespace gtsam {
28 
29  /* ************************************************************************* */
30  template <typename TERMS>
31  JacobianFactor::JacobianFactor(const TERMS& terms, const Vector& b,
32  const SharedDiagonal& model) {
33  fillTerms(terms, b, model);
34  }
35 
36  /* ************************************************************************* */
37  template <typename KEYS>
39  const VerticalBlockMatrix& augmentedMatrix,
40  const SharedDiagonal& model)
41  : Base(keys), Ab_(augmentedMatrix), model_(model) {
42  checkAb(model, augmentedMatrix);
43  }
44 
45  /* ************************************************************************* */
46  template <typename KEYS>
48  VerticalBlockMatrix&& augmentedMatrix,
49  const SharedDiagonal& model)
50  : Base(keys), Ab_(std::move(augmentedMatrix)), model_(model) {
51  checkAb(model, Ab_);
52  }
53 
54  /* ************************************************************************* */
55  template<typename TERMS>
56  void JacobianFactor::fillTerms(const TERMS& terms, const Vector& b, const SharedDiagonal& noiseModel)
57  {
58  // Check noise model dimension
59  if(noiseModel && (DenseIndex)noiseModel->dim() != b.size())
60  throw InvalidNoiseModel(b.size(), noiseModel->dim());
61 
62  // Resize base class key vector
63  Base::keys_.resize(terms.size());
64 
65  // Get dimensions of matrices
66  std::vector<size_t> dimensions;
67  dimensions.reserve(terms.size());
68  for(typename TERMS::const_iterator it = terms.begin(); it != terms.end(); ++it) {
69  const std::pair<Key, Matrix>& term = *it;
70  const Matrix& Ai = term.second;
71  dimensions.push_back(Ai.cols());
72  }
73 
74  // Construct block matrix
75  Ab_ = VerticalBlockMatrix(dimensions, b.size(), true);
76 
77  // Check and add terms
78  DenseIndex i = 0; // For block index
79  for(typename TERMS::const_iterator it = terms.begin(); it != terms.end(); ++it) {
80  const std::pair<Key, Matrix>& term = *it;
81  Key key = term.first;
82  const Matrix& Ai = term.second;
83 
84  // Check block rows
85  if(Ai.rows() != Ab_.rows())
86  throw InvalidMatrixBlock(Ab_.rows(), Ai.rows());
87 
88  // Assign key and matrix
89  Base::keys_[i] = key;
90  Ab_(i) = Ai;
91 
92  // Increment block index
93  ++ i;
94  }
95 
96  // Assign RHS vector
97  getb() = b;
98 
99  // Assign noise model
100  model_ = noiseModel;
101  }
102 
103 } // gtsam
104 
gtsam::JacobianFactor::checkAb
void checkAb(const SharedDiagonal &model, const VerticalBlockMatrix &augmentedMatrix) const
Common code between VerticalBlockMatrix constructors.
Definition: JacobianFactor.cpp:116
dimensions
const std::vector< size_t > dimensions
Definition: testVerticalBlockMatrix.cpp:27
gtsam::JacobianFactor::Ab_
VerticalBlockMatrix Ab_
Definition: JacobianFactor.h:106
keys
const KeyVector keys
Definition: testRegularImplicitSchurFactor.cpp:40
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
gtsam::InvalidNoiseModel
Definition: linearExceptions.h:106
gtsam::Factor
Definition: Factor.h:70
gtsam::Vector
Eigen::VectorXd Vector
Definition: Vector.h:39
gtsam::JacobianFactor::getb
const constBVector getb() const
Definition: JacobianFactor.h:304
gtsam::JacobianFactor::model_
noiseModel::Diagonal::shared_ptr model_
Definition: JacobianFactor.h:107
linearExceptions.h
Exceptions that may be thrown by linear solver components.
gtsam::JacobianFactor::fillTerms
void fillTerms(const TERMS &terms, const Vector &b, const SharedDiagonal &noiseModel)
Internal function to fill blocks and set dimensions.
Definition: JacobianFactor-inl.h:56
gtsam::VerticalBlockMatrix
Definition: VerticalBlockMatrix.h:44
gtsam::SharedDiagonal
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:764
model
noiseModel::Diagonal::shared_ptr model
Definition: doc/Code/Pose2SLAMExample.cpp:7
key
const gtsam::Symbol key('X', 0)
gtsam::b
const G & b
Definition: Group.h:79
move
detail::enable_if_t<!detail::move_never< T >::value, T > move(object &&obj)
Definition: cast.h:1250
gtsam
traits
Definition: SFMdata.h:40
gtsam::Factor::keys_
KeyVector keys_
The keys involved in this factor.
Definition: Factor.h:88
gtsam::DenseIndex
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:103
std
Definition: BFloat16.h:88
gtsam::VerticalBlockMatrix::rows
DenseIndex rows() const
Row size.
Definition: VerticalBlockMatrix.h:150
gtsam::JacobianFactor::JacobianFactor
JacobianFactor()
Definition: JacobianFactor.cpp:48
gtsam::InvalidMatrixBlock
Definition: linearExceptions.h:121
gtsam::Key
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:97
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9


gtsam
Author(s):
autogenerated on Wed Mar 19 2025 03:01:56