gtsam 4.2.0
gtsam
Loading...
Searching...
No Matches
utilities.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
33
34#include <exception>
35
36namespace gtsam {
37
38namespace utilities {
39
40// Create a KeyList from indices
41FastList<Key> createKeyList(const Vector& I) {
42 FastList<Key> set;
43 for (int i = 0; i < I.size(); i++)
44 set.push_back(I[i]);
45 return set;
46}
47
48// Create a KeyList from indices using symbol
49FastList<Key> createKeyList(std::string s, const Vector& I) {
50 FastList<Key> set;
51 char c = s[0];
52 for (int i = 0; i < I.size(); i++)
53 set.push_back(Symbol(c, I[i]));
54 return set;
55}
56
57// Create a KeyVector from indices
58KeyVector createKeyVector(const Vector& I) {
59 KeyVector set;
60 for (int i = 0; i < I.size(); i++)
61 set.push_back(I[i]);
62 return set;
63}
64
65// Create a KeyVector from indices using symbol
66KeyVector createKeyVector(std::string s, const Vector& I) {
67 KeyVector set;
68 char c = s[0];
69 for (int i = 0; i < I.size(); i++)
70 set.push_back(Symbol(c, I[i]));
71 return set;
72}
73
74// Create a KeySet from indices
75KeySet createKeySet(const Vector& I) {
76 KeySet set;
77 for (int i = 0; i < I.size(); i++)
78 set.insert(I[i]);
79 return set;
80}
81
82// Create a KeySet from indices using symbol
83KeySet createKeySet(std::string s, const Vector& I) {
84 KeySet set;
85 char c = s[0];
86 for (int i = 0; i < I.size(); i++)
87 set.insert(symbol(c, I[i]));
88 return set;
89}
90
92Matrix extractPoint2(const Values& values) {
93 const auto points = values.extract<gtsam::Point2>();
94 // Point2 is aliased as a gtsam::Vector in the wrapper
95 const auto points2 = values.extract<gtsam::Vector>();
96
97 Matrix result(points.size() + points2.size(), 2);
98
99 size_t j = 0;
100 for (const auto& key_value : points) {
101 result.row(j++) = key_value.second;
102 }
103 for (const auto& key_value : points2) {
104 if (key_value.second.rows() == 2) {
105 result.row(j++) = key_value.second;
106 }
107 }
108 return result;
109}
110
112Matrix extractPoint3(const Values& values) {
113 const auto points = values.extract<gtsam::Point3>();
114 // Point3 is aliased as a gtsam::Vector in the wrapper
115 const auto points2 = values.extract<gtsam::Vector>();
116
117 Matrix result(points.size() + points2.size(), 3);
118
119 size_t j = 0;
120 for (const auto& key_value : points) {
121 result.row(j++) = key_value.second;
122 }
123 for (const auto& key_value : points2) {
124 if (key_value.second.rows() == 3) {
125 result.row(j++) = key_value.second;
126 }
127 }
128 return result;
129}
130
132Values allPose2s(const Values& values) {
133 Values result;
134 for(const auto& key_value: values.extract<Pose2>())
135 result.insert(key_value.first, key_value.second);
136 return result;
137}
138
140Matrix extractPose2(const Values& values) {
141 const auto poses = values.extract<Pose2>();
142 Matrix result(poses.size(), 3);
143 size_t j = 0;
144 for(const auto& key_value: poses)
145 result.row(j++) << key_value.second.x(), key_value.second.y(), key_value.second.theta();
146 return result;
147}
148
150Values allPose3s(const Values& values) {
151 Values result;
152 for(const auto& key_value: values.extract<Pose3>())
153 result.insert(key_value.first, key_value.second);
154 return result;
155}
156
158Matrix extractPose3(const Values& values) {
159 const auto poses = values.extract<Pose3>();
160 Matrix result(poses.size(), 12);
161 size_t j = 0;
162 for(const auto& key_value: poses) {
163 result.row(j).segment(0, 3) << key_value.second.rotation().matrix().row(0);
164 result.row(j).segment(3, 3) << key_value.second.rotation().matrix().row(1);
165 result.row(j).segment(6, 3) << key_value.second.rotation().matrix().row(2);
166 result.row(j).tail(3) = key_value.second.translation();
167 j++;
168 }
169 return result;
170}
171
180Matrix extractVectors(const Values& values, char c) {
181 const auto vectors = values.extract<Vector>(Symbol::ChrTest(c));
182 if (vectors.size() == 0) {
183 return Matrix();
184 }
185 auto dim = vectors.begin()->second.size();
186 Matrix result(vectors.size(), dim);
187 Eigen::Index rowi = 0;
188 for (const auto& kv : vectors) {
189 if (kv.second.size() != dim) {
190 throw std::runtime_error(
191 "Tried to extract different-sized vectors into a single matrix");
192 }
193 result.row(rowi) = kv.second;
194 ++rowi;
195 }
196 return result;
197}
198
200void perturbPoint2(Values& values, double sigma, int32_t seed = 42u) {
201 noiseModel::Isotropic::shared_ptr model =
203 Sampler sampler(model, seed);
204 for (const auto& key_value : values.extract<Point2>()) {
205 values.update<Point2>(key_value.first,
206 key_value.second + Point2(sampler.sample()));
207 }
208 for (const auto& key_value : values.extract<gtsam::Vector>()) {
209 if (key_value.second.rows() == 2) {
210 values.update<gtsam::Vector>(key_value.first,
211 key_value.second + Point2(sampler.sample()));
212 }
213 }
214}
215
217void perturbPose2(Values& values, double sigmaT, double sigmaR, int32_t seed =
218 42u) {
219 noiseModel::Diagonal::shared_ptr model = noiseModel::Diagonal::Sigmas(
220 Vector3(sigmaT, sigmaT, sigmaR));
221 Sampler sampler(model, seed);
222 for(const auto& key_value: values.extract<Pose2>()) {
223 values.update<Pose2>(key_value.first, key_value.second.retract(sampler.sample()));
224 }
225}
226
228void perturbPoint3(Values& values, double sigma, int32_t seed = 42u) {
229 noiseModel::Isotropic::shared_ptr model =
231 Sampler sampler(model, seed);
232 for (const auto& key_value : values.extract<Point3>()) {
233 values.update<Point3>(key_value.first,
234 key_value.second + Point3(sampler.sample()));
235 }
236 for (const auto& key_value : values.extract<gtsam::Vector>()) {
237 if (key_value.second.rows() == 3) {
238 values.update<gtsam::Vector>(key_value.first,
239 key_value.second + Point3(sampler.sample()));
240 }
241 }
242}
243
254 const Vector& J, const Matrix& Z, double depth) {
255 if (Z.rows() != 2)
256 throw std::invalid_argument("insertBackProjections: Z must be 2*J");
257 if (Z.cols() != J.size())
258 throw std::invalid_argument(
259 "insertBackProjections: J and Z must have same number of entries");
260 for (int k = 0; k < Z.cols(); k++) {
261 Point2 p(Z(0, k), Z(1, k));
262 Point3 P = camera.backproject(p, depth);
263 values.insert(J(k), P);
264 }
265}
266
279 const Vector& J, const Matrix& Z, const SharedNoiseModel& model,
280 const Cal3_S2::shared_ptr K, const Pose3& body_P_sensor = Pose3()) {
281 if (Z.rows() != 2)
282 throw std::invalid_argument("addMeasurements: Z must be 2*K");
283 if (Z.cols() != J.size())
284 throw std::invalid_argument(
285 "addMeasurements: J and Z must have same number of entries");
286 for (int k = 0; k < Z.cols(); k++) {
287 graph.push_back(
288 boost::make_shared<GenericProjectionFactor<Pose3, Point3> >(
289 Point2(Z(0, k), Z(1, k)), model, i, Key(J(k)), K, body_P_sensor));
290 }
291}
292
295 const Values& values) {
296 // first count
297 size_t K = 0, k = 0;
298 for(const NonlinearFactor::shared_ptr& f: graph)
299 if (boost::dynamic_pointer_cast<const GenericProjectionFactor<Pose3, Point3> >(
300 f))
301 ++K;
302 // now fill
303 Matrix errors(2, K);
304 for(const NonlinearFactor::shared_ptr& f: graph) {
305 boost::shared_ptr<const GenericProjectionFactor<Pose3, Point3> > p =
306 boost::dynamic_pointer_cast<const GenericProjectionFactor<Pose3, Point3> >(
307 f);
308 if (p)
309 errors.col(k++) = p->unwhitenedError(values);
310 }
311 return errors;
312}
313
315Values localToWorld(const Values& local, const Pose2& base,
316 const KeyVector user_keys = KeyVector()) {
317
318 Values world;
319
320 // if no keys given, get all keys from local values
321 KeyVector keys(user_keys);
322 if (keys.size()==0)
323 keys = local.keys();
324
325 // Loop over all keys
326 for(Key key: keys) {
327 try {
328 // if value is a Pose2, compose it with base pose
329 Pose2 pose = local.at<Pose2>(key);
330 world.insert(key, base.compose(pose));
331 } catch (const std::exception& e1) {
332 try {
333 // if value is a Point2, transform it from base pose
334 Point2 point = local.at<Point2>(key);
335 world.insert(key, base.transformFrom(point));
336 } catch (const std::exception& e2) {
337 // if not Pose2 or Point2, do nothing
338 #ifndef NDEBUG
339 std::cerr << "Values[key] is neither Pose2 nor Point2, so skip" << std::endl;
340 #endif
341 }
342 }
343 }
344 return world;
345}
346
347} // namespace utilities
348
349}
350
3D Point
The most common 5DOF 3D->2D calibration.
3D Pose
2D Pose
Base class for all pinhole cameras.
2D Point
Factor Graph Values.
sampling from a NoiseModel
Non-linear factor base classes.
Factor Graph consisting of non-linear factors.
A non-templated config holding any types of Manifold-group elements.
Reprojection of a LANDMARK to a 2D point.
Global functions in a separate testing namespace.
Definition chartTesting.h:28
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition Key.h:86
Key symbol(unsigned char c, std::uint64_t j)
Create a symbol key from a character and index, i.e.
Definition Symbol.h:135
Vector2 Point2
As of GTSAM 4, in order to make GTSAM more lean, it is now possible to just typedef Point2 to Vector2...
Definition Point2.h:27
Vector3 Point3
As of GTSAM 4, in order to make GTSAM more lean, it is now possible to just typedef Point3 to Vector3...
Definition Point3.h:36
noiseModel::Base::shared_ptr SharedNoiseModel
Aliases.
Definition NoiseModel.h:724
std::uint64_t Key
Integer nonlinear key type.
Definition types.h:100
A pinhole camera class that has a Pose3 and a Calibration.
Definition PinholeCamera.h:33
Point3 backproject(const Point2 &p, double depth, OptionalJacobian< 3, 6 > Dresult_dpose=boost::none, OptionalJacobian< 3, 2 > Dresult_dp=boost::none, OptionalJacobian< 3, 1 > Dresult_ddepth=boost::none, OptionalJacobian< 3, DimK > Dresult_dcal=boost::none) const
backproject a 2-dimensional point to a 3-dimensional point at given depth
Definition PinholePose.h:139
A 2D pose (Point2,Rot2)
Definition Pose2.h:36
GTSAM_EXPORT Point2 transformFrom(const Point2 &point, OptionalJacobian< 2, 3 > Dpose=boost::none, OptionalJacobian< 2, 2 > Dpoint=boost::none) const
Return point coordinates in global frame.
Definition Pose2.cpp:226
A 3D pose (R,t) : (Rot3,Point3)
Definition Pose3.h:37
static std::function< bool(Key)> ChrTest(unsigned char c)
Return a filter function that returns true when evaluated on a Key whose character (when converted to...
Definition Symbol.cpp:65
static shared_ptr Sigmas(const Vector &sigmas, bool smart=true)
A diagonal noise model created by specifying a Vector of sigmas, i.e.
Definition NoiseModel.cpp:273
static shared_ptr Sigma(size_t dim, double sigma, bool smart=true)
An isotropic noise model created by specifying a standard devation sigma.
Definition NoiseModel.cpp:597
Sampling structure that keeps internal random number generators for diagonal distributions specified ...
Definition Sampler.h:31
Vector sample() const
sample from distribution
Definition Sampler.cpp:59
Definition NonlinearFactorGraph.h:55
A non-templated config holding any types of Manifold-group elements.
Definition Values.h:65
void update(Key j, const Value &val)
single element change of existing element
Definition Values.cpp:180
std::map< Key, ValueType > extract(const std::function< bool(Key)> &filterFcn=&_truePredicate< Key >) const
Extract a subset of values of the given type ValueType.
Definition Values-inl.h:94
const ValueType at(Key j) const
Retrieve a variable by key j.
Definition Values-inl.h:361
void insert(Key j, const Value &val)
Add a variable with the given j, throws KeyAlreadyExists<J> if j is already present.
Definition Values.cpp:157
KeyVector keys() const
Returns a vector of keys in the config.
Definition Values.cpp:230
Non-linear factor for a constraint derived from a 2D measurement.
Definition ProjectionFactor.h:40
Matrix extractPoint2(const Values &values)
Extract all Point2 values into a single matrix [x y].
Definition utilities.h:92
Values localToWorld(const Values &local, const Pose2 &base, const KeyVector user_keys=KeyVector())
Convert from local to world coordinates.
Definition utilities.h:315
Matrix extractVectors(const Values &values, char c)
Extract all Vector values with a given symbol character into an mxn matrix, where m is the number of ...
Definition utilities.h:180
Matrix extractPoint3(const Values &values)
Extract all Point3 values into a single matrix [x y z].
Definition utilities.h:112
Matrix reprojectionErrors(const NonlinearFactorGraph &graph, const Values &values)
Calculate the errors of all projection factors in a graph.
Definition utilities.h:294
void perturbPoint2(Values &values, double sigma, int32_t seed=42u)
Perturb all Point2 values using normally distributed noise.
Definition utilities.h:200
void perturbPose2(Values &values, double sigmaT, double sigmaR, int32_t seed=42u)
Perturb all Pose2 values using normally distributed noise.
Definition utilities.h:217
Values allPose2s(const Values &values)
Extract all Pose3 values.
Definition utilities.h:132
void insertProjectionFactors(NonlinearFactorGraph &graph, Key i, const Vector &J, const Matrix &Z, const SharedNoiseModel &model, const Cal3_S2::shared_ptr K, const Pose3 &body_P_sensor=Pose3())
Insert multiple projection factors for a single pose key.
Definition utilities.h:278
Values allPose3s(const Values &values)
Extract all Pose3 values.
Definition utilities.h:150
void insertBackprojections(Values &values, const PinholeCamera< Cal3_S2 > &camera, const Vector &J, const Matrix &Z, double depth)
Insert a number of initial point values by backprojecting.
Definition utilities.h:253
Matrix extractPose3(const Values &values)
Extract all Pose3 values into a single matrix [r11 r12 r13 r21 r22 r23 r31 r32 r33 x y z].
Definition utilities.h:158
Matrix extractPose2(const Values &values)
Extract all Pose2 values into a single matrix [x y theta].
Definition utilities.h:140
void perturbPoint3(Values &values, double sigma, int32_t seed=42u)
Perturb all Point3 values using normally distributed noise.
Definition utilities.h:228