g2o
Loading...
Searching...
No Matches
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
36#include "g2o/stuff/sampler.h"
37
38using namespace std;
39
41
45class VertexParams : public g2o::BaseVertex<3, Eigen::Vector3d> {
46 public:
49
50 bool read(std::istream& /*is*/) override { return false; }
51
52 bool write(std::ostream& /*os*/) const override { return false; }
53
54 void setToOriginImpl() override {}
55
56 void oplusImpl(const double* update) override {
57 Eigen::Vector3d::ConstMapType v(update);
58 _estimate += v;
59 }
60};
61
70 : public g2o::BaseUnaryEdge<1, Eigen::Vector2d, VertexParams> {
71 public:
72 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
74 bool read(std::istream& /*is*/) override {
75 cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
76 return false;
77 }
78 bool write(std::ostream& /*os*/) const override {
79 cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
80 return false;
81 }
82
83 template <typename T>
84 bool operator()(const T* params, T* error) const {
85 const T& a = params[0];
86 const T& b = params[1];
87 const T& lambda = params[2];
88 T fval = a * exp(-lambda * T(measurement()(0))) + b;
89 error[0] = fval - measurement()(1);
90 return true;
91 }
92
93 G2O_MAKE_AUTO_AD_FUNCTIONS // use autodiff
94};
95
96int main(int argc, char** argv) {
97 int numPoints;
98 int maxIterations;
99 bool verbose;
100 string dumpFilename;
102 arg.param("dump", dumpFilename, "", "dump the points into a file");
103 arg.param("numPoints", numPoints, 50,
104 "number of points sampled from the curve");
105 arg.param("i", maxIterations, 10, "perform n iterations");
106 arg.param("v", verbose, false, "verbose output of the optimization process");
107
108 arg.parseArgs(argc, argv);
109
110 // generate random data
112 double a = 2.;
113 double b = 0.4;
114 double lambda = 0.2;
115 Eigen::Vector2d* points = new Eigen::Vector2d[numPoints];
116 for (int i = 0; i < numPoints; ++i) {
117 double x = g2o::Sampler::uniformRand(0, 10);
118 double y = a * exp(-lambda * x) + b;
119 // add Gaussian noise
120 y += g2o::Sampler::gaussRand(0, 0.02);
121 points[i].x() = x;
122 points[i].y() = y;
123 }
124
125 if (dumpFilename.size() > 0) {
126 ofstream fout(dumpFilename.c_str());
127 for (int i = 0; i < numPoints; ++i) fout << points[i].transpose() << endl;
128 }
129
130 // setup the solver
131 g2o::SparseOptimizer optimizer;
132 optimizer.setVerbose(false);
133
134 // allocate the solver
136 optimizer.setAlgorithm(
137 g2o::OptimizationAlgorithmFactory::instance()->construct("lm_dense",
138 solverProperty));
139
140 // build the optimization problem given the points
141 // 1. add the parameter vertex
142 VertexParams* params = new VertexParams();
143 params->setId(0);
144 params->setEstimate(
145 Eigen::Vector3d(1, 1, 1)); // some initial value for the params
146 optimizer.addVertex(params);
147 // 2. add the points we measured to be on the curve
148 for (int i = 0; i < numPoints; ++i) {
150 e->setInformation(Eigen::Matrix<double, 1, 1>::Identity());
151 e->setVertex(0, params);
152 e->setMeasurement(points[i]);
153 optimizer.addEdge(e);
154 }
155
156 // perform the optimization
157 optimizer.initializeOptimization();
158 optimizer.setVerbose(verbose);
159 optimizer.optimize(maxIterations);
160
161 if (verbose) cout << endl;
162
163 // print out the result
164 cout << "Target curve" << endl;
165 cout << "a * exp(-lambda * x) + b" << endl;
166 cout << "Iterative least squares solution" << endl;
167 cout << "a = " << params->estimate()(0) << endl;
168 cout << "b = " << params->estimate()(1) << endl;
169 cout << "lambda = " << params->estimate()(2) << endl;
170 cout << endl;
171
172 // clean up
173 delete[] points;
174
175 return 0;
176}
#define G2O_MAKE_AUTO_AD_FUNCTIONS
measurement for a point on the curve
Definition curve_fit.cpp:70
bool read(std::istream &) override
read the vertex from a stream, i.e., the internal state of the vertex
Definition curve_fit.cpp:74
bool write(std::ostream &) const override
write the vertex to a stream
Definition curve_fit.cpp:78
bool operator()(const T *params, T *error) const
Definition curve_fit.cpp:84
EIGEN_MAKE_ALIGNED_OPERATOR_NEW EdgePointOnCurve()
Definition curve_fit.cpp:73
the params, a, b, and lambda for a * exp(-lambda * t) + b
Definition curve_fit.cpp:45
void oplusImpl(const double *update) override
Definition curve_fit.cpp:56
bool read(std::istream &) override
read the vertex from a stream, i.e., the internal state of the vertex
Definition curve_fit.cpp:50
bool write(std::ostream &) const override
write the vertex to a stream
Definition curve_fit.cpp:52
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Definition curve_fit.cpp:47
void setToOriginImpl() override
sets the node to the origin (used in the multilevel stuff)
Definition curve_fit.cpp:54
EIGEN_STRONG_INLINE const Measurement & measurement() const
accessor functions for the measurement represented by the edge
Definition base_edge.h:119
virtual void setMeasurement(const Measurement &m)
Definition base_edge.h:122
void setInformation(const InformationType &information)
Definition base_edge.h:111
const ErrorVector & error() const
Definition base_edge.h:103
Templatized BaseVertex.
Definition base_vertex.h:51
const EstimateType & estimate() const
return the current estimate of the vertex
void setEstimate(const EstimateType &et)
set the estimate for the vertex also calls updateCache()
Command line parsing of argc and argv.
bool parseArgs(int argc, char **argv, bool exitOnError=true)
void param(const std::string &name, bool &p, bool defValue, const std::string &desc)
void setVertex(size_t i, Vertex *v)
static OptimizationAlgorithmFactory * instance()
return the instance
static void seedRand()
Definition sampler.h:109
static double gaussRand(double mean, double sigma)
Definition sampler.h:89
static double uniformRand(double lowerBndr, double upperBndr)
Definition sampler.h:102
int optimize(int iterations, bool online=false)
void setVerbose(bool verbose)
virtual bool initializeOptimization(HyperGraph::EdgeSet &eset)
void setAlgorithm(OptimizationAlgorithm *algorithm)
int main()
Definition gicp_demo.cpp:44
#define __PRETTY_FUNCTION__
Definition macros.h:90
Definition jet.h:876
#define G2O_USE_OPTIMIZATION_LIBRARY(libraryname)
virtual bool addEdge(HyperGraph::Edge *e)
virtual bool addVertex(HyperGraph::Vertex *v, Data *userData)