#include <LinearContainerFactor.h>
Public Types | |
typedef std::shared_ptr< This > | shared_ptr |
Public Types inherited from gtsam::NonlinearFactor | |
typedef std::shared_ptr< This > | shared_ptr |
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 | |
NonlinearFactor::shared_ptr | clone () const override |
size_t | dim () const override |
bool | equals (const NonlinearFactor &f, double tol=1e-9) const override |
double | error (const Values &c) const override |
const GaussianFactor::shared_ptr & | factor () const |
bool | hasLinearizationPoint () const |
Casting syntactic sugar. More... | |
bool | isHessian () const |
bool | isJacobian () const |
LinearContainerFactor () | |
LinearContainerFactor (const GaussianFactor::shared_ptr &factor, const Values &linearizationPoint=Values()) | |
LinearContainerFactor (const HessianFactor &factor, const Values &linearizationPoint=Values()) | |
LinearContainerFactor (const JacobianFactor &factor, const Values &linearizationPoint=Values()) | |
const std::optional< Values > & | linearizationPoint () const |
GaussianFactor::shared_ptr | linearize (const Values &c) const override |
GaussianFactor::shared_ptr | negateToGaussian () const |
NonlinearFactor::shared_ptr | negateToNonlinear () const |
void | print (const std::string &s="", const KeyFormatter &keyFormatter=gtsam::DefaultKeyFormatter) const override |
NonlinearFactor::shared_ptr | rekey (const KeyVector &new_keys) const override |
NonlinearFactor::shared_ptr | rekey (const std::map< Key, Key > &rekey_mapping) const override |
std::shared_ptr< HessianFactor > | toHessian () const |
std::shared_ptr< JacobianFactor > | toJacobian () const |
Public Member Functions inherited from gtsam::NonlinearFactor | |
NonlinearFactor () | |
template<typename CONTAINER > | |
NonlinearFactor (const CONTAINER &keys) | |
void | print (const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override |
double | error (const HybridValues &c) const override |
virtual bool | active (const Values &) const |
virtual bool | sendable () const |
Public Member Functions inherited from gtsam::Factor | |
virtual | ~Factor ()=default |
Default destructor. More... | |
bool | empty () const |
Whether the factor is empty (involves zero variables). More... | |
Key | front () const |
First key. More... | |
Key | back () const |
Last key. More... | |
const_iterator | find (Key key) const |
find More... | |
const KeyVector & | keys () const |
Access the factor's involved variable keys. More... | |
const_iterator | begin () const |
const_iterator | end () const |
size_t | size () const |
virtual void | printKeys (const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const |
print only keys More... | |
bool | equals (const This &other, double tol=1e-9) const |
check equality More... | |
KeyVector & | keys () |
iterator | begin () |
iterator | end () |
Static Public Member Functions | |
static NonlinearFactorGraph | ConvertLinearGraph (const GaussianFactorGraph &linear_graph, const Values &linearizationPoint=Values()) |
Protected Types | |
typedef NonlinearFactor | Base |
typedef LinearContainerFactor | This |
Protected Types inherited from gtsam::NonlinearFactor | |
typedef Factor | Base |
typedef NonlinearFactor | This |
Protected Member Functions | |
void | initializeLinearizationPoint (const Values &linearizationPoint) |
LinearContainerFactor (const GaussianFactor::shared_ptr &factor, const std::optional< Values > &linearizationPoint) | |
Protected Member Functions inherited from gtsam::Factor | |
Factor () | |
template<typename CONTAINER > | |
Factor (const CONTAINER &keys) | |
template<typename ITERATOR > | |
Factor (ITERATOR first, ITERATOR last) | |
Protected Attributes | |
GaussianFactor::shared_ptr | factor_ |
std::optional< Values > | linearizationPoint_ |
Protected Attributes inherited from gtsam::Factor | |
KeyVector | keys_ |
The keys involved in this factor. More... | |
Additional Inherited Members | |
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) |
Dummy version of a generic linear factor to be injected into a nonlinear factor graph
This factor does have the ability to perform relinearization under small-angle and linearity assumptions if a linearization point is added.
Definition at line 29 of file LinearContainerFactor.h.
|
protected |
Definition at line 39 of file LinearContainerFactor.h.
typedef std::shared_ptr<This> gtsam::LinearContainerFactor::shared_ptr |
Definition at line 44 of file LinearContainerFactor.h.
|
protected |
Definition at line 40 of file LinearContainerFactor.h.
|
protected |
direct copy constructor
Definition at line 29 of file LinearContainerFactor.cpp.
|
inline |
Default constructor - necessary for serialization
Definition at line 47 of file LinearContainerFactor.h.
gtsam::LinearContainerFactor::LinearContainerFactor | ( | const JacobianFactor & | factor, |
const Values & | linearizationPoint = Values() |
||
) |
Primary constructor: store a linear factor with optional linearization point
Definition at line 35 of file LinearContainerFactor.cpp.
gtsam::LinearContainerFactor::LinearContainerFactor | ( | const HessianFactor & | factor, |
const Values & | linearizationPoint = Values() |
||
) |
Primary constructor: store a linear factor with optional linearization point
Definition at line 42 of file LinearContainerFactor.cpp.
gtsam::LinearContainerFactor::LinearContainerFactor | ( | const GaussianFactor::shared_ptr & | factor, |
const Values & | linearizationPoint = Values() |
||
) |
Constructor from shared_ptr
Definition at line 49 of file LinearContainerFactor.cpp.
|
inlineoverridevirtual |
Creates a shared_ptr clone of the factor - needs to be specialized to allow for subclasses
Clones the underlying linear factor
Reimplemented from gtsam::NonlinearFactor.
Definition at line 122 of file LinearContainerFactor.h.
|
static |
Utility function for converting linear graphs to nonlinear graphs consisting of LinearContainerFactors.
Definition at line 210 of file LinearContainerFactor.cpp.
|
overridevirtual |
get the dimension of the factor: rows of linear factor
Implements gtsam::NonlinearFactor.
Definition at line 96 of file LinearContainerFactor.cpp.
|
overridevirtual |
Check if two factors are equal
Reimplemented from gtsam::NonlinearFactor.
Definition at line 65 of file LinearContainerFactor.cpp.
|
overridevirtual |
Calculate the nonlinear error for the factor, where the error is computed by passing the delta between linearization point and c, where delta = linearizationPoint_.localCoordinates(c), into the error function of the stored linear factor.
Reimplemented from gtsam::NonlinearFactor.
Definition at line 77 of file LinearContainerFactor.cpp.
|
inline |
Definition at line 60 of file LinearContainerFactor.h.
|
inline |
Casting syntactic sugar.
Definition at line 141 of file LinearContainerFactor.h.
|
protected |
Definition at line 17 of file LinearContainerFactor.cpp.
bool gtsam::LinearContainerFactor::isHessian | ( | ) | const |
Definition at line 141 of file LinearContainerFactor.cpp.
bool gtsam::LinearContainerFactor::isJacobian | ( | ) | const |
Simple checks whether this is a Jacobian or Hessian factor
Definition at line 136 of file LinearContainerFactor.cpp.
|
inline |
Extract the linearization point used in recalculating error
Definition at line 86 of file LinearContainerFactor.h.
|
overridevirtual |
Linearize to a GaussianFactor, with method depending on the presence of a linearizationPoint
The relinearization approach used computes a linear delta between the original linearization point and the new values c, where delta = linearizationPoint_.localCoordinates(c), and substitutes this change into the system. This substitution is only really valid for linear variable manifolds, and for any variables based on a non-commutative manifold (such as Pose2, Pose3), the relinearized version will be effective for only small angles.
TODO: better approximation of relinearization TODO: switchable modes for approximation technique
Implements gtsam::NonlinearFactor.
Definition at line 104 of file LinearContainerFactor.cpp.
GaussianFactor::shared_ptr gtsam::LinearContainerFactor::negateToGaussian | ( | ) | const |
Creates an anti-factor directly
Definition at line 156 of file LinearContainerFactor.cpp.
NonlinearFactor::shared_ptr gtsam::LinearContainerFactor::negateToNonlinear | ( | ) | const |
Creates the equivalent anti-factor as another LinearContainerFactor.
Definition at line 162 of file LinearContainerFactor.cpp.
|
overridevirtual |
|
overridevirtual |
Clones a factor and fully replaces its keys
new_keys | is the full replacement set of keys |
Reimplemented from gtsam::NonlinearFactor.
Definition at line 190 of file LinearContainerFactor.cpp.
|
overridevirtual |
Creates a shared_ptr clone of the factor with different keys using a map from old->new keys
Reimplemented from gtsam::NonlinearFactor.
Definition at line 168 of file LinearContainerFactor.cpp.
HessianFactor::shared_ptr gtsam::LinearContainerFactor::toHessian | ( | ) | const |
Casts to HessianFactor
Definition at line 151 of file LinearContainerFactor.cpp.
JacobianFactor::shared_ptr gtsam::LinearContainerFactor::toJacobian | ( | ) | const |
Casts to JacobianFactor
Definition at line 146 of file LinearContainerFactor.cpp.
|
protected |
Definition at line 32 of file LinearContainerFactor.h.
|
protected |
Definition at line 33 of file LinearContainerFactor.h.