CombinedImuFactorsExample.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 
36 #include <boost/program_options.hpp>
37 
38 // GTSAM related includes.
39 #include <gtsam/inference/Symbol.h>
46 #include <gtsam/slam/dataset.h>
47 
48 #include <cstring>
49 #include <fstream>
50 #include <iostream>
51 
52 using namespace gtsam;
53 using namespace std;
54 
55 using symbol_shorthand::B; // Bias (ax,ay,az,gx,gy,gz)
56 using symbol_shorthand::V; // Vel (xdot,ydot,zdot)
57 using symbol_shorthand::X; // Pose3 (x,y,z,r,p,y)
58 
59 namespace po = boost::program_options;
60 
61 po::variables_map parseOptions(int argc, char* argv[]) {
62  po::options_description desc;
63  desc.add_options()("help,h", "produce help message") // help message
64  ("data_csv_path", po::value<string>()->default_value("imuAndGPSdata.csv"),
65  "path to the CSV file with the IMU data") // path to the data file
66  ("output_filename",
67  po::value<string>()->default_value("imuFactorExampleResults.csv"),
68  "path to the result file to use") // filename to save results to
69  ("use_isam", po::bool_switch(),
70  "use ISAM as the optimizer"); // flag for ISAM optimizer
71 
72  po::variables_map vm;
73  po::store(po::parse_command_line(argc, argv, desc), vm);
74 
75  if (vm.count("help")) {
76  cout << desc << "\n";
77  exit(1);
78  }
79 
80  return vm;
81 }
82 
83 Vector10 readInitialState(ifstream& file) {
84  string value;
85  // Format is (N,E,D,qX,qY,qZ,qW,velN,velE,velD)
86  Vector10 initial_state;
87  getline(file, value, ','); // i
88  for (int i = 0; i < 9; i++) {
89  getline(file, value, ',');
90  initial_state(i) = stof(value.c_str());
91  }
92  getline(file, value, '\n');
93  initial_state(9) = stof(value.c_str());
94 
95  return initial_state;
96 }
97 
98 std::shared_ptr<PreintegratedCombinedMeasurements::Params> imuParams() {
99  // We use the sensor specs to build the noise model for the IMU factor.
100  double accel_noise_sigma = 0.0003924;
101  double gyro_noise_sigma = 0.000205689024915;
102  double accel_bias_rw_sigma = 0.004905;
103  double gyro_bias_rw_sigma = 0.000001454441043;
104  Matrix33 measured_acc_cov = I_3x3 * pow(accel_noise_sigma, 2);
105  Matrix33 measured_omega_cov = I_3x3 * pow(gyro_noise_sigma, 2);
106  Matrix33 integration_error_cov =
107  I_3x3 * 1e-8; // error committed in integrating position from velocities
108  Matrix33 bias_acc_cov = I_3x3 * pow(accel_bias_rw_sigma, 2);
109  Matrix33 bias_omega_cov = I_3x3 * pow(gyro_bias_rw_sigma, 2);
110  Matrix66 bias_acc_omega_init =
111  I_6x6 * 1e-5; // error in the bias used for preintegration
112 
114  // PreintegrationBase params:
115  p->accelerometerCovariance =
116  measured_acc_cov; // acc white noise in continuous
117  p->integrationCovariance =
118  integration_error_cov; // integration uncertainty continuous
119  // should be using 2nd order integration
120  // PreintegratedRotation params:
121  p->gyroscopeCovariance =
122  measured_omega_cov; // gyro white noise in continuous
123  // PreintegrationCombinedMeasurements params:
124  p->biasAccCovariance = bias_acc_cov; // acc bias in continuous
125  p->biasOmegaCovariance = bias_omega_cov; // gyro bias in continuous
126  p->biasAccOmegaInt = bias_acc_omega_init;
127 
128  return p;
129 }
130 
131 int main(int argc, char* argv[]) {
132  string data_filename, output_filename;
133  po::variables_map var_map = parseOptions(argc, argv);
134 
135  data_filename = findExampleDataFile(var_map["data_csv_path"].as<string>());
136  output_filename = var_map["output_filename"].as<string>();
137 
138  // Set up output file for plotting errors
139  FILE* fp_out = fopen(output_filename.c_str(), "w+");
140  fprintf(fp_out,
141  "#time(s),x(m),y(m),z(m),qx,qy,qz,qw,gt_x(m),gt_y(m),gt_z(m),gt_qx,"
142  "gt_qy,gt_qz,gt_qw\n");
143 
144  // Begin parsing the CSV file. Input the first line for initialization.
145  // From there, we'll iterate through the file and we'll preintegrate the IMU
146  // or add in the GPS given the input.
147  ifstream file(data_filename.c_str());
148 
149  Vector10 initial_state = readInitialState(file);
150  cout << "initial state:\n" << initial_state.transpose() << "\n\n";
151 
152  // Assemble initial quaternion through GTSAM constructor
153  // ::Quaternion(w,x,y,z);
154  Rot3 prior_rotation = Rot3::Quaternion(initial_state(6), initial_state(3),
155  initial_state(4), initial_state(5));
156  Point3 prior_point(initial_state.head<3>());
157  Pose3 prior_pose(prior_rotation, prior_point);
158  Vector3 prior_velocity(initial_state.tail<3>());
159 
160  imuBias::ConstantBias prior_imu_bias; // assume zero initial bias
161 
162  int index = 0;
163 
164  Values initial_values;
165 
166  // insert pose at initialization
167  initial_values.insert(X(index), prior_pose);
168  initial_values.insert(V(index), prior_velocity);
169  initial_values.insert(B(index), prior_imu_bias);
170 
171  // Assemble prior noise model and add it the graph.`
172  auto pose_noise_model = noiseModel::Diagonal::Sigmas(
173  (Vector(6) << 0.01, 0.01, 0.01, 0.5, 0.5, 0.5)
174  .finished()); // rad,rad,rad,m, m, m
175  auto velocity_noise_model = noiseModel::Isotropic::Sigma(3, 0.1); // m/s
176  auto bias_noise_model = noiseModel::Isotropic::Sigma(6, 1e-3);
177 
178  // Add all prior factors (pose, velocity, bias) to the graph.
180  graph.addPrior<Pose3>(X(index), prior_pose, pose_noise_model);
181  graph.addPrior<Vector3>(V(index), prior_velocity, velocity_noise_model);
182  graph.addPrior<imuBias::ConstantBias>(B(index), prior_imu_bias,
183  bias_noise_model);
184 
185  auto p = imuParams();
186 
187  std::shared_ptr<PreintegrationType> preintegrated =
188  std::make_shared<PreintegratedCombinedMeasurements>(p, prior_imu_bias);
189 
190  assert(preintegrated);
191 
192  // Store previous state for imu integration and latest predicted outcome.
193  NavState prev_state(prior_pose, prior_velocity);
194  NavState prop_state = prev_state;
195  imuBias::ConstantBias prev_bias = prior_imu_bias;
196 
197  // Keep track of total error over the entire run as simple performance metric.
198  double current_position_error = 0.0, current_orientation_error = 0.0;
199 
200  double output_time = 0.0;
201  double dt = 0.005; // The real system has noise, but here, results are nearly
202  // exactly the same, so keeping this for simplicity.
203 
204  // All priors have been set up, now iterate through the data file.
205  while (file.good()) {
206  // Parse out first value
207  string value;
208  getline(file, value, ',');
209  int type = stoi(value.c_str());
210 
211  if (type == 0) { // IMU measurement
212  Vector6 imu;
213  for (int i = 0; i < 5; ++i) {
214  getline(file, value, ',');
215  imu(i) = stof(value.c_str());
216  }
217  getline(file, value, '\n');
218  imu(5) = stof(value.c_str());
219 
220  // Adding the IMU preintegration.
221  preintegrated->integrateMeasurement(imu.head<3>(), imu.tail<3>(), dt);
222 
223  } else if (type == 1) { // GPS measurement
224  Vector7 gps;
225  for (int i = 0; i < 6; ++i) {
226  getline(file, value, ',');
227  gps(i) = stof(value.c_str());
228  }
229  getline(file, value, '\n');
230  gps(6) = stof(value.c_str());
231 
232  index++;
233 
234  // Adding IMU factor and GPS factor and optimizing.
235  auto preint_imu_combined =
236  dynamic_cast<const PreintegratedCombinedMeasurements&>(
237  *preintegrated);
238  CombinedImuFactor imu_factor(X(index - 1), V(index - 1), X(index),
239  V(index), B(index - 1), B(index),
240  preint_imu_combined);
241  graph.add(imu_factor);
242 
243  auto correction_noise = noiseModel::Isotropic::Sigma(3, 1.0);
244  GPSFactor gps_factor(X(index),
245  Point3(gps(0), // N,
246  gps(1), // E,
247  gps(2)), // D,
248  correction_noise);
249  graph.add(gps_factor);
250 
251  // Now optimize and compare results.
252  prop_state = preintegrated->predict(prev_state, prev_bias);
253  initial_values.insert(X(index), prop_state.pose());
254  initial_values.insert(V(index), prop_state.v());
255  initial_values.insert(B(index), prev_bias);
256 
258  params.setVerbosityLM("SUMMARY");
259  LevenbergMarquardtOptimizer optimizer(graph, initial_values, params);
260  Values result = optimizer.optimize();
261 
262  // Overwrite the beginning of the preintegration for the next step.
263  prev_state =
264  NavState(result.at<Pose3>(X(index)), result.at<Vector3>(V(index)));
265  prev_bias = result.at<imuBias::ConstantBias>(B(index));
266 
267  // Reset the preintegration object.
268  preintegrated->resetIntegrationAndSetBias(prev_bias);
269 
270  // Print out the position and orientation error for comparison.
271  Vector3 result_position = prev_state.pose().translation();
272  Vector3 position_error = result_position - gps.head<3>();
273  current_position_error = position_error.norm();
274 
275  Quaternion result_quat = prev_state.pose().rotation().toQuaternion();
276  Quaternion gps_quat(gps(6), gps(3), gps(4), gps(5));
277  Quaternion quat_error = result_quat * gps_quat.inverse();
278  quat_error.normalize();
279  Vector3 euler_angle_error(quat_error.x() * 2, quat_error.y() * 2,
280  quat_error.z() * 2);
281  current_orientation_error = euler_angle_error.norm();
282 
283  // display statistics
284  cout << "Position error:" << current_position_error << "\t "
285  << "Angular error:" << current_orientation_error << "\n"
286  << endl;
287 
288  fprintf(fp_out, "%f,%f,%f,%f,%f,%f,%f,%f,%f,%f,%f,%f,%f,%f,%f\n",
289  output_time, result_position(0), result_position(1),
290  result_position(2), result_quat.x(), result_quat.y(),
291  result_quat.z(), result_quat.w(), gps(0), gps(1), gps(2),
292  gps_quat.x(), gps_quat.y(), gps_quat.z(), gps_quat.w());
293 
294  output_time += 1.0;
295 
296  } else {
297  cerr << "ERROR parsing file\n";
298  return 1;
299  }
300  }
301  fclose(fp_out);
302  cout << "Complete, results written to " << output_filename << "\n\n";
303 
304  return 0;
305 }
std::shared_ptr< PreintegratedCombinedMeasurements::Params > imuParams()
static std::shared_ptr< PreintegrationCombinedParams > MakeSharedD(double g=9.81)
EIGEN_DEVICE_FUNC CoeffReturnType x() const
const Point3 & translation(OptionalJacobian< 3, 6 > Hself={}) const
get translation
Definition: Pose3.cpp:308
virtual const Values & optimize()
Eigen::Vector3d Vector3
Definition: Vector.h:43
Factor Graph consisting of non-linear factors.
const ValueType at(Key j) const
Definition: Values-inl.h:204
EIGEN_DEVICE_FUNC CoeffReturnType z() const
po::variables_map parseOptions(int argc, char *argv[])
Vector10 readInitialState(ifstream &file)
int main(int argc, char *argv[])
Definition: BFloat16.h:88
EIGEN_DEVICE_FUNC CoeffReturnType w() const
Key X(std::uint64_t j)
NonlinearFactorGraph graph
Rot3 is a 3D rotation represented as a rotation matrix if the preprocessor symbol GTSAM_USE_QUATERNIO...
Definition: Rot3.h:58
IsDerived< DERIVEDFACTOR > add(std::shared_ptr< DERIVEDFACTOR > factor)
add is a synonym for push_back.
Definition: FactorGraph.h:214
static const SmartProjectionParams params
const double dt
void addPrior(Key key, const T &prior, const SharedNoiseModel &model=nullptr)
const string output_filename
static Rot3 Quaternion(double w, double x, double y, double z)
Definition: Rot3.h:204
Header file for GPS factor.
const Vector3 & v() const
Return velocity as Vector3. Computation-free.
Definition: NavState.h:107
Eigen::VectorXd Vector
Definition: Vector.h:38
Values result
EIGEN_DEVICE_FUNC Quaternion< Scalar > inverse() const
A nonlinear optimizer that uses the Levenberg-Marquardt trust-region scheme.
Array< double, 1, 3 > e(1./3., 0.5, 2.)
EIGEN_DEVICE_FUNC CoeffReturnType y() const
traits
Definition: chartTesting.h:28
gtsam::Quaternion toQuaternion() const
Definition: Rot3M.cpp:232
Key B(std::uint64_t j)
GTSAM_EXPORT std::string findExampleDataFile(const std::string &name)
Definition: dataset.cpp:70
static shared_ptr Sigmas(const Vector &sigmas, bool smart=true)
Definition: NoiseModel.cpp:270
The quaternion class used to represent 3D orientations and rotations.
Key V(std::uint64_t j)
float * p
void setVerbosityLM(const std::string &s)
void insert(Key j, const Value &val)
Definition: Values.cpp:155
const Pose3 pose() const
Definition: NavState.h:86
const Rot3 & rotation(OptionalJacobian< 3, 6 > Hself={}) const
get rotation
Definition: Pose3.cpp:315
Vector3 Point3
Definition: Point3.h:38
Jet< T, N > pow(const Jet< T, N > &f, double g)
Definition: jet.h:570
#define X
Definition: icosphere.cpp:20
utility functions for loading datasets
static shared_ptr Sigma(size_t dim, double sigma, bool smart=true)
Definition: NoiseModel.cpp:594


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