g2o
Loading...
Searching...
No Matches
optimization_algorithm_levenberg.cpp
Go to the documentation of this file.
1// g2o - General Graph Optimization
2// Copyright (C) 2011 R. Kuemmerle, G. Grisetti, W. Burgard
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
28
29#include <cassert>
30#include <iostream>
31
32#include "batch_stats.h"
33#include "g2o/stuff/logger.h"
34#include "g2o/stuff/timeutil.h"
35#include "solver.h"
36#include "sparse_optimizer.h"
37using namespace std;
38
39namespace g2o {
40
42 std::unique_ptr<Solver> solver)
43 : OptimizationAlgorithmWithHessian(*solver.get()),
44 _currentLambda(cst(-1.)),
45 _tau(cst(1e-5)),
46 _goodStepLowerScale(cst(1. / 3.)),
47 _goodStepUpperScale(cst(2. / 3.)),
48 _ni(cst(2.)),
49 _levenbergIterations(0),
50 m_solver{std::move(solver)} {
52 _properties.makeProperty<Property<double> >("initialLambda", 0.);
54 _properties.makeProperty<Property<int> >("maxTrialsAfterFailure", 10);
55}
56
58
59OptimizationAlgorithm::SolverResult OptimizationAlgorithmLevenberg::solve(
60 int iteration, bool online) {
61 assert(_optimizer && "_optimizer not set");
62 assert(_solver.optimizer() == _optimizer &&
63 "underlying linear solver operates on different graph");
64
65 if (iteration == 0 &&
66 !online) { // built up the CCS structure, here due to easy time measure
67 bool ok = _solver.buildStructure();
68 if (!ok) {
69 G2O_WARN("{}: Failure while building CCS structure", __PRETTY_FUNCTION__);
70 return OptimizationAlgorithm::Fail;
71 }
72 }
73
74 double t = get_monotonic_time();
77 if (globalStats) {
78 globalStats->timeResiduals = get_monotonic_time() - t;
80 }
81
82 double currentChi = _optimizer->activeRobustChi2();
83
85 if (globalStats) {
86 globalStats->timeQuadraticForm = get_monotonic_time() - t;
87 }
88
89 // core part of the Levenbarg algorithm
90 if (iteration == 0) {
92 _ni = 2;
93 }
94
95 double rho = 0;
96 int& qmax = _levenbergIterations;
97 qmax = 0;
98 do {
100 if (globalStats) {
101 globalStats->levenbergIterations++;
102 t = get_monotonic_time();
103 }
104 // update the diagonal of the system matrix
106 bool ok2 = _solver.solve();
107 if (globalStats) {
108 globalStats->timeLinearSolution += get_monotonic_time() - t;
109 t = get_monotonic_time();
110 }
112 if (globalStats) {
113 globalStats->timeUpdate = get_monotonic_time() - t;
114 }
115
116 // restore the diagonal
118
120 double tempChi = _optimizer->activeRobustChi2();
121
122 if (!ok2) tempChi = std::numeric_limits<double>::max();
123
124 rho = (currentChi - tempChi);
125 double scale =
126 ok2 ? computeScale() + cst(1e-3) : 1; // make sure it's non-zero :)
127 rho /= scale;
128
129 if (rho > 0 && g2o_isfinite(tempChi) && ok2) { // last step was good
130 double alpha = 1. - pow((2 * rho - 1), 3);
131 // crop lambda between minimum and maximum factors
132 alpha = (std::min)(alpha, _goodStepUpperScale);
133 double scaleFactor = (std::max)(_goodStepLowerScale, alpha);
134 _currentLambda *= scaleFactor;
135 _ni = 2;
136 currentChi = tempChi;
138 } else {
140 _ni *= 2;
141 _optimizer->pop(); // restore the last state before trying to optimize
142 if (!g2o_isfinite(_currentLambda)) break;
143 }
144 qmax++;
145 } while (rho < 0 && qmax < _maxTrialsAfterFailure->value() &&
147
148 if (qmax == _maxTrialsAfterFailure->value() || rho == 0 ||
150 return Terminate;
151 return OK;
152}
153
155 if (_userLambdaInit->value() > 0) return _userLambdaInit->value();
156 double maxDiagonal = 0;
157 for (size_t k = 0; k < _optimizer->indexMapping().size(); k++) {
159 assert(v);
160 int dim = v->dimension();
161 for (int j = 0; j < dim; ++j) {
162 maxDiagonal = std::max(fabs(v->hessian(j, j)), maxDiagonal);
163 }
164 }
165 return _tau * maxDiagonal;
166}
167
169 double scale = 0;
170 for (size_t j = 0; j < _solver.vectorSize(); j++) {
171 scale +=
172 _solver.x()[j] * (_currentLambda * _solver.x()[j] + _solver.b()[j]);
173 }
174 return scale;
175}
176
180
184
186 os << "\t schur= " << _solver.schur()
187 << "\t lambda= " << FIXED(_currentLambda)
188 << "\t levenbergIter= " << _levenbergIterations;
189}
190
191} // namespace g2o
A general case Vertex for optimization.
virtual const double & hessian(int i, int j) const =0
get the element from the hessian matrix
int dimension() const
dimension of the estimated state belonging to this node
OptimizationAlgorithmLevenberg(std::unique_ptr< Solver > solver)
virtual void printVerbose(std::ostream &os) const
virtual SolverResult solve(int iteration, bool online=false)
Base for solvers operating on the approximated Hessian, e.g., Gauss-Newton, Levenberg.
SparseOptimizer * _optimizer
the optimizer the solver is working on
P * makeProperty(const std::string &name_, const typename P::ValueType &v)
Definition property.h:116
void setValue(const T &v)
Definition property.h:58
const T & value() const
Definition property.h:59
double * b()
return b, the right hand side of the system
Definition solver.h:101
virtual void restoreDiagonal()=0
virtual bool buildStructure(bool zeroBlocks=false)=0
size_t vectorSize() const
return the size of the solution vector (x) and b
Definition solver.h:105
double * x()
return x, the solution vector
Definition solver.h:98
virtual bool setLambda(double lambda, bool backup=false)=0
virtual bool solve()=0
virtual bool buildSystem()=0
virtual bool schur()=0
should the solver perform the schur complement or not
SparseOptimizer * optimizer() const
the optimizer (graph) on which the solver works
Definition solver.h:108
void push(SparseOptimizer::VertexContainer &vlist)
push the estimate of a subset of the variables onto a stack
void pop(SparseOptimizer::VertexContainer &vlist)
pop (restore) the estimate a subset of the variables from the stack
void update(const double *update)
double activeRobustChi2() const
const VertexContainer & indexMapping() const
the index mapping of the vertices
void discardTop(SparseOptimizer::VertexContainer &vlist)
bool terminate()
if external stop flag is given, return its state. False otherwise
#define G2O_WARN(...)
Definition logger.h:88
#define g2o_isfinite(x)
Definition macros.h:102
#define __PRETTY_FUNCTION__
Definition macros.h:90
constexpr double cst(long double v)
Definition misc.h:60
double get_monotonic_time()
Definition timeutil.cpp:43
Definition jet.h:876
statistics about the optimization
Definition batch_stats.h:40
double timeResiduals
residuals
Definition batch_stats.h:49
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77
double timeUpdate
time to apply the update
Definition batch_stats.h:64
int levenbergIterations
number of iterations performed by LM
Definition batch_stats.h:52
double timeQuadraticForm
construct the quadratic form in the graph
Definition batch_stats.h:51
utility functions for handling time related stuff