gtsam 4.1.1
gtsam
|
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) . \]
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. More... | |
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. More... | |
JacobianFactor (const GaussianFactorGraph &graph) | |
Build a dense joint factor from all the factors in a factor graph. More... | |
JacobianFactor (const GaussianFactorGraph &graph, const VariableSlots &p_variableSlots) | |
Build a dense joint factor from all the factors in a factor graph. More... | |
JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering) | |
Build a dense joint factor from all the factors in a factor graph. More... | |
JacobianFactor (const GaussianFactorGraph &graph, const Ordering &ordering, const VariableSlots &p_variableSlots) | |
Build a dense joint factor from all the factors in a factor graph. More... | |
~JacobianFactor () override | |
Virtual destructor. | |
GaussianFactor::shared_ptr | clone () const override |
Clone this JacobianFactor. More... | |
void | print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override |
print More... | |
bool | equals (const GaussianFactor &lf, double tol=1e-9) const override |
Equals for testable. More... | |
Vector | unweighted_error (const VectorValues &c) const |
Vector | error_vector (const VectorValues &c) const |
(A*x-b) | |
double | error (const VectorValues &c) const override |
(A*x-b)/sigma More... | |
Matrix | augmentedInformation () const override |
0.5*(A*x-b)'D(A*x-b) More... | |
Matrix | information () const override |
Return the non-augmented information matrix represented by this GaussianFactor. More... | |
void | hessianDiagonalAdd (VectorValues &d) const override |
Add the current diagonal to a VectorValues instance. More... | |
void | hessianDiagonal (double *d) const override |
Raw memory access version of hessianDiagonal. More... | |
std::map< Key, Matrix > | hessianBlockDiagonal () const override |
Return the block diagonal of the Hessian for this factor. More... | |
std::pair< Matrix, Vector > | jacobian () const override |
Returns (dense) A,b pair associated with factor, bakes in the weights. More... | |
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. More... | |
Matrix | augmentedJacobianUnweighted () const |
Return (dense) matrix associated with factor. More... | |
const VerticalBlockMatrix & | matrixObject () const |
Return the full augmented Jacobian matrix of this factor as a VerticalBlockMatrix object. | |
VerticalBlockMatrix & | matrixObject () |
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. More... | |
bool | empty () const override |
Check if the factor is empty. More... | |
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? More... | |
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. More... | |
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. More... | |
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). More... | |
Vector | operator* (const VectorValues &x) const |
Return A*x. | |
void | transposeMultiplyAdd (double alpha, const Vector &e, VectorValues &x) const |
x += alpha * A'*e. More... | |
void | multiplyHessianAdd (double alpha, const VectorValues &x, VectorValues &y) const override |
y += alpha * A'*A*x More... | |
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. More... | |
VectorValues | gradientAtZero () const override |
A'*b for Jacobian. More... | |
void | gradientAtZero (double *d) const override |
A'*b for Jacobian (raw memory version) More... | |
Vector | gradient (Key key, const VectorValues &x) const override |
Compute the gradient wrt a key at any values. More... | |
JacobianFactor | whiten () const |
Return a whitened version of the factor, i.e. More... | |
std::pair< boost::shared_ptr< GaussianConditional >, shared_ptr > | eliminate (const Ordering &keys) |
Eliminate the requested variables. | |
void | setModel (bool anyConstrained, const Vector &sigmas) |
set noiseModel correctly | |
boost::shared_ptr< GaussianConditional > | splitConditional (size_t nrFrontals) |
splits a pre-factorized factor into a conditional, and changes the current factor to be the remaining component. More... | |
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. More... | |
virtual | ~GaussianFactor () |
Destructor. | |
void | print (const std::string &s="", const KeyFormatter &formatter=DefaultKeyFormatter) const override=0 |
print More... | |
virtual bool | equals (const GaussianFactor &lf, double tol=1e-9) const =0 |
Equals for testable. More... | |
virtual double | error (const VectorValues &c) const =0 |
Print for testable. More... | |
virtual DenseIndex | getDim (const_iterator variable) const =0 |
0.5*(A*x-b)'D(A*x-b) More... | |
virtual Matrix | augmentedJacobian () const =0 |
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. More... | |
virtual std::pair< Matrix, Vector > | jacobian () const =0 |
Return the dense Jacobian \( A \) and right-hand-side \( b \), with the noise models baked into A and b. More... | |
virtual Matrix | augmentedInformation () const =0 |
Return the augmented information matrix represented by this GaussianFactor. More... | |
virtual Matrix | information () const =0 |
Return the non-augmented information matrix represented by this GaussianFactor. More... | |
VectorValues | hessianDiagonal () const |
Return the diagonal of the Hessian for this factor. | |
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, Matrix > | hessianBlockDiagonal () const =0 |
Return the block diagonal of the Hessian for this factor. More... | |
virtual GaussianFactor::shared_ptr | clone () const =0 |
Clone a factor (make a deep copy) More... | |
virtual bool | empty () const =0 |
Test whether the factor is empty. More... | |
virtual GaussianFactor::shared_ptr | negate () const =0 |
Construct the corresponding anti-factor to negate information stored stored in this factor. More... | |
virtual void | updateHessian (const KeyVector &keys, SymmetricBlockMatrix *info) const =0 |
Update an information matrix by adding the information corresponding to this factor (used internally during elimination). More... | |
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. | |
KeyVector & | keys () |
iterator | begin () |
Iterator at beginning of involved variable keys. | |
iterator | end () |
Iterator at end of involved variable keys. | |
virtual void | printKeys (const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const |
print only keys More... | |
Key | front () const |
First key. | |
Key | back () const |
Last key. | |
const_iterator | find (Key key) const |
find | |
const KeyVector & | keys () 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 |
Public Types | |
typedef JacobianFactor | This |
Typedef to this class. | |
typedef GaussianFactor | Base |
Typedef to base class. | |
typedef boost::shared_ptr< This > | shared_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< This > | shared_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. More... | |
template<typename ITERATOR > | |
Factor (ITERATOR first, ITERATOR last) | |
Construct factor from iterator keys. More... | |
bool | equals (const This &other, double tol=1e-9) const |
check equality | |
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_ptr > | EliminateQR (const GaussianFactorGraph &factors, const Ordering &keys) |
Multiply all factors and eliminate the given keys from the resulting factor using a QR variant that handles constraints (zero sigmas). More... | |
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. More... | |
template<typename ITERATOR > | |
static Factor | FromIterators (ITERATOR first, ITERATOR last) |
Construct factor from iterator keys. More... | |
gtsam::JacobianFactor::JacobianFactor | ( | const TERMS & | terms, |
const Vector & | b, | ||
const SharedDiagonal & | model = SharedDiagonal() |
||
) |
Construct an n-ary factor.
TERMS | A container whose value type is std::pair<Key, Matrix>, specifying the collection of keys and matrices making up the factor. |
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
overridevirtual |
Return (dense) matrix associated with factor.
The returned system is an augmented matrix: [A b] weights are baked in
Implements gtsam::GaussianFactor.
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
|
inlineoverridevirtual |
Clone this JacobianFactor.
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::LinearCost, gtsam::LinearEquality, and gtsam::LinearInequality.
|
inlineoverridevirtual |
|
overridevirtual |
Equals for testable.
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::LinearCost, gtsam::LinearEquality, and gtsam::LinearInequality.
|
overridevirtual |
(A*x-b)/sigma
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::LinearCost, gtsam::LinearEquality, and gtsam::LinearInequality.
|
inline |
Get a view of the r.h.s.
vector b (non-const version)
|
inline |
Get a view of the r.h.s.
vector b, not weighted by noise
|
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.
|
overridevirtual |
Compute the gradient wrt a key at any values.
Implements gtsam::GaussianFactor.
|
overridevirtual |
A'*b for Jacobian.
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::RegularJacobianFactor< D >.
|
overridevirtual |
A'*b for Jacobian (raw memory version)
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::RegularJacobianFactor< D >.
|
overridevirtual |
Return the block diagonal of the Hessian for this factor.
Implements gtsam::GaussianFactor.
|
overridevirtual |
Raw memory access version of hessianDiagonal.
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::RegularJacobianFactor< D >, and gtsam::RegularJacobianFactor< D >.
|
overridevirtual |
Add the current diagonal to a VectorValues instance.
Implements gtsam::GaussianFactor.
|
overridevirtual |
Return the non-augmented information matrix represented by this GaussianFactor.
Implements gtsam::GaussianFactor.
|
overridevirtual |
Returns (dense) A,b pair associated with factor, bakes in the weights.
Implements gtsam::GaussianFactor.
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.
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): we should probably kill this if no longer needed
: 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
|
overridevirtual |
y += alpha * A'*A*x
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::RegularJacobianFactor< D >, and gtsam::RegularJacobianFactor< D >.
|
overridevirtual |
Construct the corresponding anti-factor to negate information stored stored in this factor.
Implements gtsam::GaussianFactor.
|
overridevirtual |
Implements gtsam::GaussianFactor.
Reimplemented in gtsam::LinearCost, gtsam::LinearEquality, and gtsam::LinearInequality.
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.
nrFrontals | number of keys to eliminate |
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.
|
overridevirtual |
Update an information matrix by adding the information corresponding to this factor (used internally during elimination).
scatter | A mapping from variable index to slot index in this HessianFactor |
info | The information matrix to be updated |
Implements gtsam::GaussianFactor.
JacobianFactor gtsam::JacobianFactor::whiten | ( | ) | const |
Return a whitened version of the factor, i.e.
with unit diagonal noise model.
|
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.