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