00001 namespace internal {
00002
00003 template <typename Scalar>
00004 void lmpar(
00005 Matrix< Scalar, Dynamic, Dynamic > &r,
00006 const VectorXi &ipvt,
00007 const Matrix< Scalar, Dynamic, 1 > &diag,
00008 const Matrix< Scalar, Dynamic, 1 > &qtb,
00009 Scalar delta,
00010 Scalar &par,
00011 Matrix< Scalar, Dynamic, 1 > &x)
00012 {
00013 typedef DenseIndex Index;
00014
00015
00016 Index i, j, l;
00017 Scalar fp;
00018 Scalar parc, parl;
00019 Index iter;
00020 Scalar temp, paru;
00021 Scalar gnorm;
00022 Scalar dxnorm;
00023
00024
00025
00026 const Scalar dwarf = std::numeric_limits<Scalar>::min();
00027 const Index n = r.cols();
00028 assert(n==diag.size());
00029 assert(n==qtb.size());
00030 assert(n==x.size());
00031
00032 Matrix< Scalar, Dynamic, 1 > wa1, wa2;
00033
00034
00035
00036 Index nsing = n-1;
00037 wa1 = qtb;
00038 for (j = 0; j < n; ++j) {
00039 if (r(j,j) == 0. && nsing == n-1)
00040 nsing = j - 1;
00041 if (nsing < n-1)
00042 wa1[j] = 0.;
00043 }
00044 for (j = nsing; j>=0; --j) {
00045 wa1[j] /= r(j,j);
00046 temp = wa1[j];
00047 for (i = 0; i < j ; ++i)
00048 wa1[i] -= r(i,j) * temp;
00049 }
00050
00051 for (j = 0; j < n; ++j)
00052 x[ipvt[j]] = wa1[j];
00053
00054
00055
00056
00057 iter = 0;
00058 wa2 = diag.cwiseProduct(x);
00059 dxnorm = wa2.blueNorm();
00060 fp = dxnorm - delta;
00061 if (fp <= Scalar(0.1) * delta) {
00062 par = 0;
00063 return;
00064 }
00065
00066
00067
00068
00069 parl = 0.;
00070 if (nsing >= n-1) {
00071 for (j = 0; j < n; ++j) {
00072 l = ipvt[j];
00073 wa1[j] = diag[l] * (wa2[l] / dxnorm);
00074 }
00075
00076
00077 for (j = 0; j < n; ++j) {
00078 Scalar sum = 0.;
00079 for (i = 0; i < j; ++i)
00080 sum += r(i,j) * wa1[i];
00081 wa1[j] = (wa1[j] - sum) / r(j,j);
00082 }
00083 temp = wa1.blueNorm();
00084 parl = fp / delta / temp / temp;
00085 }
00086
00087
00088 for (j = 0; j < n; ++j)
00089 wa1[j] = r.col(j).head(j+1).dot(qtb.head(j+1)) / diag[ipvt[j]];
00090
00091 gnorm = wa1.stableNorm();
00092 paru = gnorm / delta;
00093 if (paru == 0.)
00094 paru = dwarf / std::min(delta,Scalar(0.1));
00095
00096
00097
00098 par = std::max(par,parl);
00099 par = std::min(par,paru);
00100 if (par == 0.)
00101 par = gnorm / dxnorm;
00102
00103
00104 while (true) {
00105 ++iter;
00106
00107
00108 if (par == 0.)
00109 par = std::max(dwarf,Scalar(.001) * paru);
00110 wa1 = sqrt(par)* diag;
00111
00112 Matrix< Scalar, Dynamic, 1 > sdiag(n);
00113 qrsolv<Scalar>(r, ipvt, wa1, qtb, x, sdiag);
00114
00115 wa2 = diag.cwiseProduct(x);
00116 dxnorm = wa2.blueNorm();
00117 temp = fp;
00118 fp = dxnorm - delta;
00119
00120
00121
00122
00123 if (abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
00124 break;
00125
00126
00127 for (j = 0; j < n; ++j) {
00128 l = ipvt[j];
00129 wa1[j] = diag[l] * (wa2[l] / dxnorm);
00130 }
00131 for (j = 0; j < n; ++j) {
00132 wa1[j] /= sdiag[j];
00133 temp = wa1[j];
00134 for (i = j+1; i < n; ++i)
00135 wa1[i] -= r(i,j) * temp;
00136 }
00137 temp = wa1.blueNorm();
00138 parc = fp / delta / temp / temp;
00139
00140
00141 if (fp > 0.)
00142 parl = std::max(parl,par);
00143 if (fp < 0.)
00144 paru = std::min(paru,par);
00145
00146
00147
00148 par = std::max(parl,par+parc);
00149
00150
00151 }
00152
00153
00154 if (iter == 0)
00155 par = 0.;
00156 return;
00157 }
00158
00159 template <typename Scalar>
00160 void lmpar2(
00161 const ColPivHouseholderQR<Matrix< Scalar, Dynamic, Dynamic> > &qr,
00162 const Matrix< Scalar, Dynamic, 1 > &diag,
00163 const Matrix< Scalar, Dynamic, 1 > &qtb,
00164 Scalar delta,
00165 Scalar &par,
00166 Matrix< Scalar, Dynamic, 1 > &x)
00167
00168 {
00169 typedef DenseIndex Index;
00170
00171
00172 Index j;
00173 Scalar fp;
00174 Scalar parc, parl;
00175 Index iter;
00176 Scalar temp, paru;
00177 Scalar gnorm;
00178 Scalar dxnorm;
00179
00180
00181
00182 const Scalar dwarf = std::numeric_limits<Scalar>::min();
00183 const Index n = qr.matrixQR().cols();
00184 assert(n==diag.size());
00185 assert(n==qtb.size());
00186
00187 Matrix< Scalar, Dynamic, 1 > wa1, wa2;
00188
00189
00190
00191
00192
00193 const Index rank = qr.rank();
00194 wa1 = qtb;
00195 wa1.tail(n-rank).setZero();
00196 qr.matrixQR().topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(wa1.head(rank));
00197
00198 x = qr.colsPermutation()*wa1;
00199
00200
00201
00202
00203 iter = 0;
00204 wa2 = diag.cwiseProduct(x);
00205 dxnorm = wa2.blueNorm();
00206 fp = dxnorm - delta;
00207 if (fp <= Scalar(0.1) * delta) {
00208 par = 0;
00209 return;
00210 }
00211
00212
00213
00214
00215 parl = 0.;
00216 if (rank==n) {
00217 wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2)/dxnorm;
00218 qr.matrixQR().topLeftCorner(n, n).transpose().template triangularView<Lower>().solveInPlace(wa1);
00219 temp = wa1.blueNorm();
00220 parl = fp / delta / temp / temp;
00221 }
00222
00223
00224 for (j = 0; j < n; ++j)
00225 wa1[j] = qr.matrixQR().col(j).head(j+1).dot(qtb.head(j+1)) / diag[qr.colsPermutation().indices()(j)];
00226
00227 gnorm = wa1.stableNorm();
00228 paru = gnorm / delta;
00229 if (paru == 0.)
00230 paru = dwarf / std::min(delta,Scalar(0.1));
00231
00232
00233
00234 par = std::max(par,parl);
00235 par = std::min(par,paru);
00236 if (par == 0.)
00237 par = gnorm / dxnorm;
00238
00239
00240 Matrix< Scalar, Dynamic, Dynamic > s = qr.matrixQR();
00241 while (true) {
00242 ++iter;
00243
00244
00245 if (par == 0.)
00246 par = std::max(dwarf,Scalar(.001) * paru);
00247 wa1 = sqrt(par)* diag;
00248
00249 Matrix< Scalar, Dynamic, 1 > sdiag(n);
00250 qrsolv<Scalar>(s, qr.colsPermutation().indices(), wa1, qtb, x, sdiag);
00251
00252 wa2 = diag.cwiseProduct(x);
00253 dxnorm = wa2.blueNorm();
00254 temp = fp;
00255 fp = dxnorm - delta;
00256
00257
00258
00259
00260 if (abs(fp) <= Scalar(0.1) * delta || (parl == 0. && fp <= temp && temp < 0.) || iter == 10)
00261 break;
00262
00263
00264 wa1 = qr.colsPermutation().inverse() * diag.cwiseProduct(wa2/dxnorm);
00265
00266
00267 for (j = 0; j < n; ++j) {
00268 wa1[j] /= sdiag[j];
00269 temp = wa1[j];
00270 for (Index i = j+1; i < n; ++i)
00271 wa1[i] -= s(i,j) * temp;
00272 }
00273 temp = wa1.blueNorm();
00274 parc = fp / delta / temp / temp;
00275
00276
00277 if (fp > 0.)
00278 parl = std::max(parl,par);
00279 if (fp < 0.)
00280 paru = std::min(paru,par);
00281
00282
00283 par = std::max(parl,par+parc);
00284 }
00285 if (iter == 0)
00286 par = 0.;
00287 return;
00288 }
00289
00290 }