GTSAM  4.0.2
C++ library for smoothing and mapping (SAM)
Public Types | Public Member Functions | Protected Attributes | Friends | List of all members
gtsam::HessianFactor Class Reference

A Gaussian factor using the canonical parameters (information form) More...

#include <HessianFactor.h>

Inheritance diagram for gtsam::HessianFactor:
Inheritance graph
[legend]
Collaboration diagram for gtsam::HessianFactor:
Collaboration graph
[legend]

Public Types

typedef GaussianFactor Base
 Typedef to base class.
 
typedef HessianFactor This
 Typedef to this class.
 
typedef std::shared_ptr< Thisshared_ptr
 A shared_ptr to this class.
 
typedef SymmetricBlockMatrix::Block Block
 A block from the Hessian matrix.
 
typedef SymmetricBlockMatrix::constBlock constBlock
 A block from the Hessian matrix (const version)
 
typedef KeyVector::iterator iterator
 Iterator over keys.
 
typedef KeyVector::const_iterator const_iterator
 Const iterator over keys.
 

Public Member Functions

 HessianFactor ()
 
 HessianFactor (Key j, const Matrix &G, const Vector &g, double f)
 
 HessianFactor (Key j, const Vector &mu, const Matrix &Sigma)
 
 HessianFactor (Key j1, Key j2, const Matrix &G11, const Matrix &G12, const Vector &g1, const Matrix &G22, const Vector &g2, double f)
 
 HessianFactor (Key j1, Key j2, Key j3, const Matrix &G11, const Matrix &G12, const Matrix &G13, const Vector &g1, const Matrix &G22, const Matrix &G23, const Vector &g2, const Matrix &G33, const Vector &g3, double f)
 
 HessianFactor (const KeyVector &js, const std::vector< Matrix > &Gs, const std::vector< Vector > &gs, double f)
 
template<typename KEYS >
 HessianFactor (const KEYS &keys, const SymmetricBlockMatrix &augmentedInformation)
 
 HessianFactor (const JacobianFactor &cg)
 
 HessianFactor (const GaussianFactor &factor)
 
 HessianFactor (const GaussianFactorGraph &factors, const Scatter &scatter)
 
 HessianFactor (const GaussianFactorGraph &factors)
 
 ~HessianFactor () override
 
GaussianFactor::shared_ptr clone () const override
 
void print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override
 
bool equals (const GaussianFactor &lf, double tol=1e-9) const override
 
double error (const VectorValues &c) const override
 
DenseIndex getDim (const_iterator variable) const override
 
size_t rows () const
 
GaussianFactor::shared_ptr negate () const override
 
double constantTerm () const
 
double & constantTerm ()
 
SymmetricBlockMatrix::constBlock linearTerm (const_iterator j) const
 
SymmetricBlockMatrix::constBlock linearTerm () const
 
SymmetricBlockMatrix::Block linearTerm ()
 
const SymmetricBlockMatrixinfo () const
 Return underlying information matrix.
 
SymmetricBlockMatrixinfo ()
 
Matrix augmentedInformation () const override
 
Eigen::SelfAdjointView< SymmetricBlockMatrix::constBlock, Eigen::Upper > informationView () const
 Return self-adjoint view onto the information matrix (NOT augmented).
 
Matrix information () const override
 
void hessianDiagonalAdd (VectorValues &d) const override
 Add the current diagonal to a VectorValues instance.
 
void hessianDiagonal (double *d) const override
 Raw memory access version of hessianDiagonal.
 
std::map< Key, Matrix > hessianBlockDiagonal () const override
 Return the block diagonal of the Hessian for this factor.
 
std::pair< Matrix, Vector > jacobian () const override
 Return (dense) matrix associated with factor.
 
Matrix augmentedJacobian () const override
 
void updateHessian (const KeyVector &keys, SymmetricBlockMatrix *info) const override
 
void updateHessian (HessianFactor *other) const
 
