Public Types | List of all members
gtsam::GaussianFactor Class Referenceabstract

#include <GaussianFactor.h>

Inheritance diagram for gtsam::GaussianFactor:
Inheritance graph
[legend]

Public Types

typedef Factor Base
 Our base class. More...
 
typedef std::shared_ptr< Thisshared_ptr
 shared_ptr to this class More...
 
typedef GaussianFactor This
 This class. More...
 
- Public Types inherited from gtsam::Factor
typedef KeyVector::const_iterator const_iterator
 Const iterator over keys. More...
 
typedef KeyVector::iterator iterator
 Iterator over keys. More...
 

Public Member Functions

Standard Constructors
 GaussianFactor ()
 
template<typename CONTAINER >
 GaussianFactor (const CONTAINER &keys)
 
Testable
void print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override=0
 print with optional string More...
 
virtual bool equals (const GaussianFactor &lf, double tol=1e-9) const =0
 assert equality up to a tolerance More...
 
Standard Interface
virtual double error (const VectorValues &c) const
 
double error (const HybridValues &c) const override
 
virtual DenseIndex getDim (const_iterator variable) const =0
 
virtual Matrix augmentedJacobian () const =0
 
virtual std::pair< Matrix, Vectorjacobian () const =0
 
virtual Matrix augmentedInformation () const =0
 
virtual Matrix information () const =0
 
VectorValues hessianDiagonal () const
 Return the diagonal of the Hessian for this factor. More...
 
virtual void hessianDiagonalAdd (VectorValues &d) const =0
 Add the current diagonal to a VectorValues instance. More...
 
virtual void hessianDiagonal (double *d) const =0
 Raw memory access version of hessianDiagonal. More...
 
virtual std::map< Key, MatrixhessianBlockDiagonal () const =0
 Return the block diagonal of the Hessian for this factor. More...
 
virtual GaussianFactor::shared_ptr clone () const =0
 
virtual GaussianFactor::shared_ptr negate () const =0
 
virtual void updateHessian (const KeyVector &keys, SymmetricBlockMatrix *info) const =0
 
Operator interface
virtual void multiplyHessianAdd (double alpha, const VectorValues &x, VectorValues &y) const =0
 y += alpha * A'*A*x More...
 
virtual VectorValues gradientAtZero () const =0
 A'*b for Jacobian, eta for Hessian. More...
 
virtual void gradientAtZero (double *d) const =0
 Raw memory access version of gradientAtZero. More...
 
virtual Vector gradient (Key key, const VectorValues &x) const =0
 Gradient wrt a key at any values. More...
 
- Public Member Functions inherited from gtsam::Factor
virtual ~Factor ()=default
 Default destructor. More...
 
bool empty () const
 Whether the factor is empty (involves zero variables). More...
 
Key front () const
 First key. More...
 
Key back () const
 Last key. More...
 
const_iterator find (Key key) const
 find More...
 
const KeyVectorkeys () const
 Access the factor's involved variable keys. More...
 
const_iterator begin () const
 
const_iterator end () const
 
size_t size () const
 
