LMonestep.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 //
6 // This code initially comes from MINPACK whose original authors are:
7 // Copyright Jorge More - Argonne National Laboratory
8 // Copyright Burt Garbow - Argonne National Laboratory
9 // Copyright Ken Hillstrom - Argonne National Laboratory
10 //
11 // This Source Code Form is subject to the terms of the Minpack license
12 // (a BSD-like license) described in the campaigned CopyrightMINPACK.txt file.
13 
14 #ifndef EIGEN_LMONESTEP_H
15 #define EIGEN_LMONESTEP_H
16 
17 // IWYU pragma: private
18 #include "./InternalHeaderCheck.h"
19 
20 namespace Eigen {
21 
22 template <typename FunctorType>
24  using std::abs;
25  using std::sqrt;
26  RealScalar temp, temp1, temp2;
27  RealScalar ratio;
28  RealScalar pnorm, xnorm, fnorm1, actred, dirder, prered;
29  eigen_assert(x.size() == n); // check the caller is not cheating us
30 
31  temp = 0.0;
32  xnorm = 0.0;
33  /* calculate the jacobian matrix. */
34  Index df_ret = m_functor.df(x, m_fjac);
35  if (df_ret < 0) return LevenbergMarquardtSpace::UserAsked;
36  if (df_ret > 0)
37  // numerical diff, we evaluated the function df_ret times
38  m_nfev += df_ret;
39  else
40  m_njev++;
41 
42  /* compute the qr factorization of the jacobian. */
43  for (int j = 0; j < x.size(); ++j) m_wa2(j) = m_fjac.col(j).blueNorm();
44  QRSolver qrfac(m_fjac);
45  if (qrfac.info() != Success) {
46  m_info = NumericalIssue;
48  }
49  // Make a copy of the first factor with the associated permutation
50  m_rfactor = qrfac.matrixR();
51  m_permutation = (qrfac.colsPermutation());
52 
53  /* on the first iteration and if external scaling is not used, scale according */
54  /* to the norms of the columns of the initial jacobian. */
55  if (m_iter == 1) {
56  if (!m_useExternalScaling)
57  for (Index j = 0; j < n; ++j) m_diag[j] = (m_wa2[j] == 0.) ? 1. : m_wa2[j];
58 
59  /* on the first iteration, calculate the norm of the scaled x */
60  /* and initialize the step bound m_delta. */
61  xnorm = m_diag.cwiseProduct(x).stableNorm();
62  m_delta = m_factor * xnorm;
63  if (m_delta == 0.) m_delta = m_factor;
64  }
65 
66  /* form (q transpose)*m_fvec and store the first n components in */
67  /* m_qtf. */
68  m_wa4 = m_fvec;
69  m_wa4 = qrfac.matrixQ().adjoint() * m_fvec;
70  m_qtf = m_wa4.head(n);
71 
72  /* compute the norm of the scaled gradient. */
73  m_gnorm = 0.;
74  if (m_fnorm != 0.)
75  for (Index j = 0; j < n; ++j)
76  if (m_wa2[m_permutation.indices()[j]] != 0.)
77  m_gnorm = (std::max)(m_gnorm, abs(m_rfactor.col(j).head(j + 1).dot(m_qtf.head(j + 1) / m_fnorm) /
78  m_wa2[m_permutation.indices()[j]]));
79 
80  /* test for convergence of the gradient norm. */
81  if (m_gnorm <= m_gtol) {
82  m_info = Success;
84  }
85 
86  /* rescale if necessary. */
87  if (!m_useExternalScaling) m_diag = m_diag.cwiseMax(m_wa2);
88 
89  do {
90  /* determine the levenberg-marquardt parameter. */
91  internal::lmpar2(qrfac, m_diag, m_qtf, m_delta, m_par, m_wa1);
92 
93  /* store the direction p and x + p. calculate the norm of p. */
94  m_wa1 = -m_wa1;
95  m_wa2 = x + m_wa1;
96  pnorm = m_diag.cwiseProduct(m_wa1).stableNorm();
97 
98  /* on the first iteration, adjust the initial step bound. */
99  if (m_iter == 1) m_delta = (std::min)(m_delta, pnorm);
100 
101  /* evaluate the function at x + p and calculate its norm. */
102  if (m_functor(m_wa2, m_wa4) < 0) return LevenbergMarquardtSpace::UserAsked;
103  ++m_nfev;
104  fnorm1 = m_wa4.stableNorm();
105 
106  /* compute the scaled actual reduction. */
107  actred = -1.;
108  if (Scalar(.1) * fnorm1 < m_fnorm) actred = 1. - numext::abs2(fnorm1 / m_fnorm);
109 
110  /* compute the scaled predicted reduction and */
111  /* the scaled directional derivative. */
112  m_wa3 = m_rfactor.template triangularView<Upper>() * (m_permutation.inverse() * m_wa1);
113  temp1 = numext::abs2(m_wa3.stableNorm() / m_fnorm);
114  temp2 = numext::abs2(sqrt(m_par) * pnorm / m_fnorm);
115  prered = temp1 + temp2 / Scalar(.5);
116  dirder = -(temp1 + temp2);
117 
118  /* compute the ratio of the actual to the predicted */
119  /* reduction. */
120  ratio = 0.;
121  if (prered != 0.) ratio = actred / prered;
122 
123  /* update the step bound. */
124  if (ratio <= Scalar(.25)) {
125  if (actred >= 0.) temp = RealScalar(.5);
126  if (actred < 0.) temp = RealScalar(.5) * dirder / (dirder + RealScalar(.5) * actred);
127  if (RealScalar(.1) * fnorm1 >= m_fnorm || temp < RealScalar(.1)) temp = Scalar(.1);
128  /* Computing MIN */
129  m_delta = temp * (std::min)(m_delta, pnorm / RealScalar(.1));
130  m_par /= temp;
131  } else if (!(m_par != 0. && ratio < RealScalar(.75))) {
132  m_delta = pnorm / RealScalar(.5);
133  m_par = RealScalar(.5) * m_par;
134  }
135 
136  /* test for successful iteration. */
137  if (ratio >= RealScalar(1e-4)) {
138  /* successful iteration. update x, m_fvec, and their norms. */
139  x = m_wa2;
140  m_wa2 = m_diag.cwiseProduct(x);
141  m_fvec = m_wa4;
142  xnorm = m_wa2.stableNorm();
143  m_fnorm = fnorm1;
144  ++m_iter;
145  }
146 
147  /* tests for convergence. */
148  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1. && m_delta <= m_xtol * xnorm) {
149  m_info = Success;
151  }
152  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1.) {
153  m_info = Success;
155  }
156  if (m_delta <= m_xtol * xnorm) {
157  m_info = Success;
159  }
160 
161  /* tests for termination and stringent tolerances. */
162  if (m_nfev >= m_maxfev) {
163  m_info = NoConvergence;
165  }
166  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() &&
167  Scalar(.5) * ratio <= 1.) {
168  m_info = Success;
170  }
171  if (m_delta <= NumTraits<Scalar>::epsilon() * xnorm) {
172  m_info = Success;
174  }
175  if (m_gnorm <= NumTraits<Scalar>::epsilon()) {
176  m_info = Success;
178  }
179 
180  } while (ratio < Scalar(1e-4));
181 
183 }
184 
185 } // end namespace Eigen
186 
187 #endif // EIGEN_LMONESTEP_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define eigen_assert(x)
Definition: Macros.h:910
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
FunctorType::QRSolver QRSolver
Definition: LevenbergMarquardt/LevenbergMarquardt.h:105
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
JacobianType::RealScalar RealScalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:108
DenseIndex Index
Definition: NonLinearOptimization/LevenbergMarquardt.h:60
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
Status
Definition: LevenbergMarquardt/LevenbergMarquardt.h:27
@ RelativeReductionTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:31
@ GtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:38
@ 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
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
list x
Definition: plotDoE.py:28
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2