lmpar.h
Go to the documentation of this file.
2 
3 namespace Eigen {
4 
5 namespace internal {
6 
7 template <typename Scalar>
8 void lmpar(
10  const VectorXi &ipvt,
11  const Matrix< Scalar, Dynamic, 1 > &diag,
13  Scalar delta,
14  Scalar &par,
16 {
17  using std::abs;
18  using std::sqrt;
19  typedef DenseIndex Index;
20 
21  /* Local variables */
22  Index i, j, l;
23  Scalar fp;
24  Scalar parc, parl;
25  Index iter;
26  Scalar temp, paru;
27  Scalar gnorm;
28  Scalar dxnorm;
29 
30 
31  /* Function Body */
32  const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
33  const Index n = r.cols();
34  eigen_assert(n==diag.size());
35  eigen_assert(n==qtb.size());
36  eigen_assert(n==x.size());
37 
39 
40  /* compute and store in x the gauss-newton direction. if the */
41  /* jacobian is rank-deficient, obtain a least squares solution. */
42  Index nsing = n-1;
43  wa1 = qtb;
44  for (j = 0; j < n; ++j) {
45  if (r(j,j) == 0. && nsing == n-1)
46  nsing = j - 1;
47  if (nsing < n-1)
48  wa1[j] = 0.;
49  }
50  for (j = nsing; j>=0; --j) {
51  wa1[j] /= r(j,j);
52  temp = wa1[j];
53  for (i = 0; i < j ; ++i)
54  wa1[i] -= r(i,j) * temp;
55  }
56 
57  for (j = 0; j < n; ++j)
58  x[ipvt[j]] = wa1[j];
59 
60  /* initialize the iteration counter. */
61  /* evaluate the function at the origin, and test */
62  /* for acceptance of the gauss-newton direction. */
63  iter = 0;
64  wa2 = diag.cwiseProduct(x);
65  dxnorm = wa2.blueNorm();
66  fp = dxnorm - delta;
67  if (fp <= Scalar(0.1) * delta) {
68  par = 0;
69  return;
70  }
71 
72  /* if the jacobian is not rank deficient, the newton */
73  /* step provides a lower bound, parl, for the zero of */
74  /* the function. otherwise set this bound to zero. */
75  parl = 0.;
76  if (nsing >= n-1) {
77  for (j = 0; j < n; ++j) {
78  l = ipvt[j];
79  wa1[j] = diag[l] * (wa2[l] / dxnorm);
80  }
81  // it's actually a triangularView.solveInplace(), though in a weird
82  // way:
83  for (j = 0; j < n; ++j) {
84  Scalar sum = 0.;
85  for (i = 0; i < j; ++i)
86  sum += r(i,j) * wa1[i];
87  wa1[j] = (wa1[j] - sum) / r(j,j);
88  }
89  temp = wa1.blueNorm();
90  parl = fp / delta / temp / temp;
91  }
92 
93  /* calculate an upper bound, paru, for the zero of the function. */
94  for (j = 0; j < n; ++j)
95  wa1[j] = r.col(j).head(j+1).dot(qtb.head(j+1)) / diag[ipvt[j]];
96 
97  gnorm = wa1.stableNorm();
98  paru = gnorm / delta;
99  if (paru == 0.)
100  paru = dwarf / (std::min)(delta,Scalar(0.1));
101 
102  /* if the input par lies outside of the interval (parl,paru), */
103  /* set par to the closer endpoint. */
104  par = (std::max)(par,parl);
105  par = (std::min)(par,paru);
106  if (par == 0.)
107  par = gnorm / dxnorm;
108 
109  /* beginning of an iteration. */
110  while (true) {
111  ++iter;
112 
113  /* evaluate the function at the current value of par. */
114  if (par == 0.)
115  par = (std::max)(dwarf,Scalar(.001) * paru); /* Computing MAX */
116  wa1 = sqrt(par)* diag;
117 
119  qrsolv<Scalar>(r, ipvt, wa1, qtb, x, sdiag);
120 
121  wa2 = diag.cwiseProduct(x);
122  dxnorm = wa2.blueNorm();
123  temp = fp;
124  fp = dxnorm - delta;
125 
126  /* if the function is small enough, accept the current value */
127  /* of par. also test for the exceptional cases where parl */
128  /* is zero or the number of iterations has reached 10. */
129  if (abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
130  break;
131 
132  /* compute the newton correction. */
133  for (j = 0; j < n; ++j) {
134  l = ipvt[j];
135  wa1[j] = diag[l] * (wa2[l] / dxnorm);
136  }
137  for (j = 0; j < n; ++j) {
138  wa1[j] /= sdiag[j];
139  temp = wa1[j];
140  for (i = j+1; i < n; ++i)
141  wa1[i] -= r(i,j) * temp;
142  }
143  temp = wa1.blueNorm();
144  parc = fp / delta / temp / temp;
145 
146  /* depending on the sign of the function, update parl or paru. */
147  if (fp > 0.)
148  parl = (std::max)(parl,par);
149  if (fp < 0.)
150  paru = (std::min)(paru,par);
151 
152  /* compute an improved estimate for par. */
153  /* Computing MAX */
154  par = (std::max)(parl,par+parc);
155 
156  /* end of an iteration. */
157  }
158 
159  /* termination. */
160  if (iter == 0)
161  par = 0.;
162  return;
163 }
164 
165 template <typename Scalar>
166 void lmpar2(
168  const Matrix< Scalar, Dynamic, 1 > &diag,
169  const Matrix< Scalar, Dynamic, 1 > &qtb,
170  Scalar delta,
171  Scalar &par,
173 
174 {
175  using std::sqrt;
176  using std::abs;
177  typedef DenseIndex Index;
178 
179  /* Local variables */
180  Index j;
181  Scalar fp;
182  Scalar parc, parl;
183  Index iter;
184  Scalar temp, paru;
185  Scalar gnorm;
186  Scalar dxnorm;
187 
188 
189  /* Function Body */
190  const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
191  const Index n = qr.matrixQR().cols();
192  eigen_assert(n==diag.size());
193  eigen_assert(n==qtb.size());
194 
196 
197  /* compute and store in x the gauss-newton direction. if the */
198  /* jacobian is rank-deficient, obtain a least squares solution. */
199 
200 // const Index rank = qr.nonzeroPivots(); // exactly double(0.)
201  const Index rank = qr.rank(); // use a threshold
202  wa1 = qtb;
203  wa1.tail(n-rank).setZero();
204  qr.matrixQR().topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(wa1.head(rank));
205 
206  x = qr.colsPermutation()*wa1;
207 
208  /* initialize the iteration counter. */
209  /* evaluate the function at the origin, and test */
210  /* for acceptance of the gauss-newton direction. */
211  iter = 0;
212  wa2 = diag.cwiseProduct(x);
213  dxnorm = wa2.blueNorm();
214  fp = dxnorm - delta;
215  if (fp <= Scalar(0.1) * delta) {
216  par = 0;
217  return;
218  }
219 
220  /* if the jacobian is not rank deficient, the newton */
221  /* step provides a lower bound, parl, for the zero of */
222  /* the function. otherwise set this bound to zero. */
223  parl = 0.;
224  if (rank==n) {
225  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2)/dxnorm;
226  qr.matrixQR().topLeftCorner(n, n).transpose().template triangularView<Lower>().solveInPlace(wa1);
227  temp = wa1.blueNorm();
228  parl = fp / delta / temp / temp;
229  }
230 
231  /* calculate an upper bound, paru, for the zero of the function. */
232  for (j = 0; j < n; ++j)
233  wa1[j] = qr.matrixQR().col(j).head(j+1).dot(qtb.head(j+1)) / diag[qr.colsPermutation().indices()(j)];
234 
235  gnorm = wa1.stableNorm();
236  paru = gnorm / delta;
237  if (paru == 0.)
238  paru = dwarf / (std::min)(delta,Scalar(0.1));
239 
240  /* if the input par lies outside of the interval (parl,paru), */
241  /* set par to the closer endpoint. */
242  par = (std::max)(par,parl);
243  par = (std::min)(par,paru);
244  if (par == 0.)
245  par = gnorm / dxnorm;
246 
247  /* beginning of an iteration. */
248  Matrix< Scalar, Dynamic, Dynamic > s = qr.matrixQR();
249  while (true) {
250  ++iter;
251 
252  /* evaluate the function at the current value of par. */
253  if (par == 0.)
254  par = (std::max)(dwarf,Scalar(.001) * paru); /* Computing MAX */
255  wa1 = sqrt(par)* diag;
256 
258  qrsolv<Scalar>(s, qr.colsPermutation().indices(), wa1, qtb, x, sdiag);
259 
260  wa2 = diag.cwiseProduct(x);
261  dxnorm = wa2.blueNorm();
262  temp = fp;
263  fp = dxnorm - delta;
264 
265  /* if the function is small enough, accept the current value */
266  /* of par. also test for the exceptional cases where parl */
267  /* is zero or the number of iterations has reached 10. */
268  if (abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
269  break;
270 
271  /* compute the newton correction. */
272  wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2/dxnorm);
273  // we could almost use this here, but the diagonal is outside qr, in sdiag[]
274  // qr.matrixQR().topLeftCorner(n, n).transpose().template triangularView<Lower>().solveInPlace(wa1);
275  for (j = 0; j < n; ++j) {
276  wa1[j] /= sdiag[j];
277  temp = wa1[j];
278  for (Index i = j+1; i < n; ++i)
279  wa1[i] -= s(i,j) * temp;
280  }
281  temp = wa1.blueNorm();
282  parc = fp / delta / temp / temp;
283 
284  /* depending on the sign of the function, update parl or paru. */
285  if (fp > 0.)
286  parl = (std::max)(parl,par);
287  if (fp < 0.)
288  paru = (std::min)(paru,par);
289 
290  /* compute an improved estimate for par. */
291  par = (std::max)(parl,par+parc);
292  }
293  if (iter == 0)
294  par = 0.;
295  return;
296 }
297 
298 } // end namespace internal
299 
300 } // end namespace Eigen
int n
int i
HouseholderQR< MatrixXf > qr(A)
#define eigen_assert(x)
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Derived & setZero(Index rows, Index cols)
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
void lmpar(Matrix< Scalar, Dynamic, Dynamic > &r, const VectorXi &ipvt, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Scalar delta, Scalar &par, Matrix< Scalar, Dynamic, 1 > &x)
Definition: lmpar.h:8
: 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)
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
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