2#include "./InternalHeaderCheck.h"
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) {
13 typedef DenseIndex
Index;
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());
31 Matrix<Scalar, Dynamic, 1> wa1, wa2;
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.;
41 for (j = nsing; j >= 0; --j) {
44 for (i = 0; i < j; ++i) wa1[i] -= r(i, j) * temp;
47 for (j = 0; j < n; ++j) x[ipvt[j]] = wa1[j];
53 wa2 = diag.cwiseProduct(x);
54 dxnorm = wa2.blueNorm();
56 if (fp <= Scalar(0.1) * delta) {
66 for (j = 0; j < n; ++j) {
68 wa1[j] = diag[l] * (wa2[l] / dxnorm);
72 for (j = 0; j < n; ++j) {
74 for (i = 0; i < j; ++i) sum += r(i, j) * wa1[i];
75 wa1[j] = (wa1[j] - sum) / r(j, j);
77 temp = wa1.blueNorm();
78 parl = fp / delta / temp / temp;
82 for (j = 0; j < n; ++j) wa1[j] = r.col(j).head(j + 1).dot(qtb.head(j + 1)) / diag[ipvt[j]];
84 gnorm = wa1.stableNorm();
86 if (paru == 0.) paru = dwarf / (std::min)(delta, Scalar(0.1));
90 par = (std::max)(par, parl);
91 par = (std::min)(par, paru);
92 if (par == 0.) par = gnorm / dxnorm;
99 if (par == 0.) par = (std::max)(dwarf, Scalar(.001) * paru);
100 wa1 =
sqrt(par) * diag;
102 Matrix<Scalar, Dynamic, 1> sdiag(n);
103 qrsolv<Scalar>(r, ipvt, wa1, qtb, x, sdiag);
105 wa2 = diag.cwiseProduct(x);
106 dxnorm = wa2.blueNorm();
113 if (
abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
break;
116 for (j = 0; j < n; ++j) {
118 wa1[j] = diag[l] * (wa2[l] / dxnorm);
120 for (j = 0; j < n; ++j) {
123 for (i = j + 1; i < n; ++i) wa1[i] -= r(i, j) * temp;
125 temp = wa1.blueNorm();
126 parc = fp / delta / temp / temp;
129 if (fp > 0.) parl = (std::max)(parl, par);
130 if (fp < 0.) paru = (std::min)(paru, par);
134 par = (std::max)(parl, par + parc);
140 if (iter == 0) par = 0.;
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)
151 typedef DenseIndex
Index;
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());
168 Matrix<Scalar, Dynamic, 1> wa1, wa2;
174 const Index rank = qr.rank();
176 wa1.tail(n - rank).setZero();
177 qr.matrixQR().topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(wa1.head(rank));
179 x = qr.colsPermutation() * wa1;
185 wa2 = diag.cwiseProduct(x);
186 dxnorm = wa2.blueNorm();
188 if (fp <= Scalar(0.1) * delta) {
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;
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)];
208 gnorm = wa1.stableNorm();
209 paru = gnorm / delta;
210 if (paru == 0.) paru = dwarf / (std::min)(delta, Scalar(0.1));
214 par = (std::max)(par, parl);
215 par = (std::min)(par, paru);
216 if (par == 0.) par = gnorm / dxnorm;
219 Matrix<Scalar, Dynamic, Dynamic> s = qr.matrixQR();
224 if (par == 0.) par = (std::max)(dwarf, Scalar(.001) * paru);
225 wa1 =
sqrt(par) * diag;
227 Matrix<Scalar, Dynamic, 1> sdiag(n);
228 qrsolv<Scalar>(s, qr.colsPermutation().indices(), wa1, qtb, x, sdiag);
230 wa2 = diag.cwiseProduct(x);
231 dxnorm = wa2.blueNorm();
238 if (
abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
break;
241 wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2 / dxnorm);
244 for (j = 0; j < n; ++j) {
247 for (
Index i = j + 1; i < n; ++i) wa1[i] -= s(i, j) * temp;
249 temp = wa1.blueNorm();
250 parc = fp / delta / temp / temp;
253 if (fp > 0.) parl = (std::max)(parl, par);
254 if (fp < 0.) paru = (std::min)(paru, par);
257 par = (std::max)(parl, par + parc);
259 if (iter == 0) par = 0.;
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