5 template <
typename Scalar>
30 const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
31 const Index n = r.
cols();
42 for (j = 0; j < n; ++j) {
43 if (r(j,j) == 0. && nsing == n-1)
48 for (j = nsing; j>=0; --j) {
51 for (i = 0; i < j ; ++i)
52 wa1[i] -= r(i,j) * temp;
55 for (j = 0; j < n; ++j)
62 wa2 = diag.cwiseProduct(x);
63 dxnorm = wa2.blueNorm();
65 if (fp <= Scalar(0.1) * delta) {
75 for (j = 0; j < n; ++j) {
77 wa1[j] = diag[l] * (wa2[l] / dxnorm);
81 for (j = 0; j < n; ++j) {
83 for (i = 0; i < j; ++i)
84 sum += r(i,j) * wa1[i];
85 wa1[j] = (wa1[j] - sum) / r(j,j);
87 temp = wa1.blueNorm();
88 parl = fp / delta / temp / temp;
92 for (j = 0; j < n; ++j)
93 wa1[j] = r.col(j).head(j+1).dot(qtb.head(j+1)) / diag[ipvt[j]];
95 gnorm = wa1.stableNorm();
98 paru = dwarf / (std::min)(delta,Scalar(0.1));
102 par = (std::max)(par,parl);
103 par = (std::min)(par,paru);
105 par = gnorm / dxnorm;
113 par = (std::max)(dwarf,Scalar(.001) * paru);
114 wa1 =
sqrt(par)* diag;
117 qrsolv<Scalar>(r, ipvt, wa1, qtb, x, sdiag);
119 wa2 = diag.cwiseProduct(x);
120 dxnorm = wa2.blueNorm();
127 if (
abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
131 for (j = 0; j < n; ++j) {
133 wa1[j] = diag[l] * (wa2[l] / dxnorm);
135 for (j = 0; j < n; ++j) {
138 for (i = j+1; i < n; ++i)
139 wa1[i] -= r(i,j) * temp;
141 temp = wa1.blueNorm();
142 parc = fp / delta / temp / temp;
146 parl = (std::max)(parl,par);
148 paru = (std::min)(paru,par);
152 par = (std::max)(parl,par+parc);
163 template <
typename Scalar>
188 const Scalar dwarf = (std::numeric_limits<Scalar>::min)();
189 const Index n = qr.matrixQR().cols();
199 const Index rank = qr.rank();
202 qr.matrixQR().topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(wa1.head(rank));
204 x = qr.colsPermutation()*wa1;
210 wa2 = diag.cwiseProduct(x);
211 dxnorm = wa2.blueNorm();
213 if (fp <= Scalar(0.1) * delta) {
223 wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2)/dxnorm;
224 qr.matrixQR().topLeftCorner(n, n).transpose().template triangularView<Lower>().solveInPlace(wa1);
225 temp = wa1.blueNorm();
226 parl = fp / delta / temp / temp;
230 for (j = 0; j < n; ++j)
231 wa1[j] = qr.matrixQR().col(j).head(j+1).dot(qtb.head(j+1)) / diag[qr.colsPermutation().indices()(j)];
233 gnorm = wa1.stableNorm();
234 paru = gnorm / delta;
236 paru = dwarf / (std::min)(delta,Scalar(0.1));
240 par = (std::max)(par,parl);
241 par = (std::min)(par,paru);
243 par = gnorm / dxnorm;
252 par = (std::max)(dwarf,Scalar(.001) * paru);
253 wa1 =
sqrt(par)* diag;
256 qrsolv<Scalar>(s, qr.colsPermutation().indices(), wa1, qtb, x, sdiag);
258 wa2 = diag.cwiseProduct(x);
259 dxnorm = wa2.blueNorm();
266 if (
abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
270 wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2/dxnorm);
273 for (j = 0; j < n; ++j) {
276 for (Index i = j+1; i < n; ++i)
277 wa1[i] -= s(i,j) * temp;
279 temp = wa1.blueNorm();
280 parc = fp / delta / temp / temp;
284 parl = (std::max)(parl,par);
286 paru = (std::min)(paru,par);
289 par = (std::max)(parl,par+parc);
IntermediateState sqrt(const Expression &arg)
iterative scaling algorithm to equilibrate rows and column norms in matrices
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Derived & setZero(Index size)
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
void lmpar2(const ColPivHouseholderQR< Matrix< Scalar, Dynamic, Dynamic > > &qr, const Matrix< Scalar, Dynamic, 1 > &diag, const Matrix< Scalar, Dynamic, 1 > &qtb, Scalar delta, Scalar &par, Matrix< Scalar, Dynamic, 1 > &x)
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)
EIGEN_STRONG_INLINE Index cols() const