LevenbergMarquardt/LevenbergMarquardt.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
5 // Copyright (C) 2012 Desire Nuentsa <desire.nuentsa_wakam@inria.fr>
6 //
7 // The algorithm of this class initially comes from MINPACK whose original authors are:
8 // Copyright Jorge More - Argonne National Laboratory
9 // Copyright Burt Garbow - Argonne National Laboratory
10 // Copyright Ken Hillstrom - Argonne National Laboratory
11 //
12 // This Source Code Form is subject to the terms of the Minpack license
13 // (a BSD-like license) described in the campaigned CopyrightMINPACK.txt file.
14 //
15 // This Source Code Form is subject to the terms of the Mozilla
16 // Public License v. 2.0. If a copy of the MPL was not distributed
17 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
18 
19 #ifndef EIGEN_LEVENBERGMARQUARDT_H
20 #define EIGEN_LEVENBERGMARQUARDT_H
21 
22 // IWYU pragma: private
23 #include "./InternalHeaderCheck.h"
24 
25 namespace Eigen {
26 namespace LevenbergMarquardtSpace {
27 enum Status {
28  NotStarted = -2,
29  Running = -1,
39  UserAsked = 9
40 };
41 }
42 
43 template <typename Scalar_, int NX = Dynamic, int NY = Dynamic>
44 struct DenseFunctor {
45  typedef Scalar_ Scalar;
51  const int m_inputs, m_values;
52 
55 
56  int inputs() const { return m_inputs; }
57  int values() const { return m_values; }
58 
59  // int operator()(const InputType &x, ValueType& fvec) { }
60  // should be defined in derived classes
61 
62  // int df(const InputType &x, JacobianType& fjac) { }
63  // should be defined in derived classes
64 };
65 
66 template <typename Scalar_, typename Index_>
67 struct SparseFunctor {
68  typedef Scalar_ Scalar;
69  typedef Index_ Index;
75 
77 
78  int inputs() const { return m_inputs; }
79  int values() const { return m_values; }
80 
81  const int m_inputs, m_values;
82  // int operator()(const InputType &x, ValueType& fvec) { }
83  // to be defined in the functor
84 
85  // int df(const InputType &x, JacobianType& fjac) { }
86  // to be defined in the functor if no automatic differentiation
87 };
88 namespace internal {
89 template <typename QRSolver, typename VectorType>
90 void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta,
91  typename VectorType::Scalar &par, VectorType &x);
92 }
101 template <typename FunctorType_>
103  public:
104  typedef FunctorType_ FunctorType;
105  typedef typename FunctorType::QRSolver QRSolver;
106  typedef typename FunctorType::JacobianType JacobianType;
107  typedef typename JacobianType::Scalar Scalar;
109  typedef typename QRSolver::StorageIndex PermIndex;
112 
113  public:
115  : m_functor(functor),
116  m_nfev(0),
117  m_njev(0),
118  m_fnorm(0.0),
119  m_gnorm(0),
120  m_isInitialized(false),
122  resetParameters();
123  m_useExternalScaling = false;
124  }
125 
132 
135  using std::sqrt;
136 
137  m_factor = 100.;
138  m_maxfev = 400;
141  m_gtol = 0.;
142  m_epsfcn = 0.;
143  }
144 
147 
150 
153 
156 
158  void setEpsilon(RealScalar epsfcn) { m_epsfcn = epsfcn; }
159 
162 
165 
167  RealScalar xtol() const { return m_xtol; }
168 
170  RealScalar ftol() const { return m_ftol; }
171 
173  RealScalar gtol() const { return m_gtol; }
174 
176  RealScalar factor() const { return m_factor; }
177 
179  RealScalar epsilon() const { return m_epsfcn; }
180 
182  Index maxfev() const { return m_maxfev; }
183 
185  FVectorType &diag() { return m_diag; }
186 
188  Index iterations() { return m_iter; }
189 
191  Index nfev() { return m_nfev; }
192 
194  Index njev() { return m_njev; }
195 
197  RealScalar fnorm() { return m_fnorm; }
198 
200  RealScalar gnorm() { return m_gnorm; }
201 
203  RealScalar lm_param(void) { return m_par; }
204 
207  FVectorType &fvec() { return m_fvec; }
208 
211  JacobianType &jacobian() { return m_fjac; }
212 
217 
221 
231  ComputationInfo info() const { return m_info; }
232 
233  private:
235  JacobianType m_rfactor; // The triangular matrix R from the QR of the jacobian matrix m_fjac
242  RealScalar m_fnorm; // Norm of the current vector function
243  RealScalar m_gnorm; // Norm of the gradient of the error
245  Index m_maxfev; // Maximum number of function evaluation
246  RealScalar m_ftol; // Tolerance in the norm of the vector function
248  RealScalar m_gtol; // tolerance of the norm of the error gradient
250  Index m_iter; // Number of iterations performed
254  FVectorType m_wa1, m_wa2, m_wa3, m_wa4; // Temporary vectors
256  bool m_isInitialized; // Check whether the minimization step has been called
258 };
259 
260 template <typename FunctorType>
262  LevenbergMarquardtSpace::Status status = minimizeInit(x);
264  m_isInitialized = true;
265  return status;
266  }
267  do {
268  // std::cout << " uv " << x.transpose() << "\n";
269  status = minimizeOneStep(x);
270  } while (status == LevenbergMarquardtSpace::Running);
271  m_isInitialized = true;
272  return status;
273 }
274 
275 template <typename FunctorType>
277  n = x.size();
278  m = m_functor.values();
279 
280  m_wa1.resize(n);
281  m_wa2.resize(n);
282  m_wa3.resize(n);
283  m_wa4.resize(m);
284  m_fvec.resize(m);
285  // FIXME Sparse Case : Allocate space for the jacobian
286  m_fjac.resize(m, n);
287  // m_fjac.reserve(VectorXi::Constant(n,5)); // FIXME Find a better alternative
288  if (!m_useExternalScaling) m_diag.resize(n);
289  eigen_assert((!m_useExternalScaling || m_diag.size() == n) &&
290  "When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
291  m_qtf.resize(n);
292 
293  /* Function Body */
294  m_nfev = 0;
295  m_njev = 0;
296 
297  /* check the input parameters for errors. */
298  if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.) {
299  m_info = InvalidInput;
301  }
302 
303  if (m_useExternalScaling)
304  for (Index j = 0; j < n; ++j)
305  if (m_diag[j] <= 0.) {
306  m_info = InvalidInput;
308  }
309 
310  /* evaluate the function at the starting point */
311  /* and calculate its norm. */
312  m_nfev = 1;
313  if (m_functor(x, m_fvec) < 0) return LevenbergMarquardtSpace::UserAsked;
314  m_fnorm = m_fvec.stableNorm();
315 
316  /* initialize levenberg-marquardt parameter and iteration counter. */
317  m_par = 0.;
318  m_iter = 1;
319 
321 }
322 
323 template <typename FunctorType>
325  n = x.size();
326  m = m_functor.values();
327 
328  /* check the input parameters for errors. */
329  if (n <= 0 || m < n || tol < 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
330 
331  resetParameters();
332  m_ftol = tol;
333  m_xtol = tol;
334  m_maxfev = 100 * (n + 1);
335 
336  return minimize(x);
337 }
338 
339 template <typename FunctorType>
341  Index *nfev, const Scalar tol) {
342  Index n = x.size();
343  Index m = functor.values();
344 
345  /* check the input parameters for errors. */
346  if (n <= 0 || m < n || tol < 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
347 
348  NumericalDiff<FunctorType> numDiff(functor);
349  // embedded LevenbergMarquardt
351  lm.setFtol(tol);
352  lm.setXtol(tol);
353  lm.setMaxfev(200 * (n + 1));
354 
356  if (nfev) *nfev = lm.nfev();
357  return info;
358 }
359 
360 } // end namespace Eigen
361 
362 #endif // EIGEN_LEVENBERGMARQUARDT_H
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
HouseholderQR< MatrixXf > qr(A)
#define eigen_assert(x)
Definition: Macros.h:910
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition: ColPivHouseholderQR.h:54
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
Definition: LevenbergMarquardt/LevenbergMarquardt.h:102
RealScalar epsilon() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:179
static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
Definition: LevenbergMarquardt/LevenbergMarquardt.h:340
JacobianType::Scalar Scalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:107
FVectorType m_diag
Definition: LevenbergMarquardt/LevenbergMarquardt.h:237
FunctorType::QRSolver QRSolver
Definition: LevenbergMarquardt/LevenbergMarquardt.h:105
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
ComputationInfo info() const
Reports whether the minimization was successful.
Definition: LevenbergMarquardt/LevenbergMarquardt.h:231
JacobianType::RealScalar RealScalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:108
void setFtol(RealScalar ftol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:149
RealScalar gtol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:173
void resetParameters()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:134
FunctorType & m_functor
Definition: LevenbergMarquardt/LevenbergMarquardt.h:236
FunctorType & functor
Definition: NonLinearOptimization/LevenbergMarquardt.h:111
FVectorType m_wa2
Definition: LevenbergMarquardt/LevenbergMarquardt.h:254
FVectorType m_qtf
Definition: LevenbergMarquardt/LevenbergMarquardt.h:237
RealScalar ftol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:170
bool m_useExternalScaling
Definition: LevenbergMarquardt/LevenbergMarquardt.h:252
Index m_iter
Definition: LevenbergMarquardt/LevenbergMarquardt.h:250
RealScalar m_fnorm
Definition: LevenbergMarquardt/LevenbergMarquardt.h:242
Index m
Definition: LevenbergMarquardt/LevenbergMarquardt.h:239
Index m_maxfev
Definition: LevenbergMarquardt/LevenbergMarquardt.h:245
RealScalar m_epsfcn
Definition: LevenbergMarquardt/LevenbergMarquardt.h:249
RealScalar m_gtol
Definition: LevenbergMarquardt/LevenbergMarquardt.h:248
Index maxfev() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:182
void setExternalScaling(bool value)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:164
RealScalar xtol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:167
RealScalar m_gnorm
Definition: LevenbergMarquardt/LevenbergMarquardt.h:243
FunctorType_ FunctorType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:104
RealScalar gnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:200
LevenbergMarquardtSpace::Status minimize(FVectorType &x)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:261
FVectorType & fvec()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:207
Index iterations()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:188
JacobianType & jacobian()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:211
Matrix< Scalar, Dynamic, 1 > FVectorType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:110
JacobianType m_fjac
Definition: LevenbergMarquardt/LevenbergMarquardt.h:234
LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
Definition: LevenbergMarquardt/LevenbergMarquardt.h:324
Index nfev
Definition: NonLinearOptimization/LevenbergMarquardt.h:102
Index n
Definition: LevenbergMarquardt/LevenbergMarquardt.h:238
LevenbergMarquardt(FunctorType &functor)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:114
FVectorType & diag()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:185
FVectorType m_wa3
Definition: LevenbergMarquardt/LevenbergMarquardt.h:254
Index njev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:194
ComputationInfo m_info
Definition: LevenbergMarquardt/LevenbergMarquardt.h:257
void setGtol(RealScalar gtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:152
RealScalar lm_param(void)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:203
RealScalar fnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:197
void setEpsilon(RealScalar epsfcn)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:158
RealScalar m_ftol
Definition: LevenbergMarquardt/LevenbergMarquardt.h:246
FVectorType m_fvec
Definition: LevenbergMarquardt/LevenbergMarquardt.h:237
void setXtol(RealScalar xtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:146
RealScalar m_delta
Definition: LevenbergMarquardt/LevenbergMarquardt.h:251
PermutationType m_permutation
Definition: LevenbergMarquardt/LevenbergMarquardt.h:253
Index m_nfev
Definition: LevenbergMarquardt/LevenbergMarquardt.h:240
bool m_isInitialized
Definition: LevenbergMarquardt/LevenbergMarquardt.h:256
void setFactor(RealScalar factor)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:155
DenseIndex Index
Definition: NonLinearOptimization/LevenbergMarquardt.h:60
PermutationMatrix< Dynamic, Dynamic, int > PermutationType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:111
Index nfev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:191
RealScalar m_par
Definition: LevenbergMarquardt/LevenbergMarquardt.h:255
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:276
Index m_njev
Definition: LevenbergMarquardt/LevenbergMarquardt.h:241
RealScalar m_xtol
Definition: LevenbergMarquardt/LevenbergMarquardt.h:247
FVectorType m_wa1
Definition: LevenbergMarquardt/LevenbergMarquardt.h:254
FVectorType m_wa4
Definition: LevenbergMarquardt/LevenbergMarquardt.h:254
FunctorType::JacobianType JacobianType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:106
RealScalar m_factor
Definition: LevenbergMarquardt/LevenbergMarquardt.h:244
PermutationType permutation()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:220
RealScalar factor() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:176
QRSolver::StorageIndex PermIndex
Definition: LevenbergMarquardt/LevenbergMarquardt.h:109
JacobianType & matrixR()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:216
void setMaxfev(Index maxfev)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:161
JacobianType m_rfactor
Definition: LevenbergMarquardt/LevenbergMarquardt.h:235
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Definition: NumericalDiff.h:35
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
Sparse left-looking QR factorization with numerical column pivoting.
Definition: SparseQR.h:88
Definition: XprHelper.h:134
ComputationInfo
Definition: Constants.h:438
@ InvalidInput
Definition: Constants.h:447
int * m
Definition: level2_cplx_impl.h:294
int info
Definition: level2_cplx_impl.h:39
const char const char const char * diag
Definition: level2_impl.h:86
Status
Definition: LevenbergMarquardt/LevenbergMarquardt.h:27
@ RelativeReductionTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:31
@ GtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:38
@ NotStarted
Definition: LevenbergMarquardt/LevenbergMarquardt.h:28
@ UserAsked
Definition: LevenbergMarquardt/LevenbergMarquardt.h:39
@ Running
Definition: LevenbergMarquardt/LevenbergMarquardt.h:29
@ FtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:36
@ TooManyFunctionEvaluation
Definition: LevenbergMarquardt/LevenbergMarquardt.h:35
@ XtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:37
@ CosinusTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:34
@ RelativeErrorTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:32
@ ImproperInputParameters
Definition: LevenbergMarquardt/LevenbergMarquardt.h:30
@ RelativeErrorAndReductionTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:33
void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta, typename VectorType::Scalar &par, VectorType &x)
Definition: LMpar.h:23
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
const int Dynamic
Definition: Constants.h:25
string par
Definition: calibrate.py:135
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: LevenbergMarquardt/LevenbergMarquardt.h:44
Matrix< Scalar, InputsAtCompileTime, 1 > InputType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:47
int inputs() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:56
const int m_inputs
Definition: LevenbergMarquardt/LevenbergMarquardt.h:51
const int m_values
Definition: LevenbergMarquardt/LevenbergMarquardt.h:51
Matrix< Scalar, ValuesAtCompileTime, InputsAtCompileTime > JacobianType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:49
Scalar_ Scalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:45
ColPivHouseholderQR< JacobianType > QRSolver
Definition: LevenbergMarquardt/LevenbergMarquardt.h:50
DenseFunctor(int inputs, int values)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:54
DenseFunctor()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:53
@ InputsAtCompileTime
Definition: LevenbergMarquardt/LevenbergMarquardt.h:46
@ ValuesAtCompileTime
Definition: LevenbergMarquardt/LevenbergMarquardt.h:46
int values() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:57
Matrix< Scalar, ValuesAtCompileTime, 1 > ValueType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:48
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: LevenbergMarquardt/LevenbergMarquardt.h:67
Matrix< Scalar, Dynamic, 1 > InputType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:70
const int m_values
Definition: LevenbergMarquardt/LevenbergMarquardt.h:81
Matrix< Scalar, Dynamic, 1 > ValueType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:71
SparseMatrix< Scalar, ColMajor, Index > JacobianType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:72
@ InputsAtCompileTime
Definition: LevenbergMarquardt/LevenbergMarquardt.h:74
@ ValuesAtCompileTime
Definition: LevenbergMarquardt/LevenbergMarquardt.h:74
SparseFunctor(int inputs, int values)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:76
Index_ Index
Definition: LevenbergMarquardt/LevenbergMarquardt.h:69
const int m_inputs
Definition: LevenbergMarquardt/LevenbergMarquardt.h:81
Scalar_ Scalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:68
int inputs() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:78
int values() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:79
SparseQR< JacobianType, COLAMDOrdering< int > > QRSolver
Definition: LevenbergMarquardt/LevenbergMarquardt.h:73
Definition: fft_test_shared.h:66
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2