g2o
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
g2o::OptimizationAlgorithmLevenberg Class Reference

Implementation of the Levenberg Algorithm. More...

#include <optimization_algorithm_levenberg.h>

Inheritance diagram for g2o::OptimizationAlgorithmLevenberg:
Inheritance graph
[legend]
Collaboration diagram for g2o::OptimizationAlgorithmLevenberg:
Collaboration graph
[legend]

Public Member Functions

 OptimizationAlgorithmLevenberg (std::unique_ptr< Solver > solver)
 
virtual ~OptimizationAlgorithmLevenberg ()
 
virtual SolverResult solve (int iteration, bool online=false)
 
virtual void printVerbose (std::ostream &os) const
 
double currentLambda () const
 return the currently used damping factor
 
void setMaxTrialsAfterFailure (int max_trials)
 
int maxTrialsAfterFailure () const
 get the number of inner iterations for Levenberg-Marquardt
 
double userLambdaInit ()
 
void setUserLambdaInit (double lambda)
 
int levenbergIteration ()
 return the number of levenberg iterations performed in the last round
 
- Public Member Functions inherited from g2o::OptimizationAlgorithmWithHessian
 OptimizationAlgorithmWithHessian (Solver &solver)
 
virtual ~OptimizationAlgorithmWithHessian ()
 
virtual bool init (bool online=false)
 
virtual bool computeMarginals (SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices)
 
virtual bool buildLinearStructure ()
 
virtual void updateLinearSystem ()
 
virtual bool updateStructure (const std::vector< HyperGraph::Vertex * > &vset, const HyperGraph::EdgeSet &edges)
 
Solversolver ()
 return the underlying solver used to solve the linear system
 
virtual void setWriteDebug (bool writeDebug)
 
virtual bool writeDebug () const
 
- Public Member Functions inherited from g2o::OptimizationAlgorithm
 OptimizationAlgorithm ()
 
virtual ~OptimizationAlgorithm ()
 
const SparseOptimizeroptimizer () const
 return the optimizer operating on
 
SparseOptimizeroptimizer ()
 
void setOptimizer (SparseOptimizer *optimizer)
 
const PropertyMapproperties () const
 return the properties of the solver
 
bool updatePropertiesFromString (const std::string &propString)
 
void printProperties (std::ostream &os) const
 

Protected Member Functions

double computeLambdaInit () const
 
double computeScale () const
 

Protected Attributes

Property< int > * _maxTrialsAfterFailure
 
Property< double > * _userLambdaInit
 
double _currentLambda
 
double _tau
 
double _goodStepLowerScale
 
double _goodStepUpperScale
 
double _ni
 
int _levenbergIterations
 
- Protected Attributes inherited from g2o::OptimizationAlgorithmWithHessian
Solver_solver
 
Property< bool > * _writeDebug
 
- Protected Attributes inherited from g2o::OptimizationAlgorithm
SparseOptimizer_optimizer
 the optimizer the solver is working on
 
PropertyMap _properties
 

Private Attributes

std::unique_ptr< Solverm_solver
 

Detailed Description

Implementation of the Levenberg Algorithm.

Definition at line 40 of file optimization_algorithm_levenberg.h.

Constructor & Destructor Documentation

◆ OptimizationAlgorithmLevenberg()

g2o::OptimizationAlgorithmLevenberg::OptimizationAlgorithmLevenberg ( std::unique_ptr< Solver solver)
explicit

construct the Levenberg algorithm, which will use the given Solver for solving the linearized system.

Definition at line 41 of file optimization_algorithm_levenberg.cpp.

44 _currentLambda(cst(-1.)),
45 _tau(cst(1e-5)),
46 _goodStepLowerScale(cst(1. / 3.)),
47 _goodStepUpperScale(cst(2. / 3.)),
48 _ni(cst(2.)),
50 m_solver{std::move(solver)} {
52 _properties.makeProperty<Property<double> >("initialLambda", 0.);
54 _properties.makeProperty<Property<int> >("maxTrialsAfterFailure", 10);
55}
Solver & solver()
return the underlying solver used to solve the linear system
P * makeProperty(const std::string &name_, const typename P::ValueType &v)
Definition property.h:116
constexpr double cst(long double v)
Definition misc.h:60

References _maxTrialsAfterFailure, g2o::OptimizationAlgorithm::_properties, _userLambdaInit, and g2o::PropertyMap::makeProperty().

