gtsam 4.2.0
gtsam
Loading...
Searching...
No Matches
gtsam::JacobianFactor Class Reference

Detailed Description

A Gaussian factor in the squared-error form.

HessianFactor, which represent a Gaussian likelihood over a set of variables.

JacobianFactor implements a Gaussian, which has quadratic negative log-likelihood

\[ E(x) = \frac{1}{2} (Ax-b)^T \Sigma^{-1} (Ax-b) \]

where \( \Sigma \) is a diagonal covariance matrix. The matrix \( A \), r.h.s. vector \( b \), and diagonal noise model \( \Sigma \) are stored in this class.

This factor represents the sum-of-squares error of a linear measurement function, and is created upon linearization of a NoiseModelFactor, which in turn is a sum-of-squares factor with a nonlinear measurement function.

Here is an example of how this factor represents a sum-of-squares error:

Letting \( h(x) \) be a linear measurement prediction function, \( z \) be the actual observed measurement, the residual is

\[ f(x) = h(x) - z . \]

If we expect noise with diagonal covariance matrix \( \Sigma \) on this measurement, then the negative log-likelihood of the Gaussian induced by this measurement model is

\[ E(x) = \frac{1}{2} (h(x) - z)^T \Sigma^{-1} (h(x) - z) . \]

Because \( h(x) \) is linear, we can write it as

\[ h(x) = Ax + e \]

and thus we have

\[ E(x) = \frac{1}{2} (Ax-b)^T \Sigma^{-1} (Ax-b) \]

where \( b = z - e \).

This factor can involve an arbitrary number of variables, and in the above example \( x \) would almost always be only be a subset of the variables in the entire factor graph. There are special constructors for 1-, 2-, and 3- way JacobianFactors, and additional constructors for creating n-way JacobianFactors. The Jacobian matrix \( A \) is passed to these constructors in blocks, for example, for a 2-way factor, the constructor would accept \( A1 \) and \( A2 \), as well as the variable indices \( j1 \) and \( j2 \) and the negative log-likelihood represented by this factor would be

\[ E(x) = \frac{1}{2} (A_1 x_{j1} + A_2 x_{j2} - b)^T \Sigma^{-1} (A_1 x_{j1} + A_2 x_{j2} - b) . \]

  • Discrete factors, such as
+ Inheritance diagram for gtsam::JacobianFactor:

