testDataset.cpp
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 
19 #include <gtsam/slam/dataset.h>
21 #include <gtsam/geometry/Pose2.h>
22 #include <gtsam/inference/Symbol.h>
24 
26 
27 #include <iostream>
28 #include <sstream>
29 
30 using namespace gtsam::symbol_shorthand;
31 using namespace std;
32 using namespace gtsam;
33 
34 /* ************************************************************************* */
36  const string expected_end = "examples/Data/example.graph";
37  const string actual = findExampleDataFile("example");
38  string actual_end = actual.substr(actual.size() - expected_end.size(), expected_end.size());
39  // replace all ocurrences of \\ with / use stl
40  std::replace(actual_end.begin(), actual_end.end(), '\\', '/');
41  EXPECT(assert_equal(expected_end, actual_end));
42 }
43 
44 /* ************************************************************************* */
46 {
47  const string str = "VERTEX2 1 2.000000 3.000000 4.000000";
48  istringstream is(str);
49  string tag;
50  EXPECT(is >> tag);
51  const auto actual = parseVertexPose(is, tag);
52  EXPECT(actual);
53  if (actual) {
54  EXPECT_LONGS_EQUAL(1, actual->first);
55  EXPECT(assert_equal(Pose2(2, 3, 4), actual->second));
56  }
57 }
58 
59 /* ************************************************************************* */
61 {
62  const string str = "VERTEX_XY 1 2.000000 3.000000";
63  istringstream is(str);
64  string tag;
65  EXPECT(is >> tag);
66  const auto actual = parseVertexLandmark(is, tag);
67  EXPECT(actual);
68  if (actual) {
69  EXPECT_LONGS_EQUAL(1, actual->first);
70  EXPECT(assert_equal(Point2(2, 3), actual->second));
71  }
72 }
73 
74 /* ************************************************************************* */
75 TEST( dataSet, parseEdge)
76 {
77  const string str = "EDGE2 0 1 2.000000 3.000000 4.000000";
78  istringstream is(str);
79  string tag;
80  EXPECT(is >> tag);
81  const auto actual = parseEdge(is, tag);
82  EXPECT(actual);
83  if (actual) {
84  pair<size_t, size_t> expected(0, 1);
85  EXPECT(expected == actual->first);
86  EXPECT(assert_equal(Pose2(2, 3, 4), actual->second));
87  }
88 }
89 
90 /* ************************************************************************* */
91 TEST(dataSet, load2D) {
93  const string filename = findExampleDataFile("w100.graph");
94  const auto [graph, initial] = load2D(filename);
95  EXPECT_LONGS_EQUAL(300, graph->size());
96  EXPECT_LONGS_EQUAL(100, initial->size());
98  BetweenFactor<Pose2> expected(1, 0, Pose2(-0.99879, 0.0417574, -0.00818381),
99  model);
101  std::dynamic_pointer_cast<BetweenFactor<Pose2>>(graph->at(0));
102  EXPECT(assert_equal(expected, *actual));
103 
104  // Check binary measurements, Pose2
105  size_t maxIndex = 5;
106  auto measurements = parseMeasurements<Pose2>(filename, nullptr, maxIndex);
107  EXPECT_LONGS_EQUAL(5, measurements.size());
108 
109  // Check binary measurements, Rot2
110  auto measurements2 = parseMeasurements<Rot2>(filename);
111  EXPECT_LONGS_EQUAL(300, measurements2.size());
112 
113  // // Check factor parsing
114  const auto actualFactors = parseFactors<Pose2>(filename);
115  for (size_t i : {0, 1, 2, 3, 4, 5}) {
117  *std::dynamic_pointer_cast<BetweenFactor<Pose2>>(graph->at(i)),
118  *actualFactors[i], 1e-5));
119  }
120 
121  // Check pose parsing
122  const auto actualPoses = parseVariables<Pose2>(filename);
123  for (size_t j : {0, 1, 2, 3, 4}) {
124  EXPECT(assert_equal(initial->at<Pose2>(j), actualPoses.at(j), 1e-5));
125  }
126 
127  // Check landmark parsing
128  const auto actualLandmarks = parseVariables<Point2>(filename);
129  EXPECT_LONGS_EQUAL(0, actualLandmarks.size());
130 }
131 
132 /* ************************************************************************* */
133 TEST(dataSet, load2DVictoriaPark) {
134  const string filename = findExampleDataFile("victoria_park.txt");
135  // Load all
136  const auto [graph1, initial1] = load2D(filename);
137  EXPECT_LONGS_EQUAL(10608, graph1->size());
138  EXPECT_LONGS_EQUAL(7120, initial1->size());
139 
140  // Restrict keys
141  size_t maxIndex = 5;
142  const auto [graph2, initial2] = load2D(filename, nullptr, maxIndex);
144  EXPECT_LONGS_EQUAL(6, initial2->size()); // file has 0 as well
145  EXPECT_LONGS_EQUAL(L(5), graph2->at(4)->keys()[1]);
146 }
147 
148 /* ************************************************************************* */
149 TEST(dataSet, readG2o3D) {
150  const string g2oFile = findExampleDataFile("pose3example");
151  auto model = noiseModel::Isotropic::Precision(6, 10000);
152 
153  // Initialize 6 relative measurements with quaternion/point3 values:
154  const std::vector<Pose3> relative_poses = {
155  {{0.854230, 0.190253, 0.283162, -0.392318},
156  {1.001367, 0.015390, 0.004948}},
157  {{0.105373, 0.311512, 0.656877, -0.678505},
158  {0.523923, 0.776654, 0.326659}},
159  {{0.568551, 0.595795, -0.561677, 0.079353},
160  {0.910927, 0.055169, -0.411761}},
161  {{0.542221, -0.592077, 0.303380, -0.513226},
162  {0.775288, 0.228798, -0.596923}},
163  {{0.327419, -0.125250, -0.534379, 0.769122},
164  {-0.577841, 0.628016, -0.543592}},
165  {{0.083672, 0.104639, 0.627755, 0.766795},
166  {-0.623267, 0.086928, 0.773222}},
167  };
168 
169  // Initialize 5 poses with quaternion/point3 values:
170  const std::vector<Pose3> poses = {
171  {{1.000000, 0.000000, 0.000000, 0.000000}, {0, 0, 0}},
172  {{0.854230, 0.190253, 0.283162, -0.392318},
173  {1.001367, 0.015390, 0.004948}},
174  {{0.421446, -0.351729, -0.597838, 0.584174},
175  {1.993500, 0.023275, 0.003793}},
176  {{0.067024, 0.331798, -0.200659, 0.919323},
177  {2.004291, 1.024305, 0.018047}},
178  {{0.765488, -0.035697, -0.462490, 0.445933},
179  {0.999908, 1.055073, 0.020212}},
180  };
181 
182  // Specify connectivity:
183  using KeyPair = pair<Key, Key>;
184  std::vector<KeyPair> edges = {{0, 1}, {1, 2}, {2, 3}, {3, 4}, {1, 4}, {3, 0}};
185 
186  // Created expected factor graph
187  size_t i = 0;
189  for (const auto& keys : edges) {
190  expectedGraph.emplace_shared<BetweenFactor<Pose3>>(
191  keys.first, keys.second, relative_poses[i++], model);
192  }
193 
194  // Check factor parsing
195  const auto actualFactors = parseFactors<Pose3>(g2oFile);
196  for (size_t i : {0, 1, 2, 3, 4, 5}) {
198  *std::dynamic_pointer_cast<BetweenFactor<Pose3>>(expectedGraph[i]),
199  *actualFactors[i], 1e-5));
200  }
201 
202  // Check pose parsing
203  const auto actualPoses = parseVariables<Pose3>(g2oFile);
204  for (size_t j : {0, 1, 2, 3, 4}) {
205  EXPECT(assert_equal(poses[j], actualPoses.at(j), 1e-5));
206  }
207 
208  // Check landmark parsing
209  const auto actualLandmarks = parseVariables<Point3>(g2oFile);
210  for (size_t j : {0, 1, 2, 3, 4}) {
211  EXPECT(assert_equal(poses[j], actualPoses.at(j), 1e-5));
212  }
213 
214  // Check graph version
215  bool is3D = true;
216  const auto [actualGraph, actualValues] = readG2o(g2oFile, is3D);
217  EXPECT(assert_equal(expectedGraph, *actualGraph, 1e-5));
218  for (size_t j : {0, 1, 2, 3, 4}) {
219  EXPECT(assert_equal(poses[j], actualValues->at<Pose3>(j), 1e-5));
220  }
221 }
222 
223 /* ************************************************************************* */
224 TEST( dataSet, readG2o3DNonDiagonalNoise)
225 {
226  const string g2oFile = findExampleDataFile("pose3example-offdiagonal.txt");
227  bool is3D = true;
228  const auto [actualGraph, actualValues] = readG2o(g2oFile, is3D);
229 
230  Values expectedValues;
231  Rot3 R0 = Rot3::Quaternion(1.000000, 0.000000, 0.000000, 0.000000 );
232  Point3 p0 = Point3(0.000000, 0.000000, 0.000000);
233  expectedValues.insert(0, Pose3(R0, p0));
234 
235  Rot3 R1 = Rot3::Quaternion(0.854230, 0.190253, 0.283162, -0.392318 );
236  Point3 p1 = Point3(1.001367, 0.015390, 0.004948);
237  expectedValues.insert(1, Pose3(R1, p1));
238 
239  EXPECT(assert_equal(expectedValues,*actualValues,1e-5));
240 
241  Matrix Info = Matrix(6,6);
242  for (int i = 0; i < 6; i++){
243  for (int j = i; j < 6; j++){
244  if(i==j)
245  Info(i, j) = 10000;
246  else{
247  Info(i, j) = i+1; // arbitrary nonzero number
248  Info(j, i) = i+1;
249  }
250  }
251  }
254  Point3 p01 = Point3(1.001367, 0.015390, 0.004948);
255  Rot3 R01 = Rot3::Quaternion(0.854230, 0.190253, 0.283162, -0.392318 );
256  expectedGraph.emplace_shared<BetweenFactor<Pose3> >(0, 1, Pose3(R01,p01), model);
257 
258  EXPECT(assert_equal(expectedGraph,*actualGraph,1e-2));
259 }
260 
261 /* ************************************************************************* */
262 TEST(dataSet, readG2oCheckDeterminants) {
263  const string g2oFile = findExampleDataFile("toyExample.g2o");
264 
265  // Check determinants in factors
266  auto factors = parseFactors<Pose3>(g2oFile);
268  for (const auto& factor : factors) {
269  const Rot3 R = factor->measured().rotation();
270  EXPECT_DOUBLES_EQUAL(1.0, R.matrix().determinant(), 1e-9);
271  }
272 
273  // Check determinants in initial values
274  const map<size_t, Pose3> poses = parseVariables<Pose3>(g2oFile);
275  EXPECT_LONGS_EQUAL(5, poses.size());
276  for (const auto& key_value : poses) {
277  const Rot3 R = key_value.second.rotation();
278  EXPECT_DOUBLES_EQUAL(1.0, R.matrix().determinant(), 1e-9);
279  }
280  const map<size_t, Point3> landmarks = parseVariables<Point3>(g2oFile);
281  EXPECT_LONGS_EQUAL(0, landmarks.size());
282 }
283 
284 /* ************************************************************************* */
285 TEST(dataSet, readG2oLandmarks) {
286  const string g2oFile = findExampleDataFile("example_with_vertices.g2o");
287 
288  // Check number of poses and landmarks. Should be 8 each.
289  const map<size_t, Pose3> poses = parseVariables<Pose3>(g2oFile);
290  EXPECT_LONGS_EQUAL(8, poses.size());
291  const map<size_t, Point3> landmarks = parseVariables<Point3>(g2oFile);
292  EXPECT_LONGS_EQUAL(8, landmarks.size());
293 
294  auto graphAndValues = load3D(g2oFile);
295  EXPECT(graphAndValues.second->exists(L(0)));
296 }
297 
298 /* ************************************************************************* */
302  g.emplace_shared<Factor>(0, 1, Pose2(1.030390, 0.011350, -0.081596), model);
303  g.emplace_shared<Factor>(1, 2, Pose2(1.013900, -0.058639, -0.220291), model);
304  g.emplace_shared<Factor>(2, 3, Pose2(1.027650, -0.007456, -0.043627), model);
305  g.emplace_shared<Factor>(3, 4, Pose2(-0.012016, 1.004360, 1.560229), model);
306  g.emplace_shared<Factor>(4, 5, Pose2(1.016030, 0.014565, -0.030930), model);
307  g.emplace_shared<Factor>(5, 6, Pose2(1.023890, 0.006808, -0.007452), model);
308  g.emplace_shared<Factor>(6, 7, Pose2(0.957734, 0.003159, 0.082836), model);
309  g.emplace_shared<Factor>(7, 8, Pose2(-1.023820, -0.013668, -3.084560), model);
310  g.emplace_shared<Factor>(8, 9, Pose2(1.023440, 0.013984, -0.127624), model);
311  g.emplace_shared<Factor>(9, 10, Pose2(1.003350, 0.022250, -0.195918), model);
312  g.emplace_shared<Factor>(5, 9, Pose2(0.033943, 0.032439, 3.073637), model);
313  g.emplace_shared<Factor>(3, 10, Pose2(0.044020, 0.988477, -1.553511), model);
314  return g;
315 }
316 
317 /* ************************************************************************* */
318 TEST(dataSet, readG2o) {
319  const string g2oFile = findExampleDataFile("pose2example");
320  const auto [actualGraph, actualValues] = readG2o(g2oFile);
321 
323  Vector3(44.721360, 44.721360, 30.901699));
324  EXPECT(assert_equal(expectedGraph(model), *actualGraph, 1e-5));
325 
326  Values expectedValues;
327  expectedValues.insert(0, Pose2(0.000000, 0.000000, 0.000000));
328  expectedValues.insert(1, Pose2(1.030390, 0.011350, -0.081596));
329  expectedValues.insert(2, Pose2(2.036137, -0.129733, -0.301887));
330  expectedValues.insert(3, Pose2(3.015097, -0.442395, -0.345514));
331  expectedValues.insert(4, Pose2(3.343949, 0.506678, 1.214715));
332  expectedValues.insert(5, Pose2(3.684491, 1.464049, 1.183785));
333  expectedValues.insert(6, Pose2(4.064626, 2.414783, 1.176333));
334  expectedValues.insert(7, Pose2(4.429778, 3.300180, 1.259169));
335  expectedValues.insert(8, Pose2(4.128877, 2.321481, -1.825391));
336  expectedValues.insert(9, Pose2(3.884653, 1.327509, -1.953016));
337  expectedValues.insert(10, Pose2(3.531067, 0.388263, -2.148934));
338  EXPECT(assert_equal(expectedValues, *actualValues, 1e-5));
339 }
340 
341 /* ************************************************************************* */
342 TEST(dataSet, readG2oHuber) {
343  const string g2oFile = findExampleDataFile("pose2example");
344  bool is3D = false;
345  const auto [actualGraph, actualValues] =
346  readG2o(g2oFile, is3D, KernelFunctionTypeHUBER);
347 
348  auto baseModel = noiseModel::Diagonal::Precisions(
349  Vector3(44.721360, 44.721360, 30.901699));
351  noiseModel::mEstimator::Huber::Create(1.345), baseModel);
352 
353  EXPECT(assert_equal(expectedGraph(model), *actualGraph, 1e-5));
354 }
355 
356 /* ************************************************************************* */
357 TEST(dataSet, readG2oTukey) {
358  const string g2oFile = findExampleDataFile("pose2example");
359  bool is3D = false;
360  const auto [actualGraph, actualValues] =
361  readG2o(g2oFile, is3D, KernelFunctionTypeTUKEY);
362 
363  auto baseModel = noiseModel::Diagonal::Precisions(
364  Vector3(44.721360, 44.721360, 30.901699));
366  noiseModel::mEstimator::Tukey::Create(4.6851), baseModel);
367 
368  EXPECT(assert_equal(expectedGraph(model), *actualGraph, 1e-5));
369 }
370 
371 /* ************************************************************************* */
372 TEST( dataSet, writeG2o)
373 {
374  const string g2oFile = findExampleDataFile("pose2example");
375  const auto [expectedGraph, expectedValues] = readG2o(g2oFile);
376 
377  const string filenameToWrite = createRewrittenFileName(g2oFile);
378  writeG2o(*expectedGraph, *expectedValues, filenameToWrite);
379 
380  const auto [actualGraph, actualValues] = readG2o(filenameToWrite);
381  EXPECT(assert_equal(*expectedValues,*actualValues,1e-5));
382  EXPECT(assert_equal(*expectedGraph,*actualGraph,1e-5));
383 }
384 
385 /* ************************************************************************* */
386 TEST( dataSet, writeG2o3D)
387 {
388  const string g2oFile = findExampleDataFile("pose3example");
389  bool is3D = true;
390  const auto [expectedGraph, expectedValues] = readG2o(g2oFile, is3D);
391 
392  const string filenameToWrite = createRewrittenFileName(g2oFile);
393  writeG2o(*expectedGraph, *expectedValues, filenameToWrite);
394 
395  const auto [actualGraph, actualValues] = readG2o(filenameToWrite, is3D);
396  EXPECT(assert_equal(*expectedValues,*actualValues,1e-4));
397  EXPECT(assert_equal(*expectedGraph,*actualGraph,1e-4));
398 }
399 
400 /* ************************************************************************* */
401 TEST( dataSet, writeG2o3DNonDiagonalNoise)
402 {
403  const string g2oFile = findExampleDataFile("pose3example-offdiagonal");
404  bool is3D = true;
405  const auto [expectedGraph, expectedValues] = readG2o(g2oFile, is3D);
406 
407  const string filenameToWrite = createRewrittenFileName(g2oFile);
408  writeG2o(*expectedGraph, *expectedValues, filenameToWrite);
409 
410  const auto [actualGraph, actualValues] = readG2o(filenameToWrite, is3D);
411  EXPECT(assert_equal(*expectedValues,*actualValues,1e-4));
412  EXPECT(assert_equal(*expectedGraph,*actualGraph,1e-4));
413 }
414 
415 /* ************************************************************************* */
416 int main() { TestResult tr; return TestRegistry::runAllTests(tr); }
417 /* ************************************************************************* */
vector< MFAS::KeyPair > edges
Definition: testMFAS.cpp:25
Provides additional testing facilities for common data structures.
void writeG2o(const NonlinearFactorGraph &graph, const Values &estimate, const std::string &filename)
This function writes a g2o file from NonlinearFactorGraph and a Values structure. ...
Definition: dataset.cpp:636
static shared_ptr Covariance(const Matrix &covariance, bool smart=true)
Definition: NoiseModel.cpp:114
IsDerived< DERIVEDFACTOR > emplace_shared(Args &&... args)
Emplace a shared pointer to factor of given type.
Definition: FactorGraph.h:196
static int runAllTests(TestResult &result)
Eigen::Vector3d Vector3
Definition: Vector.h:43
Vector3f p1
noiseModel::Diagonal::shared_ptr model
Matrix expected
Definition: testMatrix.cpp:971
Vector2 Point2
Definition: Point2.h:32
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:40
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
const GaussianFactorGraph factors
Values initial
Definition: BFloat16.h:88
Vector3f p0
GTSAM_EXPORT std::map< size_t, Pose3 > parseVariables< Pose3 >(const std::string &filename, size_t maxIndex)
Definition: dataset.cpp:775
static shared_ptr Create(const RobustModel::shared_ptr &robust, const NoiseModel::shared_ptr noise)
Definition: NoiseModel.cpp:704
NonlinearFactorGraph graph
Rot3 is a 3D rotation represented as a rotation matrix if the preprocessor symbol GTSAM_USE_QUATERNIO...
Definition: Rot3.h:58
#define EXPECT_DOUBLES_EQUAL(expected, actual, threshold)
Definition: Test.h:161
void g(const string &key, int i)
Definition: testBTree.cpp:41
static shared_ptr Create(size_t dim)
Definition: NoiseModel.h:610
static NonlinearFactorGraph expectedGraph(const SharedNoiseModel &model)
static Rot3 Quaternion(double w, double x, double y, double z)
Definition: Rot3.h:204
size_t size() const
Definition: FactorGraph.h:334
GTSAM_EXPORT std::map< size_t, Point2 > parseVariables< Point2 >(const std::string &filename, size_t maxIndex)
Definition: dataset.cpp:209
std::shared_ptr< BetweenFactor > shared_ptr
Definition: BetweenFactor.h:63
std::optional< IndexedEdge > parseEdge(std::istream &is, const std::string &tag)
Definition: dataset.cpp:299
Definition: pytypes.h:1403
static shared_ptr Create(double k, const ReweightScheme reweight=Block)
#define EXPECT(condition)
Definition: Test.h:150
GraphAndValues load3D(const std::string &filename)
Load TORO 3D Graph.
Definition: dataset.cpp:922
GraphAndValues load2D(const std::string &filename, SharedNoiseModel model, size_t maxIndex, bool addNoise, bool smart, NoiseFormat noiseFormat, KernelFunctionType kernelFunctionType)
Definition: dataset.cpp:505
int main()
Array< double, 1, 3 > e(1./3., 0.5, 2.)
GraphAndValues readG2o(const std::string &g2oFile, const bool is3D, KernelFunctionType kernelFunctionType)
This function parses a g2o file and stores the measurements into a NonlinearFactorGraph and the initi...
Definition: dataset.cpp:621
static shared_ptr Precisions(const Vector &precisions, bool smart=true)
Definition: NoiseModel.cpp:287
NonlinearFactorGraph graph2()
GTSAM_EXPORT std::map< size_t, Pose2 > parseVariables< Pose2 >(const std::string &filename, size_t maxIndex)
Definition: dataset.cpp:187
traits
Definition: chartTesting.h:28
GTSAM_EXPORT std::string findExampleDataFile(const std::string &name)
Definition: dataset.cpp:70
#define EXPECT_LONGS_EQUAL(expected, actual)
Definition: Test.h:154
Key R(std::uint64_t j)
const sharedFactor at(size_t i) const
Definition: FactorGraph.h:343
Matrix3 matrix() const
Definition: Rot3M.cpp:218
Key L(std::uint64_t j)
std::optional< IndexedLandmark > parseVertexLandmark(std::istream &is, const std::string &tag)
Definition: dataset.cpp:193
GTSAM_EXPORT std::vector< BetweenFactor< Pose2 >::shared_ptr > parseFactors< Pose2 >(const std::string &filename, const noiseModel::Diagonal::shared_ptr &model, size_t maxIndex)
Definition: dataset.cpp:440
void insert(Key j, const Value &val)
Definition: Values.cpp:155
static Rot3 R0
Vector3 Point3
Definition: Point3.h:38
TEST(SmartFactorBase, Pinhole)
std::shared_ptr< Gaussian > shared_ptr
Definition: NoiseModel.h:189
const KeyVector keys
2D Pose
GTSAM_EXPORT std::string createRewrittenFileName(const std::string &name)
Definition: dataset.cpp:105
std::optional< IndexedPose > parseVertexPose(std::istream &is, const std::string &tag)
Definition: dataset.cpp:173
GTSAM_EXPORT std::map< size_t, Point3 > parseVariables< Point3 >(const std::string &filename, size_t maxIndex)
Definition: dataset.cpp:793
static shared_ptr Precision(size_t dim, double precision, bool smart=true)
Definition: NoiseModel.h:560
utility functions for loading datasets
std::ptrdiff_t j
static shared_ptr Create(double k, const ReweightScheme reweight=Block)
noiseModel::Base::shared_ptr SharedNoiseModel
Definition: NoiseModel.h:741
GTSAM_EXPORT std::vector< BetweenFactor< Pose3 >::shared_ptr > parseFactors< Pose3 >(const std::string &filename, const noiseModel::Diagonal::shared_ptr &model, size_t maxIndex)
Definition: dataset.cpp:914


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:38:01