◆ ~OptimizationAlgorithmLevenberg()

g2o::OptimizationAlgorithmLevenberg::~OptimizationAlgorithmLevenberg ( )
virtual

Definition at line 57 of file optimization_algorithm_levenberg.cpp.

57{}

Member Function Documentation

◆ computeLambdaInit()

double g2o::OptimizationAlgorithmLevenberg::computeLambdaInit ( ) const
protected

helper for Levenberg, this function computes the initial damping factor, if the user did not specify an own value, see setUserLambdaInit()

Definition at line 154 of file optimization_algorithm_levenberg.cpp.

154 {
155 if (_userLambdaInit->value() > 0) return _userLambdaInit->value();
156 double maxDiagonal = 0;
157 for (size_t k = 0; k < _optimizer->indexMapping().size(); k++) {
158 OptimizableGraph::Vertex* v = _optimizer->indexMapping()[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}
SparseOptimizer * _optimizer
the optimizer the solver is working on
const T & value() const
Definition property.h:59
const VertexContainer & indexMapping() const
the index mapping of the vertices

References g2o::OptimizationAlgorithm::_optimizer, _tau, _userLambdaInit, g2o::OptimizableGraph::Vertex::dimension(), g2o::OptimizableGraph::Vertex::hessian(), g2o::SparseOptimizer::indexMapping(), and g2o::Property< T >::value().

Referenced by solve().

◆ computeScale()

double g2o::OptimizationAlgorithmLevenberg::computeScale ( ) const
protected

Definition at line 168 of file optimization_algorithm_levenberg.cpp.

168 {
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}
double * b()
return b, the right hand side of the system
Definition solver.h:101
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

References _currentLambda, g2o::OptimizationAlgorithmWithHessian::_solver, g2o::Solver::b(), g2o::Solver::vectorSize(), and g2o::Solver::x().

Referenced by solve().

◆ currentLambda()

double g2o::OptimizationAlgorithmLevenberg::currentLambda ( ) const
inline

return the currently used damping factor

Definition at line 55 of file optimization_algorithm_levenberg.h.

55{ return _currentLambda; }

◆ levenbergIteration()

int g2o::OptimizationAlgorithmLevenberg::levenbergIteration ( )
inline

return the number of levenberg iterations performed in the last round

Definition at line 72 of file optimization_algorithm_levenberg.h.

72{ return _levenbergIterations; }

◆ maxTrialsAfterFailure()

int g2o::OptimizationAlgorithmLevenberg::maxTrialsAfterFailure ( ) const
inline

get the number of inner iterations for Levenberg-Marquardt

Definition at line 62 of file optimization_algorithm_levenberg.h.

62{ return _maxTrialsAfterFailure->value(); }

◆ printVerbose()

void g2o::OptimizationAlgorithmLevenberg::printVerbose ( std::ostream &  os) const
virtual

called by the optimizer if verbose. re-implement, if you want to print something

Reimplemented from g2o::OptimizationAlgorithm.

Definition at line 185 of file optimization_algorithm_levenberg.cpp.

185 {
186 os << "\t schur= " << _solver.schur()
187 << "\t lambda= " << FIXED(_currentLambda)
188 << "\t levenbergIter= " << _levenbergIterations;
189}
virtual bool schur()=0
should the solver perform the schur complement or not

References _currentLambda, _levenbergIterations, g2o::OptimizationAlgorithmWithHessian::_solver, and g2o::Solver::schur().

◆ setMaxTrialsAfterFailure()

void g2o::OptimizationAlgorithmLevenberg::setMaxTrialsAfterFailure ( int  max_trials)

the number of internal iteration if an update step increases chi^2 within Levenberg-Marquardt

Definition at line 177 of file optimization_algorithm_levenberg.cpp.

177 {
178 _maxTrialsAfterFailure->setValue(max_trials);
179}
void setValue(const T &v)
Definition property.h:58

References _maxTrialsAfterFailure, and g2o::Property< T >::setValue().

◆ setUserLambdaInit()

void g2o::OptimizationAlgorithmLevenberg::setUserLambdaInit ( double  lambda)

specify the initial lambda used for the first iteraion, if not given the SparseOptimizer tries to compute a suitable value

Definition at line 181 of file optimization_algorithm_levenberg.cpp.

181 {
182 _userLambdaInit->setValue(lambda);
183}

References _userLambdaInit, and g2o::Property< T >::setValue().

◆ solve()

OptimizationAlgorithm::SolverResult g2o::OptimizationAlgorithmLevenberg::solve ( int  iteration,
bool  online = false 
)
virtual

Solve one iteration. The SparseOptimizer running on-top will call this for the given number of iterations.

Parameters
iterationindicates the current iteration

Implements g2o::OptimizationAlgorithm.

Definition at line 59 of file optimization_algorithm_levenberg.cpp.

60 {
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();
76 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
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}
virtual void restoreDiagonal()=0
virtual bool buildStructure(bool zeroBlocks=false)=0
virtual bool setLambda(double lambda, bool backup=false)=0
virtual bool solve()=0
virtual bool buildSystem()=0
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
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
Jet< T, N > pow(const Jet< T, N > &f, double g)
Definition jet.h:743
double get_monotonic_time()
Definition timeutil.cpp:43
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77

References __PRETTY_FUNCTION__, _currentLambda, _goodStepLowerScale, _goodStepUpperScale, _levenbergIterations, _maxTrialsAfterFailure, _ni, g2o::OptimizationAlgorithm::_optimizer, g2o::OptimizationAlgorithmWithHessian::_solver, g2o::SparseOptimizer::activeRobustChi2(), g2o::Solver::buildStructure(), g2o::Solver::buildSystem(), g2o::SparseOptimizer::computeActiveErrors(), computeLambdaInit(), computeScale(), g2o::cst(), g2o::SparseOptimizer::discardTop(), g2o_isfinite, G2O_WARN, g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), g2o::G2OBatchStatistics::levenbergIterations, OK, g2o::Solver::optimizer(), g2o::SparseOptimizer::pop(), g2o::SparseOptimizer::push(), g2o::Solver::restoreDiagonal(), g2o::Solver::setLambda(), g2o::Solver::solve(), Terminate, g2o::SparseOptimizer::terminate(), g2o::G2OBatchStatistics::timeLinearSolution, g2o::G2OBatchStatistics::timeQuadraticForm, g2o::G2OBatchStatistics::timeResiduals, g2o::G2OBatchStatistics::timeUpdate, g2o::SparseOptimizer::update(), g2o::Property< T >::value(), and g2o::Solver::x().

◆ userLambdaInit()

double g2o::OptimizationAlgorithmLevenberg::userLambdaInit ( )
inline

return the lambda set by the user, if < 0 the SparseOptimizer will compute the initial lambda

Definition at line 66 of file optimization_algorithm_levenberg.h.

66{ return _userLambdaInit->value(); }

Member Data Documentation

◆ _currentLambda

double g2o::OptimizationAlgorithmLevenberg::_currentLambda
protected

Definition at line 78 of file optimization_algorithm_levenberg.h.

Referenced by computeScale(), printVerbose(), and solve().

◆ _goodStepLowerScale

double g2o::OptimizationAlgorithmLevenberg::_goodStepLowerScale
protected

lower bound for lambda decrease if a good LM step

Definition at line 80 of file optimization_algorithm_levenberg.h.

Referenced by solve().

◆ _goodStepUpperScale

double g2o::OptimizationAlgorithmLevenberg::_goodStepUpperScale
protected

upper bound for lambda decrease if a good LM step

Definition at line 82 of file optimization_algorithm_levenberg.h.

Referenced by solve().

◆ _levenbergIterations

int g2o::OptimizationAlgorithmLevenberg::_levenbergIterations
protected

the number of levenberg iterations performed to accept the last step

Definition at line 85 of file optimization_algorithm_levenberg.h.

Referenced by printVerbose(), and solve().

◆ _maxTrialsAfterFailure

Property<int>* g2o::OptimizationAlgorithmLevenberg::_maxTrialsAfterFailure
protected

◆ _ni

double g2o::OptimizationAlgorithmLevenberg::_ni
protected

Definition at line 84 of file optimization_algorithm_levenberg.h.

Referenced by solve().

◆ _tau

double g2o::OptimizationAlgorithmLevenberg::_tau
protected

Definition at line 79 of file optimization_algorithm_levenberg.h.

Referenced by computeLambdaInit().

◆ _userLambdaInit

Property<double>* g2o::OptimizationAlgorithmLevenberg::_userLambdaInit
protected

◆ m_solver

std::unique_ptr<Solver> g2o::OptimizationAlgorithmLevenberg::m_solver
private

Definition at line 96 of file optimization_algorithm_levenberg.h.


The documentation for this class was generated from the following files: