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 #include "./InternalHeaderCheck.h"
18 
19 namespace Eigen {
20 
21 template<typename FunctorType>
24 {
25  using std::abs;
26  using std::sqrt;
27  RealScalar temp, temp1,temp2;
28  RealScalar ratio;
29  RealScalar pnorm, xnorm, fnorm1, actred, dirder, prered;
30  eigen_assert(x.size()==n); // check the caller is not cheating us
31 
32  temp = 0.0; xnorm = 0.0;
33  /* calculate the jacobian matrix. */
34  Index df_ret = m_functor.df(x, m_fjac);
35  if (df_ret<0)
37  if (df_ret>0)
38  // numerical diff, we evaluated the function df_ret times
39  m_nfev += df_ret;
40  else m_njev++;
41 
42  /* compute the qr factorization of the jacobian. */
43  for (int j = 0; j < x.size(); ++j)
44  m_wa2(j) = m_fjac.col(j).blueNorm();
45  QRSolver qrfac(m_fjac);
46  if(qrfac.info() != Success) {
47  m_info = NumericalIssue;
49  }
50  // Make a copy of the first factor with the associated permutation
51  m_rfactor = qrfac.matrixR();
52  m_permutation = (qrfac.colsPermutation());
53 
54  /* on the first iteration and if external scaling is not used, scale according */
55  /* to the norms of the columns of the initial jacobian. */
56  if (m_iter == 1) {
57  if (!m_useExternalScaling)
58  for (Index j = 0; j < n; ++j)
59  m_diag[j] = (m_wa2[j]==0.)? 1. : m_wa2[j];
60 
61  /* on the first iteration, calculate the norm of the scaled x */
62  /* and initialize the step bound m_delta. */
63  xnorm = m_diag.cwiseProduct(x).stableNorm();
64  m_delta = m_factor * xnorm;
65  if (m_delta == 0.)
66  m_delta = m_factor;
67  }
68 
69  /* form (q transpose)*m_fvec and store the first n components in */
70  /* m_qtf. */
71  m_wa4 = m_fvec;
72  m_wa4 = qrfac.matrixQ().adjoint() * m_fvec;
73  m_qtf = m_wa4.head(n);
74 
75  /* compute the norm of the scaled gradient. */
76  m_gnorm = 0.;
77  if (m_fnorm != 0.)
78  for (Index j = 0; j < n; ++j)
79  if (m_wa2[m_permutation.indices()[j]] != 0.)
80  m_gnorm = (std::max)(m_gnorm, abs( m_rfactor.col(j).head(j+1).dot(m_qtf.head(j+1)/m_fnorm) / m_wa2[m_permutation.indices()[j]]));
81 
82  /* test for convergence of the gradient norm. */
83  if (m_gnorm <= m_gtol) {
84  m_info = Success;
86  }
87 
88  /* rescale if necessary. */
89  if (!m_useExternalScaling)
90  m_diag = m_diag.cwiseMax(m_wa2);
91 
92  do {
93  /* determine the levenberg-marquardt parameter. */
94  internal::lmpar2(qrfac, m_diag, m_qtf, m_delta, m_par, m_wa1);
95 
96  /* store the direction p and x + p. calculate the norm of p. */
97  m_wa1 = -m_wa1;
98  m_wa2 = x + m_wa1;
99  pnorm = m_diag.cwiseProduct(m_wa1).stableNorm();
100 
101  /* on the first iteration, adjust the initial step bound. */
102  if (m_iter == 1)
103  m_delta = (std::min)(m_delta,pnorm);
104 
105  /* evaluate the function at x + p and calculate its norm. */
106  if ( m_functor(m_wa2, m_wa4) < 0)
108  ++m_nfev;
109  fnorm1 = m_wa4.stableNorm();
110 
111  /* compute the scaled actual reduction. */
112  actred = -1.;
113  if (Scalar(.1) * fnorm1 < m_fnorm)
114  actred = 1. - numext::abs2(fnorm1 / m_fnorm);
115 
116  /* compute the scaled predicted reduction and */
117  /* the scaled directional derivative. */
118  m_wa3 = m_rfactor.template triangularView<Upper>() * (m_permutation.inverse() *m_wa1);
119  temp1 = numext::abs2(m_wa3.stableNorm() / m_fnorm);
120  temp2 = numext::abs2(sqrt(m_par) * pnorm / m_fnorm);
121  prered = temp1 + temp2 / Scalar(.5);
122  dirder = -(temp1 + temp2);
123 
124  /* compute the ratio of the actual to the predicted */
125  /* reduction. */
126  ratio = 0.;
127  if (prered != 0.)
128  ratio = actred / prered;
129 
130  /* update the step bound. */
131  if (ratio <= Scalar(.25)) {
132  if (actred >= 0.)
133  temp = RealScalar(.5);
134  if (actred < 0.)
135  temp = RealScalar(.5) * dirder / (dirder + RealScalar(.5) * actred);
136  if (RealScalar(.1) * fnorm1 >= m_fnorm || temp < RealScalar(.1))
137  temp = Scalar(.1);
138  /* Computing MIN */
139  m_delta = temp * (std::min)(m_delta, pnorm / RealScalar(.1));
140  m_par /= temp;
141  } else if (!(m_par != 0. && ratio < RealScalar(.75))) {
142  m_delta = pnorm / RealScalar(.5);
143  m_par = RealScalar(.5) * m_par;
144  }
145 
146  /* test for successful iteration. */
147  if (ratio >= RealScalar(1e-4)) {
148  /* successful iteration. update x, m_fvec, and their norms. */
149  x = m_wa2;
150  m_wa2 = m_diag.cwiseProduct(x);
151  m_fvec = m_wa4;
152  xnorm = m_wa2.stableNorm();
153  m_fnorm = fnorm1;
154  ++m_iter;
155  }
156 
157  /* tests for convergence. */
158  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1. && m_delta <= m_xtol * xnorm)
159  {
160  m_info = Success;
162  }
163  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1.)
164  {
165  m_info = Success;
167  }
168  if (m_delta <= m_xtol * xnorm)
169  {
170  m_info = Success;
172  }
173 
174  /* tests for termination and stringent tolerances. */
175  if (m_nfev >= m_maxfev)
176  {
177  m_info = NoConvergence;
179  }
180  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
181  {
182  m_info = Success;
184  }
185  if (m_delta <= NumTraits<Scalar>::epsilon() * xnorm)
186  {
187  m_info = Success;
189  }
190  if (m_gnorm <= NumTraits<Scalar>::epsilon())
191  {
192  m_info = Success;
194  }
195 
196  } while (ratio < Scalar(1e-4));
197 
199 }
200 
201 
202 } // end namespace Eigen
203 
204 #endif // EIGEN_LMONESTEP_H
int n
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define eigen_assert(x)
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
NumericalIssue
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
bool abs2(bool x)
: 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)
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
CleanedUpDerType< DerType >::type() max(const AutoDiffScalar< DerType > &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)
Definition: AdolcForward:74
std::ptrdiff_t j