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 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 
21  template <typename QRSolver, typename VectorType>
22  void lmpar2(
23  const QRSolver &qr,
24  const VectorType &diag,
25  const VectorType &qtb,
26  typename VectorType::Scalar m_delta,
27  typename VectorType::Scalar &par,
28  VectorType &x)
29 
30  {
31  using std::sqrt;
32  using std::abs;
33  typedef typename QRSolver::MatrixType MatrixType;
34  typedef typename QRSolver::Scalar Scalar;
35 // typedef typename QRSolver::StorageIndex StorageIndex;
36 
37  /* Local variables */
38  Index j;
39  Scalar fp;
40  Scalar parc, parl;
41  Index iter;
42  Scalar temp, paru;
43  Scalar gnorm;
44  Scalar dxnorm;
45 
46  // Make a copy of the triangular factor.
47  // This copy is modified during call the qrsolv
48  MatrixType s;
49  s = qr.matrixR();
50 
51  /* Function Body */
52  const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
53  const Index n = qr.matrixR().cols();
54  eigen_assert(n==diag.size());
55  eigen_assert(n==qtb.size());
56 
57  VectorType wa1, wa2;
58 
59  /* compute and store in x the gauss-newton direction. if the */
60  /* jacobian is rank-deficient, obtain a least squares solution. */
61 
62  // const Index rank = qr.nonzeroPivots(); // exactly double(0.)
63  const Index rank = qr.rank(); // use a threshold
64  wa1 = qtb;
65  wa1.tail(n-rank).setZero();
66  //FIXME There is no solve in place for sparse triangularView
67  wa1.head(rank) = s.topLeftCorner(rank,rank).template triangularView<Upper>().solve(qtb.head(rank));
68 
69  x = qr.colsPermutation()*wa1;
70 
71  /* initialize the iteration counter. */
72  /* evaluate the function at the origin, and test */
73  /* for acceptance of the gauss-newton direction. */
74  iter = 0;
75  wa2 = diag.cwiseProduct(x);
76  dxnorm = wa2.blueNorm();
77  fp = dxnorm - m_delta;
78  if (fp <= Scalar(0.1) * m_delta) {
79  par = 0;
80  return;
81  }
82 
83  /* if the jacobian is not rank deficient, the newton */
84  /* step provides a lower bound, parl, for the zero of */
85  /* the function. otherwise set this bound to zero. */
86  parl = 0.;
87  if (rank==n) {
88  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2)/dxnorm;
89  s.topLeftCorner(n,n).transpose().template triangularView<Lower>().solveInPlace(wa1);
90  temp = wa1.blueNorm();
91  parl = fp / m_delta / temp / temp;
92  }
93 
94  /* calculate an upper bound, paru, for the zero of the function. */
95  for (j = 0; j < n; ++j)
96  wa1[j] = s.col(j).head(j+1).dot(qtb.head(j+1)) / diag[qr.colsPermutation().indices()(j)];
97 
98  gnorm = wa1.stableNorm();
99  paru = gnorm / m_delta;
100  if (paru == 0.)
101  paru = dwarf / (std::min)(m_delta,Scalar(0.1));
102 
103  /* if the input par lies outside of the interval (parl,paru), */
104  /* set par to the closer endpoint. */
105  par = (std::max)(par,parl);
106  par = (std::min)(par,paru);
107  if (par == 0.)
108  par = gnorm / dxnorm;
109 
110  /* beginning of an iteration. */
111  while (true) {
112  ++iter;
113 
114  /* evaluate the function at the current value of par. */
115  if (par == 0.)
116  par = (std::max)(dwarf,Scalar(.001) * paru); /* Computing MAX */
117  wa1 = sqrt(par)* diag;
118 
119  VectorType sdiag(n);
120  lmqrsolv(s, qr.colsPermutation(), wa1, qtb, x, sdiag);
121 
122  wa2 = diag.cwiseProduct(x);
123  dxnorm = wa2.blueNorm();
124  temp = fp;
125  fp = dxnorm - m_delta;
126 
127  /* if the function is small enough, accept the current value */
128  /* of par. also test for the exceptional cases where parl */
129  /* is zero or the number of iterations has reached 10. */
130  if (abs(fp) <= Scalar(0.1) * m_delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
131  break;
132 
133  /* compute the newton correction. */
134  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2/dxnorm);
135  // we could almost use this here, but the diagonal is outside qr, in sdiag[]
136  for (j = 0; j < n; ++j) {
137  wa1[j] /= sdiag[j];
138  temp = wa1[j];
139  for (Index i = j+1; i < n; ++i)
140  wa1[i] -= s.coeff(i,j) * temp;
141  }
142  temp = wa1.blueNorm();
143  parc = fp / m_delta / temp / temp;
144 
145  /* depending on the sign of the function, update parl or paru. */
146  if (fp > 0.)
147  parl = (std::max)(parl,par);
148  if (fp < 0.)
149  paru = (std::min)(paru,par);
150 
151  /* compute an improved estimate for par. */
152  par = (std::max)(parl,par+parc);
153  }
154  if (iter == 0)
155  par = 0.;
156  return;
157  }
158 } // end namespace internal
159 
160 } // end namespace Eigen
161 
162 #endif // EIGEN_LMPAR_H
int n
int i
HouseholderQR< MatrixXf > qr(A)
#define eigen_assert(x)
Matrix< float, 1, Dynamic > MatrixType
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:25
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
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