g2o
Loading...
Searching...
No Matches
circle_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 <Eigen/Geometry>
29#include <iostream>
30
37#include "g2o/stuff/sampler.h"
38
39using namespace std;
40
42
43double errorOfSolution(int numPoints, Eigen::Vector2d* points,
44 const Eigen::Vector3d& circle) {
45 Eigen::Vector2d center = circle.head<2>();
46 double radius = circle(2);
47 double error = 0.;
48 for (int i = 0; i < numPoints; ++i) {
49 double d = (points[i] - center).norm() - radius;
50 error += d * d;
51 }
52 return error;
53}
54
58class VertexCircle : public g2o::BaseVertex<3, Eigen::Vector3d> {
59 public:
62
63 bool read(std::istream& /*is*/) override { return false; }
64
65 bool write(std::ostream& /*os*/) const override { return false; }
66
67 void setToOriginImpl() override {
68 cerr << __PRETTY_FUNCTION__ << " not implemented yet" << endl;
69 }
70
71 void oplusImpl(const double* update) override {
72 Eigen::Vector3d::ConstMapType v(update);
73 _estimate += v;
74 }
75};
76
85 : public g2o::BaseUnaryEdge<1, Eigen::Vector2d, VertexCircle> {
86 public:
87 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
89 bool read(std::istream& /*is*/) override { return false; }
90 bool write(std::ostream& /*os*/) const override { return false; }
91
92 template <typename T>
93 bool operator()(const T* circle, T* error) const {
94 typename g2o::VectorN<2, T>::ConstMapType center(circle);
95 const T& radius = circle[2];
96
97 error[0] = (measurement().cast<T>() - center).norm() - radius;
98 return true;
99 }
100
101 G2O_MAKE_AUTO_AD_FUNCTIONS // use autodiff
102};
103
104int main(int argc, char** argv) {
105 int numPoints;
106 int maxIterations;
107 bool verbose;
109 arg.param("numPoints", numPoints, 100,
110 "number of points sampled from the circle");
111 arg.param("i", maxIterations, 10, "perform n iterations");
112 arg.param("v", verbose, false, "verbose output of the optimization process");
113 arg.parseArgs(argc, argv);
114
115 // generate random data
116 Eigen::Vector2d center(4.0, 2.0);
117 double radius = 2.0;
118 Eigen::Vector2d* points = new Eigen::Vector2d[numPoints];
119
121 for (int i = 0; i < numPoints; ++i) {
122 double r = g2o::Sampler::gaussRand(radius, 0.05);
123 double angle = g2o::Sampler::uniformRand(0.0, 2.0 * M_PI);
124 points[i].x() = center.x() + r * cos(angle);
125 points[i].y() = center.y() + r * sin(angle);
126 }
127
128 // setup the solver
129 g2o::SparseOptimizer optimizer;
130 optimizer.setVerbose(false);
131
132 // allocate the solver
134 optimizer.setAlgorithm(
135 g2o::OptimizationAlgorithmFactory::instance()->construct("lm_dense",
136 solverProperty));
137
138 // build the optimization problem given the points
139 // 1. add the circle vertex
140 VertexCircle* circle = new VertexCircle();
141 circle->setId(0);
142 circle->setEstimate(
143 Eigen::Vector3d(3.0, 3.0, 3.0)); // some initial value for the circle
144 optimizer.addVertex(circle);
145 // 2. add the points we measured
146 for (int i = 0; i < numPoints; ++i) {
148 e->setInformation(Eigen::Matrix<double, 1, 1>::Identity());
149 e->setVertex(0, circle);
150 e->setMeasurement(points[i]);
151 optimizer.addEdge(e);
152 }
153
154 // perform the optimization
155 optimizer.initializeOptimization();
156 optimizer.setVerbose(verbose);
157 optimizer.optimize(maxIterations);
158
159 if (verbose) cout << endl;
160
161 // print out the result
162 cout << "Iterative least squares solution" << endl;
163 cout << "center of the circle " << circle->estimate().head<2>().transpose()
164 << endl;
165 cout << "radius of the cirlce " << circle->estimate()(2) << endl;
166 cout << "error " << errorOfSolution(numPoints, points, circle->estimate())
167 << endl;
168 cout << endl;
169
170 // solve by linear least squares
171 // Let (a, b) be the center of the circle and r the radius of the circle.
172 // For a point (x, y) on the circle we have:
173 // (x - a)^2 + (y - b)^2 = r^2
174 // This leads to
175 // (-2x -2y 1)^T * (a b c) = -x^2 - y^2 (1)
176 // where c = a^2 + b^2 - r^2.
177 // Since we have a bunch of points, we accumulate Eqn (1) in a matrix and
178 // compute the normal equation to obtain a solution for (a b c).
179 // Afterwards the radius r is recovered.
180 Eigen::MatrixXd A(numPoints, 3);
181 Eigen::VectorXd b(numPoints);
182 for (int i = 0; i < numPoints; ++i) {
183 A(i, 0) = -2 * points[i].x();
184 A(i, 1) = -2 * points[i].y();
185 A(i, 2) = 1;
186 b(i) = -pow(points[i].x(), 2) - pow(points[i].y(), 2);
187 }
188 Eigen::Vector3d solution =
189 (A.transpose() * A).ldlt().solve(A.transpose() * b);
190 // calculate the radius of the circle given the solution so far
191 solution(2) = sqrt(pow(solution(0), 2) + pow(solution(1), 2) - solution(2));
192 cout << "Linear least squares solution" << endl;
193 cout << "center of the circle " << solution.head<2>().transpose() << endl;
194 cout << "radius of the cirlce " << solution(2) << endl;
195 cout << "error " << errorOfSolution(numPoints, points, solution) << endl;
196
197 // clean up
198 delete[] points;
199
200 return 0;
201}
#define G2O_MAKE_AUTO_AD_FUNCTIONS
double errorOfSolution(int numPoints, Eigen::Vector2d *points, const Eigen::Vector3d &circle)
measurement for a point on the circle
bool operator()(const T *circle, T *error) const
bool write(std::ostream &) const override
write the vertex to a stream
EIGEN_MAKE_ALIGNED_OPERATOR_NEW EdgePointOnCircle()
bool read(std::istream &) override
read the vertex from a stream, i.e., the internal state of the vertex
a circle located at x,y with radius r
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
void oplusImpl(const double *update) override
void setToOriginImpl() override
sets the node to the origin (used in the multilevel stuff)
bool write(std::ostream &) const override
write the vertex to a stream
bool read(std::istream &) override
read the vertex from a stream, i.e., the internal state of the vertex
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
Eigen::Matrix< T, N, 1, Eigen::ColMajor > VectorN
Definition eigen_types.h:49
Definition jet.h:876
#define G2O_USE_OPTIMIZATION_LIBRARY(libraryname)
virtual bool addEdge(HyperGraph::Edge *e)
virtual bool addVertex(HyperGraph::Vertex *v, Data *userData)