Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
gtsam::JacobianFactor Class Reference

#include <JacobianFactor.h>

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

Public Types

typedef VerticalBlockMatrix::Block ABlock
 
typedef GaussianFactor Base
 Typedef to base class. More...
 
typedef ABlock::ColXpr BVector
 
typedef VerticalBlockMatrix::constBlock constABlock
 
typedef constABlock::ConstColXpr constBVector
 
typedef boost::shared_ptr< Thisshared_ptr
 shared_ptr to this class More...
 
typedef JacobianFactor This
 Typedef to this class. More...
 
- Public Types inherited from gtsam::GaussianFactor
typedef Factor Base
 Our base class. More...
 
typedef boost::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

Matrix augmentedInformation () const override
 
Matrix augmentedJacobian () const override
 
Matrix augmentedJacobianUnweighted () const
 
GaussianFactor::shared_ptr clone () const override
 
size_t cols () const
 
std::pair< boost::shared_ptr< GaussianConditional >, shared_ptreliminate (const Ordering &keys)
 
bool empty () const override
 
bool equals (const GaussianFactor &lf, double tol=1e-9) const override
 
double error (const VectorValues &c) const override
 
Vector error_vector (const VectorValues &c) const
 
const SharedDiagonalget_model () const
 
SharedDiagonalget_model ()
 
constABlock getA (const_iterator variable) const
 
constABlock getA () const
 
ABlock getA (iterator variable)
 
ABlock getA ()
 
const constBVector getb () const
 
BVector getb ()
 
DenseIndex getDim (const_iterator variable) const override
 
Vector gradient (Key key, const VectorValues &x) const override
 Compute the gradient wrt a key at any values. More...
 
VectorValues gradientAtZero () const override
 A'*b for Jacobian. More...
 
void gradientAtZero (double *d) const override
 A'*b for Jacobian (raw memory version) More...
 
std::map< Key, MatrixhessianBlockDiagonal () const override
 Return the block diagonal of the Hessian for this factor. More...
 
void hessianDiagonal (double *d) const override
 Raw memory access version of hessianDiagonal. More...
 
void hessianDiagonalAdd (VectorValues &d) const override
 Add the current diagonal to a VectorValues instance. More...
 
Matrix information () const override
 
bool isConstrained () const
 
std::pair< Matrix, Vectorjacobian () const override
 Returns (dense) A,b pair associated with factor, bakes in the weights. More...
 
 JacobianFactor (const GaussianFactor &gf)
 
 JacobianFactor (const JacobianFactor &jf)
 
 JacobianFactor (const HessianFactor &hf)
 
 JacobianFactor ()
 
 JacobianFactor (const Vector &b_in)
 
 JacobianFactor (Key i1, const Matrix &A1, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 
 JacobianFactor (Key i1, const Matrix &A1, Key i2, const Matrix &A2, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 
 JacobianFactor (Key i1, const Matrix &A1, Key i2, const Matrix &A2, Key i3, const Matrix &A3, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 
template<typename TERMS >
 JacobianFactor (const TERMS &terms, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
 
template<typename KEYS >
 JacobianFactor (const KEYS &keys, const VerticalBlockMatrix &augmentedMatrix, const SharedDiagonal &sigmas=SharedDiagonal())
 
 JacobianFactor (const GaussianFactorGraph &graph)
 
 JacobianFactor (const GaussianFactorGraph &graph, const VariableSlots &p_variableSlots)
 
 JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering)
 
 JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering, const VariableSlots &p_variableSlots)
 
std::pair< Matrix, VectorjacobianUnweighted () const
 Returns (dense) A,b pair associated with factor, does not bake in weights. More...
 
const VerticalBlockMatrixmatrixObject () const
 
VerticalBlockMatrixmatrixObject ()
 
void multiplyHessianAdd (double alpha, const VectorValues &x, VectorValues &y) const override
 
void multiplyHessianAdd (double alpha, const double *x, double *y, const std::vector< size_t > &accumulatedDims) const
 
GaussianFactor::shared_ptr negate () const override
 
Vector operator* (const VectorValues &x) const
 
void print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override
 print More...
 
size_t rows () const
 
void setModel (bool anyConstrained, const Vector &sigmas)
 
boost::shared_ptr< GaussianConditionalsplitConditional (size_t nrFrontals)
 
void transposeMultiplyAdd (double alpha, const Vector &e, VectorValues &x) const
 
Vector unweighted_error (const VectorValues &c) const
 
void updateHessian (const KeyVector &keys, SymmetricBlockMatrix *info) const override
 
JacobianFactor whiten () const
 
 ~JacobianFactor () override
 
- Public Member Functions inherited from gtsam::GaussianFactor
 GaussianFactor ()
 
template<typename CONTAINER >
 GaussianFactor (const CONTAINER &keys)
 
VectorValues hessianDiagonal () const
 Return the diagonal of the Hessian for this factor. More...
 
virtual ~GaussianFactor ()
 
- Public Member Functions inherited from gtsam::Factor
virtual ~Factor ()=default
 Default destructor. 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...
 
KeyVectorkeys ()
 
iterator begin ()
 
iterator end ()
 

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. More...
 
- Protected Member Functions inherited from gtsam::Factor
 Factor ()
 
template<typename CONTAINER >
 Factor (const CONTAINER &keys)
 
template<typename ITERATOR >
 Factor (ITERATOR first, ITERATOR last)
 
bool equals (const This &other, double tol=1e-9) const
 check equality More...
 

Protected Attributes

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

Private Member Functions

template<class KEYS , class DIMENSIONS >
 JacobianFactor (const KEYS &keys, const DIMENSIONS &dims, DenseIndex m, const SharedDiagonal &model=SharedDiagonal())
 
void JacobianFactorHelper (const GaussianFactorGraph &graph, const FastVector< VariableSlots::const_iterator > &orderedSlots)
 
template<class ARCHIVE >
void load (ARCHIVE &ar, const unsigned int version)
 
template<class ARCHIVE >
void save (ARCHIVE &ar, const unsigned int version) const
 

Friends

class boost::serialization::access
 
GTSAM_EXPORT std::pair< boost::shared_ptr< GaussianConditional >, shared_ptrEliminateQR (const GaussianFactorGraph &factors, const Ordering &keys)
 
template<typename T >
class ExpressionFactor
 

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)
 
template<typename ITERATOR >
static Factor FromIterators (ITERATOR first, ITERATOR last)
 

Detailed Description

A Gaussian factor in the squared-error form.

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) . \]

