gtsam/examples/SFMExample_bal_COLAMD_METIS...

139 lines
4.6 KiB
C++
Raw Normal View History

2015-02-19 12:36:31 +08:00
/* ----------------------------------------------------------------------------
* GTSAM Copyright 2010, Georgia Tech Research Corporation,
* Atlanta, Georgia 30332-0415
* All Rights Reserved
* Authors: Frank Dellaert, et al. (see THANKS for the full author list)
* See LICENSE for the license information
* -------------------------------------------------------------------------- */
/**
2020-05-10 07:08:31 +08:00
* @file SFMExample_bal_COLAMD_METIS.cpp
2015-02-19 12:36:31 +08:00
* @brief This file is to compare the ordering performance for COLAMD vs METIS.
* Example problem is to solve a structure-from-motion problem from a "Bundle Adjustment in the Large" file.
* @author Frank Dellaert, Zhaoyang Lv
*/
// For an explanation of headers, see SFMExample.cpp
#include <gtsam/slam/GeneralSFMFactor.h>
2022-02-01 01:29:13 +08:00
#include <gtsam/sfm/SfmData.h> // for loading BAL datasets !
#include <gtsam/slam/dataset.h>
#include <gtsam/nonlinear/NonlinearFactorGraph.h>
#include <gtsam/nonlinear/LevenbergMarquardtOptimizer.h>
#include <gtsam/inference/Symbol.h>
#include <gtsam/inference/Ordering.h>
2015-02-23 23:24:34 +08:00
#include <gtsam/base/timing.h>
2015-02-19 12:36:31 +08:00
#include <vector>
using namespace std;
using namespace gtsam;
using symbol_shorthand::C;
using symbol_shorthand::P;
// We will be using a projection factor that ties a SFM_Camera to a 3D point.
2020-05-10 07:08:31 +08:00
// An SFM_Camera is defined in datase.h as a camera with unknown Cal3Bundler
// calibration and has a total of 9 free parameters
typedef GeneralSFMFactor<SfmCamera, Point3> MyFactor;
2015-02-19 12:36:31 +08:00
2020-05-10 07:08:31 +08:00
int main(int argc, char* argv[]) {
2015-02-19 12:36:31 +08:00
// Find default file, but if an argument is given, try loading a file
string filename = findExampleDataFile("dubrovnik-3-7-pre");
2020-05-10 07:08:31 +08:00
if (argc > 1) filename = string(argv[1]);
2015-02-19 12:36:31 +08:00
// Load the SfM data from file
2022-02-01 01:29:13 +08:00
SfmData mydata = SfmData::FromBalFile(filename);
2023-01-21 05:04:12 +08:00
cout << "read " << mydata.numberTracks() << " tracks on " << mydata.numberCameras() << " cameras" << endl;
2015-02-19 12:36:31 +08:00
// Create a factor graph
2015-02-20 05:00:21 +08:00
NonlinearFactorGraph graph;
2015-02-19 12:36:31 +08:00
// We share *one* noiseModel between all projection factors
2020-05-10 07:08:31 +08:00
auto noise = noiseModel::Isotropic::Sigma(2, 1.0); // one pixel in u and v
2015-02-19 12:36:31 +08:00
// Add measurements to the factor graph
size_t j = 0;
2020-05-10 07:08:31 +08:00
for (const SfmTrack& track : mydata.tracks) {
2023-02-05 02:35:42 +08:00
for (const auto& [i, uv] : track.measurements) {
2020-05-10 07:08:31 +08:00
graph.emplace_shared<MyFactor>(
uv, noise, C(i), P(j)); // note use of shorthand symbols C and P
2015-02-19 12:36:31 +08:00
}
j += 1;
}
// Add a prior on pose x1. This indirectly specifies where the origin is.
// and a prior on the position of the first landmark to fix the scale
2020-05-10 07:08:31 +08:00
graph.addPrior(C(0), mydata.cameras[0], noiseModel::Isotropic::Sigma(9, 0.1));
graph.addPrior(P(0), mydata.tracks[0].p,
noiseModel::Isotropic::Sigma(3, 0.1));
2015-02-19 12:36:31 +08:00
// Create initial estimate
Values initial;
2020-05-10 07:08:31 +08:00
size_t i = 0;
j = 0;
for (const SfmCamera& camera : mydata.cameras) initial.insert(C(i++), camera);
for (const SfmTrack& track : mydata.tracks) initial.insert(P(j++), track.p);
2015-02-19 12:36:31 +08:00
/** --------------- COMPARISON -----------------------**/
/** ----------------------------------------------------**/
2015-02-19 12:36:31 +08:00
2015-02-23 23:24:34 +08:00
LevenbergMarquardtParams params_using_COLAMD, params_using_METIS;
2015-02-19 12:36:31 +08:00
try {
params_using_METIS.setVerbosity("ERROR");
2015-02-23 23:24:34 +08:00
gttic_(METIS_ORDERING);
params_using_METIS.ordering = Ordering::Create(Ordering::METIS, graph);
2015-02-23 23:24:34 +08:00
gttoc_(METIS_ORDERING);
2015-02-19 12:36:31 +08:00
params_using_COLAMD.setVerbosity("ERROR");
2015-02-23 23:24:34 +08:00
gttic_(COLAMD_ORDERING);
2015-02-20 05:00:21 +08:00
params_using_COLAMD.ordering = Ordering::Create(Ordering::COLAMD, graph);
2015-02-23 23:24:34 +08:00
gttoc_(COLAMD_ORDERING);
2015-02-19 12:36:31 +08:00
} catch (exception& e) {
cout << e.what();
}
// expect they have different ordering results
2020-05-10 07:08:31 +08:00
if (params_using_COLAMD.ordering == params_using_METIS.ordering) {
cout << "COLAMD and METIS produce the same ordering. "
<< "Problem here!!!" << endl;
}
2015-02-19 12:36:31 +08:00
2015-02-23 23:24:34 +08:00
/* Optimize the graph with METIS and COLAMD and time the results */
2015-02-19 12:36:31 +08:00
2015-02-23 23:24:34 +08:00
Values result_METIS, result_COLAMD;
2015-02-19 12:36:31 +08:00
try {
2015-02-23 23:24:34 +08:00
gttic_(OPTIMIZE_WITH_METIS);
LevenbergMarquardtOptimizer lm_METIS(graph, initial, params_using_METIS);
result_METIS = lm_METIS.optimize();
2015-02-23 23:24:34 +08:00
gttoc_(OPTIMIZE_WITH_METIS);
2015-02-19 12:36:31 +08:00
2015-02-23 23:24:34 +08:00
gttic_(OPTIMIZE_WITH_COLAMD);
LevenbergMarquardtOptimizer lm_COLAMD(graph, initial, params_using_COLAMD);
result_COLAMD = lm_COLAMD.optimize();
2015-02-23 23:24:34 +08:00
gttoc_(OPTIMIZE_WITH_COLAMD);
2015-02-19 12:36:31 +08:00
} catch (exception& e) {
cout << e.what();
}
2020-05-10 07:08:31 +08:00
{ // printing the result
2015-02-19 12:36:31 +08:00
2015-02-23 23:24:34 +08:00
cout << "COLAMD final error: " << graph.error(result_COLAMD) << endl;
cout << "METIS final error: " << graph.error(result_METIS) << endl;
2015-02-23 23:31:25 +08:00
cout << endl << endl;
2015-02-20 05:00:21 +08:00
2015-02-23 23:31:25 +08:00
cout << "Time comparison by solving " << filename << " results:" << endl;
2023-01-21 05:04:12 +08:00
cout << mydata.numberTracks() << " point tracks and " << mydata.numberCameras()
<< " cameras" << endl;
2015-02-23 23:31:25 +08:00
tictoc_print_();
2015-02-20 05:00:21 +08:00
}
2015-02-19 12:36:31 +08:00
return 0;
}