gtsam 4.1.1
gtsam
LevenbergMarquardtParams.h
Go to the documentation of this file.
1/* ----------------------------------------------------------------------------
2
3 * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4 * Atlanta, Georgia 30332-0415
5 * All Rights Reserved
6 * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7
8 * See LICENSE for the license information
9
10 * -------------------------------------------------------------------------- */
11
21#pragma once
22
25
26namespace gtsam {
27
28class LevenbergMarquardtOptimizer;
29
36
37public:
40 SILENT = 0, SUMMARY, TERMINATION, LAMBDA, TRYLAMBDA, TRYCONFIG, DAMPED, TRYDELTA
41 };
42
43 static VerbosityLM verbosityLMTranslator(const std::string &s);
44 static std::string verbosityLMTranslator(VerbosityLM value);
45 using OptimizerType = LevenbergMarquardtOptimizer;
46
47public:
48
50 double lambdaFactor;
55 std::string logFile;
58 double minDiagonal;
59 double maxDiagonal;
60
62 : verbosityLM(SILENT),
63 diagonalDamping(false),
64 minDiagonal(1e-6),
65 maxDiagonal(1e32) {
66 SetLegacyDefaults(this);
67 }
68
69 static void SetLegacyDefaults(LevenbergMarquardtParams* p) {
70 // Relevant NonlinearOptimizerParams:
71 p->maxIterations = 100;
72 p->relativeErrorTol = 1e-5;
73 p->absoluteErrorTol = 1e-5;
74 // LM-specific:
75 p->lambdaInitial = 1e-5;
76 p->lambdaFactor = 10.0;
77 p->lambdaUpperBound = 1e5;
78 p->lambdaLowerBound = 0.0;
79 p->minModelFidelity = 1e-3;
80 p->diagonalDamping = false;
81 p->useFixedLambdaFactor = true;
82 }
83
84 // these do seem to work better for SFM
85 static void SetCeresDefaults(LevenbergMarquardtParams* p) {
86 // Relevant NonlinearOptimizerParams:
87 p->maxIterations = 50;
88 p->absoluteErrorTol = 0; // No corresponding option in CERES
89 p->relativeErrorTol = 1e-6; // This is function_tolerance
90 // LM-specific:
91 p->lambdaUpperBound = 1e32;
92 p->lambdaLowerBound = 1e-16;
93 p->lambdaInitial = 1e-04;
94 p->lambdaFactor = 2.0;
95 p->minModelFidelity = 1e-3; // options.min_relative_decrease in CERES
96 p->diagonalDamping = true;
97 p->useFixedLambdaFactor = false; // This is important
98 }
99
100 static LevenbergMarquardtParams LegacyDefaults() {
101 LevenbergMarquardtParams p;
102 SetLegacyDefaults(&p);
103 return p;
104 }
105
106 static LevenbergMarquardtParams CeresDefaults() {
107 LevenbergMarquardtParams p;
108 SetCeresDefaults(&p);
109 return p;
110 }
111
112 static LevenbergMarquardtParams EnsureHasOrdering(LevenbergMarquardtParams params,
113 const NonlinearFactorGraph& graph) {
114 if (!params.ordering)
115 params.ordering = Ordering::Create(params.orderingType, graph);
116 return params;
117 }
118
119 static LevenbergMarquardtParams ReplaceOrdering(LevenbergMarquardtParams params,
120 const Ordering& ordering) {
121 params.ordering = ordering;
122 return params;
123 }
124
125 ~LevenbergMarquardtParams() override {}
126 void print(const std::string& str = "") const override;
127
130 bool getDiagonalDamping() const { return diagonalDamping; }
131 double getlambdaFactor() const { return lambdaFactor; }
132 double getlambdaInitial() const { return lambdaInitial; }
133 double getlambdaLowerBound() const { return lambdaLowerBound; }
134 double getlambdaUpperBound() const { return lambdaUpperBound; }
135 bool getUseFixedLambdaFactor() { return useFixedLambdaFactor; }
136 std::string getLogFile() const { return logFile; }
137 std::string getVerbosityLM() const { return verbosityLMTranslator(verbosityLM);}
138
139 void setDiagonalDamping(bool flag) { diagonalDamping = flag; }
140 void setlambdaFactor(double value) { lambdaFactor = value; }
141 void setlambdaInitial(double value) { lambdaInitial = value; }
142 void setlambdaLowerBound(double value) { lambdaLowerBound = value; }
143 void setlambdaUpperBound(double value) { lambdaUpperBound = value; }
144 void setUseFixedLambdaFactor(bool flag) { useFixedLambdaFactor = flag;}
145 void setLogFile(const std::string& s) { logFile = s; }
146 void setVerbosityLM(const std::string& s) { verbosityLM = verbosityLMTranslator(s);}
147 // @}
150
152 boost::shared_ptr<NonlinearOptimizerParams> clone() const {
153 return boost::shared_ptr<NonlinearOptimizerParams>(new LevenbergMarquardtParams(*this));
154 }
155
157};
158
159}
Parameters for nonlinear optimization.
Factor Graph consisting of non-linear factors.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
void print(const Matrix &A, const string &s, ostream &stream)
print without optional string, must specify cout yourself
Definition: Matrix.cpp:155
This class performs Levenberg-Marquardt nonlinear optimization.
Definition: LevenbergMarquardtOptimizer.h:35
Parameters for Levenberg-Marquardt optimization.
Definition: LevenbergMarquardtParams.h:35
double lambdaFactor
The amount by which to multiply or divide lambda when adjusting lambda (default: 10....
Definition: LevenbergMarquardtParams.h:50
double minDiagonal
when using diagonal damping saturates the minimum diagonal entries (default: 1e-6)
Definition: LevenbergMarquardtParams.h:58
double lambdaUpperBound
The maximum lambda to try before assuming the optimization has failed (default: 1e5)
Definition: LevenbergMarquardtParams.h:51
double lambdaInitial
The initial Levenberg-Marquardt damping term (default: 1e-5)
Definition: LevenbergMarquardtParams.h:49
double maxDiagonal
when using diagonal damping saturates the maximum diagonal entries (default: 1e32)
Definition: LevenbergMarquardtParams.h:59
double minModelFidelity
Lower bound for the modelFidelity to accept the result of an LM iteration.
Definition: LevenbergMarquardtParams.h:54
double lambdaLowerBound
The minimum lambda used in LM (default: 0)
Definition: LevenbergMarquardtParams.h:52
boost::shared_ptr< NonlinearOptimizerParams > clone() const
Definition: LevenbergMarquardtParams.h:152
bool diagonalDamping
if true, use diagonal of Hessian
Definition: LevenbergMarquardtParams.h:56
bool useFixedLambdaFactor
if true applies constant increase (or decrease) to lambda according to lambdaFactor
Definition: LevenbergMarquardtParams.h:57
std::string logFile
an optional CSV log file, with [iteration, time, error, lambda]
Definition: LevenbergMarquardtParams.h:55
VerbosityLM
See LevenbergMarquardtParams::lmVerbosity.
Definition: LevenbergMarquardtParams.h:39
VerbosityLM verbosityLM
The verbosity level for Levenberg-Marquardt (default: SILENT), see also NonlinearOptimizerParams::ver...
Definition: LevenbergMarquardtParams.h:53
The common parameters for Nonlinear optimizers.
Definition: NonlinearOptimizerParams.h:34
double absoluteErrorTol
The maximum absolute error decrease to stop iterating (default 1e-5)
Definition: NonlinearOptimizerParams.h:43
size_t maxIterations
The maximum iterations to stop iterating (default 100)
Definition: NonlinearOptimizerParams.h:41
double relativeErrorTol
The maximum relative error decrease to stop iterating (default 1e-5)
Definition: NonlinearOptimizerParams.h:42