Definition at line 90 of file JacobianFactor.h.

Member Typedef Documentation

Definition at line 98 of file JacobianFactor.h.

Typedef to base class.

Definition at line 95 of file JacobianFactor.h.

Definition at line 100 of file JacobianFactor.h.

Definition at line 99 of file JacobianFactor.h.

Definition at line 101 of file JacobianFactor.h.

typedef boost::shared_ptr<This> gtsam::JacobianFactor::shared_ptr

shared_ptr to this class

Definition at line 96 of file JacobianFactor.h.

Typedef to this class.

Definition at line 94 of file JacobianFactor.h.

Constructor & Destructor Documentation

gtsam::JacobianFactor::JacobianFactor ( const GaussianFactor gf)
explicit

Convert from other GaussianFactor

Definition at line 66 of file JacobianFactor.cpp.

gtsam::JacobianFactor::JacobianFactor ( const JacobianFactor jf)
inline

Copy constructor

Definition at line 114 of file JacobianFactor.h.

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

Conversion from HessianFactor (does Cholesky to obtain Jacobian matrix)

Definition at line 105 of file JacobianFactor.cpp.

gtsam::JacobianFactor::JacobianFactor ( )

default constructor for I/O

Definition at line 60 of file JacobianFactor.cpp.

gtsam::JacobianFactor::JacobianFactor ( const Vector b_in)
explicit

Construct Null factor

Definition at line 78 of file JacobianFactor.cpp.

gtsam::JacobianFactor::JacobianFactor ( Key  i1,
const Matrix A1,
const Vector b,
const SharedDiagonal model = SharedDiagonal() 
)

Construct unary factor

Definition at line 84 of file JacobianFactor.cpp.

gtsam::JacobianFactor::JacobianFactor ( Key  i1,
const Matrix A1,
Key  i2,
const Matrix A2,
const Vector b,
const SharedDiagonal model = SharedDiagonal() 
)

Construct binary factor

Definition at line 90 of file JacobianFactor.cpp.

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

Definition at line 96 of file JacobianFactor.cpp.

template<typename TERMS >
gtsam::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.

Definition at line 27 of file JacobianFactor-inl.h.

template<typename KEYS >
gtsam::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.

Definition at line 34 of file JacobianFactor-inl.h.

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

Definition at line 343 of file JacobianFactor.cpp.

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

Definition at line 368 of file JacobianFactor.cpp.

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

Definition at line 393 of file JacobianFactor.cpp.

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

Definition at line 409 of file JacobianFactor.cpp.

gtsam::JacobianFactor::~JacobianFactor ( )
inlineoverride

Virtual destructor

Definition at line 186 of file JacobianFactor.h.

template<class KEYS , class DIMENSIONS >
gtsam::JacobianFactor::JacobianFactor ( const KEYS &  keys,
const DIMENSIONS &  dims,
DenseIndex  m,
const SharedDiagonal model = SharedDiagonal() 
)
inlineprivate

Unsafe Constructor that creates an uninitialized Jacobian of right size

Parameters
keysin some order
diemnsionsof the variables in same order
moutput dimension
modelnoise model (default nullptr)

Definition at line 407 of file JacobianFactor.h.

Member Function Documentation

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

0.5*(A*x-b)'D(A*x-b) 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.

Definition at line 513 of file JacobianFactor.cpp.

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

Definition at line 749 of file JacobianFactor.cpp.

Matrix gtsam::JacobianFactor::augmentedJacobianUnweighted ( ) const

Return (dense) matrix associated with factor. The returned system is an augmented matrix: [A b] weights are not baked in

Definition at line 757 of file JacobianFactor.cpp.

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

Clone this JacobianFactor

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::LinearInequality, gtsam::LinearEquality, and gtsam::LinearCost.

Definition at line 189 of file JacobianFactor.h.

size_t gtsam::JacobianFactor::cols ( void  ) const
inline

return the number of columns in the corresponding linear system

Definition at line 286 of file JacobianFactor.h.

std::pair< GaussianConditional::shared_ptr, JacobianFactor::shared_ptr > gtsam::JacobianFactor::eliminate ( const Ordering keys)

Eliminate the requested variables.

Definition at line 779 of file JacobianFactor.cpp.

bool gtsam::JacobianFactor::empty ( ) const
inlineoverridevirtual

Check if the factor is empty. TODO: How should this be defined?

Implements gtsam::GaussianFactor.

Definition at line 264 of file JacobianFactor.h.

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

Equals for testable

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::LinearInequality, gtsam::LinearEquality, and gtsam::LinearCost.

Definition at line 439 of file JacobianFactor.cpp.

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

(A*x-b)/sigma

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::LinearInequality, gtsam::LinearEquality, and gtsam::LinearCost.

Definition at line 505 of file JacobianFactor.cpp.

Vector gtsam::JacobianFactor::error_vector ( const VectorValues c) const

(A*x-b)

Definition at line 498 of file JacobianFactor.cpp.

template<typename TERMS >
void gtsam::JacobianFactor::fillTerms ( const TERMS &  terms,
const Vector b,
const SharedDiagonal noiseModel 
)
protected

Internal function to fill blocks and set dimensions.

Definition at line 60 of file JacobianFactor-inl.h.

const SharedDiagonal& gtsam::JacobianFactor::get_model ( ) const
inline

get a copy of model

Definition at line 289 of file JacobianFactor.h.

SharedDiagonal& gtsam::JacobianFactor::get_model ( )
inline

get a copy of model (non-const version)

Definition at line 292 of file JacobianFactor.h.

constABlock gtsam::JacobianFactor::getA ( const_iterator  variable) const
inline

Get a view of the A matrix for the variable pointed to by the given key iterator

Definition at line 298 of file JacobianFactor.h.

constABlock gtsam::JacobianFactor::getA ( ) const
inline

Get a view of the A matrix, not weighted by noise

Definition at line 301 of file JacobianFactor.h.

ABlock gtsam::JacobianFactor::getA ( iterator  variable)
inline

Get a view of the A matrix for the variable pointed to by the given key iterator (non-const version)

Definition at line 307 of file JacobianFactor.h.

ABlock gtsam::JacobianFactor::getA ( )
inline

Get a view of the A matrix

Definition at line 310 of file JacobianFactor.h.

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

Get a view of the r.h.s. vector b, not weighted by noise

Definition at line 295 of file JacobianFactor.h.

BVector gtsam::JacobianFactor::getb ( )
inline

Get a view of the r.h.s. vector b (non-const version)

Definition at line 304 of file JacobianFactor.h.

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.

Definition at line 274 of file JacobianFactor.h.

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

Compute the gradient wrt a key at any values.

Implements gtsam::GaussianFactor.

Definition at line 726 of file JacobianFactor.cpp.

VectorValues gtsam::JacobianFactor::gradientAtZero ( ) const
overridevirtual

A'*b for Jacobian.

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

Definition at line 709 of file JacobianFactor.cpp.

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

A'*b for Jacobian (raw memory version)

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

Definition at line 721 of file JacobianFactor.cpp.

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

Return the block diagonal of the Hessian for this factor.

Implements gtsam::GaussianFactor.

Definition at line 569 of file JacobianFactor.cpp.

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

Raw memory access version of hessianDiagonal.

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

Definition at line 564 of file JacobianFactor.cpp.

void gtsam::JacobianFactor::hessianDiagonalAdd ( VectorValues d) const
overridevirtual

Add the current diagonal to a VectorValues instance.

Implements gtsam::GaussianFactor.

Definition at line 535 of file JacobianFactor.cpp.

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

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

Implements gtsam::GaussianFactor.

Definition at line 524 of file JacobianFactor.cpp.

bool gtsam::JacobianFactor::isConstrained ( ) const
inline

is noise model constrained ?

Definition at line 267 of file JacobianFactor.h.

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

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

Implements gtsam::GaussianFactor.

Definition at line 733 of file JacobianFactor.cpp.

void gtsam::JacobianFactor::JacobianFactorHelper ( const GaussianFactorGraph graph,
const FastVector< VariableSlots::const_iterator > &  orderedSlots 
)
private

Helper function for public constructors: Build a dense joint factor from all the factors in a factor graph. Takes in ordered variable slots

Definition at line 223 of file JacobianFactor.cpp.

pair< Matrix, Vector > gtsam::JacobianFactor::jacobianUnweighted ( ) const

Returns (dense) A,b pair associated with factor, does not bake in weights.

Definition at line 742 of file JacobianFactor.cpp.

template<class ARCHIVE >
void gtsam::JacobianFactor::load ( ARCHIVE &  ar,
const unsigned int  version 
)
inlineprivate

Definition at line 436 of file JacobianFactor.h.

const VerticalBlockMatrix& gtsam::JacobianFactor::matrixObject ( ) const
inline

Return the full augmented Jacobian matrix of this factor as a VerticalBlockMatrix object.

Definition at line 251 of file JacobianFactor.h.

VerticalBlockMatrix& gtsam::JacobianFactor::matrixObject ( )
inline

Mutable access to the full augmented Jacobian matrix of this factor as a VerticalBlockMatrix object.

Definition at line 254 of file JacobianFactor.h.

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

y += alpha * A'*A*x

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::RegularJacobianFactor< D >.

Definition at line 657 of file JacobianFactor.cpp.

void gtsam::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

Definition at line 672 of file JacobianFactor.cpp.

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

Definition at line 772 of file JacobianFactor.cpp.

Vector gtsam::JacobianFactor::operator* ( const VectorValues x) const

Return A*x

Definition at line 620 of file JacobianFactor.cpp.

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

print

Implements gtsam::GaussianFactor.

Reimplemented in gtsam::LinearInequality, gtsam::LinearEquality, and gtsam::LinearCost.

Definition at line 422 of file JacobianFactor.cpp.

size_t gtsam::JacobianFactor::rows ( void  ) const
inline

return the number of rows in the corresponding linear system

Definition at line 281 of file JacobianFactor.h.

template<class ARCHIVE >
void gtsam::JacobianFactor::save ( ARCHIVE &  ar,
const unsigned int  version 
) const
inlineprivate

Definition at line 418 of file JacobianFactor.h.

void gtsam::JacobianFactor::setModel ( bool  anyConstrained,
const Vector sigmas 
)

set noiseModel correctly

Definition at line 787 of file JacobianFactor.cpp.

GaussianConditional::shared_ptr gtsam::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

Definition at line 838 of file JacobianFactor.cpp.

void gtsam::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.

Definition at line 637 of file JacobianFactor.cpp.

Vector gtsam::JacobianFactor::unweighted_error ( const VectorValues c) const

Definition at line 490 of file JacobianFactor.cpp.

void gtsam::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.

Definition at line 582 of file JacobianFactor.cpp.

JacobianFactor gtsam::JacobianFactor::whiten ( ) const

Return a whitened version of the factor, i.e. with unit diagonal noise model.

Definition at line 762 of file JacobianFactor.cpp.

Friends And Related Function Documentation

friend class boost::serialization::access
friend

Serialization function

Definition at line 416 of file JacobianFactor.h.

GTSAM_EXPORT std::pair<boost::shared_ptr<GaussianConditional>, shared_ptr> EliminateQR ( const GaussianFactorGraph factors,
const Ordering keys 
)
friend

Multiply all factors and eliminate the given keys from the resulting factor using a QR variant that handles constraints (zero sigmas). Computation happens in noiseModel::Gaussian::QR Returns a conditional on those keys, and a new factor on the separator.

Definition at line 797 of file JacobianFactor.cpp.

template<typename T >
friend class ExpressionFactor
friend

Definition at line 413 of file JacobianFactor.h.

Member Data Documentation

VerticalBlockMatrix gtsam::JacobianFactor::Ab_
protected

Definition at line 105 of file JacobianFactor.h.

noiseModel::Diagonal::shared_ptr gtsam::JacobianFactor::model_
protected

Definition at line 106 of file JacobianFactor.h.


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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:58:17