26#include <boost/make_shared.hpp>
34 class GaussianConditional;
35 class GaussianBayesNet;
36 class GaussianFactorGraph;
111 typedef SymmetricBlockMatrix::Block
Block;
146 const Matrix& G11,
const Matrix& G12,
const Vector& g1,
147 const Matrix& G22,
const Vector& g2,
double f);
154 const Matrix& G11,
const Matrix& G12,
const Matrix& G13,
const Vector& g1,
155 const Matrix& G22,
const Matrix& G23,
const Vector& g2,
156 const Matrix& G33,
const Vector& g3,
double f);
163 const std::vector<Vector>& gs,
double f);
167 template<
typename KEYS>
190 return boost::make_shared<HessianFactor>(*
this); }
193 void print(
const std::string& s =
"",
194 const KeyFormatter& formatter = DefaultKeyFormatter)
const override;
211 return info_.
getDim(std::distance(begin(), variable));
225 bool empty()
const override {
return size() == 0 ; }
286 Matrix augmentedInformation()
const override;
289 Eigen::SelfAdjointView<SymmetricBlockMatrix::constBlock, Eigen::Upper> informationView()
const;
294 Matrix information()
const override;
297 void hessianDiagonalAdd(
VectorValues& d)
const override;
300 using Base::hessianDiagonal;
303 void hessianDiagonal(
double* d)
const override;
306 std::map<Key,Matrix> hessianBlockDiagonal()
const override;
309 std::pair<Matrix, Vector> jacobian()
const override;
316 Matrix augmentedJacobian()
const override;
340 void gradientAtZero(
double* d)
const override;
352 boost::shared_ptr<GaussianConditional> eliminateCholesky(
const Ordering& keys);
359 void Allocate(
const Scatter& scatter);
368 friend class boost::serialization::access;
369 template<
class ARCHIVE>
370 void serialize(ARCHIVE & ar,
const unsigned int ) {
372 ar & BOOST_SERIALIZATION_NVP(info_);
392GTSAM_EXPORT std::pair<boost::shared_ptr<GaussianConditional>, boost::shared_ptr<HessianFactor> >
393 EliminateCholesky(
const GaussianFactorGraph& factors,
const Ordering& keys);
410GTSAM_EXPORT std::pair<boost::shared_ptr<GaussianConditional>, boost::shared_ptr<GaussianFactor> >
411 EliminatePreferCholesky(
const GaussianFactorGraph& factors,
const Ordering& keys);
A thin wrapper around std::vector that uses a custom allocator.
Access to matrices via blocks of pre-defined sizes.
Maps global variable indices to slot indices.
A factor with a quadratic error function - a Gaussian.
Contains the HessianFactor class, a general quadratic factor.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:86
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:75
std::string serialize(const T &input)
serializes to a string
Definition: serialization.h:112
void print(const Matrix &A, const string &s, ostream &stream)
print without optional string, must specify cout yourself
Definition: Matrix.cpp:155
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
Definition: SymmetricBlockMatrix.h:52
DenseIndex getDim(DenseIndex block) const
Number of dimensions for variable on this diagonal block.
Definition: SymmetricBlockMatrix.h:128
constBlock aboveDiagonalBlock(DenseIndex I, DenseIndex J) const
Get block above the diagonal (I, J).
Definition: SymmetricBlockMatrix.h:155
DenseIndex rows() const
Row size.
Definition: SymmetricBlockMatrix.h:119
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:175
Eigen::SelfAdjointView< Block, Eigen::Upper > diagonalBlock(DenseIndex J)
Return the J'th diagonal block as a self adjoint view.
Definition: SymmetricBlockMatrix.h:140
Template to create a binary predicate.
Definition: Testable.h:111
A helper that implements the traits interface for GTSAM types.
Definition: Testable.h:151
KeyVector keys_
The keys involved in this factor.
Definition: Factor.h:73
KeyVector::const_iterator const_iterator
Const iterator over keys.
Definition: Factor.h:68
Definition: Ordering.h:34
An abstract virtual base class for JacobianFactor and HessianFactor.
Definition: GaussianFactor.h:39
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactor.h:42
A Linear Factor Graph is a factor graph where all factors are Gaussian, i.e.
Definition: GaussianFactorGraph.h:69
A Gaussian factor using the canonical parameters (information form)
Definition: HessianFactor.h:101
HessianFactor This
Typedef to this class.
Definition: HessianFactor.h:109
void updateHessian(HessianFactor *other) const
Update another Hessian factor.
Definition: HessianFactor.h:328
boost::shared_ptr< This > shared_ptr
A shared_ptr to this class.
Definition: HessianFactor.h:110
GaussianFactor::shared_ptr clone() const override
Clone this HessianFactor.
Definition: HessianFactor.h:189
const SymmetricBlockMatrix & info() const
Return underlying information matrix.
Definition: HessianFactor.h:265
SymmetricBlockMatrix::Block linearTerm()
Return the complete linear term as described above.
Definition: HessianFactor.h:259
double constantTerm() const
Return the constant term as described above.
Definition: HessianFactor.h:230
SymmetricBlockMatrix::constBlock linearTerm() const
Return the complete linear term as described above.
Definition: HessianFactor.h:251
~HessianFactor() override
Destructor.
Definition: HessianFactor.h:186
HessianFactor(const GaussianFactorGraph &factors)
Combine a set of factors into a single dense HessianFactor.
Definition: HessianFactor.h:182
DenseIndex getDim(const_iterator variable) const override
Return the dimension of the variable pointed to by the given key iterator todo: Remove this in favor ...
Definition: HessianFactor.h:210
SymmetricBlockMatrix::constBlock linearTerm(const_iterator j) const
Return the part of linear term as described above corresponding to the requested variable.
Definition: HessianFactor.h:244
SymmetricBlockMatrix info_
The full augmented information matrix, s.t. the quadratic error is 0.5*[x -1]'H[x -1].
Definition: HessianFactor.h:104
GaussianFactor Base
Typedef to base class.
Definition: HessianFactor.h:108
size_t rows() const
Return the number of columns and rows of the Hessian matrix, including the information vector.
Definition: HessianFactor.h:215
SymmetricBlockMatrix::constBlock constBlock
A block from the Hessian matrix (const version)
Definition: HessianFactor.h:112
SymmetricBlockMatrix & info()
Return non-const information matrix.
Definition: HessianFactor.h:269
SymmetricBlockMatrix::Block Block
A block from the Hessian matrix.
Definition: HessianFactor.h:111
double & constantTerm()
Return the constant term as described above.
Definition: HessianFactor.h:238
bool empty() const override
Check if the factor is empty.
Definition: HessianFactor.h:225
A Gaussian factor in the squared-error form.
Definition: JacobianFactor.h:91
Scatter is an intermediate data structure used when building a HessianFactor incrementally,...
Definition: Scatter.h:49
This class represents a collection of vector-valued variables associated each with a unique integer i...
Definition: VectorValues.h:74
A non-linear factor graph is a graph of non-Gaussian, i.e.
Definition: NonlinearFactorGraph.h:78
Definition: NonlinearClusterTree.h:14