g2o
curve_fit.cpp
Go to the documentation of this file.
1 // g2o - General Graph Optimization
2 // Copyright (C) 2012 R. Kümmerle
3 // All rights reserved.
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are
7 // met:
8 //
9 // * Redistributions of source code must retain the above copyright notice,
10 // this list of conditions and the following disclaimer.
11 // * Redistributions in binary form must reproduce the above copyright
12 // notice, this list of conditions and the following disclaimer in the
13 // documentation and/or other materials provided with the distribution.
14 //
15 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
16 // IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
17 // TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
18 // PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
19 // HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
20 // SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
21 // TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
22 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
23 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
24 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
25 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26 
27 #include <Eigen/Core>
28 #include <iostream>
29 
30 #include "g2o/stuff/sampler.h"
31 #include "g2o/stuff/command_args.h"
33 #include "g2o/core/block_solver.h"
34 #include "g2o/core/solver.h"
36 #include "g2o/core/base_vertex.h"
39 
40 using namespace std;
41 
45 class VertexParams : public g2o::BaseVertex<3, Eigen::Vector3d>
46 {
47  public:
50  {
51  }
52 
53  virtual bool read(std::istream& /*is*/)
54  {
55  cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
56  return false;
57  }
58 
59  virtual bool write(std::ostream& /*os*/) const
60  {
61  cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
62  return false;
63  }
64 
65  virtual void setToOriginImpl()
66  {
67  cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
68  }
69 
70  virtual void oplusImpl(const double* update)
71  {
72  Eigen::Vector3d::ConstMapType v(update);
73  _estimate += v;
74  }
75 };
76 
84 class EdgePointOnCurve : public g2o::BaseUnaryEdge<1, Eigen::Vector2d, VertexParams>
85 {
86  public:
87  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
89  {
90  }
91  virtual bool read(std::istream& /*is*/)
92  {
93  cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
94  return false;
95  }
96  virtual bool write(std::ostream& /*os*/) const
97  {
98  cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
99  return false;
100  }
101 
103  {
104  const VertexParams* params = static_cast<const VertexParams*>(vertex(0));
105  const double& a = params->estimate()(0);
106  const double& b = params->estimate()(1);
107  const double& lambda = params->estimate()(2);
108  double fval = a * exp(-lambda * measurement()(0)) + b;
109  _error(0) = fval - measurement()(1);
110  }
111 };
112 
113 int main(int argc, char** argv)
114 {
115  int numPoints;
116  int maxIterations;
117  bool verbose;
118  std::vector<int> gaugeList;
119  string dumpFilename;
120  g2o::CommandArgs arg;
121  arg.param("dump", dumpFilename, "", "dump the points into a file");
122  arg.param("numPoints", numPoints, 50, "number of points sampled from the curve");
123  arg.param("i", maxIterations, 10, "perform n iterations");
124  arg.param("v", verbose, false, "verbose output of the optimization process");
125 
126  arg.parseArgs(argc, argv);
127 
128  // generate random data
129  double a = 2.;
130  double b = 0.4;
131  double lambda = 0.2;
132  Eigen::Vector2d* points = new Eigen::Vector2d[numPoints];
133  for (int i = 0; i < numPoints; ++i) {
134  double x = g2o::Sampler::uniformRand(0, 10);
135  double y = a * exp(-lambda * x) + b;
136  // add Gaussian noise
137  y += g2o::Sampler::gaussRand(0, 0.02);
138  points[i].x() = x;
139  points[i].y() = y;
140  }
141 
142  if (dumpFilename.size() > 0) {
143  ofstream fout(dumpFilename.c_str());
144  for (int i = 0; i < numPoints; ++i)
145  fout << points[i].transpose() << endl;
146  }
147 
148  // some handy typedefs
151 
152  // setup the solver
153  g2o::SparseOptimizer optimizer;
154  optimizer.setVerbose(false);
155  MyLinearSolver* linearSolver = new MyLinearSolver();
156  MyBlockSolver* solver_ptr = new MyBlockSolver(linearSolver);
158  optimizer.setAlgorithm(solver);
159 
160  // build the optimization problem given the points
161  // 1. add the parameter vertex
162  VertexParams* params = new VertexParams();
163  params->setId(0);
164  params->setEstimate(Eigen::Vector3d(1,1,1)); // some initial value for the params
165  optimizer.addVertex(params);
166  // 2. add the points we measured to be on the curve
167  for (int i = 0; i < numPoints; ++i) {
169  e->setInformation(Eigen::Matrix<double, 1, 1>::Identity());
170  e->setVertex(0, params);
171  e->setMeasurement(points[i]);
172  optimizer.addEdge(e);
173  }
174 
175  // perform the optimization
176  optimizer.initializeOptimization();
177  optimizer.setVerbose(verbose);
178  optimizer.optimize(maxIterations);
179 
180  if (verbose)
181  cout << endl;
182 
183  // print out the result
184  cout << "Target curve" << endl;
185  cout << "a * exp(-lambda * x) + b" << endl;
186  cout << "Iterative least squares solution" << endl;
187  cout << "a = " << params->estimate()(0) << endl;
188  cout << "b = " << params->estimate()(1) << endl;
189  cout << "lambda = " << params->estimate()(2) << endl;
190  cout << endl;
191 
192  // clean up
193  delete[] points;
194 
195  return 0;
196 }
virtual bool read(std::istream &)
read the vertex from a stream, i.e., the internal state of the vertex
Definition: curve_fit.cpp:91
#define __PRETTY_FUNCTION__
Definition: macros.h:89
Command line parsing of argc and argv.
Definition: command_args.h:46
measurement for a point on the curve
Definition: curve_fit.cpp:84
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Definition: curve_fit.cpp:48
static double gaussRand(double mean, double sigma)
Definition: sampler.h:86
Implementation of the Levenberg Algorithm.
int optimize(int iterations, bool online=false)
Templatized BaseVertex.
Definition: base_vertex.h:50
virtual void setMeasurement(const Measurement &m)
Definition: base_edge.h:76
linear solver using dense cholesky decomposition
int main(int argc, char **argv)
Definition: curve_fit.cpp:113
void setVertex(size_t i, Vertex *v)
Definition: hyper_graph.h:194
virtual bool read(std::istream &)
read the vertex from a stream, i.e., the internal state of the vertex
Definition: curve_fit.cpp:53
virtual void setId(int id)
sets the id of the node in the graph be sure that the graph keeps consistent after changing the id ...
static double uniformRand(double lowerBndr, double upperBndr)
Definition: sampler.h:100
bool parseArgs(int argc, char **argv, bool exitOnError=true)
void setVerbose(bool verbose)
virtual bool write(std::ostream &) const
write the vertex to a stream
Definition: curve_fit.cpp:59
void setAlgorithm(OptimizationAlgorithm *algorithm)
EIGEN_STRONG_INLINE void setInformation(const InformationType &information)
Definition: base_edge.h:69
virtual void oplusImpl(const double *update)
Definition: curve_fit.cpp:70
EIGEN_MAKE_ALIGNED_OPERATOR_NEW EdgePointOnCurve()
Definition: curve_fit.cpp:88
void param(const std::string &name, bool &p, bool defValue, const std::string &desc)
void setEstimate(const EstimateType &et)
set the estimate for the vertex also calls updateCache()
Definition: base_vertex.h:101
the params, a, b, and lambda for a * exp(-lambda * t) + b
Definition: curve_fit.cpp:45
virtual bool write(std::ostream &) const
write the vertex to a stream
Definition: curve_fit.cpp:96
const EstimateType & estimate() const
return the current estimate of the vertex
Definition: base_vertex.h:99
Implementation of a solver operating on the blocks of the Hessian.
Definition: block_solver.h:96
virtual void setToOriginImpl()
sets the node to the origin (used in the multilevel stuff)
Definition: curve_fit.cpp:65
virtual bool initializeOptimization(HyperGraph::EdgeSet &eset)
virtual bool addEdge(HyperGraph::Edge *e)
virtual bool addVertex(HyperGraph::Vertex *v, Data *userData)