gtsam 4.1.1
gtsam
PartialPriorFactor.h
Go to the documentation of this file.
1/* ----------------------------------------------------------------------------
2
3 * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4 * Atlanta, Georgia 30332-0415
5 * All Rights Reserved
6 * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7
8 * See LICENSE for the license information
9
10 * -------------------------------------------------------------------------- */
11
18#pragma once
19
21#include <gtsam/base/Lie.h>
22
23namespace gtsam {
24
37 template<class VALUE>
39
40 public:
41 typedef VALUE T;
42
43 protected:
44 // Concept checks on the variable type - currently requires Lie
45 GTSAM_CONCEPT_LIE_TYPE(VALUE)
46
49
50 Vector prior_;
51 std::vector<size_t> indices_;
52
55
61 : Base(model, key) {}
62
63 public:
64
65 ~PartialPriorFactor() override {}
66
68 PartialPriorFactor(Key key, size_t idx, double prior, const SharedNoiseModel& model) :
69 Base(model, key),
70 prior_((Vector(1) << prior).finished()),
71 indices_(1, idx) {
72 assert(model->dim() == 1);
73 }
74
76 PartialPriorFactor(Key key, const std::vector<size_t>& indices, const Vector& prior,
77 const SharedNoiseModel& model) :
78 Base(model, key),
79 prior_(prior),
80 indices_(indices) {
81 assert((size_t)prior_.size() == indices_.size());
82 assert(model->dim() == (size_t)prior.size());
83 }
84
86 gtsam::NonlinearFactor::shared_ptr clone() const override {
87 return boost::static_pointer_cast<gtsam::NonlinearFactor>(
88 gtsam::NonlinearFactor::shared_ptr(new This(*this))); }
89
93 void print(const std::string& s, const KeyFormatter& keyFormatter = DefaultKeyFormatter) const override {
94 Base::print(s, keyFormatter);
95 gtsam::print(prior_, "prior");
96 }
97
99 bool equals(const NonlinearFactor& expected, double tol=1e-9) const override {
100 const This *e = dynamic_cast<const This*> (&expected);
101 return e != nullptr && Base::equals(*e, tol) &&
102 gtsam::equal_with_abs_tol(this->prior_, e->prior_, tol) &&
103 this->indices_ == e->indices_;
104 }
105
109 Vector evaluateError(const T& p, boost::optional<Matrix&> H = boost::none) const override {
110 Eigen::Matrix<double, T::dimension, T::dimension> H_local;
111
112 // If the Rot3 Cayley map is used, Rot3::LocalCoordinates will throw a runtime error
113 // when asked to compute the Jacobian matrix (see Rot3M.cpp).
114 #ifdef GTSAM_ROT3_EXPMAP
115 const Vector full_tangent = T::LocalCoordinates(p, H ? &H_local : nullptr);
116 #else
117 const Vector full_tangent = T::Logmap(p, H ? &H_local : nullptr);
118 #endif
119
120 if (H) {
121 (*H) = Matrix::Zero(indices_.size(), T::dimension);
122 for (size_t i = 0; i < indices_.size(); ++i) {
123 (*H).row(i) = H_local.row(indices_.at(i));
124 }
125 }
126 // Select relevant parameters from the tangent vector.
127 Vector partial_tangent = Vector::Zero(indices_.size());
128 for (size_t i = 0; i < indices_.size(); ++i) {
129 partial_tangent(i) = full_tangent(indices_.at(i));
130 }
131
132 return partial_tangent - prior_;
133 }
134
135 // access
136 const Vector& prior() const { return prior_; }
137 const std::vector<size_t>& indices() const { return indices_; }
138
139 private:
142 template<class ARCHIVE>
143 void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
144 ar & boost::serialization::make_nvp("NoiseModelFactor1",
145 boost::serialization::base_object<Base>(*this));
146 ar & BOOST_SERIALIZATION_NVP(prior_);
147 ar & BOOST_SERIALIZATION_NVP(indices_);
148 // ar & BOOST_SERIALIZATION_NVP(H_);
149 }
150 }; // \class PartialPriorFactor
151
152}
Base class and basic functions for Lie types.
Non-linear factor base classes.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
void print(const Matrix &A, const string &s, ostream &stream)
print without optional string, must specify cout yourself
Definition: Matrix.cpp:155
noiseModel::Base::shared_ptr SharedNoiseModel
Note, deliberately not in noiseModel namespace.
Definition: NoiseModel.h:736
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
equals with a tolerance
Definition: Matrix.h:84
This is the base class for all factor types.
Definition: Factor.h:56
Nonlinear factor base class.
Definition: NonlinearFactor.h:43
bool equals(const NonlinearFactor &f, double tol=1e-9) const override
Check if two factors are equal.
Definition: NonlinearFactor.cpp:71
void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
Print.
Definition: NonlinearFactor.cpp:63
A convenient base class for creating your own NoiseModelFactor with 1 variable.
Definition: NonlinearFactor.h:285
A class for a soft partial prior on any Lie type, with a mask over Expmap parameters.
Definition: PartialPriorFactor.h:38
PartialPriorFactor(Key key, const SharedNoiseModel &model)
constructor with just minimum requirements for a factor - allows more computation in the constructor.
Definition: PartialPriorFactor.h:60
gtsam::NonlinearFactor::shared_ptr clone() const override
Definition: PartialPriorFactor.h:86
PartialPriorFactor()
default constructor - only use for serialization
Definition: PartialPriorFactor.h:54
PartialPriorFactor(Key key, const std::vector< size_t > &indices, const Vector &prior, const SharedNoiseModel &model)
Indices Constructor: Specify the relevant measured indices in the tangent vector.
Definition: PartialPriorFactor.h:76
Vector prior_
Measurement on tangent space parameters, in compressed form.
Definition: PartialPriorFactor.h:50
Vector evaluateError(const T &p, boost::optional< Matrix & > H=boost::none) const override
implement functions needed to derive from Factor
Definition: PartialPriorFactor.h:109
bool equals(const NonlinearFactor &expected, double tol=1e-9) const override
equals
Definition: PartialPriorFactor.h:99
void print(const std::string &s, const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
implement functions needed for Testable
Definition: PartialPriorFactor.h:93
friend class boost::serialization::access
Serialization function.
Definition: PartialPriorFactor.h:141
PartialPriorFactor(Key key, size_t idx, double prior, const SharedNoiseModel &model)
Single Element Constructor: Prior on a single parameter at index 'idx' in the tangent vector.
Definition: PartialPriorFactor.h:68
std::vector< size_t > indices_
Indices of the measured tangent space parameters.
Definition: PartialPriorFactor.h:51