g2o
Loading...
Searching...
No Matches
Classes | Functions
circle_fit.cpp File Reference
#include <Eigen/Core>
#include <Eigen/Geometry>
#include <iostream>
#include "g2o/core/auto_differentiation.h"
#include "g2o/core/base_unary_edge.h"
#include "g2o/core/base_vertex.h"
#include "g2o/core/optimization_algorithm_factory.h"
#include "g2o/core/sparse_optimizer.h"
#include "g2o/stuff/command_args.h"
#include "g2o/stuff/sampler.h"
Include dependency graph for circle_fit.cpp:

Go to the source code of this file.

Classes

class  VertexCircle
 a circle located at x,y with radius r More...
 
class  EdgePointOnCircle
 measurement for a point on the circle More...
 

Functions

 G2O_USE_OPTIMIZATION_LIBRARY (dense)
 
double errorOfSolution (int numPoints, Eigen::Vector2d *points, const Eigen::Vector3d &circle)
 
int main (int argc, char **argv)
 

Function Documentation

◆ errorOfSolution()

double errorOfSolution ( int  numPoints,
Eigen::Vector2d *  points,
const Eigen::Vector3d &  circle 
)

Definition at line 43 of file circle_fit.cpp.

44 {
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}

Referenced by main().

◆ G2O_USE_OPTIMIZATION_LIBRARY()

G2O_USE_OPTIMIZATION_LIBRARY ( dense  )

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 104 of file circle_fit.cpp.

104 {
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}
double errorOfSolution(int numPoints, Eigen::Vector2d *points, const Eigen::Vector3d &circle)
measurement for a point on the circle
a circle located at x,y with radius r
virtual void setMeasurement(const Measurement &m)
Definition base_edge.h:122
void setInformation(const InformationType &information)
Definition base_edge.h:111
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)
Jet< T, N > pow(const Jet< T, N > &f, double g)
Definition jet.h:743
Jet< T, N > cos(const Jet< T, N > &f)
Definition jet.h:452
Jet< T, N > sin(const Jet< T, N > &f)
Definition jet.h:465
Jet< T, N > sqrt(const Jet< T, N > &f)
Definition jet.h:444
virtual bool addEdge(HyperGraph::Edge *e)
virtual bool addVertex(HyperGraph::Vertex *v, Data *userData)

References g2o::OptimizableGraph::addEdge(), g2o::OptimizableGraph::addVertex(), errorOfSolution(), g2o::BaseVertex< D, T >::estimate(), g2o::Sampler::gaussRand(), g2o::SparseOptimizer::initializeOptimization(), g2o::OptimizationAlgorithmFactory::instance(), g2o::SparseOptimizer::optimize(), g2o::CommandArgs::param(), g2o::CommandArgs::parseArgs(), g2o::Sampler::seedRand(), g2o::SparseOptimizer::setAlgorithm(), g2o::BaseVertex< D, T >::setEstimate(), g2o::OptimizableGraph::Vertex::setId(), g2o::BaseEdge< D, E >::setInformation(), g2o::BaseEdge< D, E >::setMeasurement(), g2o::SparseOptimizer::setVerbose(), g2o::HyperGraph::Edge::setVertex(), and g2o::Sampler::uniformRand().