virtual void printKeys (const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const
 print only keys More...
 
bool equals (const This &other, double tol=1e-9) const
 check equality More...
 
KeyVectorkeys ()
 
iterator begin ()
 
iterator end ()
 

Static Public Member Functions

Advanced Interface
template<typename CONTAINER >
static DenseIndex Slot (const CONTAINER &keys, Key key)
 

Additional Inherited Members

- Protected Member Functions inherited from gtsam::Factor
 Factor ()
 
template<typename CONTAINER >
 Factor (const CONTAINER &keys)
 
template<typename ITERATOR >
 Factor (ITERATOR first, ITERATOR last)
 
- Static Protected Member Functions inherited from gtsam::Factor
template<typename CONTAINER >
static Factor FromKeys (const CONTAINER &keys)
 
template<typename ITERATOR >
static Factor FromIterators (ITERATOR first, ITERATOR last)
 
- Protected Attributes inherited from gtsam::Factor
KeyVector keys_
 The keys involved in this factor. More...
 

Detailed Description

An abstract virtual base class for JacobianFactor and HessianFactor. A GaussianFactor has a quadratic error function. GaussianFactor is non-mutable (all methods const!). The factor value is exp(-0.5*||Ax-b||^2)

Definition at line 38 of file GaussianFactor.h.

Member Typedef Documentation

◆ Base

Our base class.

Definition at line 43 of file GaussianFactor.h.

◆ shared_ptr

typedef std::shared_ptr<This> gtsam::GaussianFactor::shared_ptr

shared_ptr to this class

Definition at line 42 of file GaussianFactor.h.

◆ This

This class.

Definition at line 41 of file GaussianFactor.h.

Constructor & Destructor Documentation

◆ GaussianFactor() [1/2]

gtsam::GaussianFactor::GaussianFactor ( )
inline

Default constructor creates empty factor

Definition at line 49 of file GaussianFactor.h.

◆ GaussianFactor() [2/2]

template<typename CONTAINER >
gtsam::GaussianFactor::GaussianFactor ( const CONTAINER &  keys)
inline

Construct from container of keys. This constructor is used internally from derived factor constructors, either from a container of keys or from a boost::assign::list_of.

Definition at line 54 of file GaussianFactor.h.

Member Function Documentation

◆ augmentedInformation()

virtual Matrix gtsam::GaussianFactor::augmentedInformation ( ) const
pure virtual

Return the augmented information matrix represented by this GaussianFactor. The augmented information matrix contains the information matrix with an additional column holding the information vector, and an additional row holding the transpose of the information vector. The lower-right entry contains the constant error term (when $ \delta x = 0 $). The augmented information matrix is described in more detail in HessianFactor, which in fact stores an augmented information matrix.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ augmentedJacobian()

virtual Matrix gtsam::GaussianFactor::augmentedJacobian ( ) const
pure virtual

Return a dense $ [ \;A\;b\; ] \in \mathbb{R}^{m \times n+1} $ Jacobian matrix, augmented with b with the noise models baked into A and b. The negative log-likelihood is $ \frac{1}{2} \Vert Ax-b \Vert^2 $. See also GaussianFactorGraph::jacobian and GaussianFactorGraph::sparseJacobian.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ clone()

virtual GaussianFactor::shared_ptr gtsam::GaussianFactor::clone ( ) const
pure virtual

◆ equals()

virtual bool gtsam::GaussianFactor::equals ( const GaussianFactor lf,
double  tol = 1e-9 
) const
pure virtual

◆ error() [1/2]

double GaussianFactor::error ( const HybridValues c) const
overridevirtual

All factor types need to implement an error function. In factor graphs, this is the negative log-likelihood.

Reimplemented from gtsam::Factor.

Definition at line 31 of file GaussianFactor.cpp.

◆ error() [2/2]

double GaussianFactor::error ( const VectorValues c) const
virtual

◆ getDim()

virtual DenseIndex gtsam::GaussianFactor::getDim ( const_iterator  variable) const
pure virtual

Return the dimension of the variable pointed to by the given key iterator

Implemented in gtsam::JacobianFactor, gtsam::HessianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ gradient()

virtual Vector gtsam::GaussianFactor::gradient ( Key  key,
const VectorValues x 
) const
pure virtual

Gradient wrt a key at any values.

Implemented in gtsam::RegularImplicitSchurFactor< CAMERA >, gtsam::JacobianFactor, and gtsam::HessianFactor.

◆ gradientAtZero() [1/2]

virtual VectorValues gtsam::GaussianFactor::gradientAtZero ( ) const
pure virtual

◆ gradientAtZero() [2/2]

virtual void gtsam::GaussianFactor::gradientAtZero ( double *  d) const
pure virtual

◆ hessianBlockDiagonal()

virtual std::map<Key,Matrix> gtsam::GaussianFactor::hessianBlockDiagonal ( ) const
pure virtual

Return the block diagonal of the Hessian for this factor.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ hessianDiagonal() [1/2]

VectorValues GaussianFactor::hessianDiagonal ( ) const

Return the diagonal of the Hessian for this factor.

Definition at line 35 of file GaussianFactor.cpp.

◆ hessianDiagonal() [2/2]

virtual void gtsam::GaussianFactor::hessianDiagonal ( double *  d) const
pure virtual

◆ hessianDiagonalAdd()

virtual void gtsam::GaussianFactor::hessianDiagonalAdd ( VectorValues d) const
pure virtual

Add the current diagonal to a VectorValues instance.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ information()

virtual Matrix gtsam::GaussianFactor::information ( ) const
pure virtual

Return the non-augmented information matrix represented by this GaussianFactor.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ jacobian()

virtual std::pair<Matrix,Vector> gtsam::GaussianFactor::jacobian ( ) const
pure virtual

Return the dense Jacobian $ A $ and right-hand-side $ b $, with the noise models baked into A and b. The negative log-likelihood is $ \frac{1}{2} \Vert Ax-b \Vert^2 $. See also GaussianFactorGraph::augmentedJacobian and GaussianFactorGraph::sparseJacobian.

Implemented in gtsam::HessianFactor, gtsam::JacobianFactor, and gtsam::RegularImplicitSchurFactor< CAMERA >.

◆ multiplyHessianAdd()

virtual void gtsam::GaussianFactor::multiplyHessianAdd ( double  alpha,
const VectorValues x,
VectorValues y 
) const
pure virtual

◆ negate()

virtual GaussianFactor::shared_ptr gtsam::GaussianFactor::negate ( ) const
pure virtual

Construct the corresponding anti-factor to negate information stored stored in this factor.

Returns
a HessianFactor with negated Hessian matrices

Implemented in gtsam::JacobianFactor, gtsam::RegularImplicitSchurFactor< CAMERA >, and gtsam::HessianFactor.

◆ print()

void gtsam::GaussianFactor::print ( const std::string &  s = "",
const KeyFormatter formatter = DefaultKeyFormatter 
) const
overridepure virtual

◆ Slot()

template<typename CONTAINER >
static DenseIndex gtsam::GaussianFactor::Slot ( const CONTAINER &  keys,
Key  key 
)
inlinestatic

Definition at line 175 of file GaussianFactor.h.

◆ updateHessian()

virtual void gtsam::GaussianFactor::updateHessian ( const KeyVector keys,
SymmetricBlockMatrix info 
) const
pure virtual

Update an information matrix by adding the information corresponding to this factor (used internally during elimination).

Parameters
scatterA mapping from variable index to slot index in this HessianFactor
infoThe information matrix to be updated

Implemented in gtsam::JacobianFactor, gtsam::HessianFactor, gtsam::RegularImplicitSchurFactor< CAMERA >, and gtsam::BinaryJacobianFactor< M, N1, N2 >.


The documentation for this class was generated from the following files:


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:15:33