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 
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>
45 {
46  typedef Scalar_ Scalar;
47  enum {
50  };
55  const int m_inputs, m_values;
56 
59 
60  int inputs() const { return m_inputs; }
61  int values() const { return m_values; }
62 
63  //int operator()(const InputType &x, ValueType& fvec) { }
64  // should be defined in derived classes
65 
66  //int df(const InputType &x, JacobianType& fjac) { }
67  // should be defined in derived classes
68 };
69 
70 template <typename Scalar_, typename Index_>
72 {
73  typedef Scalar_ Scalar;
74  typedef Index_ Index;
79  enum {
82  };
83 
85 
86  int inputs() const { return m_inputs; }
87  int values() const { return m_values; }
88 
89  const int m_inputs, m_values;
90  //int operator()(const InputType &x, ValueType& fvec) { }
91  // to be defined in the functor
92 
93  //int df(const InputType &x, JacobianType& fjac) { }
94  // to be defined in the functor if no automatic differentiation
95 
96 };
97 namespace internal {
98 template <typename QRSolver, typename VectorType>
99 void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb,
100  typename VectorType::Scalar m_delta, typename VectorType::Scalar &par,
101  VectorType &x);
102  }
111 template<typename FunctorType_>
112 class LevenbergMarquardt : internal::no_assignment_operator
113 {
114  public:
115  typedef FunctorType_ FunctorType;
116  typedef typename FunctorType::QRSolver QRSolver;
117  typedef typename FunctorType::JacobianType JacobianType;
118  typedef typename JacobianType::Scalar Scalar;
119  typedef typename JacobianType::RealScalar RealScalar;
120  typedef typename QRSolver::StorageIndex PermIndex;
123  public:
125  : m_functor(functor),m_nfev(0),m_njev(0),m_fnorm(0.0),m_gnorm(0),
127  {
128  resetParameters();
129  m_useExternalScaling=false;
130  }
131 
136  FVectorType &x,
138  );
141  FVectorType &x,
142  Index *nfev,
144  );
145 
148  {
149  using std::sqrt;
150 
151  m_factor = 100.;
152  m_maxfev = 400;
155  m_gtol = 0. ;
156  m_epsfcn = 0. ;
157  }
158 
161 
164 
167 
170 
172  void setEpsilon (RealScalar epsfcn) { m_epsfcn = epsfcn; }
173 
176 
178  void setExternalScaling(bool value) {m_useExternalScaling = value; }
179 
181  RealScalar xtol() const {return m_xtol; }
182 
184  RealScalar ftol() const {return m_ftol; }
185 
187  RealScalar gtol() const {return m_gtol; }
188 
190  RealScalar factor() const {return m_factor; }
191 
193  RealScalar epsilon() const {return m_epsfcn; }
194 
196  Index maxfev() const {return m_maxfev; }
197 
199  FVectorType& diag() {return m_diag; }
200 
202  Index iterations() { return m_iter; }
203 
205  Index nfev() { return m_nfev; }
206 
208  Index njev() { return m_njev; }
209 
211  RealScalar fnorm() {return m_fnorm; }
212 
214  RealScalar gnorm() {return m_gnorm; }
215 
217  RealScalar lm_param(void) { return m_par; }
218 
221  FVectorType& fvec() {return m_fvec; }
222 
226 
231 
235 
246  {
247 
248  return m_info;
249  }
250  private:
252  JacobianType m_rfactor; // The triangular matrix R from the QR of the jacobian matrix m_fjac
259  RealScalar m_fnorm; // Norm of the current vector function
260  RealScalar m_gnorm; //Norm of the gradient of the error
262  Index m_maxfev; // Maximum number of function evaluation
263  RealScalar m_ftol; //Tolerance in the norm of the vector function
265  RealScalar m_gtol; //tolerance of the norm of the error gradient
267  Index m_iter; // Number of iterations performed
271  FVectorType m_wa1, m_wa2, m_wa3, m_wa4; //Temporary vectors
273  bool m_isInitialized; // Check whether the minimization step has been called
275 };
276 
277 template<typename FunctorType>
280 {
281  LevenbergMarquardtSpace::Status status = minimizeInit(x);
283  m_isInitialized = true;
284  return status;
285  }
286  do {
287 // std::cout << " uv " << x.transpose() << "\n";
288  status = minimizeOneStep(x);
289  } while (status==LevenbergMarquardtSpace::Running);
290  m_isInitialized = true;
291  return status;
292 }
293 
294 template<typename FunctorType>
297 {
298  n = x.size();
299  m = m_functor.values();
300 
301  m_wa1.resize(n); m_wa2.resize(n); m_wa3.resize(n);
302  m_wa4.resize(m);
303  m_fvec.resize(m);
304  //FIXME Sparse Case : Allocate space for the jacobian
305  m_fjac.resize(m, n);
306 // m_fjac.reserve(VectorXi::Constant(n,5)); // FIXME Find a better alternative
307  if (!m_useExternalScaling)
308  m_diag.resize(n);
309  eigen_assert( (!m_useExternalScaling || m_diag.size()==n) && "When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
310  m_qtf.resize(n);
311 
312  /* Function Body */
313  m_nfev = 0;
314  m_njev = 0;
315 
316  /* check the input parameters for errors. */
317  if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.){
318  m_info = InvalidInput;
320  }
321 
322  if (m_useExternalScaling)
323  for (Index j = 0; j < n; ++j)
324  if (m_diag[j] <= 0.)
325  {
326  m_info = InvalidInput;
328  }
329 
330  /* evaluate the function at the starting point */
331  /* and calculate its norm. */
332  m_nfev = 1;
333  if ( m_functor(x, m_fvec) < 0)
335  m_fnorm = m_fvec.stableNorm();
336 
337  /* initialize levenberg-marquardt parameter and iteration counter. */
338  m_par = 0.;
339  m_iter = 1;
340 
342 }
343 
344 template<typename FunctorType>
347  FVectorType &x,
348  const Scalar tol
349  )
350 {
351  n = x.size();
352  m = m_functor.values();
353 
354  /* check the input parameters for errors. */
355  if (n <= 0 || m < n || tol < 0.)
357 
358  resetParameters();
359  m_ftol = tol;
360  m_xtol = tol;
361  m_maxfev = 100*(n+1);
362 
363  return minimize(x);
364 }
365 
366 
367 template<typename FunctorType>
370  FunctorType &functor,
371  FVectorType &x,
372  Index *nfev,
373  const Scalar tol
374  )
375 {
376  Index n = x.size();
377  Index m = functor.values();
378 
379  /* check the input parameters for errors. */
380  if (n <= 0 || m < n || tol < 0.)
382 
383  NumericalDiff<FunctorType> numDiff(functor);
384  // embedded LevenbergMarquardt
386  lm.setFtol(tol);
387  lm.setXtol(tol);
388  lm.setMaxfev(200*(n+1));
389 
391  if (nfev)
392  * nfev = lm.nfev();
393  return info;
394 }
395 
396 } // end namespace Eigen
397 
398 #endif // EIGEN_LEVENBERGMARQUARDT_H
Matrix3f m
int n
HouseholderQR< MatrixXf > qr(A)
#define eigen_assert(x)
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
ComputationInfo info() const
Reports whether the minimization was successful.
LevenbergMarquardtSpace::Status minimize(FVectorType &x)
LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)
PermutationMatrix< Dynamic, Dynamic, int > PermutationType
ComputationInfo
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:22
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Eigen::internal::remove_all_t< DerType >, typename Eigen::internal::traits< Eigen::internal::remove_all_t< DerType >>::Scalar, product) > sqrt(const Eigen::AutoDiffScalar< DerType > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const int Dynamic
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Matrix< Scalar, ValuesAtCompileTime, 1 > ValueType
ColPivHouseholderQR< JacobianType > QRSolver
Matrix< Scalar, ValuesAtCompileTime, InputsAtCompileTime > JacobianType
Matrix< Scalar, InputsAtCompileTime, 1 > InputType
Matrix< Scalar, Dynamic, 1 > ValueType
SparseMatrix< Scalar, ColMajor, Index > JacobianType
Matrix< Scalar, Dynamic, 1 > InputType
SparseQR< JacobianType, COLAMDOrdering< int > > QRSolver
std::ptrdiff_t j