00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include <stdio.h>
00013
00014 #include "main.h"
00015 #include <unsupported/Eigen/LevenbergMarquardt>
00016
00017
00018
00019 #include <Eigen/src/Core/util/DisableStupidWarnings.h>
00020
00021 using std::sqrt;
00022
00023 struct lmder_functor : DenseFunctor<double>
00024 {
00025 lmder_functor(void): DenseFunctor<double>(3,15) {}
00026 int operator()(const VectorXd &x, VectorXd &fvec) const
00027 {
00028 double tmp1, tmp2, tmp3;
00029 static const double y[15] = {1.4e-1, 1.8e-1, 2.2e-1, 2.5e-1, 2.9e-1, 3.2e-1, 3.5e-1,
00030 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
00031
00032 for (int i = 0; i < values(); i++)
00033 {
00034 tmp1 = i+1;
00035 tmp2 = 16 - i - 1;
00036 tmp3 = (i>=8)? tmp2 : tmp1;
00037 fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
00038 }
00039 return 0;
00040 }
00041
00042 int df(const VectorXd &x, MatrixXd &fjac) const
00043 {
00044 double tmp1, tmp2, tmp3, tmp4;
00045 for (int i = 0; i < values(); i++)
00046 {
00047 tmp1 = i+1;
00048 tmp2 = 16 - i - 1;
00049 tmp3 = (i>=8)? tmp2 : tmp1;
00050 tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
00051 fjac(i,0) = -1;
00052 fjac(i,1) = tmp1*tmp2/tmp4;
00053 fjac(i,2) = tmp1*tmp3/tmp4;
00054 }
00055 return 0;
00056 }
00057 };
00058
00059 void testLmder1()
00060 {
00061 int n=3, info;
00062
00063 VectorXd x;
00064
00065
00066 x.setConstant(n, 1.);
00067
00068
00069 lmder_functor functor;
00070 LevenbergMarquardt<lmder_functor> lm(functor);
00071 info = lm.lmder1(x);
00072
00073
00074 VERIFY_IS_EQUAL(info, 1);
00075 VERIFY_IS_EQUAL(lm.nfev(), 6);
00076 VERIFY_IS_EQUAL(lm.njev(), 5);
00077
00078
00079 VERIFY_IS_APPROX(lm.fvec().blueNorm(), 0.09063596);
00080
00081
00082 VectorXd x_ref(n);
00083 x_ref << 0.08241058, 1.133037, 2.343695;
00084 VERIFY_IS_APPROX(x, x_ref);
00085 }
00086
00087 void testLmder()
00088 {
00089 const int m=15, n=3;
00090 int info;
00091 double fnorm, covfac;
00092 VectorXd x;
00093
00094
00095 x.setConstant(n, 1.);
00096
00097
00098 lmder_functor functor;
00099 LevenbergMarquardt<lmder_functor> lm(functor);
00100 info = lm.minimize(x);
00101
00102
00103 VERIFY_IS_EQUAL(info, 1);
00104 VERIFY_IS_EQUAL(lm.nfev(), 6);
00105 VERIFY_IS_EQUAL(lm.njev(), 5);
00106
00107
00108 fnorm = lm.fvec().blueNorm();
00109 VERIFY_IS_APPROX(fnorm, 0.09063596);
00110
00111
00112 VectorXd x_ref(n);
00113 x_ref << 0.08241058, 1.133037, 2.343695;
00114 VERIFY_IS_APPROX(x, x_ref);
00115
00116
00117 covfac = fnorm*fnorm/(m-n);
00118 internal::covar(lm.matrixR(), lm.permutation().indices());
00119
00120 MatrixXd cov_ref(n,n);
00121 cov_ref <<
00122 0.0001531202, 0.002869941, -0.002656662,
00123 0.002869941, 0.09480935, -0.09098995,
00124 -0.002656662, -0.09098995, 0.08778727;
00125
00126
00127
00128 MatrixXd cov;
00129 cov = covfac*lm.matrixR().topLeftCorner<n,n>();
00130 VERIFY_IS_APPROX( cov, cov_ref);
00131
00132
00133 }
00134
00135 struct lmdif_functor : DenseFunctor<double>
00136 {
00137 lmdif_functor(void) : DenseFunctor<double>(3,15) {}
00138 int operator()(const VectorXd &x, VectorXd &fvec) const
00139 {
00140 int i;
00141 double tmp1,tmp2,tmp3;
00142 static const double y[15]={1.4e-1,1.8e-1,2.2e-1,2.5e-1,2.9e-1,3.2e-1,3.5e-1,3.9e-1,
00143 3.7e-1,5.8e-1,7.3e-1,9.6e-1,1.34e0,2.1e0,4.39e0};
00144
00145 assert(x.size()==3);
00146 assert(fvec.size()==15);
00147 for (i=0; i<15; i++)
00148 {
00149 tmp1 = i+1;
00150 tmp2 = 15 - i;
00151 tmp3 = tmp1;
00152
00153 if (i >= 8) tmp3 = tmp2;
00154 fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
00155 }
00156 return 0;
00157 }
00158 };
00159
00160 void testLmdif1()
00161 {
00162 const int n=3;
00163 int info;
00164
00165 VectorXd x(n), fvec(15);
00166
00167
00168 x.setConstant(n, 1.);
00169
00170
00171 lmdif_functor functor;
00172 DenseIndex nfev;
00173 info = LevenbergMarquardt<lmdif_functor>::lmdif1(functor, x, &nfev);
00174
00175
00176 VERIFY_IS_EQUAL(info, 1);
00177
00178
00179
00180 functor(x, fvec);
00181 VERIFY_IS_APPROX(fvec.blueNorm(), 0.09063596);
00182
00183
00184 VectorXd x_ref(n);
00185 x_ref << 0.0824106, 1.1330366, 2.3436947;
00186 VERIFY_IS_APPROX(x, x_ref);
00187
00188 }
00189
00190 void testLmdif()
00191 {
00192 const int m=15, n=3;
00193 int info;
00194 double fnorm, covfac;
00195 VectorXd x(n);
00196
00197
00198 x.setConstant(n, 1.);
00199
00200
00201 lmdif_functor functor;
00202 NumericalDiff<lmdif_functor> numDiff(functor);
00203 LevenbergMarquardt<NumericalDiff<lmdif_functor> > lm(numDiff);
00204 info = lm.minimize(x);
00205
00206
00207 VERIFY_IS_EQUAL(info, 1);
00208
00209
00210
00211 fnorm = lm.fvec().blueNorm();
00212 VERIFY_IS_APPROX(fnorm, 0.09063596);
00213
00214
00215 VectorXd x_ref(n);
00216 x_ref << 0.08241058, 1.133037, 2.343695;
00217 VERIFY_IS_APPROX(x, x_ref);
00218
00219
00220 covfac = fnorm*fnorm/(m-n);
00221 internal::covar(lm.matrixR(), lm.permutation().indices());
00222
00223 MatrixXd cov_ref(n,n);
00224 cov_ref <<
00225 0.0001531202, 0.002869942, -0.002656662,
00226 0.002869942, 0.09480937, -0.09098997,
00227 -0.002656662, -0.09098997, 0.08778729;
00228
00229
00230
00231 MatrixXd cov;
00232 cov = covfac*lm.matrixR().topLeftCorner<n,n>();
00233 VERIFY_IS_APPROX( cov, cov_ref);
00234
00235
00236 }
00237
00238 struct chwirut2_functor : DenseFunctor<double>
00239 {
00240 chwirut2_functor(void) : DenseFunctor<double>(3,54) {}
00241 static const double m_x[54];
00242 static const double m_y[54];
00243 int operator()(const VectorXd &b, VectorXd &fvec)
00244 {
00245 int i;
00246
00247 assert(b.size()==3);
00248 assert(fvec.size()==54);
00249 for(i=0; i<54; i++) {
00250 double x = m_x[i];
00251 fvec[i] = exp(-b[0]*x)/(b[1]+b[2]*x) - m_y[i];
00252 }
00253 return 0;
00254 }
00255 int df(const VectorXd &b, MatrixXd &fjac)
00256 {
00257 assert(b.size()==3);
00258 assert(fjac.rows()==54);
00259 assert(fjac.cols()==3);
00260 for(int i=0; i<54; i++) {
00261 double x = m_x[i];
00262 double factor = 1./(b[1]+b[2]*x);
00263 double e = exp(-b[0]*x);
00264 fjac(i,0) = -x*e*factor;
00265 fjac(i,1) = -e*factor*factor;
00266 fjac(i,2) = -x*e*factor*factor;
00267 }
00268 return 0;
00269 }
00270 };
00271 const double chwirut2_functor::m_x[54] = { 0.500E0, 1.000E0, 1.750E0, 3.750E0, 5.750E0, 0.875E0, 2.250E0, 3.250E0, 5.250E0, 0.750E0, 1.750E0, 2.750E0, 4.750E0, 0.625E0, 1.250E0, 2.250E0, 4.250E0, .500E0, 3.000E0, .750E0, 3.000E0, 1.500E0, 6.000E0, 3.000E0, 6.000E0, 1.500E0, 3.000E0, .500E0, 2.000E0, 4.000E0, .750E0, 2.000E0, 5.000E0, .750E0, 2.250E0, 3.750E0, 5.750E0, 3.000E0, .750E0, 2.500E0, 4.000E0, .750E0, 2.500E0, 4.000E0, .750E0, 2.500E0, 4.000E0, .500E0, 6.000E0, 3.000E0, .500E0, 2.750E0, .500E0, 1.750E0};
00272 const double chwirut2_functor::m_y[54] = { 92.9000E0 ,57.1000E0 ,31.0500E0 ,11.5875E0 ,8.0250E0 ,63.6000E0 ,21.4000E0 ,14.2500E0 ,8.4750E0 ,63.8000E0 ,26.8000E0 ,16.4625E0 ,7.1250E0 ,67.3000E0 ,41.0000E0 ,21.1500E0 ,8.1750E0 ,81.5000E0 ,13.1200E0 ,59.9000E0 ,14.6200E0 ,32.9000E0 ,5.4400E0 ,12.5600E0 ,5.4400E0 ,32.0000E0 ,13.9500E0 ,75.8000E0 ,20.0000E0 ,10.4200E0 ,59.5000E0 ,21.6700E0 ,8.5500E0 ,62.0000E0 ,20.2000E0 ,7.7600E0 ,3.7500E0 ,11.8100E0 ,54.7000E0 ,23.7000E0 ,11.5500E0 ,61.3000E0 ,17.7000E0 ,8.7400E0 ,59.2000E0 ,16.3000E0 ,8.6200E0 ,81.0000E0 ,4.8700E0 ,14.6200E0 ,81.7000E0 ,17.1700E0 ,81.3000E0 ,28.9000E0 };
00273
00274
00275 void testNistChwirut2(void)
00276 {
00277 const int n=3;
00278 int info;
00279
00280 VectorXd x(n);
00281
00282
00283
00284
00285 x<< 0.1, 0.01, 0.02;
00286
00287 chwirut2_functor functor;
00288 LevenbergMarquardt<chwirut2_functor> lm(functor);
00289 info = lm.minimize(x);
00290
00291
00292 VERIFY_IS_EQUAL(info, 1);
00293
00294 VERIFY_IS_EQUAL(lm.njev(), 8);
00295
00296 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.1304802941E+02);
00297
00298 VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
00299 VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
00300 VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
00301
00302
00303
00304
00305 x<< 0.15, 0.008, 0.010;
00306
00307 lm.resetParameters();
00308 lm.setFtol(1.E6*NumTraits<double>::epsilon());
00309 lm.setXtol(1.E6*NumTraits<double>::epsilon());
00310 info = lm.minimize(x);
00311
00312
00313 VERIFY_IS_EQUAL(info, 1);
00314
00315 VERIFY_IS_EQUAL(lm.njev(), 6);
00316
00317 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.1304802941E+02);
00318
00319 VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
00320 VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
00321 VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
00322 }
00323
00324
00325 struct misra1a_functor : DenseFunctor<double>
00326 {
00327 misra1a_functor(void) : DenseFunctor<double>(2,14) {}
00328 static const double m_x[14];
00329 static const double m_y[14];
00330 int operator()(const VectorXd &b, VectorXd &fvec)
00331 {
00332 assert(b.size()==2);
00333 assert(fvec.size()==14);
00334 for(int i=0; i<14; i++) {
00335 fvec[i] = b[0]*(1.-exp(-b[1]*m_x[i])) - m_y[i] ;
00336 }
00337 return 0;
00338 }
00339 int df(const VectorXd &b, MatrixXd &fjac)
00340 {
00341 assert(b.size()==2);
00342 assert(fjac.rows()==14);
00343 assert(fjac.cols()==2);
00344 for(int i=0; i<14; i++) {
00345 fjac(i,0) = (1.-exp(-b[1]*m_x[i]));
00346 fjac(i,1) = (b[0]*m_x[i]*exp(-b[1]*m_x[i]));
00347 }
00348 return 0;
00349 }
00350 };
00351 const double misra1a_functor::m_x[14] = { 77.6E0, 114.9E0, 141.1E0, 190.8E0, 239.9E0, 289.0E0, 332.8E0, 378.4E0, 434.8E0, 477.3E0, 536.8E0, 593.1E0, 689.1E0, 760.0E0};
00352 const double misra1a_functor::m_y[14] = { 10.07E0, 14.73E0, 17.94E0, 23.93E0, 29.61E0, 35.18E0, 40.02E0, 44.82E0, 50.76E0, 55.05E0, 61.01E0, 66.40E0, 75.47E0, 81.78E0};
00353
00354
00355 void testNistMisra1a(void)
00356 {
00357 const int n=2;
00358 int info;
00359
00360 VectorXd x(n);
00361
00362
00363
00364
00365 x<< 500., 0.0001;
00366
00367 misra1a_functor functor;
00368 LevenbergMarquardt<misra1a_functor> lm(functor);
00369 info = lm.minimize(x);
00370
00371
00372 VERIFY_IS_EQUAL(info, 1);
00373 VERIFY_IS_EQUAL(lm.nfev(), 19);
00374 VERIFY_IS_EQUAL(lm.njev(), 15);
00375
00376 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.2455138894E-01);
00377
00378 VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
00379 VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
00380
00381
00382
00383
00384 x<< 250., 0.0005;
00385
00386 info = lm.minimize(x);
00387
00388
00389 VERIFY_IS_EQUAL(info, 1);
00390 VERIFY_IS_EQUAL(lm.nfev(), 5);
00391 VERIFY_IS_EQUAL(lm.njev(), 4);
00392
00393 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.2455138894E-01);
00394
00395 VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
00396 VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
00397 }
00398
00399 struct hahn1_functor : DenseFunctor<double>
00400 {
00401 hahn1_functor(void) : DenseFunctor<double>(7,236) {}
00402 static const double m_x[236];
00403 int operator()(const VectorXd &b, VectorXd &fvec)
00404 {
00405 static const double m_y[236] = { .591E0 , 1.547E0 , 2.902E0 , 2.894E0 , 4.703E0 , 6.307E0 , 7.03E0 , 7.898E0 , 9.470E0 , 9.484E0 , 10.072E0 , 10.163E0 , 11.615E0 , 12.005E0 , 12.478E0 , 12.982E0 , 12.970E0 , 13.926E0 , 14.452E0 , 14.404E0 , 15.190E0 , 15.550E0 , 15.528E0 , 15.499E0 , 16.131E0 , 16.438E0 , 16.387E0 , 16.549E0 , 16.872E0 , 16.830E0 , 16.926E0 , 16.907E0 , 16.966E0 , 17.060E0 , 17.122E0 , 17.311E0 , 17.355E0 , 17.668E0 , 17.767E0 , 17.803E0 , 17.765E0 , 17.768E0 , 17.736E0 , 17.858E0 , 17.877E0 , 17.912E0 , 18.046E0 , 18.085E0 , 18.291E0 , 18.357E0 , 18.426E0 , 18.584E0 , 18.610E0 , 18.870E0 , 18.795E0 , 19.111E0 , .367E0 , .796E0 , 0.892E0 , 1.903E0 , 2.150E0 , 3.697E0 , 5.870E0 , 6.421E0 , 7.422E0 , 9.944E0 , 11.023E0 , 11.87E0 , 12.786E0 , 14.067E0 , 13.974E0 , 14.462E0 , 14.464E0 , 15.381E0 , 15.483E0 , 15.59E0 , 16.075E0 , 16.347E0 , 16.181E0 , 16.915E0 , 17.003E0 , 16.978E0 , 17.756E0 , 17.808E0 , 17.868E0 , 18.481E0 , 18.486E0 , 19.090E0 , 16.062E0 , 16.337E0 , 16.345E0 ,
00406 16.388E0 , 17.159E0 , 17.116E0 , 17.164E0 , 17.123E0 , 17.979E0 , 17.974E0 , 18.007E0 , 17.993E0 , 18.523E0 , 18.669E0 , 18.617E0 , 19.371E0 , 19.330E0 , 0.080E0 , 0.248E0 , 1.089E0 , 1.418E0 , 2.278E0 , 3.624E0 , 4.574E0 , 5.556E0 , 7.267E0 , 7.695E0 , 9.136E0 , 9.959E0 , 9.957E0 , 11.600E0 , 13.138E0 , 13.564E0 , 13.871E0 , 13.994E0 , 14.947E0 , 15.473E0 , 15.379E0 , 15.455E0 , 15.908E0 , 16.114E0 , 17.071E0 , 17.135E0 , 17.282E0 , 17.368E0 , 17.483E0 , 17.764E0 , 18.185E0 , 18.271E0 , 18.236E0 , 18.237E0 , 18.523E0 , 18.627E0 , 18.665E0 , 19.086E0 , 0.214E0 , 0.943E0 , 1.429E0 , 2.241E0 , 2.951E0 , 3.782E0 , 4.757E0 , 5.602E0 , 7.169E0 , 8.920E0 , 10.055E0 , 12.035E0 , 12.861E0 , 13.436E0 , 14.167E0 , 14.755E0 , 15.168E0 , 15.651E0 , 15.746E0 , 16.216E0 , 16.445E0 , 16.965E0 , 17.121E0 , 17.206E0 , 17.250E0 , 17.339E0 , 17.793E0 , 18.123E0 , 18.49E0 , 18.566E0 , 18.645E0 , 18.706E0 , 18.924E0 , 19.1E0 , 0.375E0 , 0.471E0 , 1.504E0 , 2.204E0 , 2.813E0 , 4.765E0 , 9.835E0 , 10.040E0 , 11.946E0 ,
00407 12.596E0 ,
00408 13.303E0 , 13.922E0 , 14.440E0 , 14.951E0 , 15.627E0 , 15.639E0 , 15.814E0 , 16.315E0 , 16.334E0 , 16.430E0 , 16.423E0 , 17.024E0 , 17.009E0 , 17.165E0 , 17.134E0 , 17.349E0 , 17.576E0 , 17.848E0 , 18.090E0 , 18.276E0 , 18.404E0 , 18.519E0 , 19.133E0 , 19.074E0 , 19.239E0 , 19.280E0 , 19.101E0 , 19.398E0 , 19.252E0 , 19.89E0 , 20.007E0 , 19.929E0 , 19.268E0 , 19.324E0 , 20.049E0 , 20.107E0 , 20.062E0 , 20.065E0 , 19.286E0 , 19.972E0 , 20.088E0 , 20.743E0 , 20.83E0 , 20.935E0 , 21.035E0 , 20.93E0 , 21.074E0 , 21.085E0 , 20.935E0 };
00409
00410
00411
00412 assert(b.size()==7);
00413 assert(fvec.size()==236);
00414 for(int i=0; i<236; i++) {
00415 double x=m_x[i], xx=x*x, xxx=xx*x;
00416 fvec[i] = (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) / (1.+b[4]*x+b[5]*xx+b[6]*xxx) - m_y[i];
00417 }
00418 return 0;
00419 }
00420
00421 int df(const VectorXd &b, MatrixXd &fjac)
00422 {
00423 assert(b.size()==7);
00424 assert(fjac.rows()==236);
00425 assert(fjac.cols()==7);
00426 for(int i=0; i<236; i++) {
00427 double x=m_x[i], xx=x*x, xxx=xx*x;
00428 double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
00429 fjac(i,0) = 1.*fact;
00430 fjac(i,1) = x*fact;
00431 fjac(i,2) = xx*fact;
00432 fjac(i,3) = xxx*fact;
00433 fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
00434 fjac(i,4) = x*fact;
00435 fjac(i,5) = xx*fact;
00436 fjac(i,6) = xxx*fact;
00437 }
00438 return 0;
00439 }
00440 };
00441 const double hahn1_functor::m_x[236] = { 24.41E0 , 34.82E0 , 44.09E0 , 45.07E0 , 54.98E0 , 65.51E0 , 70.53E0 , 75.70E0 , 89.57E0 , 91.14E0 , 96.40E0 , 97.19E0 , 114.26E0 , 120.25E0 , 127.08E0 , 133.55E0 , 133.61E0 , 158.67E0 , 172.74E0 , 171.31E0 , 202.14E0 , 220.55E0 , 221.05E0 , 221.39E0 , 250.99E0 , 268.99E0 , 271.80E0 , 271.97E0 , 321.31E0 , 321.69E0 , 330.14E0 , 333.03E0 , 333.47E0 , 340.77E0 , 345.65E0 , 373.11E0 , 373.79E0 , 411.82E0 , 419.51E0 , 421.59E0 , 422.02E0 , 422.47E0 , 422.61E0 , 441.75E0 , 447.41E0 , 448.7E0 , 472.89E0 , 476.69E0 , 522.47E0 , 522.62E0 , 524.43E0 , 546.75E0 , 549.53E0 , 575.29E0 , 576.00E0 , 625.55E0 , 20.15E0 , 28.78E0 , 29.57E0 , 37.41E0 , 39.12E0 , 50.24E0 , 61.38E0 , 66.25E0 , 73.42E0 , 95.52E0 , 107.32E0 , 122.04E0 , 134.03E0 , 163.19E0 , 163.48E0 , 175.70E0 , 179.86E0 , 211.27E0 , 217.78E0 , 219.14E0 , 262.52E0 , 268.01E0 , 268.62E0 , 336.25E0 , 337.23E0 , 339.33E0 , 427.38E0 , 428.58E0 , 432.68E0 , 528.99E0 , 531.08E0 , 628.34E0 , 253.24E0 , 273.13E0 , 273.66E0 ,
00442 282.10E0 , 346.62E0 , 347.19E0 , 348.78E0 , 351.18E0 , 450.10E0 , 450.35E0 , 451.92E0 , 455.56E0 , 552.22E0 , 553.56E0 , 555.74E0 , 652.59E0 , 656.20E0 , 14.13E0 , 20.41E0 , 31.30E0 , 33.84E0 , 39.70E0 , 48.83E0 , 54.50E0 , 60.41E0 , 72.77E0 , 75.25E0 , 86.84E0 , 94.88E0 , 96.40E0 , 117.37E0 , 139.08E0 , 147.73E0 , 158.63E0 , 161.84E0 , 192.11E0 , 206.76E0 , 209.07E0 , 213.32E0 , 226.44E0 , 237.12E0 , 330.90E0 , 358.72E0 , 370.77E0 , 372.72E0 , 396.24E0 , 416.59E0 , 484.02E0 , 495.47E0 , 514.78E0 , 515.65E0 , 519.47E0 , 544.47E0 , 560.11E0 , 620.77E0 , 18.97E0 , 28.93E0 , 33.91E0 , 40.03E0 , 44.66E0 , 49.87E0 , 55.16E0 , 60.90E0 , 72.08E0 , 85.15E0 , 97.06E0 , 119.63E0 , 133.27E0 , 143.84E0 , 161.91E0 , 180.67E0 , 198.44E0 , 226.86E0 , 229.65E0 , 258.27E0 , 273.77E0 , 339.15E0 , 350.13E0 , 362.75E0 , 371.03E0 , 393.32E0 , 448.53E0 , 473.78E0 , 511.12E0 , 524.70E0 , 548.75E0 , 551.64E0 , 574.02E0 , 623.86E0 , 21.46E0 , 24.33E0 , 33.43E0 , 39.22E0 , 44.18E0 , 55.02E0 , 94.33E0 , 96.44E0 , 118.82E0 , 128.48E0 ,
00443 141.94E0 , 156.92E0 , 171.65E0 , 190.00E0 , 223.26E0 , 223.88E0 , 231.50E0 , 265.05E0 , 269.44E0 , 271.78E0 , 273.46E0 , 334.61E0 , 339.79E0 , 349.52E0 , 358.18E0 , 377.98E0 , 394.77E0 , 429.66E0 , 468.22E0 , 487.27E0 , 519.54E0 , 523.03E0 , 612.99E0 , 638.59E0 , 641.36E0 , 622.05E0 , 631.50E0 , 663.97E0 , 646.9E0 , 748.29E0 , 749.21E0 , 750.14E0 , 647.04E0 , 646.89E0 , 746.9E0 , 748.43E0 , 747.35E0 , 749.27E0 , 647.61E0 , 747.78E0 , 750.51E0 , 851.37E0 , 845.97E0 , 847.54E0 , 849.93E0 , 851.61E0 , 849.75E0 , 850.98E0 , 848.23E0};
00444
00445
00446 void testNistHahn1(void)
00447 {
00448 const int n=7;
00449 int info;
00450
00451 VectorXd x(n);
00452
00453
00454
00455
00456 x<< 10., -1., .05, -.00001, -.05, .001, -.000001;
00457
00458 hahn1_functor functor;
00459 LevenbergMarquardt<hahn1_functor> lm(functor);
00460 info = lm.minimize(x);
00461
00462
00463 VERIFY_IS_EQUAL(info, 1);
00464 VERIFY_IS_EQUAL(lm.nfev(), 11);
00465 VERIFY_IS_EQUAL(lm.njev(), 10);
00466
00467 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.5324382854E+00);
00468
00469 VERIFY_IS_APPROX(x[0], 1.0776351733E+00);
00470 VERIFY_IS_APPROX(x[1],-1.2269296921E-01);
00471 VERIFY_IS_APPROX(x[2], 4.0863750610E-03);
00472 VERIFY_IS_APPROX(x[3],-1.426264e-06);
00473 VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
00474 VERIFY_IS_APPROX(x[5], 2.4053735503E-04);
00475 VERIFY_IS_APPROX(x[6],-1.2314450199E-07);
00476
00477
00478
00479
00480 x<< .1, -.1, .005, -.000001, -.005, .0001, -.0000001;
00481
00482 info = lm.minimize(x);
00483
00484
00485 VERIFY_IS_EQUAL(info, 1);
00486
00487 VERIFY_IS_EQUAL(lm.njev(), 10);
00488
00489 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.5324382854E+00);
00490
00491 VERIFY_IS_APPROX(x[0], 1.077640);
00492 VERIFY_IS_APPROX(x[1], -0.1226933);
00493 VERIFY_IS_APPROX(x[2], 0.004086383);
00494 VERIFY_IS_APPROX(x[3], -1.426277e-06);
00495 VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
00496 VERIFY_IS_APPROX(x[5], 0.00024053772);
00497 VERIFY_IS_APPROX(x[6], -1.231450e-07);
00498
00499 }
00500
00501 struct misra1d_functor : DenseFunctor<double>
00502 {
00503 misra1d_functor(void) : DenseFunctor<double>(2,14) {}
00504 static const double x[14];
00505 static const double y[14];
00506 int operator()(const VectorXd &b, VectorXd &fvec)
00507 {
00508 assert(b.size()==2);
00509 assert(fvec.size()==14);
00510 for(int i=0; i<14; i++) {
00511 fvec[i] = b[0]*b[1]*x[i]/(1.+b[1]*x[i]) - y[i];
00512 }
00513 return 0;
00514 }
00515 int df(const VectorXd &b, MatrixXd &fjac)
00516 {
00517 assert(b.size()==2);
00518 assert(fjac.rows()==14);
00519 assert(fjac.cols()==2);
00520 for(int i=0; i<14; i++) {
00521 double den = 1.+b[1]*x[i];
00522 fjac(i,0) = b[1]*x[i] / den;
00523 fjac(i,1) = b[0]*x[i]*(den-b[1]*x[i])/den/den;
00524 }
00525 return 0;
00526 }
00527 };
00528 const double misra1d_functor::x[14] = { 77.6E0, 114.9E0, 141.1E0, 190.8E0, 239.9E0, 289.0E0, 332.8E0, 378.4E0, 434.8E0, 477.3E0, 536.8E0, 593.1E0, 689.1E0, 760.0E0};
00529 const double misra1d_functor::y[14] = { 10.07E0, 14.73E0, 17.94E0, 23.93E0, 29.61E0, 35.18E0, 40.02E0, 44.82E0, 50.76E0, 55.05E0, 61.01E0, 66.40E0, 75.47E0, 81.78E0};
00530
00531
00532 void testNistMisra1d(void)
00533 {
00534 const int n=2;
00535 int info;
00536
00537 VectorXd x(n);
00538
00539
00540
00541
00542 x<< 500., 0.0001;
00543
00544 misra1d_functor functor;
00545 LevenbergMarquardt<misra1d_functor> lm(functor);
00546 info = lm.minimize(x);
00547
00548
00549 VERIFY_IS_EQUAL(info, 1);
00550 VERIFY_IS_EQUAL(lm.nfev(), 9);
00551 VERIFY_IS_EQUAL(lm.njev(), 7);
00552
00553 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6419295283E-02);
00554
00555 VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
00556 VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
00557
00558
00559
00560
00561 x<< 450., 0.0003;
00562
00563 info = lm.minimize(x);
00564
00565
00566 VERIFY_IS_EQUAL(info, 1);
00567 VERIFY_IS_EQUAL(lm.nfev(), 4);
00568 VERIFY_IS_EQUAL(lm.njev(), 3);
00569
00570 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6419295283E-02);
00571
00572 VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
00573 VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
00574 }
00575
00576
00577 struct lanczos1_functor : DenseFunctor<double>
00578 {
00579 lanczos1_functor(void) : DenseFunctor<double>(6,24) {}
00580 static const double x[24];
00581 static const double y[24];
00582 int operator()(const VectorXd &b, VectorXd &fvec)
00583 {
00584 assert(b.size()==6);
00585 assert(fvec.size()==24);
00586 for(int i=0; i<24; i++)
00587 fvec[i] = b[0]*exp(-b[1]*x[i]) + b[2]*exp(-b[3]*x[i]) + b[4]*exp(-b[5]*x[i]) - y[i];
00588 return 0;
00589 }
00590 int df(const VectorXd &b, MatrixXd &fjac)
00591 {
00592 assert(b.size()==6);
00593 assert(fjac.rows()==24);
00594 assert(fjac.cols()==6);
00595 for(int i=0; i<24; i++) {
00596 fjac(i,0) = exp(-b[1]*x[i]);
00597 fjac(i,1) = -b[0]*x[i]*exp(-b[1]*x[i]);
00598 fjac(i,2) = exp(-b[3]*x[i]);
00599 fjac(i,3) = -b[2]*x[i]*exp(-b[3]*x[i]);
00600 fjac(i,4) = exp(-b[5]*x[i]);
00601 fjac(i,5) = -b[4]*x[i]*exp(-b[5]*x[i]);
00602 }
00603 return 0;
00604 }
00605 };
00606 const double lanczos1_functor::x[24] = { 0.000000000000E+00, 5.000000000000E-02, 1.000000000000E-01, 1.500000000000E-01, 2.000000000000E-01, 2.500000000000E-01, 3.000000000000E-01, 3.500000000000E-01, 4.000000000000E-01, 4.500000000000E-01, 5.000000000000E-01, 5.500000000000E-01, 6.000000000000E-01, 6.500000000000E-01, 7.000000000000E-01, 7.500000000000E-01, 8.000000000000E-01, 8.500000000000E-01, 9.000000000000E-01, 9.500000000000E-01, 1.000000000000E+00, 1.050000000000E+00, 1.100000000000E+00, 1.150000000000E+00 };
00607 const double lanczos1_functor::y[24] = { 2.513400000000E+00 ,2.044333373291E+00 ,1.668404436564E+00 ,1.366418021208E+00 ,1.123232487372E+00 ,9.268897180037E-01 ,7.679338563728E-01 ,6.388775523106E-01 ,5.337835317402E-01 ,4.479363617347E-01 ,3.775847884350E-01 ,3.197393199326E-01 ,2.720130773746E-01 ,2.324965529032E-01 ,1.996589546065E-01 ,1.722704126914E-01 ,1.493405660168E-01 ,1.300700206922E-01 ,1.138119324644E-01 ,1.000415587559E-01 ,8.833209084540E-02 ,7.833544019350E-02 ,6.976693743449E-02 ,6.239312536719E-02 };
00608
00609
00610 void testNistLanczos1(void)
00611 {
00612 const int n=6;
00613 int info;
00614
00615 VectorXd x(n);
00616
00617
00618
00619
00620 x<< 1.2, 0.3, 5.6, 5.5, 6.5, 7.6;
00621
00622 lanczos1_functor functor;
00623 LevenbergMarquardt<lanczos1_functor> lm(functor);
00624 info = lm.minimize(x);
00625
00626
00627 VERIFY_IS_EQUAL(info, 2);
00628 VERIFY_IS_EQUAL(lm.nfev(), 79);
00629 VERIFY_IS_EQUAL(lm.njev(), 72);
00630
00631
00632
00633 VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
00634 VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
00635 VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
00636 VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
00637 VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
00638 VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
00639
00640
00641
00642
00643 x<< 0.5, 0.7, 3.6, 4.2, 4., 6.3;
00644
00645 info = lm.minimize(x);
00646
00647
00648 VERIFY_IS_EQUAL(info, 2);
00649 VERIFY_IS_EQUAL(lm.nfev(), 9);
00650 VERIFY_IS_EQUAL(lm.njev(), 8);
00651
00652
00653
00654 VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
00655 VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
00656 VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
00657 VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
00658 VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
00659 VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
00660
00661 }
00662
00663 struct rat42_functor : DenseFunctor<double>
00664 {
00665 rat42_functor(void) : DenseFunctor<double>(3,9) {}
00666 static const double x[9];
00667 static const double y[9];
00668 int operator()(const VectorXd &b, VectorXd &fvec)
00669 {
00670 assert(b.size()==3);
00671 assert(fvec.size()==9);
00672 for(int i=0; i<9; i++) {
00673 fvec[i] = b[0] / (1.+exp(b[1]-b[2]*x[i])) - y[i];
00674 }
00675 return 0;
00676 }
00677
00678 int df(const VectorXd &b, MatrixXd &fjac)
00679 {
00680 assert(b.size()==3);
00681 assert(fjac.rows()==9);
00682 assert(fjac.cols()==3);
00683 for(int i=0; i<9; i++) {
00684 double e = exp(b[1]-b[2]*x[i]);
00685 fjac(i,0) = 1./(1.+e);
00686 fjac(i,1) = -b[0]*e/(1.+e)/(1.+e);
00687 fjac(i,2) = +b[0]*e*x[i]/(1.+e)/(1.+e);
00688 }
00689 return 0;
00690 }
00691 };
00692 const double rat42_functor::x[9] = { 9.000E0, 14.000E0, 21.000E0, 28.000E0, 42.000E0, 57.000E0, 63.000E0, 70.000E0, 79.000E0 };
00693 const double rat42_functor::y[9] = { 8.930E0 ,10.800E0 ,18.590E0 ,22.330E0 ,39.350E0 ,56.110E0 ,61.730E0 ,64.620E0 ,67.080E0 };
00694
00695
00696 void testNistRat42(void)
00697 {
00698 const int n=3;
00699 int info;
00700
00701 VectorXd x(n);
00702
00703
00704
00705
00706 x<< 100., 1., 0.1;
00707
00708 rat42_functor functor;
00709 LevenbergMarquardt<rat42_functor> lm(functor);
00710 info = lm.minimize(x);
00711
00712
00713 VERIFY_IS_EQUAL(info, 1);
00714 VERIFY_IS_EQUAL(lm.nfev(), 10);
00715 VERIFY_IS_EQUAL(lm.njev(), 8);
00716
00717 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.0565229338E+00);
00718
00719 VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
00720 VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
00721 VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
00722
00723
00724
00725
00726 x<< 75., 2.5, 0.07;
00727
00728 info = lm.minimize(x);
00729
00730
00731 VERIFY_IS_EQUAL(info, 1);
00732 VERIFY_IS_EQUAL(lm.nfev(), 6);
00733 VERIFY_IS_EQUAL(lm.njev(), 5);
00734
00735 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.0565229338E+00);
00736
00737 VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
00738 VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
00739 VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
00740 }
00741
00742 struct MGH10_functor : DenseFunctor<double>
00743 {
00744 MGH10_functor(void) : DenseFunctor<double>(3,16) {}
00745 static const double x[16];
00746 static const double y[16];
00747 int operator()(const VectorXd &b, VectorXd &fvec)
00748 {
00749 assert(b.size()==3);
00750 assert(fvec.size()==16);
00751 for(int i=0; i<16; i++)
00752 fvec[i] = b[0] * exp(b[1]/(x[i]+b[2])) - y[i];
00753 return 0;
00754 }
00755 int df(const VectorXd &b, MatrixXd &fjac)
00756 {
00757 assert(b.size()==3);
00758 assert(fjac.rows()==16);
00759 assert(fjac.cols()==3);
00760 for(int i=0; i<16; i++) {
00761 double factor = 1./(x[i]+b[2]);
00762 double e = exp(b[1]*factor);
00763 fjac(i,0) = e;
00764 fjac(i,1) = b[0]*factor*e;
00765 fjac(i,2) = -b[1]*b[0]*factor*factor*e;
00766 }
00767 return 0;
00768 }
00769 };
00770 const double MGH10_functor::x[16] = { 5.000000E+01, 5.500000E+01, 6.000000E+01, 6.500000E+01, 7.000000E+01, 7.500000E+01, 8.000000E+01, 8.500000E+01, 9.000000E+01, 9.500000E+01, 1.000000E+02, 1.050000E+02, 1.100000E+02, 1.150000E+02, 1.200000E+02, 1.250000E+02 };
00771 const double MGH10_functor::y[16] = { 3.478000E+04, 2.861000E+04, 2.365000E+04, 1.963000E+04, 1.637000E+04, 1.372000E+04, 1.154000E+04, 9.744000E+03, 8.261000E+03, 7.030000E+03, 6.005000E+03, 5.147000E+03, 4.427000E+03, 3.820000E+03, 3.307000E+03, 2.872000E+03 };
00772
00773
00774 void testNistMGH10(void)
00775 {
00776 const int n=3;
00777 int info;
00778
00779 VectorXd x(n);
00780
00781
00782
00783
00784 x<< 2., 400000., 25000.;
00785
00786 MGH10_functor functor;
00787 LevenbergMarquardt<MGH10_functor> lm(functor);
00788 info = lm.minimize(x);
00789
00790
00791 VERIFY_IS_EQUAL(info, 1);
00792 VERIFY_IS_EQUAL(lm.nfev(), 284 );
00793 VERIFY_IS_EQUAL(lm.njev(), 249 );
00794
00795 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7945855171E+01);
00796
00797 VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
00798 VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
00799 VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
00800
00801
00802
00803
00804 x<< 0.02, 4000., 250.;
00805
00806 info = lm.minimize(x);
00807
00808
00809 VERIFY_IS_EQUAL(info, 1);
00810 VERIFY_IS_EQUAL(lm.nfev(), 126);
00811 VERIFY_IS_EQUAL(lm.njev(), 116);
00812
00813 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7945855171E+01);
00814
00815 VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
00816 VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
00817 VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
00818 }
00819
00820
00821 struct BoxBOD_functor : DenseFunctor<double>
00822 {
00823 BoxBOD_functor(void) : DenseFunctor<double>(2,6) {}
00824 static const double x[6];
00825 int operator()(const VectorXd &b, VectorXd &fvec)
00826 {
00827 static const double y[6] = { 109., 149., 149., 191., 213., 224. };
00828 assert(b.size()==2);
00829 assert(fvec.size()==6);
00830 for(int i=0; i<6; i++)
00831 fvec[i] = b[0]*(1.-exp(-b[1]*x[i])) - y[i];
00832 return 0;
00833 }
00834 int df(const VectorXd &b, MatrixXd &fjac)
00835 {
00836 assert(b.size()==2);
00837 assert(fjac.rows()==6);
00838 assert(fjac.cols()==2);
00839 for(int i=0; i<6; i++) {
00840 double e = exp(-b[1]*x[i]);
00841 fjac(i,0) = 1.-e;
00842 fjac(i,1) = b[0]*x[i]*e;
00843 }
00844 return 0;
00845 }
00846 };
00847 const double BoxBOD_functor::x[6] = { 1., 2., 3., 5., 7., 10. };
00848
00849
00850 void testNistBoxBOD(void)
00851 {
00852 const int n=2;
00853 int info;
00854
00855 VectorXd x(n);
00856
00857
00858
00859
00860 x<< 1., 1.;
00861
00862 BoxBOD_functor functor;
00863 LevenbergMarquardt<BoxBOD_functor> lm(functor);
00864 lm.setFtol(1.E6*NumTraits<double>::epsilon());
00865 lm.setXtol(1.E6*NumTraits<double>::epsilon());
00866 lm.setFactor(10);
00867 info = lm.minimize(x);
00868
00869
00870 VERIFY_IS_EQUAL(info, 1);
00871 VERIFY_IS_EQUAL(lm.nfev(), 31);
00872 VERIFY_IS_EQUAL(lm.njev(), 25);
00873
00874 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.1680088766E+03);
00875
00876 VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
00877 VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
00878
00879
00880
00881
00882 x<< 100., 0.75;
00883
00884 lm.resetParameters();
00885 lm.setFtol(NumTraits<double>::epsilon());
00886 lm.setXtol( NumTraits<double>::epsilon());
00887 info = lm.minimize(x);
00888
00889
00890 VERIFY_IS_EQUAL(info, 1);
00891 VERIFY_IS_EQUAL(lm.nfev(), 15 );
00892 VERIFY_IS_EQUAL(lm.njev(), 14 );
00893
00894 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.1680088766E+03);
00895
00896 VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
00897 VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
00898 }
00899
00900 struct MGH17_functor : DenseFunctor<double>
00901 {
00902 MGH17_functor(void) : DenseFunctor<double>(5,33) {}
00903 static const double x[33];
00904 static const double y[33];
00905 int operator()(const VectorXd &b, VectorXd &fvec)
00906 {
00907 assert(b.size()==5);
00908 assert(fvec.size()==33);
00909 for(int i=0; i<33; i++)
00910 fvec[i] = b[0] + b[1]*exp(-b[3]*x[i]) + b[2]*exp(-b[4]*x[i]) - y[i];
00911 return 0;
00912 }
00913 int df(const VectorXd &b, MatrixXd &fjac)
00914 {
00915 assert(b.size()==5);
00916 assert(fjac.rows()==33);
00917 assert(fjac.cols()==5);
00918 for(int i=0; i<33; i++) {
00919 fjac(i,0) = 1.;
00920 fjac(i,1) = exp(-b[3]*x[i]);
00921 fjac(i,2) = exp(-b[4]*x[i]);
00922 fjac(i,3) = -x[i]*b[1]*exp(-b[3]*x[i]);
00923 fjac(i,4) = -x[i]*b[2]*exp(-b[4]*x[i]);
00924 }
00925 return 0;
00926 }
00927 };
00928 const double MGH17_functor::x[33] = { 0.000000E+00, 1.000000E+01, 2.000000E+01, 3.000000E+01, 4.000000E+01, 5.000000E+01, 6.000000E+01, 7.000000E+01, 8.000000E+01, 9.000000E+01, 1.000000E+02, 1.100000E+02, 1.200000E+02, 1.300000E+02, 1.400000E+02, 1.500000E+02, 1.600000E+02, 1.700000E+02, 1.800000E+02, 1.900000E+02, 2.000000E+02, 2.100000E+02, 2.200000E+02, 2.300000E+02, 2.400000E+02, 2.500000E+02, 2.600000E+02, 2.700000E+02, 2.800000E+02, 2.900000E+02, 3.000000E+02, 3.100000E+02, 3.200000E+02 };
00929 const double MGH17_functor::y[33] = { 8.440000E-01, 9.080000E-01, 9.320000E-01, 9.360000E-01, 9.250000E-01, 9.080000E-01, 8.810000E-01, 8.500000E-01, 8.180000E-01, 7.840000E-01, 7.510000E-01, 7.180000E-01, 6.850000E-01, 6.580000E-01, 6.280000E-01, 6.030000E-01, 5.800000E-01, 5.580000E-01, 5.380000E-01, 5.220000E-01, 5.060000E-01, 4.900000E-01, 4.780000E-01, 4.670000E-01, 4.570000E-01, 4.480000E-01, 4.380000E-01, 4.310000E-01, 4.240000E-01, 4.200000E-01, 4.140000E-01, 4.110000E-01, 4.060000E-01 };
00930
00931
00932 void testNistMGH17(void)
00933 {
00934 const int n=5;
00935 int info;
00936
00937 VectorXd x(n);
00938
00939
00940
00941
00942 x<< 50., 150., -100., 1., 2.;
00943
00944 MGH17_functor functor;
00945 LevenbergMarquardt<MGH17_functor> lm(functor);
00946 lm.setFtol(NumTraits<double>::epsilon());
00947 lm.setXtol(NumTraits<double>::epsilon());
00948 lm.setMaxfev(1000);
00949 info = lm.minimize(x);
00950
00951
00952
00953
00954 VERIFY_IS_EQUAL(lm.njev(), 545 );
00955
00956 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.4648946975E-05);
00957
00958 VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
00959 VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
00960 VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
00961 VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
00962 VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
00963
00964
00965
00966
00967 x<< 0.5 ,1.5 ,-1 ,0.01 ,0.02;
00968
00969 lm.resetParameters();
00970 info = lm.minimize(x);
00971
00972
00973 VERIFY_IS_EQUAL(info, 1);
00974 VERIFY_IS_EQUAL(lm.nfev(), 18);
00975 VERIFY_IS_EQUAL(lm.njev(), 15);
00976
00977 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.4648946975E-05);
00978
00979 VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
00980 VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
00981 VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
00982 VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
00983 VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
00984 }
00985
00986 struct MGH09_functor : DenseFunctor<double>
00987 {
00988 MGH09_functor(void) : DenseFunctor<double>(4,11) {}
00989 static const double _x[11];
00990 static const double y[11];
00991 int operator()(const VectorXd &b, VectorXd &fvec)
00992 {
00993 assert(b.size()==4);
00994 assert(fvec.size()==11);
00995 for(int i=0; i<11; i++) {
00996 double x = _x[i], xx=x*x;
00997 fvec[i] = b[0]*(xx+x*b[1])/(xx+x*b[2]+b[3]) - y[i];
00998 }
00999 return 0;
01000 }
01001 int df(const VectorXd &b, MatrixXd &fjac)
01002 {
01003 assert(b.size()==4);
01004 assert(fjac.rows()==11);
01005 assert(fjac.cols()==4);
01006 for(int i=0; i<11; i++) {
01007 double x = _x[i], xx=x*x;
01008 double factor = 1./(xx+x*b[2]+b[3]);
01009 fjac(i,0) = (xx+x*b[1]) * factor;
01010 fjac(i,1) = b[0]*x* factor;
01011 fjac(i,2) = - b[0]*(xx+x*b[1]) * x * factor * factor;
01012 fjac(i,3) = - b[0]*(xx+x*b[1]) * factor * factor;
01013 }
01014 return 0;
01015 }
01016 };
01017 const double MGH09_functor::_x[11] = { 4., 2., 1., 5.E-1 , 2.5E-01, 1.670000E-01, 1.250000E-01, 1.E-01, 8.330000E-02, 7.140000E-02, 6.250000E-02 };
01018 const double MGH09_functor::y[11] = { 1.957000E-01, 1.947000E-01, 1.735000E-01, 1.600000E-01, 8.440000E-02, 6.270000E-02, 4.560000E-02, 3.420000E-02, 3.230000E-02, 2.350000E-02, 2.460000E-02 };
01019
01020
01021 void testNistMGH09(void)
01022 {
01023 const int n=4;
01024 int info;
01025
01026 VectorXd x(n);
01027
01028
01029
01030
01031 x<< 25., 39, 41.5, 39.;
01032
01033 MGH09_functor functor;
01034 LevenbergMarquardt<MGH09_functor> lm(functor);
01035 lm.setMaxfev(1000);
01036 info = lm.minimize(x);
01037
01038
01039 VERIFY_IS_EQUAL(info, 1);
01040 VERIFY_IS_EQUAL(lm.nfev(), 490 );
01041 VERIFY_IS_EQUAL(lm.njev(), 376 );
01042
01043 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 3.0750560385E-04);
01044
01045 VERIFY_IS_APPROX(x[0], 0.1928077089);
01046 VERIFY_IS_APPROX(x[1], 0.19126423573);
01047 VERIFY_IS_APPROX(x[2], 0.12305309914);
01048 VERIFY_IS_APPROX(x[3], 0.13605395375);
01049
01050
01051
01052
01053 x<< 0.25, 0.39, 0.415, 0.39;
01054
01055 lm.resetParameters();
01056 info = lm.minimize(x);
01057
01058
01059 VERIFY_IS_EQUAL(info, 1);
01060 VERIFY_IS_EQUAL(lm.nfev(), 18);
01061 VERIFY_IS_EQUAL(lm.njev(), 16);
01062
01063 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 3.0750560385E-04);
01064
01065 VERIFY_IS_APPROX(x[0], 0.19280781);
01066 VERIFY_IS_APPROX(x[1], 0.19126265);
01067 VERIFY_IS_APPROX(x[2], 0.12305280);
01068 VERIFY_IS_APPROX(x[3], 0.13605322);
01069 }
01070
01071
01072
01073 struct Bennett5_functor : DenseFunctor<double>
01074 {
01075 Bennett5_functor(void) : DenseFunctor<double>(3,154) {}
01076 static const double x[154];
01077 static const double y[154];
01078 int operator()(const VectorXd &b, VectorXd &fvec)
01079 {
01080 assert(b.size()==3);
01081 assert(fvec.size()==154);
01082 for(int i=0; i<154; i++)
01083 fvec[i] = b[0]* pow(b[1]+x[i],-1./b[2]) - y[i];
01084 return 0;
01085 }
01086 int df(const VectorXd &b, MatrixXd &fjac)
01087 {
01088 assert(b.size()==3);
01089 assert(fjac.rows()==154);
01090 assert(fjac.cols()==3);
01091 for(int i=0; i<154; i++) {
01092 double e = pow(b[1]+x[i],-1./b[2]);
01093 fjac(i,0) = e;
01094 fjac(i,1) = - b[0]*e/b[2]/(b[1]+x[i]);
01095 fjac(i,2) = b[0]*e*log(b[1]+x[i])/b[2]/b[2];
01096 }
01097 return 0;
01098 }
01099 };
01100 const double Bennett5_functor::x[154] = { 7.447168E0, 8.102586E0, 8.452547E0, 8.711278E0, 8.916774E0, 9.087155E0, 9.232590E0, 9.359535E0, 9.472166E0, 9.573384E0, 9.665293E0, 9.749461E0, 9.827092E0, 9.899128E0, 9.966321E0, 10.029280E0, 10.088510E0, 10.144430E0, 10.197380E0, 10.247670E0, 10.295560E0, 10.341250E0, 10.384950E0, 10.426820E0, 10.467000E0, 10.505640E0, 10.542830E0, 10.578690E0, 10.613310E0, 10.646780E0, 10.679150E0, 10.710520E0, 10.740920E0, 10.770440E0, 10.799100E0, 10.826970E0, 10.854080E0, 10.880470E0, 10.906190E0, 10.931260E0, 10.955720E0, 10.979590E0, 11.002910E0, 11.025700E0, 11.047980E0, 11.069770E0, 11.091100E0, 11.111980E0, 11.132440E0, 11.152480E0, 11.172130E0, 11.191410E0, 11.210310E0, 11.228870E0, 11.247090E0, 11.264980E0, 11.282560E0, 11.299840E0, 11.316820E0, 11.333520E0, 11.349940E0, 11.366100E0, 11.382000E0, 11.397660E0, 11.413070E0, 11.428240E0, 11.443200E0, 11.457930E0, 11.472440E0, 11.486750E0, 11.500860E0, 11.514770E0, 11.528490E0, 11.542020E0, 11.555380E0, 11.568550E0,
01101 11.581560E0, 11.594420E0, 11.607121E0, 11.619640E0, 11.632000E0, 11.644210E0, 11.656280E0, 11.668200E0, 11.679980E0, 11.691620E0, 11.703130E0, 11.714510E0, 11.725760E0, 11.736880E0, 11.747890E0, 11.758780E0, 11.769550E0, 11.780200E0, 11.790730E0, 11.801160E0, 11.811480E0, 11.821700E0, 11.831810E0, 11.841820E0, 11.851730E0, 11.861550E0, 11.871270E0, 11.880890E0, 11.890420E0, 11.899870E0, 11.909220E0, 11.918490E0, 11.927680E0, 11.936780E0, 11.945790E0, 11.954730E0, 11.963590E0, 11.972370E0, 11.981070E0, 11.989700E0, 11.998260E0, 12.006740E0, 12.015150E0, 12.023490E0, 12.031760E0, 12.039970E0, 12.048100E0, 12.056170E0, 12.064180E0, 12.072120E0, 12.080010E0, 12.087820E0, 12.095580E0, 12.103280E0, 12.110920E0, 12.118500E0, 12.126030E0, 12.133500E0, 12.140910E0, 12.148270E0, 12.155570E0, 12.162830E0, 12.170030E0, 12.177170E0, 12.184270E0, 12.191320E0, 12.198320E0, 12.205270E0, 12.212170E0, 12.219030E0, 12.225840E0, 12.232600E0, 12.239320E0, 12.245990E0, 12.252620E0, 12.259200E0, 12.265750E0, 12.272240E0 };
01102 const double Bennett5_functor::y[154] = { -34.834702E0 ,-34.393200E0 ,-34.152901E0 ,-33.979099E0 ,-33.845901E0 ,-33.732899E0 ,-33.640301E0 ,-33.559200E0 ,-33.486801E0 ,-33.423100E0 ,-33.365101E0 ,-33.313000E0 ,-33.260899E0 ,-33.217400E0 ,-33.176899E0 ,-33.139198E0 ,-33.101601E0 ,-33.066799E0 ,-33.035000E0 ,-33.003101E0 ,-32.971298E0 ,-32.942299E0 ,-32.916302E0 ,-32.890202E0 ,-32.864101E0 ,-32.841000E0 ,-32.817799E0 ,-32.797501E0 ,-32.774300E0 ,-32.757000E0 ,-32.733799E0 ,-32.716400E0 ,-32.699100E0 ,-32.678799E0 ,-32.661400E0 ,-32.644001E0 ,-32.626701E0 ,-32.612202E0 ,-32.597698E0 ,-32.583199E0 ,-32.568699E0 ,-32.554298E0 ,-32.539799E0 ,-32.525299E0 ,-32.510799E0 ,-32.499199E0 ,-32.487598E0 ,-32.473202E0 ,-32.461601E0 ,-32.435501E0 ,-32.435501E0 ,-32.426800E0 ,-32.412300E0 ,-32.400799E0 ,-32.392101E0 ,-32.380501E0 ,-32.366001E0 ,-32.357300E0 ,-32.348598E0 ,-32.339901E0 ,-32.328400E0 ,-32.319698E0 ,-32.311001E0 ,-32.299400E0 ,-32.290699E0 ,-32.282001E0 ,-32.273300E0 ,-32.264599E0 ,-32.256001E0 ,-32.247299E0
01103 ,-32.238602E0 ,-32.229900E0 ,-32.224098E0 ,-32.215401E0 ,-32.203800E0 ,-32.198002E0 ,-32.189400E0 ,-32.183601E0 ,-32.174900E0 ,-32.169102E0 ,-32.163300E0 ,-32.154598E0 ,-32.145901E0 ,-32.140099E0 ,-32.131401E0 ,-32.125599E0 ,-32.119801E0 ,-32.111198E0 ,-32.105400E0 ,-32.096699E0 ,-32.090900E0 ,-32.088001E0 ,-32.079300E0 ,-32.073502E0 ,-32.067699E0 ,-32.061901E0 ,-32.056099E0 ,-32.050301E0 ,-32.044498E0 ,-32.038799E0 ,-32.033001E0 ,-32.027199E0 ,-32.024300E0 ,-32.018501E0 ,-32.012699E0 ,-32.004002E0 ,-32.001099E0 ,-31.995300E0 ,-31.989500E0 ,-31.983700E0 ,-31.977900E0 ,-31.972099E0 ,-31.969299E0 ,-31.963501E0 ,-31.957701E0 ,-31.951900E0 ,-31.946100E0 ,-31.940300E0 ,-31.937401E0 ,-31.931601E0 ,-31.925800E0 ,-31.922899E0 ,-31.917101E0 ,-31.911301E0 ,-31.908400E0 ,-31.902599E0 ,-31.896900E0 ,-31.893999E0 ,-31.888201E0 ,-31.885300E0 ,-31.882401E0 ,-31.876600E0 ,-31.873699E0 ,-31.867901E0 ,-31.862101E0 ,-31.859200E0 ,-31.856300E0 ,-31.850500E0 ,-31.844700E0 ,-31.841801E0 ,-31.838900E0 ,-31.833099E0 ,-31.830200E0 ,
01104 -31.827299E0 ,-31.821600E0 ,-31.818701E0 ,-31.812901E0 ,-31.809999E0 ,-31.807100E0 ,-31.801300E0 ,-31.798401E0 ,-31.795500E0 ,-31.789700E0 ,-31.786800E0 };
01105
01106
01107 void testNistBennett5(void)
01108 {
01109 const int n=3;
01110 int info;
01111
01112 VectorXd x(n);
01113
01114
01115
01116
01117 x<< -2000., 50., 0.8;
01118
01119 Bennett5_functor functor;
01120 LevenbergMarquardt<Bennett5_functor> lm(functor);
01121 lm.setMaxfev(1000);
01122 info = lm.minimize(x);
01123
01124
01125 VERIFY_IS_EQUAL(info, 1);
01126 VERIFY_IS_EQUAL(lm.nfev(), 758);
01127 VERIFY_IS_EQUAL(lm.njev(), 744);
01128
01129 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.2404744073E-04);
01130
01131 VERIFY_IS_APPROX(x[0], -2.5235058043E+03);
01132 VERIFY_IS_APPROX(x[1], 4.6736564644E+01);
01133 VERIFY_IS_APPROX(x[2], 9.3218483193E-01);
01134
01135
01136
01137 x<< -1500., 45., 0.85;
01138
01139 lm.resetParameters();
01140 info = lm.minimize(x);
01141
01142
01143 VERIFY_IS_EQUAL(info, 1);
01144 VERIFY_IS_EQUAL(lm.nfev(), 203);
01145 VERIFY_IS_EQUAL(lm.njev(), 192);
01146
01147 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.2404744073E-04);
01148
01149 VERIFY_IS_APPROX(x[0], -2523.3007865);
01150 VERIFY_IS_APPROX(x[1], 46.735705771);
01151 VERIFY_IS_APPROX(x[2], 0.93219881891);
01152 }
01153
01154 struct thurber_functor : DenseFunctor<double>
01155 {
01156 thurber_functor(void) : DenseFunctor<double>(7,37) {}
01157 static const double _x[37];
01158 static const double _y[37];
01159 int operator()(const VectorXd &b, VectorXd &fvec)
01160 {
01161
01162 assert(b.size()==7);
01163 assert(fvec.size()==37);
01164 for(int i=0; i<37; i++) {
01165 double x=_x[i], xx=x*x, xxx=xx*x;
01166 fvec[i] = (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) / (1.+b[4]*x+b[5]*xx+b[6]*xxx) - _y[i];
01167 }
01168 return 0;
01169 }
01170 int df(const VectorXd &b, MatrixXd &fjac)
01171 {
01172 assert(b.size()==7);
01173 assert(fjac.rows()==37);
01174 assert(fjac.cols()==7);
01175 for(int i=0; i<37; i++) {
01176 double x=_x[i], xx=x*x, xxx=xx*x;
01177 double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
01178 fjac(i,0) = 1.*fact;
01179 fjac(i,1) = x*fact;
01180 fjac(i,2) = xx*fact;
01181 fjac(i,3) = xxx*fact;
01182 fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
01183 fjac(i,4) = x*fact;
01184 fjac(i,5) = xx*fact;
01185 fjac(i,6) = xxx*fact;
01186 }
01187 return 0;
01188 }
01189 };
01190 const double thurber_functor::_x[37] = { -3.067E0, -2.981E0, -2.921E0, -2.912E0, -2.840E0, -2.797E0, -2.702E0, -2.699E0, -2.633E0, -2.481E0, -2.363E0, -2.322E0, -1.501E0, -1.460E0, -1.274E0, -1.212E0, -1.100E0, -1.046E0, -0.915E0, -0.714E0, -0.566E0, -0.545E0, -0.400E0, -0.309E0, -0.109E0, -0.103E0, 0.010E0, 0.119E0, 0.377E0, 0.790E0, 0.963E0, 1.006E0, 1.115E0, 1.572E0, 1.841E0, 2.047E0, 2.200E0 };
01191 const double thurber_functor::_y[37] = { 80.574E0, 84.248E0, 87.264E0, 87.195E0, 89.076E0, 89.608E0, 89.868E0, 90.101E0, 92.405E0, 95.854E0, 100.696E0, 101.060E0, 401.672E0, 390.724E0, 567.534E0, 635.316E0, 733.054E0, 759.087E0, 894.206E0, 990.785E0, 1090.109E0, 1080.914E0, 1122.643E0, 1178.351E0, 1260.531E0, 1273.514E0, 1288.339E0, 1327.543E0, 1353.863E0, 1414.509E0, 1425.208E0, 1421.384E0, 1442.962E0, 1464.350E0, 1468.705E0, 1447.894E0, 1457.628E0};
01192
01193
01194 void testNistThurber(void)
01195 {
01196 const int n=7;
01197 int info;
01198
01199 VectorXd x(n);
01200
01201
01202
01203
01204 x<< 1000 ,1000 ,400 ,40 ,0.7,0.3,0.0 ;
01205
01206 thurber_functor functor;
01207 LevenbergMarquardt<thurber_functor> lm(functor);
01208 lm.setFtol(1.E4*NumTraits<double>::epsilon());
01209 lm.setXtol(1.E4*NumTraits<double>::epsilon());
01210 info = lm.minimize(x);
01211
01212
01213 VERIFY_IS_EQUAL(info, 1);
01214 VERIFY_IS_EQUAL(lm.nfev(), 39);
01215 VERIFY_IS_EQUAL(lm.njev(), 36);
01216
01217 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6427082397E+03);
01218
01219 VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
01220 VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
01221 VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
01222 VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
01223 VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
01224 VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
01225 VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
01226
01227
01228
01229
01230 x<< 1300 ,1500 ,500 ,75 ,1 ,0.4 ,0.05 ;
01231
01232 lm.resetParameters();
01233 lm.setFtol(1.E4*NumTraits<double>::epsilon());
01234 lm.setXtol(1.E4*NumTraits<double>::epsilon());
01235 info = lm.minimize(x);
01236
01237
01238 VERIFY_IS_EQUAL(info, 1);
01239 VERIFY_IS_EQUAL(lm.nfev(), 29);
01240 VERIFY_IS_EQUAL(lm.njev(), 28);
01241
01242 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6427082397E+03);
01243
01244 VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
01245 VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
01246 VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
01247 VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
01248 VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
01249 VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
01250 VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
01251 }
01252
01253 struct rat43_functor : DenseFunctor<double>
01254 {
01255 rat43_functor(void) : DenseFunctor<double>(4,15) {}
01256 static const double x[15];
01257 static const double y[15];
01258 int operator()(const VectorXd &b, VectorXd &fvec)
01259 {
01260 assert(b.size()==4);
01261 assert(fvec.size()==15);
01262 for(int i=0; i<15; i++)
01263 fvec[i] = b[0] * pow(1.+exp(b[1]-b[2]*x[i]),-1./b[3]) - y[i];
01264 return 0;
01265 }
01266 int df(const VectorXd &b, MatrixXd &fjac)
01267 {
01268 assert(b.size()==4);
01269 assert(fjac.rows()==15);
01270 assert(fjac.cols()==4);
01271 for(int i=0; i<15; i++) {
01272 double e = exp(b[1]-b[2]*x[i]);
01273 double power = -1./b[3];
01274 fjac(i,0) = pow(1.+e, power);
01275 fjac(i,1) = power*b[0]*e*pow(1.+e, power-1.);
01276 fjac(i,2) = -power*b[0]*e*x[i]*pow(1.+e, power-1.);
01277 fjac(i,3) = b[0]*power*power*log(1.+e)*pow(1.+e, power);
01278 }
01279 return 0;
01280 }
01281 };
01282 const double rat43_functor::x[15] = { 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14., 15. };
01283 const double rat43_functor::y[15] = { 16.08, 33.83, 65.80, 97.20, 191.55, 326.20, 386.87, 520.53, 590.03, 651.92, 724.93, 699.56, 689.96, 637.56, 717.41 };
01284
01285
01286 void testNistRat43(void)
01287 {
01288 const int n=4;
01289 int info;
01290
01291 VectorXd x(n);
01292
01293
01294
01295
01296 x<< 100., 10., 1., 1.;
01297
01298 rat43_functor functor;
01299 LevenbergMarquardt<rat43_functor> lm(functor);
01300 lm.setFtol(1.E6*NumTraits<double>::epsilon());
01301 lm.setXtol(1.E6*NumTraits<double>::epsilon());
01302 info = lm.minimize(x);
01303
01304
01305 VERIFY_IS_EQUAL(info, 1);
01306 VERIFY_IS_EQUAL(lm.nfev(), 27);
01307 VERIFY_IS_EQUAL(lm.njev(), 20);
01308
01309 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7864049080E+03);
01310
01311 VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
01312 VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
01313 VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
01314 VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
01315
01316
01317
01318
01319 x<< 700., 5., 0.75, 1.3;
01320
01321 lm.resetParameters();
01322 lm.setFtol(1.E5*NumTraits<double>::epsilon());
01323 lm.setXtol(1.E5*NumTraits<double>::epsilon());
01324 info = lm.minimize(x);
01325
01326
01327 VERIFY_IS_EQUAL(info, 1);
01328 VERIFY_IS_EQUAL(lm.nfev(), 9);
01329 VERIFY_IS_EQUAL(lm.njev(), 8);
01330
01331 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7864049080E+03);
01332
01333 VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
01334 VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
01335 VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
01336 VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
01337 }
01338
01339
01340
01341 struct eckerle4_functor : DenseFunctor<double>
01342 {
01343 eckerle4_functor(void) : DenseFunctor<double>(3,35) {}
01344 static const double x[35];
01345 static const double y[35];
01346 int operator()(const VectorXd &b, VectorXd &fvec)
01347 {
01348 assert(b.size()==3);
01349 assert(fvec.size()==35);
01350 for(int i=0; i<35; i++)
01351 fvec[i] = b[0]/b[1] * exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/(b[1]*b[1])) - y[i];
01352 return 0;
01353 }
01354 int df(const VectorXd &b, MatrixXd &fjac)
01355 {
01356 assert(b.size()==3);
01357 assert(fjac.rows()==35);
01358 assert(fjac.cols()==3);
01359 for(int i=0; i<35; i++) {
01360 double b12 = b[1]*b[1];
01361 double e = exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/b12);
01362 fjac(i,0) = e / b[1];
01363 fjac(i,1) = ((x[i]-b[2])*(x[i]-b[2])/b12-1.) * b[0]*e/b12;
01364 fjac(i,2) = (x[i]-b[2])*e*b[0]/b[1]/b12;
01365 }
01366 return 0;
01367 }
01368 };
01369 const double eckerle4_functor::x[35] = { 400.0, 405.0, 410.0, 415.0, 420.0, 425.0, 430.0, 435.0, 436.5, 438.0, 439.5, 441.0, 442.5, 444.0, 445.5, 447.0, 448.5, 450.0, 451.5, 453.0, 454.5, 456.0, 457.5, 459.0, 460.5, 462.0, 463.5, 465.0, 470.0, 475.0, 480.0, 485.0, 490.0, 495.0, 500.0};
01370 const double eckerle4_functor::y[35] = { 0.0001575, 0.0001699, 0.0002350, 0.0003102, 0.0004917, 0.0008710, 0.0017418, 0.0046400, 0.0065895, 0.0097302, 0.0149002, 0.0237310, 0.0401683, 0.0712559, 0.1264458, 0.2073413, 0.2902366, 0.3445623, 0.3698049, 0.3668534, 0.3106727, 0.2078154, 0.1164354, 0.0616764, 0.0337200, 0.0194023, 0.0117831, 0.0074357, 0.0022732, 0.0008800, 0.0004579, 0.0002345, 0.0001586, 0.0001143, 0.0000710 };
01371
01372
01373 void testNistEckerle4(void)
01374 {
01375 const int n=3;
01376 int info;
01377
01378 VectorXd x(n);
01379
01380
01381
01382
01383 x<< 1., 10., 500.;
01384
01385 eckerle4_functor functor;
01386 LevenbergMarquardt<eckerle4_functor> lm(functor);
01387 info = lm.minimize(x);
01388
01389
01390 VERIFY_IS_EQUAL(info, 1);
01391 VERIFY_IS_EQUAL(lm.nfev(), 18);
01392 VERIFY_IS_EQUAL(lm.njev(), 15);
01393
01394 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.4635887487E-03);
01395
01396 VERIFY_IS_APPROX(x[0], 1.5543827178);
01397 VERIFY_IS_APPROX(x[1], 4.0888321754);
01398 VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
01399
01400
01401
01402
01403 x<< 1.5, 5., 450.;
01404
01405 info = lm.minimize(x);
01406
01407
01408 VERIFY_IS_EQUAL(info, 1);
01409 VERIFY_IS_EQUAL(lm.nfev(), 7);
01410 VERIFY_IS_EQUAL(lm.njev(), 6);
01411
01412 VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.4635887487E-03);
01413
01414 VERIFY_IS_APPROX(x[0], 1.5543827178);
01415 VERIFY_IS_APPROX(x[1], 4.0888321754);
01416 VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
01417 }
01418
01419 void test_levenberg_marquardt()
01420 {
01421
01422 CALL_SUBTEST(testLmder1());
01423 CALL_SUBTEST(testLmder());
01424 CALL_SUBTEST(testLmdif1());
01425
01426
01427 CALL_SUBTEST(testLmdif());
01428
01429
01430 CALL_SUBTEST(testNistMisra1a());
01431 CALL_SUBTEST(testNistChwirut2());
01432
01433
01434 CALL_SUBTEST(testNistHahn1());
01435 CALL_SUBTEST(testNistMisra1d());
01436 CALL_SUBTEST(testNistMGH17());
01437 CALL_SUBTEST(testNistLanczos1());
01438
01439
01440 CALL_SUBTEST(testNistRat42());
01441 CALL_SUBTEST(testNistMGH10());
01442 CALL_SUBTEST(testNistBoxBOD());
01443
01444 CALL_SUBTEST(testNistBennett5());
01445 CALL_SUBTEST(testNistThurber());
01446 CALL_SUBTEST(testNistRat43());
01447 CALL_SUBTEST(testNistEckerle4());
01448 }