LMpar.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 // This code initially comes from MINPACK whose original authors are:
5 // Copyright Jorge More - Argonne National Laboratory
6 // Copyright Burt Garbow - Argonne National Laboratory
7 // Copyright Ken Hillstrom - Argonne National Laboratory
8 //
9 // This Source Code Form is subject to the terms of the Minpack license
10 // (a BSD-like license) described in the campaigned CopyrightMINPACK.txt file.
11 
12 #ifndef EIGEN_LMPAR_H
13 #define EIGEN_LMPAR_H
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21 
22 template <typename QRSolver, typename VectorType>
23 void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta,
24  typename VectorType::Scalar &par, VectorType &x)
25 
26 {
27  using std::abs;
28  using std::sqrt;
29  typedef typename QRSolver::MatrixType MatrixType;
30  typedef typename QRSolver::Scalar Scalar;
31  // typedef typename QRSolver::StorageIndex StorageIndex;
32 
33  /* Local variables */
34  Index j;
35  Scalar fp;
36  Scalar parc, parl;
37  Index iter;
38  Scalar temp, paru;
39  Scalar gnorm;
40  Scalar dxnorm;
41 
42  // Make a copy of the triangular factor.
43  // This copy is modified during call the qrsolv
44  MatrixType s;
45  s = qr.matrixR();
46 
47  /* Function Body */
48  const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
49  const Index n = qr.matrixR().cols();
50  eigen_assert(n == diag.size());
51  eigen_assert(n == qtb.size());
52 
53  VectorType wa1, wa2;
54 
55  /* compute and store in x the gauss-newton direction. if the */
56  /* jacobian is rank-deficient, obtain a least squares solution. */
57 
58  // const Index rank = qr.nonzeroPivots(); // exactly double(0.)
59  const Index rank = qr.rank(); // use a threshold
60  wa1 = qtb;
61  wa1.tail(n - rank).setZero();
62  // FIXME There is no solve in place for sparse triangularView
63  wa1.head(rank) = s.topLeftCorner(rank, rank).template triangularView<Upper>().solve(qtb.head(rank));
64 
65  x = qr.colsPermutation() * wa1;
66 
67  /* initialize the iteration counter. */
68  /* evaluate the function at the origin, and test */
69  /* for acceptance of the gauss-newton direction. */
70  iter = 0;
71  wa2 = diag.cwiseProduct(x);
72  dxnorm = wa2.blueNorm();
73  fp = dxnorm - m_delta;
74  if (fp <= Scalar(0.1) * m_delta) {
75  par = 0;
76  return;
77  }
78 
79  /* if the jacobian is not rank deficient, the newton */
80  /* step provides a lower bound, parl, for the zero of */
81  /* the function. otherwise set this bound to zero. */
82  parl = 0.;
83  if (rank == n) {
84  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2) / dxnorm;
85  s.topLeftCorner(n, n).transpose().template triangularView<Lower>().solveInPlace(wa1);
86  temp = wa1.blueNorm();
87  parl = fp / m_delta / temp / temp;
88  }
89 
90  /* calculate an upper bound, paru, for the zero of the function. */
91  for (j = 0; j < n; ++j) wa1[j] = s.col(j).head(j + 1).dot(qtb.head(j + 1)) / diag[qr.colsPermutation().indices()(j)];
92 
93  gnorm = wa1.stableNorm();
94  paru = gnorm / m_delta;
95  if (paru == 0.) paru = dwarf / (std::min)(m_delta, Scalar(0.1));
96 
97  /* if the input par lies outside of the interval (parl,paru), */
98  /* set par to the closer endpoint. */
99  par = (std::max)(par, parl);
100  par = (std::min)(par, paru);
101  if (par == 0.) par = gnorm / dxnorm;
102 
103  /* beginning of an iteration. */
104  while (true) {
105  ++iter;
106 
107  /* evaluate the function at the current value of par. */
108  if (par == 0.) par = (std::max)(dwarf, Scalar(.001) * paru); /* Computing MAX */
109  wa1 = sqrt(par) * diag;
110 
111  VectorType sdiag(n);
112  lmqrsolv(s, qr.colsPermutation(), wa1, qtb, x, sdiag);
113 
114  wa2 = diag.cwiseProduct(x);
115  dxnorm = wa2.blueNorm();
116  temp = fp;
117  fp = dxnorm - m_delta;
118 
119  /* if the function is small enough, accept the current value */
120  /* of par. also test for the exceptional cases where parl */
121  /* is zero or the number of iterations has reached 10. */
122  if (abs(fp) <= Scalar(0.1) * m_delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10) break;
123 
124  /* compute the newton correction. */
125  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2 / dxnorm);
126  // we could almost use this here, but the diagonal is outside qr, in sdiag[]
127  for (j = 0; j < n; ++j) {
128  wa1[j] /= sdiag[j];
129  temp = wa1[j];
130  for (Index i = j + 1; i < n; ++i) wa1[i] -= s.coeff(i, j) * temp;
131  }
132  temp = wa1.blueNorm();
133  parc = fp / m_delta / temp / temp;
134 
135  /* depending on the sign of the function, update parl or paru. */
136  if (fp > 0.) parl = (std::max)(parl, par);
137  if (fp < 0.) paru = (std::min)(paru, par);
138 
139  /* compute an improved estimate for par. */
140  par = (std::max)(parl, par + parc);
141  }
142  if (iter == 0) par = 0.;
143  return;
144 }
145 } // end namespace internal
146 
147 } // end namespace Eigen
148 
149 #endif // EIGEN_LMPAR_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
RealScalar s
Definition: level1_cplx_impl.h:130
const char const char const char * diag
Definition: level2_impl.h:86
void lmqrsolv(Matrix< Scalar, Rows, Cols > &s, const PermutationMatrix< Dynamic, Dynamic, PermIndex > &iPerm, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &sdiag)
Definition: LMqrsolv.h:26
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
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
string par
Definition: calibrate.py:135
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: fft_test_shared.h:66
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2