void multiplyHessianAdd (double alpha, const VectorValues &x, VectorValues &y) const override
 
VectorValues gradientAtZero () const override
 eta for Hessian
 
void gradientAtZero (double *d) const override
 Raw memory access version of gradientAtZero.
 
Vector gradient (Key key, const VectorValues &x) const override
 
std::shared_ptr< GaussianConditionaleliminateCholesky (const Ordering &keys)
 
VectorValues solve ()
 Solve the system A'*A delta = A'*b in-place, return delta as VectorValues.
 
Testable
bool equals (const This &other, double tol=1e-9) const
 check equality
 
virtual void printKeys (const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const
 print only keys
 
Standard Interface
double error (const HybridValues &c) const override
 
VectorValues hessianDiagonal () const
 Return the diagonal of the Hessian for this factor.
 
Standard Interface
bool empty () const
 Whether the factor is empty (involves zero variables).
 
Key front () const
 First key.
 
Key back () const
 Last key.
 
const_iterator find (Key key) const
 find
 
const KeyVectorkeys () const
 Access the factor's involved variable keys.
 
const_iterator begin () const
 
const_iterator end () const
 
size_t size () const
 
Advanced Interface
KeyVectorkeys ()
 
iterator begin ()
 
iterator end ()
 

Static Public Member Functions

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

Static Protected Member Functions

Standard Constructors
template<typename CONTAINER >
static Factor FromKeys (const CONTAINER &keys)
 
template<typename ITERATOR >
static Factor FromIterators (ITERATOR first, ITERATOR last)
 

Protected Attributes

SymmetricBlockMatrix info_
 The full augmented information matrix, s.t. the quadratic error is 0.5*[x -1]'H[x -1].
 
KeyVector keys_
 The keys involved in this factor.
 

Friends

class NonlinearFactorGraph
 
class NonlinearClusterTree
 

Detailed Description

A Gaussian factor using the canonical parameters (information form)

HessianFactor implements a general quadratic factor of the form

\[ E(x) = 0.5 x^T G x - x^T g + 0.5 f \]

that stores the matrix \( G \), the vector \( g \), and the constant term \( f \).

When \( G \) is positive semidefinite, this factor represents a Gaussian, in which case \( G \) is the information matrix \( \Lambda \), \( g \) is the information vector \( \eta \), and \( f \) is the residual sum-square-error at the mean, when \( x = \mu \).

Indeed, the negative log-likelihood of a Gaussian is (up to a constant) \( E(x) = 0.5(x-\mu)^T P^{-1} (x-\mu) \) with \( \mu \) the mean and \( P \) the covariance matrix. Expanding the product we get

\[ E(x) = 0.5 x^T P^{-1} x - x^T P^{-1} \mu + 0.5 \mu^T P^{-1} \mu \]

We define the Information matrix (or Hessian) \( \Lambda = P^{-1} \) and the information vector \( \eta = P^{-1} \mu = \Lambda \mu \) to arrive at the canonical form of the Gaussian:

\[ E(x) = 0.5 x^T \Lambda x - x^T \eta + 0.5 \mu^T \Lambda \mu \]

This factor is one of the factors that can be in a GaussianFactorGraph. It may be returned from NonlinearFactor::linearize(), but is also used internally to store the Hessian during Cholesky elimination.

This can represent a quadratic factor with characteristics that cannot be represented using a JacobianFactor (which has the form \( E(x) = \Vert Ax - b \Vert^2 \) and stores the Jacobian \( A \) and error vector \( b \), i.e. is a sum-of-squares factor). For example, a HessianFactor need not be positive semidefinite, it can be indefinite or even negative semidefinite.

If a HessianFactor is indefinite or negative semi-definite, then in order for solving the linear system to be possible, the Hessian of the full system must be positive definite (i.e. when all small Hessians are combined, the result must be positive definite). If this is not the case, an error will occur during elimination.

This class stores G, g, and f as an augmented matrix HessianFactor::matrix_. The upper-left n x n blocks of HessianFactor::matrix_ store the upper-right triangle of G, the upper-right-most column of length n of HessianFactor::matrix_ stores g, and the lower-right entry of HessianFactor::matrix_ stores f, i.e.

HessianFactor::matrix_ = [ G11 G12 G13 ... g1
0 G22 G23 ... g2
0 0 G33 ... g3
: : : :
0 0 0 ... f ]

Blocks can be accessed as follows:

G11 = info(begin(), begin());
G12 = info(begin(), begin()+1);
G23 = info(begin()+1, begin()+2);
g2 = linearTerm(begin()+1);
.......

Constructor & Destructor Documentation

◆ HessianFactor() [1/11]

gtsam::HessianFactor::HessianFactor ( )

default constructor for I/O

◆ HessianFactor() [2/11]

gtsam::HessianFactor::HessianFactor ( Key  j,
const Matrix &  G,
const Vector &  g,
double  f 
)

Construct a unary factor. G is the quadratic term (Hessian matrix), g the linear term (a vector), and f the constant term. The quadratic error is: 0.5*(f - 2*x'*g + x'*G*x)

◆ HessianFactor() [3/11]

gtsam::HessianFactor::HessianFactor ( Key  j,
const Vector &  mu,
const Matrix &  Sigma 
)

Construct a unary factor, given a mean and covariance matrix. error is 0.5*(x-mu)'inv(Sigma)(x-mu)

◆ HessianFactor() [4/11]

gtsam::HessianFactor::HessianFactor ( Key  j1,
Key  j2,
const Matrix &  G11,
const Matrix &  G12,
const Vector &  g1,
const Matrix &  G22,
const Vector &  g2,
double  f 
)

Construct a binary factor. Gxx are the upper-triangle blocks of the quadratic term (the Hessian matrix), gx the pieces of the linear vector term, and f the constant term. JacobianFactor error is

\[ 0.5* (Ax-b)' M (Ax-b) = 0.5*x'A'MAx - x'A'Mb + 0.5*b'Mb \]

HessianFactor error is

\[ 0.5*(x'Gx - 2x'g + f) = 0.5*x'Gx - x'*g + 0.5*f \]

So, with \( A = [A1 A2] \) and \( G=A*'M*A = [A1';A2']*M*[A1 A2] \) we have

n1*n1 G11 = A1'*M*A1
n1*n2 G12 = A1'*M*A2
n2*n2 G22 = A2'*M*A2
n1*1 g1 = A1'*M*b
n2*1 g2 = A2'*M*b
1*1 f = b'*M*b

◆ HessianFactor() [5/11]

gtsam::HessianFactor::HessianFactor ( Key  j1,
Key  j2,
Key  j3,
const Matrix &  G11,
const Matrix &  G12,
const Matrix &  G13,
const Vector &  g1,
const Matrix &  G22,
const Matrix &  G23,
const Vector &  g2,
const Matrix &  G33,
const Vector &  g3,
double  f 
)

Construct a ternary factor. Gxx are the upper-triangle blocks of the quadratic term (the Hessian matrix), gx the pieces of the linear vector term, and f the constant term.

◆ HessianFactor() [6/11]

gtsam::HessianFactor::HessianFactor ( const KeyVector js,
const std::vector< Matrix > &  Gs,
const std::vector< Vector > &  gs,
double  f 
)

Construct an n-way factor. Gs contains the upper-triangle blocks of the quadratic term (the Hessian matrix) provided in row-order, gs the pieces of the linear vector term, and f the constant term.

◆ HessianFactor() [7/11]

template<typename KEYS >
gtsam::HessianFactor::HessianFactor ( const KEYS &  keys,
const SymmetricBlockMatrix augmentedInformation 
)

Constructor with an arbitrary number of keys and with the augmented information matrix specified as a block matrix.

◆ HessianFactor() [8/11]

gtsam::HessianFactor::HessianFactor ( const JacobianFactor cg)
explicit

Construct from a JacobianFactor (or from a GaussianConditional since it derives from it)

◆ HessianFactor() [9/11]

gtsam::HessianFactor::HessianFactor ( const GaussianFactor factor)
explicit

Attempt to construct from any GaussianFactor - currently supports JacobianFactor, HessianFactor, GaussianConditional, or any derived classes.

◆ HessianFactor() [10/11]

gtsam::HessianFactor::HessianFactor ( const GaussianFactorGraph factors,
const Scatter scatter 
)
explicit

Combine a set of factors into a single dense HessianFactor

◆ HessianFactor() [11/11]

gtsam::HessianFactor::HessianFactor ( const GaussianFactorGraph factors)
inlineexplicit

Combine a set of factors into a single dense HessianFactor

◆ ~HessianFactor()

gtsam::HessianFactor::~HessianFactor ( )
inlineoverride

Destructor

Member Function Documentation

◆ augmentedInformation()

Matrix gtsam::HessianFactor::augmentedInformation ( ) const
overridevirtual

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.

For HessianFactor, this is the same as info() except that this function returns a complete symmetric matrix whereas info() returns a matrix where only the upper triangle is valid, but should be interpreted as symmetric. This is because info() returns only a reference to the internal representation of the augmented information matrix, which stores only the upper triangle.

Implements gtsam::GaussianFactor.

◆ augmentedJacobian()

Matrix gtsam::HessianFactor::augmentedJacobian ( ) const
overridevirtual

Return (dense) matrix associated with factor The returned system is an augmented matrix: [A b]

Parameters
setweight to use whitening to bake in weights

Implements gtsam::GaussianFactor.

◆ begin() [1/2]

const_iterator gtsam::Factor::begin ( ) const
inlineinherited

Iterator at beginning of involved variable keys

◆ begin() [2/2]

iterator gtsam::Factor::begin ( )
inlineinherited

Iterator at beginning of involved variable keys

◆ clone()

GaussianFactor::shared_ptr gtsam::HessianFactor::clone ( ) const
inlineoverridevirtual

Clone this HessianFactor

Implements gtsam::GaussianFactor.

◆ constantTerm() [1/2]

double gtsam::HessianFactor::constantTerm ( ) const
inline

Return the constant term \( f \) as described above

Returns
The constant term \( f \)

◆ constantTerm() [2/2]

double& gtsam::HessianFactor::constantTerm ( )
inline

Return the constant term \( f \) as described above

Returns
The constant term \( f \)

◆ eliminateCholesky()

std::shared_ptr<GaussianConditional> gtsam::HessianFactor::eliminateCholesky ( const Ordering keys)

In-place elimination that returns a conditional on (ordered) keys specified, and leaves this factor to be on the remaining keys (separator) only. Does dense partial Cholesky.

◆ end() [1/2]

const_iterator gtsam::Factor::end ( ) const
inlineinherited

Iterator at end of involved variable keys

◆ end() [2/2]

iterator gtsam::Factor::end ( )
inlineinherited

Iterator at end of involved variable keys

◆ equals()

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

Compare to another factor for testing (implementing Testable)

Implements gtsam::GaussianFactor.

◆ error() [1/2]

double gtsam::GaussianFactor::error ( const HybridValues c) const
overridevirtualinherited

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

Reimplemented from gtsam::Factor.

◆ error() [2/2]

double gtsam::HessianFactor::error ( const VectorValues c) const
overridevirtual

Evaluate the factor error f(x). returns 0.5*[x -1]'H[x -1] (also see constructor documentation)

Reimplemented from gtsam::GaussianFactor.

◆ FromIterators()

template<typename ITERATOR >
static Factor gtsam::Factor::FromIterators ( ITERATOR  first,
ITERATOR  last 
)
inlinestaticprotectedinherited

Construct factor from iterator keys. This is called internally from derived factor static factor methods, as a workaround for not being able to call the protected constructors above.

◆ FromKeys()

template<typename CONTAINER >
static Factor gtsam::Factor::FromKeys ( const CONTAINER &  keys)
inlinestaticprotectedinherited

Construct factor from container of keys. This is called internally from derived factor static factor methods, as a workaround for not being able to call the protected constructors above.

◆ getDim()

DenseIndex gtsam::HessianFactor::getDim ( const_iterator  variable) const
inlineoverridevirtual

Return the dimension of the variable pointed to by the given key iterator todo: Remove this in favor of keeping track of dimensions with variables?

Parameters
variableAn iterator pointing to the slot in this factor. You can use, for example, begin() + 2 to get the 3rd variable in this factor.

Implements gtsam::GaussianFactor.

◆ gradient()

Vector gtsam::HessianFactor::gradient ( Key  key,
const VectorValues x 
) const
overridevirtual

Compute the gradient at a key: \( \grad f(x_i) = \sum_j G_ij*x_j - g_i \)

Implements gtsam::GaussianFactor.

◆ info()

SymmetricBlockMatrix& gtsam::HessianFactor::info ( )
inline

Return non-const information matrix. TODO(gareth): Review the sanity of having non-const access to this.

◆ information()

Matrix gtsam::HessianFactor::information ( ) const
overridevirtual

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

Implements gtsam::GaussianFactor.

◆ keys()

KeyVector& gtsam::Factor::keys ( )
inlineinherited
Returns
keys involved in this factor

◆ linearTerm() [1/3]

SymmetricBlockMatrix::constBlock gtsam::HessianFactor::linearTerm ( const_iterator  j) const
inline

Return the part of linear term \( g \) as described above corresponding to the requested variable.

Parameters
jWhich block row to get, as an iterator pointing to the slot in this factor. You can use, for example, begin() + 2 to get the 3rd variable in this factor.
Returns
The linear term \( g \)

◆ linearTerm() [2/3]

SymmetricBlockMatrix::constBlock gtsam::HessianFactor::linearTerm ( ) const
inline

Return the complete linear term \( g \) as described above.

Returns
The linear term \( g \)

◆ linearTerm() [3/3]

SymmetricBlockMatrix::Block gtsam::HessianFactor::linearTerm ( )
inline

Return the complete linear term \( g \) as described above.

Returns
The linear term \( g \)

◆ multiplyHessianAdd()

void gtsam::HessianFactor::multiplyHessianAdd ( double  alpha,
const VectorValues x,
VectorValues y 
) const
overridevirtual

y += alpha * A'*A*x

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularHessianFactor< D >.

◆ negate()

GaussianFactor::shared_ptr gtsam::HessianFactor::negate ( ) const
overridevirtual

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

Returns
a HessianFactor with negated Hessian matrices

Implements gtsam::GaussianFactor.

◆ print()

void gtsam::HessianFactor::print ( const std::string &  s = "",
const KeyFormatter formatter = DefaultKeyFormatter 
) const
overridevirtual

Print the factor for debugging and testing (implementing Testable)

Implements gtsam::GaussianFactor.

◆ rows()

size_t gtsam::HessianFactor::rows ( ) const
inline

Return the number of columns and rows of the Hessian matrix, including the information vector.

◆ size()

size_t gtsam::Factor::size ( ) const
inlineinherited
Returns
the number of variables involved in this factor

◆ updateHessian() [1/2]

void gtsam::HessianFactor::updateHessian ( const KeyVector keys,
SymmetricBlockMatrix info 
) const
overridevirtual

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

Parameters
keysTHe ordered vector of keys for the information matrix to be updated
infoThe information matrix to be updated

Implements gtsam::GaussianFactor.

◆ updateHessian() [2/2]

void gtsam::HessianFactor::updateHessian ( HessianFactor other) const
inline

Update another Hessian factor

Parameters
otherthe HessianFactor to be updated

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