Public Member Functions

 JacobianFactor (const GaussianFactor &gf)
 Convert from other GaussianFactor.
 
 JacobianFactor (const JacobianFactor &jf)
 Copy constructor.
 
 JacobianFactor (const HessianFactor &hf)
 Conversion from HessianFactor (does Cholesky to obtain Jacobian matrix)
 
 JacobianFactor ()
 default constructor for I/O
 
 JacobianFactor (const Vector &b_in)
 Construct Null factor.
 
 JacobianFactor (Key i1, const Matrix &A1, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 Construct unary factor.
 
 JacobianFactor (Key i1, const Matrix &A1, Key i2, const Matrix &A2, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 Construct binary factor.
 
 JacobianFactor (Key i1, const Matrix &A1, Key i2, const Matrix &A2, Key i3, const Matrix &A3, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 Construct ternary factor.
 
template<typename TERMS >
 JacobianFactor (const TERMS &terms, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 Construct an n-ary factor.
 
template<typename KEYS >
 JacobianFactor (const KEYS &keys, const VerticalBlockMatrix &augmentedMatrix, const SharedDiagonal &sigmas=SharedDiagonal())
 Constructor with arbitrary number keys, and where the augmented matrix is given all together instead of in block terms.
 
 JacobianFactor (const GaussianFactorGraph &graph)
 Build a dense joint factor from all the factors in a factor graph.
 
 JacobianFactor (const GaussianFactorGraph &graph, const VariableSlots &p_variableSlots)
 Build a dense joint factor from all the factors in a factor graph.
 
 JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering)
 Build a dense joint factor from all the factors in a factor graph.
 
 JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering, const VariableSlots &p_variableSlots)
 Build a dense joint factor from all the factors in a factor graph.
 
 ~JacobianFactor () override
 Virtual destructor.
 
GaussianFactor::shared_ptr clone () const override
 Clone this JacobianFactor.
 
void print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override
 print
 
bool equals (const GaussianFactor &lf, double tol=1e-9) const override
 Equals for testable.
 
Vector unweighted_error (const VectorValues &c) const
 
Vector error_vector (const VectorValues &c) const
 (A*x-b)
 
double error (const VectorValues &c) const override
 
Matrix augmentedInformation () const override
 Return the augmented information matrix represented by this GaussianFactor.
 
Matrix information () const override
 Return the non-augmented information matrix represented by this GaussianFactor.
 
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
 Returns (dense) A,b pair associated with factor, bakes in the weights.
 
std::pair< Matrix, Vector > jacobianUnweighted () const
 Returns (dense) A,b pair associated with factor, does not bake in weights.
 
Matrix augmentedJacobian () const override
 Return (dense) matrix associated with factor.
 
Matrix augmentedJacobianUnweighted () const
 Return (dense) matrix associated with factor.
 
const VerticalBlockMatrixmatrixObject () const
 Return the full augmented Jacobian matrix of this factor as a VerticalBlockMatrix object.
 
VerticalBlockMatrixmatrixObject ()
 Mutable access to the full augmented Jacobian matrix of this factor as a VerticalBlockMatrix object.
 
GaussianFactor::shared_ptr negate () const override
 Construct the corresponding anti-factor to negate information stored stored in this factor.
 
bool isConstrained () const
 is noise model constrained ?
 
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 of keeping track of dimensions with variables?
 
size_t rows () const
 return the number of rows in the corresponding linear system
 
size_t cols () const
 return the number of columns in the corresponding linear system
 
const SharedDiagonal & get_model () const
 get a copy of model
 
SharedDiagonal & get_model ()
 get a copy of model (non-const version)
 
const constBVector getb () const
 Get a view of the r.h.s.
 
constABlock getA (const_iterator variable) const
 Get a view of the A matrix for the variable pointed to by the given key iterator.
 
constABlock getA () const
 Get a view of the A matrix, not weighted by noise.
 
BVector getb ()
 Get a view of the r.h.s.
 
ABlock getA (iterator variable)
 Get a view of the A matrix for the variable pointed to by the given key iterator (non-const version)
 
ABlock getA ()
 Get a view of the A matrix.
 
void updateHessian (const KeyVector &keys, SymmetricBlockMatrix *info) const override
 Update an information matrix by adding the information corresponding to this factor (used internally during elimination).
 
Vector operator* (const VectorValues &x) const
 Return A*x.
 
void transposeMultiplyAdd (double alpha, const Vector &e, VectorValues &x) const
 x += alpha * A'*e.
 
void multiplyHessianAdd (double alpha, const VectorValues &x, VectorValues &y) const override
 y += alpha * A'*A*x
 
void multiplyHessianAdd (double alpha, const double *x, double *y, const std::vector< size_t > &accumulatedDims) const
 Raw memory access version of multiplyHessianAdd y += alpha * A'*A*x Requires the vector accumulatedDims to tell the dimension of each variable: e.g.: x0 has dim 3, x2 has dim 6, x3 has dim 2, then accumulatedDims is [0 3 9 11 13] NOTE: size of accumulatedDims is size of keys + 1!! TODO(frank): we should probably kill this if no longer needed.
 
VectorValues gradientAtZero () const override
 A'*b for Jacobian.
 
void gradientAtZero (double *d) const override
 A'*b for Jacobian (raw memory version)
 
Vector gradient (Key key, const VectorValues &x) const override
 Compute the gradient wrt a key at any values.
 
JacobianFactor whiten () const
 Return a whitened version of the factor, i.e.
 
std::pair< boost::shared_ptr< GaussianConditional >, shared_ptreliminate (const Ordering &keys)
 Eliminate the requested variables.
 
void setModel (bool anyConstrained, const Vector &sigmas)
 set noiseModel correctly
 
boost::shared_ptr< GaussianConditionalsplitConditional (size_t nrFrontals)
 splits a pre-factorized factor into a conditional, and changes the current factor to be the remaining component.
 
virtual double error (const VectorValues &c) const
 
double error (const HybridValues &c) const override
 All factor types need to implement an error function.
 
VectorValues hessianDiagonal () const
 Using the base method.
 
virtual void hessianDiagonal (double *d) const=0
 Using the base method.
 
- Public Member Functions inherited from gtsam::GaussianFactor
 GaussianFactor ()
 Default constructor creates empty factor.
 
template<typename CONTAINER >
 GaussianFactor (const CONTAINER &keys)
 Construct from container of keys.
 
virtual ~GaussianFactor ()
 Destructor.
 
VectorValues hessianDiagonal () const
 Return the diagonal of the Hessian for this factor.
 
- Public Member Functions inherited from gtsam::Factor
virtual ~Factor ()=default
 Default destructor.
 
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
 Iterator at beginning of involved variable keys.
 
const_iterator end () const
 Iterator at end of involved variable keys.
 
size_t size () const
 
virtual void printKeys (const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const
 print only keys
 
bool equals (const This &other, double tol=1e-9) const
 check equality
 
KeyVectorkeys ()
 
iterator begin ()
 Iterator at beginning of involved variable keys.
 
iterator end ()
 Iterator at end of involved variable keys.
 

Public Types

typedef JacobianFactor This
 Typedef to this class.
 
typedef GaussianFactor Base
 Typedef to base class.
 
typedef boost::shared_ptr< Thisshared_ptr
 shared_ptr to this class
 
typedef VerticalBlockMatrix::Block ABlock
 
typedef VerticalBlockMatrix::constBlock constABlock
 
typedef ABlock::ColXpr BVector
 
typedef constABlock::ConstColXpr constBVector
 
- Public Types inherited from gtsam::GaussianFactor
typedef GaussianFactor This
 This class.
 
typedef boost::shared_ptr< Thisshared_ptr
 shared_ptr to this class
 
typedef Factor Base
 Our base class.
 
- Public Types inherited from gtsam::Factor
typedef KeyVector::iterator iterator
 Iterator over keys.
 
typedef KeyVector::const_iterator const_iterator
 Const iterator over keys.
 

Protected Member Functions

template<typename TERMS >
void fillTerms (const TERMS &terms, const Vector &b, const SharedDiagonal &noiseModel)
 Internal function to fill blocks and set dimensions.
 
- Protected Member Functions inherited from gtsam::Factor
 Factor ()
 Default constructor for I/O.
 
template<typename CONTAINER >
 Factor (const CONTAINER &keys)
 Construct factor from container of keys.
 
template<typename ITERATOR >
 Factor (ITERATOR first, ITERATOR last)
 Construct factor from iterator keys.
 

Protected Attributes

VerticalBlockMatrix Ab_
 
noiseModel::Diagonal::shared_ptr model_
 
- Protected Attributes inherited from gtsam::Factor
KeyVector keys_
 The keys involved in this factor.
 

Friends

template<typename T >
class ExpressionFactor
 
class boost::serialization::access
 Serialization function.
 
GTSAM_EXPORT std::pair< boost::shared_ptr< GaussianConditional >, shared_ptrEliminateQR (const GaussianFactorGraph &factors, const Ordering &keys)
 Densely partially eliminate with QR factorization, this is usually provided as an argument to one of the factor graph elimination functions (see EliminateableFactorGraph).
 

Additional Inherited Members

- Static Public Member Functions inherited from gtsam::GaussianFactor
template<typename CONTAINER >
static DenseIndex Slot (const CONTAINER &keys, Key key)
 
- Static Protected Member Functions inherited from gtsam::Factor
template<typename CONTAINER >
static Factor FromKeys (const CONTAINER &keys)
 Construct factor from container of keys.
 
template<typename ITERATOR >
static Factor FromIterators (ITERATOR first, ITERATOR last)
 Construct factor from iterator keys.
 

Constructor & Destructor Documentation

◆ JacobianFactor() [1/6]

template<typename TERMS >
JacobianFactor::JacobianFactor ( const TERMS &  terms,
const Vector &  b,
const SharedDiagonal &  model = SharedDiagonal() 
)

Construct an n-ary factor.

Template Parameters
TERMSA container whose value type is std::pair<Key, Matrix>, specifying the collection of keys and matrices making up the factor.

◆ JacobianFactor() [2/6]

template<typename KEYS >
JacobianFactor::JacobianFactor ( const KEYS &  keys,
const VerticalBlockMatrix augmentedMatrix,
const SharedDiagonal &  sigmas = SharedDiagonal() 
)

Constructor with arbitrary number keys, and where the augmented matrix is given all together instead of in block terms.

Note that only the active view of the provided augmented matrix is used, and that the matrix data is copied into a newly-allocated matrix in the constructed factor.

◆ JacobianFactor() [3/6]

JacobianFactor::JacobianFactor ( const GaussianFactorGraph graph)
explicit

Build a dense joint factor from all the factors in a factor graph.

If a VariableSlots structure computed for graph is already available, providing it will reduce the amount of computation performed.

◆ JacobianFactor() [4/6]

JacobianFactor::JacobianFactor ( const GaussianFactorGraph graph,
const VariableSlots p_variableSlots 
)
explicit

Build a dense joint factor from all the factors in a factor graph.

If a VariableSlots structure computed for graph is already available, providing it will reduce the amount of computation performed.

◆ JacobianFactor() [5/6]

JacobianFactor::JacobianFactor ( const GaussianFactorGraph graph,
const Ordering ordering 
)
explicit

Build a dense joint factor from all the factors in a factor graph.

If a VariableSlots structure computed for graph is already available, providing it will reduce the amount of computation performed.

◆ JacobianFactor() [6/6]

JacobianFactor::JacobianFactor ( const GaussianFactorGraph graph,
const Ordering ordering,
const VariableSlots p_variableSlots 
)
explicit

Build a dense joint factor from all the factors in a factor graph.

If a VariableSlots structure computed for graph is already available, providing it will reduce the amount of computation performed.

Member Function Documentation

◆ augmentedInformation()

Matrix JacobianFactor::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.

Implements gtsam::GaussianFactor.

◆ augmentedJacobian()

Matrix JacobianFactor::augmentedJacobian ( ) const
overridevirtual

Return (dense) matrix associated with factor.

The returned system is an augmented matrix: [A b] weights are baked in

Implements gtsam::GaussianFactor.

◆ augmentedJacobianUnweighted()

Matrix JacobianFactor::augmentedJacobianUnweighted ( ) const

Return (dense) matrix associated with factor.

The returned system is an augmented matrix: [A b] weights are not baked in

◆ clone()

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

Clone this JacobianFactor.

Implements gtsam::GaussianFactor.

◆ equals()

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

Equals for testable.

Implements gtsam::GaussianFactor.

◆ error() [1/3]

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::GaussianFactor.

◆ error() [2/3]

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

Reimplemented from gtsam::GaussianFactor.

◆ error() [3/3]

double JacobianFactor::error ( const VectorValues c) const
overridevirtual

Reimplemented from gtsam::GaussianFactor.

◆ getb() [1/2]

BVector gtsam::JacobianFactor::getb ( )
inline

Get a view of the r.h.s.

vector b (non-const version)

◆ getb() [2/2]

const constBVector gtsam::JacobianFactor::getb ( ) const
inline

Get a view of the r.h.s.

vector b, not weighted by noise

◆ getDim()

DenseIndex gtsam::JacobianFactor::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?

Implements gtsam::GaussianFactor.

◆ gradient()

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

Compute the gradient wrt a key at any values.

Implements gtsam::GaussianFactor.

◆ gradientAtZero() [1/2]

VectorValues JacobianFactor::gradientAtZero ( ) const
overridevirtual

A'*b for Jacobian.

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

◆ gradientAtZero() [2/2]

void JacobianFactor::gradientAtZero ( double *  d) const
overridevirtual

A'*b for Jacobian (raw memory version)

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

◆ hessianBlockDiagonal()

map< Key, Matrix > JacobianFactor::hessianBlockDiagonal ( ) const
overridevirtual

Return the block diagonal of the Hessian for this factor.

Implements gtsam::GaussianFactor.

◆ hessianDiagonal() [1/2]

void JacobianFactor::hessianDiagonal ( double *  d) const
overridevirtual

Raw memory access version of hessianDiagonal.

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >, and gtsam::RegularJacobianFactor< D >.

◆ hessianDiagonal() [2/2]

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

Using the base method.

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >, and gtsam::RegularJacobianFactor< D >.

◆ hessianDiagonalAdd()

void JacobianFactor::hessianDiagonalAdd ( VectorValues d) const
overridevirtual

Add the current diagonal to a VectorValues instance.

Implements gtsam::GaussianFactor.

◆ information()

Matrix JacobianFactor::information ( ) const
overridevirtual

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

Implements gtsam::GaussianFactor.

◆ jacobian()

pair< Matrix, Vector > JacobianFactor::jacobian ( ) const
overridevirtual

Returns (dense) A,b pair associated with factor, bakes in the weights.

Implements gtsam::GaussianFactor.

◆ multiplyHessianAdd() [1/2]

void JacobianFactor::multiplyHessianAdd ( double  alpha,
const double *  x,
double *  y,
const std::vector< size_t > &  accumulatedDims 
) const

Raw memory access version of multiplyHessianAdd y += alpha * A'*A*x Requires the vector accumulatedDims to tell the dimension of each variable: e.g.: x0 has dim 3, x2 has dim 6, x3 has dim 2, then accumulatedDims is [0 3 9 11 13] NOTE: size of accumulatedDims is size of keys + 1!! TODO(frank): we should probably kill this if no longer needed.

Raw memory access version of multiplyHessianAdd y += alpha * A'*A*x Note: this is not assuming a fixed dimension for the variables, but requires the vector accumulatedDims to tell the dimension of each variable: e.g.: x0 has dim 3, x2 has dim 6, x3 has dim 2, then accumulatedDims is [0 3 9 11 13] NOTE: size of accumulatedDims is size of keys + 1!! TODO Frank asks: why is this here if not regular ????

Use Eigen magic to access raw memory

Just iterate over all A matrices and multiply in correct config part (looping over keys) E.g.: Jacobian A = [A0 A1 A2] multiplies x = [x0 x1 x2]' Hence: Ax = A0 x0 + A1 x1 + A2 x2 (hence we loop over the keys and accumulate)

Deal with noise properly, need to Double* whiten as we are dividing by variance

multiply with alpha

Again iterate over all A matrices and insert Ai^T into y

◆ multiplyHessianAdd() [2/2]

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

y += alpha * A'*A*x

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >, and gtsam::RegularJacobianFactor< D >.

◆ negate()

GaussianFactor::shared_ptr JacobianFactor::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 JacobianFactor::print ( const std::string &  s = "",
const KeyFormatter formatter = DefaultKeyFormatter 
) const
overridevirtual

print

Implements gtsam::GaussianFactor.

◆ splitConditional()

GaussianConditional::shared_ptr JacobianFactor::splitConditional ( size_t  nrFrontals)

splits a pre-factorized factor into a conditional, and changes the current factor to be the remaining component.

Performs same operation as eliminate(), but without running QR. NOTE: looks at dimension of noise model to determine how many rows to keep.

Parameters
nrFrontalsnumber of keys to eliminate

◆ transposeMultiplyAdd()

void JacobianFactor::transposeMultiplyAdd ( double  alpha,
const Vector &  e,
VectorValues x 
) const

x += alpha * A'*e.

If x is initially missing any values, they are created and assumed to start as zero vectors.

◆ updateHessian()

void JacobianFactor::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
scatterA mapping from variable index to slot index in this HessianFactor
infoThe information matrix to be updated

Implements gtsam::GaussianFactor.

◆ whiten()

JacobianFactor JacobianFactor::whiten ( ) const

Return a whitened version of the factor, i.e.

with unit diagonal noise model.


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