$darkmode
Eigen-unsupported  5.0.1-dev
LevenbergMarquardt.h
1 // -*- coding: utf-8
2 // vim: set fileencoding=utf-8
3 
4 // This file is part of Eigen, a lightweight C++ template library
5 // for linear algebra.
6 //
7 // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
8 //
9 // This Source Code Form is subject to the terms of the Mozilla
10 // Public License v. 2.0. If a copy of the MPL was not distributed
11 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12 
13 #ifndef EIGEN_LEVENBERGMARQUARDT__H
14 #define EIGEN_LEVENBERGMARQUARDT__H
15 
16 // IWYU pragma: private
17 #include "./InternalHeaderCheck.h"
18 
19 namespace Eigen {
20 
21 namespace LevenbergMarquardtSpace {
22 enum Status {
23  NotStarted = -2,
24  Running = -1,
25  ImproperInputParameters = 0,
26  RelativeReductionTooSmall = 1,
27  RelativeErrorTooSmall = 2,
28  RelativeErrorAndReductionTooSmall = 3,
29  CosinusTooSmall = 4,
30  TooManyFunctionEvaluation = 5,
31  FtolTooSmall = 6,
32  XtolTooSmall = 7,
33  GtolTooSmall = 8,
34  UserAsked = 9
35 };
36 }
37 
46 template <typename FunctorType, typename Scalar = double>
47 class LevenbergMarquardt {
48  static Scalar sqrt_epsilon() {
49  using std::sqrt;
50  return sqrt(NumTraits<Scalar>::epsilon());
51  }
52 
53  public:
54  LevenbergMarquardt(FunctorType &_functor) : functor(_functor) {
55  nfev = njev = iter = 0;
56  fnorm = gnorm = 0.;
57  useExternalScaling = false;
58  }
59 
60  typedef DenseIndex Index;
61 
62  struct Parameters {
63  Parameters()
64  : factor(Scalar(100.)),
65  maxfev(400),
66  ftol(sqrt_epsilon()),
67  xtol(sqrt_epsilon()),
68  gtol(Scalar(0.)),
69  epsfcn(Scalar(0.)) {}
70  Scalar factor;
71  Index maxfev; // maximum number of function evaluation
72  Scalar ftol;
73  Scalar xtol;
74  Scalar gtol;
75  Scalar epsfcn;
76  };
77 
78  typedef Matrix<Scalar, Dynamic, 1> FVectorType;
79  typedef Matrix<Scalar, Dynamic, Dynamic> JacobianType;
80 
81  LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol = sqrt_epsilon());
82 
83  LevenbergMarquardtSpace::Status minimize(FVectorType &x);
84  LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x);
85  LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x);
86 
87  static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev,
88  const Scalar tol = sqrt_epsilon());
89 
90  LevenbergMarquardtSpace::Status lmstr1(FVectorType &x, const Scalar tol = sqrt_epsilon());
91 
92  LevenbergMarquardtSpace::Status minimizeOptimumStorage(FVectorType &x);
93  LevenbergMarquardtSpace::Status minimizeOptimumStorageInit(FVectorType &x);
94  LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep(FVectorType &x);
95 
96  void resetParameters(void) { parameters = Parameters(); }
97 
98  Parameters parameters;
99  FVectorType fvec, qtf, diag;
100  JacobianType fjac;
101  PermutationMatrix<Dynamic, Dynamic> permutation;
102  Index nfev;
103  Index njev;
104  Index iter;
105  Scalar fnorm, gnorm;
106  bool useExternalScaling;
107 
108  Scalar lm_param(void) { return par; }
109 
110  private:
111  FunctorType &functor;
112  Index n;
113  Index m;
114  FVectorType wa1, wa2, wa3, wa4;
115 
116  Scalar par, sum;
117  Scalar temp, temp1, temp2;
118  Scalar delta;
119  Scalar ratio;
120  Scalar pnorm, xnorm, fnorm1, actred, dirder, prered;
121 
122  LevenbergMarquardt &operator=(const LevenbergMarquardt &);
123 };
124 
125 template <typename FunctorType, typename Scalar>
126 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::lmder1(FVectorType &x, const Scalar tol) {
127  n = x.size();
128  m = functor.values();
129 
130  /* check the input parameters for errors. */
131  if (n <= 0 || m < n || tol < 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
132 
133  resetParameters();
134  parameters.ftol = tol;
135  parameters.xtol = tol;
136  parameters.maxfev = 100 * (n + 1);
137 
138  return minimize(x);
139 }
140 
141 template <typename FunctorType, typename Scalar>
142 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimize(FVectorType &x) {
143  LevenbergMarquardtSpace::Status status = minimizeInit(x);
144  if (status == LevenbergMarquardtSpace::ImproperInputParameters) return status;
145  do {
146  status = minimizeOneStep(x);
147  } while (status == LevenbergMarquardtSpace::Running);
148  return status;
149 }
150 
151 template <typename FunctorType, typename Scalar>
152 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimizeInit(FVectorType &x) {
153  n = x.size();
154  m = functor.values();
155 
156  wa1.resize(n);
157  wa2.resize(n);
158  wa3.resize(n);
159  wa4.resize(m);
160  fvec.resize(m);
161  fjac.resize(m, n);
162  if (!useExternalScaling) diag.resize(n);
163  eigen_assert((!useExternalScaling || diag.size() == n) &&
164  "When useExternalScaling is set, the caller must provide a valid 'diag'");
165  qtf.resize(n);
166 
167  /* Function Body */
168  nfev = 0;
169  njev = 0;
170 
171  /* check the input parameters for errors. */
172  if (n <= 0 || m < n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. ||
173  parameters.maxfev <= 0 || parameters.factor <= 0.)
174  return LevenbergMarquardtSpace::ImproperInputParameters;
175 
176  if (useExternalScaling)
177  for (Index j = 0; j < n; ++j)
178  if (diag[j] <= 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
179 
180  /* evaluate the function at the starting point */
181  /* and calculate its norm. */
182  nfev = 1;
183  if (functor(x, fvec) < 0) return LevenbergMarquardtSpace::UserAsked;
184  fnorm = fvec.stableNorm();
185 
186  /* initialize levenberg-marquardt parameter and iteration counter. */
187  par = 0.;
188  iter = 1;
189 
190  return LevenbergMarquardtSpace::NotStarted;
191 }
192 
193 template <typename FunctorType, typename Scalar>
194 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimizeOneStep(FVectorType &x) {
195  using std::abs;
196  using std::sqrt;
197 
198  eigen_assert(x.size() == n); // check the caller is not cheating us
199 
200  /* calculate the jacobian matrix. */
201  Index df_ret = functor.df(x, fjac);
202  if (df_ret < 0) return LevenbergMarquardtSpace::UserAsked;
203  if (df_ret > 0)
204  // numerical diff, we evaluated the function df_ret times
205  nfev += df_ret;
206  else
207  njev++;
208 
209  /* compute the qr factorization of the jacobian. */
210  wa2 = fjac.colwise().blueNorm();
211  ColPivHouseholderQR<JacobianType> qrfac(fjac);
212  fjac = qrfac.matrixQR();
213  permutation = qrfac.colsPermutation();
214 
215  /* on the first iteration and if external scaling is not used, scale according */
216  /* to the norms of the columns of the initial jacobian. */
217  if (iter == 1) {
218  if (!useExternalScaling)
219  for (Index j = 0; j < n; ++j) diag[j] = (wa2[j] == 0.) ? 1. : wa2[j];
220 
221  /* on the first iteration, calculate the norm of the scaled x */
222  /* and initialize the step bound delta. */
223  xnorm = diag.cwiseProduct(x).stableNorm();
224  delta = parameters.factor * xnorm;
225  if (delta == 0.) delta = parameters.factor;
226  }
227 
228  /* form (q transpose)*fvec and store the first n components in */
229  /* qtf. */
230  wa4 = fvec;
231  wa4.applyOnTheLeft(qrfac.householderQ().adjoint());
232  qtf = wa4.head(n);
233 
234  /* compute the norm of the scaled gradient. */
235  gnorm = 0.;
236  if (fnorm != 0.)
237  for (Index j = 0; j < n; ++j)
238  if (wa2[permutation.indices()[j]] != 0.)
239  gnorm = (std::max)(gnorm,
240  abs(fjac.col(j).head(j + 1).dot(qtf.head(j + 1) / fnorm) / wa2[permutation.indices()[j]]));
241 
242  /* test for convergence of the gradient norm. */
243  if (gnorm <= parameters.gtol) return LevenbergMarquardtSpace::CosinusTooSmall;
244 
245  /* rescale if necessary. */
246  if (!useExternalScaling) diag = diag.cwiseMax(wa2);
247 
248  do {
249  /* determine the levenberg-marquardt parameter. */
250  internal::lmpar2<Scalar>(qrfac, diag, qtf, delta, par, wa1);
251 
252  /* store the direction p and x + p. calculate the norm of p. */
253  wa1 = -wa1;
254  wa2 = x + wa1;
255  pnorm = diag.cwiseProduct(wa1).stableNorm();
256 
257  /* on the first iteration, adjust the initial step bound. */
258  if (iter == 1) delta = (std::min)(delta, pnorm);
259 
260  /* evaluate the function at x + p and calculate its norm. */
261  if (functor(wa2, wa4) < 0) return LevenbergMarquardtSpace::UserAsked;
262  ++nfev;
263  fnorm1 = wa4.stableNorm();
264 
265  /* compute the scaled actual reduction. */
266  actred = -1.;
267  if (Scalar(.1) * fnorm1 < fnorm) actred = 1. - numext::abs2(fnorm1 / fnorm);
268 
269  /* compute the scaled predicted reduction and */
270  /* the scaled directional derivative. */
271  wa3.noalias() = fjac.template triangularView<Upper>() * (qrfac.colsPermutation().inverse() * wa1);
272  temp1 = numext::abs2(wa3.stableNorm() / fnorm);
273  temp2 = numext::abs2(sqrt(par) * pnorm / fnorm);
274  prered = temp1 + temp2 / Scalar(.5);
275  dirder = -(temp1 + temp2);
276 
277  /* compute the ratio of the actual to the predicted */
278  /* reduction. */
279  ratio = 0.;
280  if (prered != 0.) ratio = actred / prered;
281 
282  /* update the step bound. */
283  if (ratio <= Scalar(.25)) {
284  if (actred >= 0.) temp = Scalar(.5);
285  if (actred < 0.) temp = Scalar(.5) * dirder / (dirder + Scalar(.5) * actred);
286  if (Scalar(.1) * fnorm1 >= fnorm || temp < Scalar(.1)) temp = Scalar(.1);
287  /* Computing MIN */
288  delta = temp * (std::min)(delta, pnorm / Scalar(.1));
289  par /= temp;
290  } else if (!(par != 0. && ratio < Scalar(.75))) {
291  delta = pnorm / Scalar(.5);
292  par = Scalar(.5) * par;
293  }
294 
295  /* test for successful iteration. */
296  if (ratio >= Scalar(1e-4)) {
297  /* successful iteration. update x, fvec, and their norms. */
298  x = wa2;
299  wa2 = diag.cwiseProduct(x);
300  fvec = wa4;
301  xnorm = wa2.stableNorm();
302  fnorm = fnorm1;
303  ++iter;
304  }
305 
306  /* tests for convergence. */
307  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1. &&
308  delta <= parameters.xtol * xnorm)
309  return LevenbergMarquardtSpace::RelativeErrorAndReductionTooSmall;
310  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1.)
311  return LevenbergMarquardtSpace::RelativeReductionTooSmall;
312  if (delta <= parameters.xtol * xnorm) return LevenbergMarquardtSpace::RelativeErrorTooSmall;
313 
314  /* tests for termination and stringent tolerances. */
315  if (nfev >= parameters.maxfev) return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
316  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() &&
317  Scalar(.5) * ratio <= 1.)
318  return LevenbergMarquardtSpace::FtolTooSmall;
319  if (delta <= NumTraits<Scalar>::epsilon() * xnorm) return LevenbergMarquardtSpace::XtolTooSmall;
320  if (gnorm <= NumTraits<Scalar>::epsilon()) return LevenbergMarquardtSpace::GtolTooSmall;
321 
322  } while (ratio < Scalar(1e-4));
323 
324  return LevenbergMarquardtSpace::Running;
325 }
326 
327 template <typename FunctorType, typename Scalar>
328 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::lmstr1(FVectorType &x, const Scalar tol) {
329  n = x.size();
330  m = functor.values();
331 
332  /* check the input parameters for errors. */
333  if (n <= 0 || m < n || tol < 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
334 
335  resetParameters();
336  parameters.ftol = tol;
337  parameters.xtol = tol;
338  parameters.maxfev = 100 * (n + 1);
339 
340  return minimizeOptimumStorage(x);
341 }
342 
343 template <typename FunctorType, typename Scalar>
344 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimizeOptimumStorageInit(FVectorType &x) {
345  n = x.size();
346  m = functor.values();
347 
348  wa1.resize(n);
349  wa2.resize(n);
350  wa3.resize(n);
351  wa4.resize(m);
352  fvec.resize(m);
353  // Only R is stored in fjac. Q is only used to compute 'qtf', which is
354  // Q.transpose()*rhs. qtf will be updated using givens rotation,
355  // instead of storing them in Q.
356  // The purpose it to only use a nxn matrix, instead of mxn here, so
357  // that we can handle cases where m>>n :
358  fjac.resize(n, n);
359  if (!useExternalScaling) diag.resize(n);
360  eigen_assert((!useExternalScaling || diag.size() == n) &&
361  "When useExternalScaling is set, the caller must provide a valid 'diag'");
362  qtf.resize(n);
363 
364  /* Function Body */
365  nfev = 0;
366  njev = 0;
367 
368  /* check the input parameters for errors. */
369  if (n <= 0 || m < n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. ||
370  parameters.maxfev <= 0 || parameters.factor <= 0.)
371  return LevenbergMarquardtSpace::ImproperInputParameters;
372 
373  if (useExternalScaling)
374  for (Index j = 0; j < n; ++j)
375  if (diag[j] <= 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
376 
377  /* evaluate the function at the starting point */
378  /* and calculate its norm. */
379  nfev = 1;
380  if (functor(x, fvec) < 0) return LevenbergMarquardtSpace::UserAsked;
381  fnorm = fvec.stableNorm();
382 
383  /* initialize levenberg-marquardt parameter and iteration counter. */
384  par = 0.;
385  iter = 1;
386 
387  return LevenbergMarquardtSpace::NotStarted;
388 }
389 
390 template <typename FunctorType, typename Scalar>
391 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimizeOptimumStorageOneStep(FVectorType &x) {
392  using std::abs;
393  using std::sqrt;
394 
395  eigen_assert(x.size() == n); // check the caller is not cheating us
396 
397  Index i, j;
398  bool sing;
399 
400  /* compute the qr factorization of the jacobian matrix */
401  /* calculated one row at a time, while simultaneously */
402  /* forming (q transpose)*fvec and storing the first */
403  /* n components in qtf. */
404  qtf.fill(0.);
405  fjac.fill(0.);
406  Index rownb = 2;
407  for (i = 0; i < m; ++i) {
408  if (functor.df(x, wa3, rownb) < 0) return LevenbergMarquardtSpace::UserAsked;
409  internal::rwupdt<Scalar>(fjac, wa3, qtf, fvec[i]);
410  ++rownb;
411  }
412  ++njev;
413 
414  /* if the jacobian is rank deficient, call qrfac to */
415  /* reorder its columns and update the components of qtf. */
416  sing = false;
417  for (j = 0; j < n; ++j) {
418  if (fjac(j, j) == 0.) sing = true;
419  wa2[j] = fjac.col(j).head(j).stableNorm();
420  }
421  permutation.setIdentity(n);
422  if (sing) {
423  wa2 = fjac.colwise().blueNorm();
424  // TODO We have no unit test covering this code path, do not modify
425  // until it is carefully tested
426  ColPivHouseholderQR<JacobianType> qrfac(fjac);
427  fjac = qrfac.matrixQR();
428  wa1 = fjac.diagonal();
429  fjac.diagonal() = qrfac.hCoeffs();
430  permutation = qrfac.colsPermutation();
431  // TODO : avoid this:
432  for (Index ii = 0; ii < fjac.cols(); ii++)
433  fjac.col(ii).segment(ii + 1, fjac.rows() - ii - 1) *= fjac(ii, ii); // rescale vectors
434 
435  for (j = 0; j < n; ++j) {
436  if (fjac(j, j) != 0.) {
437  sum = 0.;
438  for (i = j; i < n; ++i) sum += fjac(i, j) * qtf[i];
439  temp = -sum / fjac(j, j);
440  for (i = j; i < n; ++i) qtf[i] += fjac(i, j) * temp;
441  }
442  fjac(j, j) = wa1[j];
443  }
444  }
445 
446  /* on the first iteration and if external scaling is not used, scale according */
447  /* to the norms of the columns of the initial jacobian. */
448  if (iter == 1) {
449  if (!useExternalScaling)
450  for (j = 0; j < n; ++j) diag[j] = (wa2[j] == 0.) ? 1. : wa2[j];
451 
452  /* on the first iteration, calculate the norm of the scaled x */
453  /* and initialize the step bound delta. */
454  xnorm = diag.cwiseProduct(x).stableNorm();
455  delta = parameters.factor * xnorm;
456  if (delta == 0.) delta = parameters.factor;
457  }
458 
459  /* compute the norm of the scaled gradient. */
460  gnorm = 0.;
461  if (fnorm != 0.)
462  for (j = 0; j < n; ++j)
463  if (wa2[permutation.indices()[j]] != 0.)
464  gnorm = (std::max)(gnorm,
465  abs(fjac.col(j).head(j + 1).dot(qtf.head(j + 1) / fnorm) / wa2[permutation.indices()[j]]));
466 
467  /* test for convergence of the gradient norm. */
468  if (gnorm <= parameters.gtol) return LevenbergMarquardtSpace::CosinusTooSmall;
469 
470  /* rescale if necessary. */
471  if (!useExternalScaling) diag = diag.cwiseMax(wa2);
472 
473  do {
474  /* determine the levenberg-marquardt parameter. */
475  internal::lmpar<Scalar>(fjac, permutation.indices(), diag, qtf, delta, par, wa1);
476 
477  /* store the direction p and x + p. calculate the norm of p. */
478  wa1 = -wa1;
479  wa2 = x + wa1;
480  pnorm = diag.cwiseProduct(wa1).stableNorm();
481 
482  /* on the first iteration, adjust the initial step bound. */
483  if (iter == 1) delta = (std::min)(delta, pnorm);
484 
485  /* evaluate the function at x + p and calculate its norm. */
486  if (functor(wa2, wa4) < 0) return LevenbergMarquardtSpace::UserAsked;
487  ++nfev;
488  fnorm1 = wa4.stableNorm();
489 
490  /* compute the scaled actual reduction. */
491  actred = -1.;
492  if (Scalar(.1) * fnorm1 < fnorm) actred = 1. - numext::abs2(fnorm1 / fnorm);
493 
494  /* compute the scaled predicted reduction and */
495  /* the scaled directional derivative. */
496  wa3 = fjac.topLeftCorner(n, n).template triangularView<Upper>() * (permutation.inverse() * wa1);
497  temp1 = numext::abs2(wa3.stableNorm() / fnorm);
498  temp2 = numext::abs2(sqrt(par) * pnorm / fnorm);
499  prered = temp1 + temp2 / Scalar(.5);
500  dirder = -(temp1 + temp2);
501 
502  /* compute the ratio of the actual to the predicted */
503  /* reduction. */
504  ratio = 0.;
505  if (prered != 0.) ratio = actred / prered;
506 
507  /* update the step bound. */
508  if (ratio <= Scalar(.25)) {
509  if (actred >= 0.) temp = Scalar(.5);
510  if (actred < 0.) temp = Scalar(.5) * dirder / (dirder + Scalar(.5) * actred);
511  if (Scalar(.1) * fnorm1 >= fnorm || temp < Scalar(.1)) temp = Scalar(.1);
512  /* Computing MIN */
513  delta = temp * (std::min)(delta, pnorm / Scalar(.1));
514  par /= temp;
515  } else if (!(par != 0. && ratio < Scalar(.75))) {
516  delta = pnorm / Scalar(.5);
517  par = Scalar(.5) * par;
518  }
519 
520  /* test for successful iteration. */
521  if (ratio >= Scalar(1e-4)) {
522  /* successful iteration. update x, fvec, and their norms. */
523  x = wa2;
524  wa2 = diag.cwiseProduct(x);
525  fvec = wa4;
526  xnorm = wa2.stableNorm();
527  fnorm = fnorm1;
528  ++iter;
529  }
530 
531  /* tests for convergence. */
532  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1. &&
533  delta <= parameters.xtol * xnorm)
534  return LevenbergMarquardtSpace::RelativeErrorAndReductionTooSmall;
535  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1.)
536  return LevenbergMarquardtSpace::RelativeReductionTooSmall;
537  if (delta <= parameters.xtol * xnorm) return LevenbergMarquardtSpace::RelativeErrorTooSmall;
538 
539  /* tests for termination and stringent tolerances. */
540  if (nfev >= parameters.maxfev) return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
541  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() &&
542  Scalar(.5) * ratio <= 1.)
543  return LevenbergMarquardtSpace::FtolTooSmall;
544  if (delta <= NumTraits<Scalar>::epsilon() * xnorm) return LevenbergMarquardtSpace::XtolTooSmall;
545  if (gnorm <= NumTraits<Scalar>::epsilon()) return LevenbergMarquardtSpace::GtolTooSmall;
546 
547  } while (ratio < Scalar(1e-4));
548 
549  return LevenbergMarquardtSpace::Running;
550 }
551 
552 template <typename FunctorType, typename Scalar>
553 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::minimizeOptimumStorage(FVectorType &x) {
554  LevenbergMarquardtSpace::Status status = minimizeOptimumStorageInit(x);
555  if (status == LevenbergMarquardtSpace::ImproperInputParameters) return status;
556  do {
557  status = minimizeOptimumStorageOneStep(x);
558  } while (status == LevenbergMarquardtSpace::Running);
559  return status;
560 }
561 
562 template <typename FunctorType, typename Scalar>
563 LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType, Scalar>::lmdif1(FunctorType &functor, FVectorType &x,
564  Index *nfev, const Scalar tol) {
565  Index n = x.size();
566  Index m = functor.values();
567 
568  /* check the input parameters for errors. */
569  if (n <= 0 || m < n || tol < 0.) return LevenbergMarquardtSpace::ImproperInputParameters;
570 
571  NumericalDiff<FunctorType> numDiff(functor);
572  // embedded LevenbergMarquardt
573  LevenbergMarquardt<NumericalDiff<FunctorType>, Scalar> lm(numDiff);
574  lm.parameters.ftol = tol;
575  lm.parameters.xtol = tol;
576  lm.parameters.maxfev = 200 * (n + 1);
577 
578  LevenbergMarquardtSpace::Status info = LevenbergMarquardtSpace::Status(lm.minimize(x));
579  if (nfev) *nfev = lm.nfev;
580  return info;
581 }
582 
583 } // end namespace Eigen
584 
585 #endif // EIGEN_LEVENBERGMARQUARDT__H
586 
587 // vim: ai ts=4 sts=4 et sw=4
RealScalar xtol() const
Definition: LevenbergMarquardt.h:167
RealScalar fnorm()
Definition: LevenbergMarquardt.h:197
Index nfev()
Definition: LevenbergMarquardt.h:191
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
RealScalar factor() const
Definition: LevenbergMarquardt.h:176
Index maxfev() const
Definition: LevenbergMarquardt.h:182
void resetParameters()
Definition: LevenbergMarquardt.h:134
RealScalar ftol() const
Definition: LevenbergMarquardt.h:170
Index njev()
Definition: LevenbergMarquardt.h:194
RealScalar gnorm()
Definition: LevenbergMarquardt.h:200
FVectorType & fvec()
Definition: LevenbergMarquardt.h:207
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
PermutationType permutation()
Definition: LevenbergMarquardt.h:220
FVectorType & diag()
Definition: LevenbergMarquardt.h:185
RealScalar lm_param(void)
Definition: LevenbergMarquardt.h:203
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
RealScalar gtol() const
Definition: LevenbergMarquardt.h:173