g2o
Loading...
Searching...
No Matches
static_dynamic_function_fit.cpp
Go to the documentation of this file.
1// This example illustrates how to mix static and dynamic vertices in
2// a graph.
3
4// The goal is to fit a function of the form y(x) = f(x) + x^3 * p(x)
5// to a set of data:
6// * f(x) is a quadratic
7// * p(x) is a polynomial whose dimensions are established at runtime.
8
9// The ith observation consists of m_i pairs of values
10// Z_i={(x_1,z_1),(x_2,z_2),(x_m_i,z_m_i)}, where z_i=y(x_i)+w_i, where
11// w_i is additive white noise with information matrix Omega.
12
13// In this example both the measurement edges and one vertex can
14// changed dynamically.
15
16#include <random>
17#include <unsupported/Eigen/Polynomials>
18
27#include "g2o/stuff/sampler.h"
28
29// Declare the custom types used in the graph
30
31// This vertex stores the coefficients of the f(x) polynomial. This is
32// quadratic, and always has a degree of three.
33
35 : public g2o::BaseVertex<3, Eigen::Vector3d> {
36 public:
38
39 // Create the vertex
41
42 // Read the vertex
43 virtual bool read(std::istream& is) {
44 // Read the state
46 }
47
48 // Write the vertex
49 virtual bool write(std::ostream& os) const {
51 }
52
53 // Reset to zero
54 virtual void setToOriginImpl() { _estimate.setZero(); }
55
56 // Direct linear add
57 virtual void oplusImpl(const double* update) {
58 Eigen::Vector3d::ConstMapType v(update, 3);
59 _estimate += v;
60 }
61};
62
63// This vertex stores the coefficients of the p(x) polynomial. It is dynamic
64// because we can change it at runtime.
65
67 : public g2o::BaseDynamicVertex<Eigen::VectorXd> {
68 public:
70
71 // Create the vertex
73
74 // Read the vertex
75 virtual bool read(std::istream& is) {
76 // Read the dimension
77 int dimension;
78 is >> dimension;
79 if (is.good() == false) {
80 return false;
81 }
82
83 // Set the dimension; we call the method here to ensure stuff like
84 // cache and the workspace is setup
86
87 // Read the state
89 }
90
91 // Write the vertex
92 virtual bool write(std::ostream& os) const {
93 os << _estimate.size() << " ";
95 }
96
97 // Reset to zero
98 virtual void setToOriginImpl() { _estimate.setZero(); }
99
100 // Direct linear add
101 virtual void oplusImpl(const double* update) {
102 Eigen::VectorXd::ConstMapType v(update, _dimension);
103 _estimate += v;
104 }
105
106 // Resize the vertex state. In this case, we simply trash whatever
107 // was there before.
108 virtual bool setDimensionImpl(int newDimension) {
109 _estimate.resize(newDimension);
110 _estimate.setZero();
111 return true;
112 }
113};
114
115// Helper structure
116
118 Eigen::VectorXd x;
119 Eigen::VectorXd z;
120};
121
122// The edge which encodes the observations
123
125 : public g2o::BaseBinaryEdge<Eigen::Dynamic, Eigen::VectorXd,
126 FPolynomialCoefficientVertex,
127 PPolynomialCoefficientVertex> {
128 public:
130
131 MultipleValueEdge(const FunctionObservation& obs, double omega) : _x(obs.x) {
132 setDimension(obs.z.size());
133 setMeasurement(obs.z);
134 InformationType I = Eigen::MatrixXd::Identity(_x.size(), _x.size()) * omega;
136 }
137
138 virtual bool read(std::istream& is) {
139 Eigen::VectorXd z;
143
144 return readInformationMatrix(is);
145 }
146
147 virtual bool write(std::ostream& os) const {
150 return writeInformationMatrix(os);
151 }
152
153 // Compute the measurement from the eigen polynomial module
154 virtual void computeError() {
155 const FPolynomialCoefficientVertex* fvertex =
156 dynamic_cast<const FPolynomialCoefficientVertex*>(_vertices[0]);
157 const PPolynomialCoefficientVertex* pvertex =
158 dynamic_cast<const PPolynomialCoefficientVertex*>(_vertices[1]);
159 for (int i = 0; i < _measurement.size(); ++i) {
160 double x3 = pow(_x[i], 3);
161 _error[i] = _measurement[i] -
162 Eigen::poly_eval(fvertex->estimate(), _x[i]) -
163 x3 * (Eigen::poly_eval(pvertex->estimate(), _x[i]));
164 }
165 }
166
167 private:
168 // The points that the polynomial is computed at
169 Eigen::VectorXd _x;
170};
171
172int main(int argc, const char* argv[]) {
173 // Random number generator
174 std::default_random_engine generator;
175
176 // Create the coefficients for the f-polynomial (all drawn randomly)
177 Eigen::Vector3d f;
178 for (int i = 0; i < 3; ++i) {
179 f[i] = g2o::sampleUniform(-1, 1);
180 }
181
182 // Number of dimensions of the polynomial; the default is 4
183 int polynomialDimension = 4;
184 if (argc > 1) {
185 polynomialDimension = atoi(argv[1]);
186 }
187
188 // Create the coefficients for the polynomial (all drawn randomly)
189 Eigen::VectorXd p(polynomialDimension);
190 for (int i = 0; i < polynomialDimension; ++i) {
191 p[i] = g2o::sampleUniform(-1, 1);
192 }
193
194 std::cout << "Ground truth vectors f=" << f.transpose()
195 << "; p=" << p.transpose() << std::endl;
196
197 // The number of observations in the polynomial; the defaultis 6
198 int obs = 6;
199 if (argc > 2) {
200 obs = atoi(argv[2]);
201 }
202
203 // The number of observations will be randomly sampled each time.
204
205 // Sample the observations. This is a set of function
206 // observations. The cardinality of each observation set is random.
207 double sigmaZ = 0.1;
208 std::vector<FunctionObservation> observations(obs);
209 std::uniform_int_distribution<int> cardinalitySampler(1, 5);
210
211 for (int i = 0; i < obs; ++i) {
212 FunctionObservation& fo = observations[i];
213 int numObs = cardinalitySampler(generator);
214 fo.x.resize(numObs);
215 fo.z.resize(numObs);
216 for (int o = 0; o < numObs; ++o) {
217 fo.x[o] = g2o::sampleUniform(-5, 5);
218 double x3 = pow(fo.x[o], 3);
219 fo.z[o] = Eigen::poly_eval(f, fo.x[o]) +
220 x3 * (Eigen::poly_eval(p, fo.x[o])) +
221 sigmaZ * g2o::sampleGaussian();
222 }
223 }
224
225 // Construct the graph and set up the solver and optimiser
226 std::unique_ptr<g2o::BlockSolverX::LinearSolverType> linearSolver =
227 std::make_unique<
229
230 // Set up the solver
231 std::unique_ptr<g2o::BlockSolverX> blockSolver =
232 std::make_unique<g2o::BlockSolverX>(std::move(linearSolver));
233
234 // Set up the optimisation algorithm
235 g2o::OptimizationAlgorithm* optimisationAlgorithm =
236 new g2o::OptimizationAlgorithmLevenberg(std::move(blockSolver));
237
238 // Create the graph and configure it
239 std::unique_ptr<g2o::SparseOptimizer> optimizer =
240 std::make_unique<g2o::SparseOptimizer>();
241 optimizer->setVerbose(true);
242 optimizer->setAlgorithm(optimisationAlgorithm);
243
244 // Create the f vertex; its dimensions are known
246 pf->setId(0);
247 optimizer->addVertex(pf);
248
249 // Create the vertex; note its dimension is currently is undefined
251 pv->setId(1);
252 optimizer->addVertex(pv);
253
254 // Create the information matrix
255 double omega = 1 / (sigmaZ * sigmaZ);
256
257 // Create the edges
258 for (int i = 0; i < obs; ++i) {
259 MultipleValueEdge* mve = new MultipleValueEdge(observations[i], omega);
260 mve->setVertex(0, pf);
261 mve->setVertex(1, pv);
262 optimizer->addEdge(mve);
263 }
264
265 // Now run the same optimization problem for different choices of
266 // dimension of the polynomial vertex. This shows how we can
267 // dynamically change the vertex dimensions in an alreacy
268 // constructed graph. Note that you must call initializeOptimization
269 // before you can optimize after a state dimension has changed.
270 for (int testDimension = 1; testDimension <= polynomialDimension;
271 ++testDimension) {
272 pv->setDimension(testDimension);
273 optimizer->initializeOptimization();
274 optimizer->optimize(10);
275 std::cout << "Computed parameters: f=" << pf->estimate().transpose()
276 << "; p=" << pv->estimate().transpose() << std::endl;
277 }
278 for (int testDimension = polynomialDimension - 1; testDimension >= 1;
279 --testDimension) {
280 pv->setDimension(testDimension);
281 optimizer->initializeOptimization();
282 optimizer->optimize(10);
283 std::cout << "Computed parameters: f= " << pf->estimate().transpose()
284 << "; p=" << pv->estimate().transpose() << std::endl;
285 }
286}
virtual void oplusImpl(const double *update)
virtual void setToOriginImpl()
sets the node to the origin (used in the multilevel stuff)
virtual bool write(std::ostream &os) const
write the vertex to a stream
virtual bool read(std::istream &is)
read the vertex from a stream, i.e., the internal state of the vertex
virtual bool read(std::istream &is)
read the vertex from a stream, i.e., the internal state of the vertex
MultipleValueEdge(const FunctionObservation &obs, double omega)
virtual bool write(std::ostream &os) const
write the vertex to a stream
virtual void oplusImpl(const double *update)
virtual bool write(std::ostream &os) const
write the vertex to a stream
virtual bool read(std::istream &is)
read the vertex from a stream, i.e., the internal state of the vertex
virtual bool setDimensionImpl(int newDimension)
virtual void setToOriginImpl()
sets the node to the origin (used in the multilevel stuff)
virtual bool setDimension(int newDimension)
bool writeInformationMatrix(std::ostream &os) const
write the upper trinagular part of the information matrix into the stream
Definition base_edge.h:165
bool readInformationMatrix(std::istream &is)
Definition base_edge.h:173
virtual void setMeasurement(const Measurement &m)
Definition base_edge.h:122
void setInformation(const InformationType &information)
Definition base_edge.h:111
std::enable_if< Dim==-1, void >::type setDimension(int dim)
Definition base_edge.h:139
Measurement _measurement
the measurement of the edge
Definition base_edge.h:146
ErrorVector _error
Definition base_edge.h:149
BaseEdge< D, E >::InformationType InformationType
Templatized BaseVertex.
Definition base_vertex.h:51
const EstimateType & estimate() const
return the current estimate of the vertex
void setVertex(size_t i, Vertex *v)
VertexContainer _vertices
linear solver which uses the sparse Cholesky solver from Eigen
int dimension() const
dimension of the estimated state belonging to this node
void setToOrigin()
sets the node to the origin (used in the multilevel stuff)
Implementation of the Levenberg Algorithm.
Generic interface for a non-linear solver operating on a graph.
int main()
Definition gicp_demo.cpp:44
bool writeVector(std::ostream &os, const Eigen::DenseBase< Derived > &b)
Definition io_helper.h:36
bool readVector(std::istream &is, Eigen::DenseBase< Derived > &b)
Definition io_helper.h:42
double sampleUniform(double min, double max, std::mt19937 *generator)
Definition sampler.cpp:35
double sampleGaussian(std::mt19937 *generator)
Definition sampler.cpp:40