You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

1457 lines
54 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
  5. // Copyright (C) 2012 desire Nuentsa <desire.nuentsa_wakam@inria.fr
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. // FIXME: These tests all check for hard-coded values. Ideally, parameters and start estimates should be randomized.
  11. #include <stdio.h>
  12. #include "main.h"
  13. #include <unsupported/Eigen/LevenbergMarquardt>
  14. // This disables some useless Warnings on MSVC.
  15. // It is intended to be done for this test only.
  16. #include <Eigen/src/Core/util/DisableStupidWarnings.h>
  17. using std::sqrt;
  18. struct lmder_functor : DenseFunctor<double>
  19. {
  20. lmder_functor(void): DenseFunctor<double>(3,15) {}
  21. int operator()(const VectorXd &x, VectorXd &fvec) const
  22. {
  23. double tmp1, tmp2, tmp3;
  24. 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,
  25. 3.9e-1, 3.7e-1, 5.8e-1, 7.3e-1, 9.6e-1, 1.34, 2.1, 4.39};
  26. for (int i = 0; i < values(); i++)
  27. {
  28. tmp1 = i+1;
  29. tmp2 = 16 - i - 1;
  30. tmp3 = (i>=8)? tmp2 : tmp1;
  31. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  32. }
  33. return 0;
  34. }
  35. int df(const VectorXd &x, MatrixXd &fjac) const
  36. {
  37. double tmp1, tmp2, tmp3, tmp4;
  38. for (int i = 0; i < values(); i++)
  39. {
  40. tmp1 = i+1;
  41. tmp2 = 16 - i - 1;
  42. tmp3 = (i>=8)? tmp2 : tmp1;
  43. tmp4 = (x[1]*tmp2 + x[2]*tmp3); tmp4 = tmp4*tmp4;
  44. fjac(i,0) = -1;
  45. fjac(i,1) = tmp1*tmp2/tmp4;
  46. fjac(i,2) = tmp1*tmp3/tmp4;
  47. }
  48. return 0;
  49. }
  50. };
  51. void testLmder1()
  52. {
  53. int n=3, info;
  54. VectorXd x;
  55. /* the following starting values provide a rough fit. */
  56. x.setConstant(n, 1.);
  57. // do the computation
  58. lmder_functor functor;
  59. LevenbergMarquardt<lmder_functor> lm(functor);
  60. info = lm.lmder1(x);
  61. // check return value
  62. VERIFY_IS_EQUAL(info, 1);
  63. VERIFY_IS_EQUAL(lm.nfev(), 6);
  64. VERIFY_IS_EQUAL(lm.njev(), 5);
  65. // check norm
  66. VERIFY_IS_APPROX(lm.fvec().blueNorm(), 0.09063596);
  67. // check x
  68. VectorXd x_ref(n);
  69. x_ref << 0.08241058, 1.133037, 2.343695;
  70. VERIFY_IS_APPROX(x, x_ref);
  71. }
  72. void testLmder()
  73. {
  74. const int m=15, n=3;
  75. int info;
  76. double fnorm, covfac;
  77. VectorXd x;
  78. /* the following starting values provide a rough fit. */
  79. x.setConstant(n, 1.);
  80. // do the computation
  81. lmder_functor functor;
  82. LevenbergMarquardt<lmder_functor> lm(functor);
  83. info = lm.minimize(x);
  84. // check return values
  85. VERIFY_IS_EQUAL(info, 1);
  86. VERIFY_IS_EQUAL(lm.nfev(), 6);
  87. VERIFY_IS_EQUAL(lm.njev(), 5);
  88. // check norm
  89. fnorm = lm.fvec().blueNorm();
  90. VERIFY_IS_APPROX(fnorm, 0.09063596);
  91. // check x
  92. VectorXd x_ref(n);
  93. x_ref << 0.08241058, 1.133037, 2.343695;
  94. VERIFY_IS_APPROX(x, x_ref);
  95. // check covariance
  96. covfac = fnorm*fnorm/(m-n);
  97. internal::covar(lm.matrixR(), lm.permutation().indices()); // TODO : move this as a function of lm
  98. MatrixXd cov_ref(n,n);
  99. cov_ref <<
  100. 0.0001531202, 0.002869941, -0.002656662,
  101. 0.002869941, 0.09480935, -0.09098995,
  102. -0.002656662, -0.09098995, 0.08778727;
  103. // std::cout << fjac*covfac << std::endl;
  104. MatrixXd cov;
  105. cov = covfac*lm.matrixR().topLeftCorner<n,n>();
  106. VERIFY_IS_APPROX( cov, cov_ref);
  107. // TODO: why isn't this allowed ? :
  108. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  109. }
  110. struct lmdif_functor : DenseFunctor<double>
  111. {
  112. lmdif_functor(void) : DenseFunctor<double>(3,15) {}
  113. int operator()(const VectorXd &x, VectorXd &fvec) const
  114. {
  115. int i;
  116. double tmp1,tmp2,tmp3;
  117. 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,
  118. 3.7e-1,5.8e-1,7.3e-1,9.6e-1,1.34e0,2.1e0,4.39e0};
  119. assert(x.size()==3);
  120. assert(fvec.size()==15);
  121. for (i=0; i<15; i++)
  122. {
  123. tmp1 = i+1;
  124. tmp2 = 15 - i;
  125. tmp3 = tmp1;
  126. if (i >= 8) tmp3 = tmp2;
  127. fvec[i] = y[i] - (x[0] + tmp1/(x[1]*tmp2 + x[2]*tmp3));
  128. }
  129. return 0;
  130. }
  131. };
  132. void testLmdif1()
  133. {
  134. const int n=3;
  135. int info;
  136. VectorXd x(n), fvec(15);
  137. /* the following starting values provide a rough fit. */
  138. x.setConstant(n, 1.);
  139. // do the computation
  140. lmdif_functor functor;
  141. DenseIndex nfev;
  142. info = LevenbergMarquardt<lmdif_functor>::lmdif1(functor, x, &nfev);
  143. // check return value
  144. VERIFY_IS_EQUAL(info, 1);
  145. // VERIFY_IS_EQUAL(nfev, 26);
  146. // check norm
  147. functor(x, fvec);
  148. VERIFY_IS_APPROX(fvec.blueNorm(), 0.09063596);
  149. // check x
  150. VectorXd x_ref(n);
  151. x_ref << 0.0824106, 1.1330366, 2.3436947;
  152. VERIFY_IS_APPROX(x, x_ref);
  153. }
  154. void testLmdif()
  155. {
  156. const int m=15, n=3;
  157. int info;
  158. double fnorm, covfac;
  159. VectorXd x(n);
  160. /* the following starting values provide a rough fit. */
  161. x.setConstant(n, 1.);
  162. // do the computation
  163. lmdif_functor functor;
  164. NumericalDiff<lmdif_functor> numDiff(functor);
  165. LevenbergMarquardt<NumericalDiff<lmdif_functor> > lm(numDiff);
  166. info = lm.minimize(x);
  167. // check return values
  168. VERIFY_IS_EQUAL(info, 1);
  169. // VERIFY_IS_EQUAL(lm.nfev(), 26);
  170. // check norm
  171. fnorm = lm.fvec().blueNorm();
  172. VERIFY_IS_APPROX(fnorm, 0.09063596);
  173. // check x
  174. VectorXd x_ref(n);
  175. x_ref << 0.08241058, 1.133037, 2.343695;
  176. VERIFY_IS_APPROX(x, x_ref);
  177. // check covariance
  178. covfac = fnorm*fnorm/(m-n);
  179. internal::covar(lm.matrixR(), lm.permutation().indices()); // TODO : move this as a function of lm
  180. MatrixXd cov_ref(n,n);
  181. cov_ref <<
  182. 0.0001531202, 0.002869942, -0.002656662,
  183. 0.002869942, 0.09480937, -0.09098997,
  184. -0.002656662, -0.09098997, 0.08778729;
  185. // std::cout << fjac*covfac << std::endl;
  186. MatrixXd cov;
  187. cov = covfac*lm.matrixR().topLeftCorner<n,n>();
  188. VERIFY_IS_APPROX( cov, cov_ref);
  189. // TODO: why isn't this allowed ? :
  190. // VERIFY_IS_APPROX( covfac*fjac.topLeftCorner<n,n>() , cov_ref);
  191. }
  192. struct chwirut2_functor : DenseFunctor<double>
  193. {
  194. chwirut2_functor(void) : DenseFunctor<double>(3,54) {}
  195. static const double m_x[54];
  196. static const double m_y[54];
  197. int operator()(const VectorXd &b, VectorXd &fvec)
  198. {
  199. int i;
  200. assert(b.size()==3);
  201. assert(fvec.size()==54);
  202. for(i=0; i<54; i++) {
  203. double x = m_x[i];
  204. fvec[i] = exp(-b[0]*x)/(b[1]+b[2]*x) - m_y[i];
  205. }
  206. return 0;
  207. }
  208. int df(const VectorXd &b, MatrixXd &fjac)
  209. {
  210. assert(b.size()==3);
  211. assert(fjac.rows()==54);
  212. assert(fjac.cols()==3);
  213. for(int i=0; i<54; i++) {
  214. double x = m_x[i];
  215. double factor = 1./(b[1]+b[2]*x);
  216. double e = exp(-b[0]*x);
  217. fjac(i,0) = -x*e*factor;
  218. fjac(i,1) = -e*factor*factor;
  219. fjac(i,2) = -x*e*factor*factor;
  220. }
  221. return 0;
  222. }
  223. };
  224. 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};
  225. 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 };
  226. // http://www.itl.nist.gov/div898/strd/nls/data/chwirut2.shtml
  227. void testNistChwirut2(void)
  228. {
  229. const int n=3;
  230. LevenbergMarquardtSpace::Status info;
  231. VectorXd x(n);
  232. /*
  233. * First try
  234. */
  235. x<< 0.1, 0.01, 0.02;
  236. // do the computation
  237. chwirut2_functor functor;
  238. LevenbergMarquardt<chwirut2_functor> lm(functor);
  239. info = lm.minimize(x);
  240. // check return value
  241. VERIFY_IS_EQUAL(info, 1);
  242. // VERIFY_IS_EQUAL(lm.nfev(), 10);
  243. VERIFY_IS_EQUAL(lm.njev(), 8);
  244. // check norm^2
  245. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.1304802941E+02);
  246. // check x
  247. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  248. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  249. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  250. /*
  251. * Second try
  252. */
  253. x<< 0.15, 0.008, 0.010;
  254. // do the computation
  255. lm.resetParameters();
  256. lm.setFtol(1.E6*NumTraits<double>::epsilon());
  257. lm.setXtol(1.E6*NumTraits<double>::epsilon());
  258. info = lm.minimize(x);
  259. // check return value
  260. VERIFY_IS_EQUAL(info, 1);
  261. // VERIFY_IS_EQUAL(lm.nfev(), 7);
  262. VERIFY_IS_EQUAL(lm.njev(), 6);
  263. // check norm^2
  264. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.1304802941E+02);
  265. // check x
  266. VERIFY_IS_APPROX(x[0], 1.6657666537E-01);
  267. VERIFY_IS_APPROX(x[1], 5.1653291286E-03);
  268. VERIFY_IS_APPROX(x[2], 1.2150007096E-02);
  269. }
  270. struct misra1a_functor : DenseFunctor<double>
  271. {
  272. misra1a_functor(void) : DenseFunctor<double>(2,14) {}
  273. static const double m_x[14];
  274. static const double m_y[14];
  275. int operator()(const VectorXd &b, VectorXd &fvec)
  276. {
  277. assert(b.size()==2);
  278. assert(fvec.size()==14);
  279. for(int i=0; i<14; i++) {
  280. fvec[i] = b[0]*(1.-exp(-b[1]*m_x[i])) - m_y[i] ;
  281. }
  282. return 0;
  283. }
  284. int df(const VectorXd &b, MatrixXd &fjac)
  285. {
  286. assert(b.size()==2);
  287. assert(fjac.rows()==14);
  288. assert(fjac.cols()==2);
  289. for(int i=0; i<14; i++) {
  290. fjac(i,0) = (1.-exp(-b[1]*m_x[i]));
  291. fjac(i,1) = (b[0]*m_x[i]*exp(-b[1]*m_x[i]));
  292. }
  293. return 0;
  294. }
  295. };
  296. 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};
  297. 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};
  298. // http://www.itl.nist.gov/div898/strd/nls/data/misra1a.shtml
  299. void testNistMisra1a(void)
  300. {
  301. const int n=2;
  302. int info;
  303. VectorXd x(n);
  304. /*
  305. * First try
  306. */
  307. x<< 500., 0.0001;
  308. // do the computation
  309. misra1a_functor functor;
  310. LevenbergMarquardt<misra1a_functor> lm(functor);
  311. info = lm.minimize(x);
  312. // check return value
  313. VERIFY_IS_EQUAL(info, 1);
  314. VERIFY_IS_EQUAL(lm.nfev(), 19);
  315. VERIFY_IS_EQUAL(lm.njev(), 15);
  316. // check norm^2
  317. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.2455138894E-01);
  318. // check x
  319. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  320. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  321. /*
  322. * Second try
  323. */
  324. x<< 250., 0.0005;
  325. // do the computation
  326. info = lm.minimize(x);
  327. // check return value
  328. VERIFY_IS_EQUAL(info, 1);
  329. VERIFY_IS_EQUAL(lm.nfev(), 5);
  330. VERIFY_IS_EQUAL(lm.njev(), 4);
  331. // check norm^2
  332. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.2455138894E-01);
  333. // check x
  334. VERIFY_IS_APPROX(x[0], 2.3894212918E+02);
  335. VERIFY_IS_APPROX(x[1], 5.5015643181E-04);
  336. }
  337. struct hahn1_functor : DenseFunctor<double>
  338. {
  339. hahn1_functor(void) : DenseFunctor<double>(7,236) {}
  340. static const double m_x[236];
  341. int operator()(const VectorXd &b, VectorXd &fvec)
  342. {
  343. 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 ,
  344. 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 ,
  345. 12.596E0 ,
  346. 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 };
  347. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  348. assert(b.size()==7);
  349. assert(fvec.size()==236);
  350. for(int i=0; i<236; i++) {
  351. double x=m_x[i], xx=x*x, xxx=xx*x;
  352. 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];
  353. }
  354. return 0;
  355. }
  356. int df(const VectorXd &b, MatrixXd &fjac)
  357. {
  358. assert(b.size()==7);
  359. assert(fjac.rows()==236);
  360. assert(fjac.cols()==7);
  361. for(int i=0; i<236; i++) {
  362. double x=m_x[i], xx=x*x, xxx=xx*x;
  363. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  364. fjac(i,0) = 1.*fact;
  365. fjac(i,1) = x*fact;
  366. fjac(i,2) = xx*fact;
  367. fjac(i,3) = xxx*fact;
  368. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  369. fjac(i,4) = x*fact;
  370. fjac(i,5) = xx*fact;
  371. fjac(i,6) = xxx*fact;
  372. }
  373. return 0;
  374. }
  375. };
  376. 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 ,
  377. 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 ,
  378. 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};
  379. // http://www.itl.nist.gov/div898/strd/nls/data/hahn1.shtml
  380. void testNistHahn1(void)
  381. {
  382. const int n=7;
  383. int info;
  384. VectorXd x(n);
  385. /*
  386. * First try
  387. */
  388. x<< 10., -1., .05, -.00001, -.05, .001, -.000001;
  389. // do the computation
  390. hahn1_functor functor;
  391. LevenbergMarquardt<hahn1_functor> lm(functor);
  392. info = lm.minimize(x);
  393. // check return value
  394. VERIFY_IS_EQUAL(info, 1);
  395. VERIFY_IS_EQUAL(lm.nfev(), 11);
  396. VERIFY_IS_EQUAL(lm.njev(), 10);
  397. // check norm^2
  398. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.5324382854E+00);
  399. // check x
  400. VERIFY_IS_APPROX(x[0], 1.0776351733E+00);
  401. VERIFY_IS_APPROX(x[1],-1.2269296921E-01);
  402. VERIFY_IS_APPROX(x[2], 4.0863750610E-03);
  403. VERIFY_IS_APPROX(x[3],-1.426264e-06); // shoulde be : -1.4262662514E-06
  404. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  405. VERIFY_IS_APPROX(x[5], 2.4053735503E-04);
  406. VERIFY_IS_APPROX(x[6],-1.2314450199E-07);
  407. /*
  408. * Second try
  409. */
  410. x<< .1, -.1, .005, -.000001, -.005, .0001, -.0000001;
  411. // do the computation
  412. info = lm.minimize(x);
  413. // check return value
  414. VERIFY_IS_EQUAL(info, 1);
  415. // VERIFY_IS_EQUAL(lm.nfev(), 11);
  416. VERIFY_IS_EQUAL(lm.njev(), 10);
  417. // check norm^2
  418. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.5324382854E+00);
  419. // check x
  420. VERIFY_IS_APPROX(x[0], 1.077640); // should be : 1.0776351733E+00
  421. VERIFY_IS_APPROX(x[1], -0.1226933); // should be : -1.2269296921E-01
  422. VERIFY_IS_APPROX(x[2], 0.004086383); // should be : 4.0863750610E-03
  423. VERIFY_IS_APPROX(x[3], -1.426277e-06); // shoulde be : -1.4262662514E-06
  424. VERIFY_IS_APPROX(x[4],-5.7609940901E-03);
  425. VERIFY_IS_APPROX(x[5], 0.00024053772); // should be : 2.4053735503E-04
  426. VERIFY_IS_APPROX(x[6], -1.231450e-07); // should be : -1.2314450199E-07
  427. }
  428. struct misra1d_functor : DenseFunctor<double>
  429. {
  430. misra1d_functor(void) : DenseFunctor<double>(2,14) {}
  431. static const double x[14];
  432. static const double y[14];
  433. int operator()(const VectorXd &b, VectorXd &fvec)
  434. {
  435. assert(b.size()==2);
  436. assert(fvec.size()==14);
  437. for(int i=0; i<14; i++) {
  438. fvec[i] = b[0]*b[1]*x[i]/(1.+b[1]*x[i]) - y[i];
  439. }
  440. return 0;
  441. }
  442. int df(const VectorXd &b, MatrixXd &fjac)
  443. {
  444. assert(b.size()==2);
  445. assert(fjac.rows()==14);
  446. assert(fjac.cols()==2);
  447. for(int i=0; i<14; i++) {
  448. double den = 1.+b[1]*x[i];
  449. fjac(i,0) = b[1]*x[i] / den;
  450. fjac(i,1) = b[0]*x[i]*(den-b[1]*x[i])/den/den;
  451. }
  452. return 0;
  453. }
  454. };
  455. 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};
  456. 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};
  457. // http://www.itl.nist.gov/div898/strd/nls/data/misra1d.shtml
  458. void testNistMisra1d(void)
  459. {
  460. const int n=2;
  461. int info;
  462. VectorXd x(n);
  463. /*
  464. * First try
  465. */
  466. x<< 500., 0.0001;
  467. // do the computation
  468. misra1d_functor functor;
  469. LevenbergMarquardt<misra1d_functor> lm(functor);
  470. info = lm.minimize(x);
  471. // check return value
  472. VERIFY_IS_EQUAL(info, 1);
  473. VERIFY_IS_EQUAL(lm.nfev(), 9);
  474. VERIFY_IS_EQUAL(lm.njev(), 7);
  475. // check norm^2
  476. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6419295283E-02);
  477. // check x
  478. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  479. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  480. /*
  481. * Second try
  482. */
  483. x<< 450., 0.0003;
  484. // do the computation
  485. info = lm.minimize(x);
  486. // check return value
  487. VERIFY_IS_EQUAL(info, 1);
  488. VERIFY_IS_EQUAL(lm.nfev(), 4);
  489. VERIFY_IS_EQUAL(lm.njev(), 3);
  490. // check norm^2
  491. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6419295283E-02);
  492. // check x
  493. VERIFY_IS_APPROX(x[0], 4.3736970754E+02);
  494. VERIFY_IS_APPROX(x[1], 3.0227324449E-04);
  495. }
  496. struct lanczos1_functor : DenseFunctor<double>
  497. {
  498. lanczos1_functor(void) : DenseFunctor<double>(6,24) {}
  499. static const double x[24];
  500. static const double y[24];
  501. int operator()(const VectorXd &b, VectorXd &fvec)
  502. {
  503. assert(b.size()==6);
  504. assert(fvec.size()==24);
  505. for(int i=0; i<24; i++)
  506. 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];
  507. return 0;
  508. }
  509. int df(const VectorXd &b, MatrixXd &fjac)
  510. {
  511. assert(b.size()==6);
  512. assert(fjac.rows()==24);
  513. assert(fjac.cols()==6);
  514. for(int i=0; i<24; i++) {
  515. fjac(i,0) = exp(-b[1]*x[i]);
  516. fjac(i,1) = -b[0]*x[i]*exp(-b[1]*x[i]);
  517. fjac(i,2) = exp(-b[3]*x[i]);
  518. fjac(i,3) = -b[2]*x[i]*exp(-b[3]*x[i]);
  519. fjac(i,4) = exp(-b[5]*x[i]);
  520. fjac(i,5) = -b[4]*x[i]*exp(-b[5]*x[i]);
  521. }
  522. return 0;
  523. }
  524. };
  525. 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 };
  526. 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 };
  527. // http://www.itl.nist.gov/div898/strd/nls/data/lanczos1.shtml
  528. void testNistLanczos1(void)
  529. {
  530. const int n=6;
  531. LevenbergMarquardtSpace::Status info;
  532. VectorXd x(n);
  533. /*
  534. * First try
  535. */
  536. x<< 1.2, 0.3, 5.6, 5.5, 6.5, 7.6;
  537. // do the computation
  538. lanczos1_functor functor;
  539. LevenbergMarquardt<lanczos1_functor> lm(functor);
  540. info = lm.minimize(x);
  541. // check return value
  542. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeErrorTooSmall);
  543. VERIFY_IS_EQUAL(lm.nfev(), 79);
  544. VERIFY_IS_EQUAL(lm.njev(), 72);
  545. // check norm^2
  546. // VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.430899764097e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  547. // check x
  548. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  549. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  550. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  551. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  552. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  553. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  554. /*
  555. * Second try
  556. */
  557. x<< 0.5, 0.7, 3.6, 4.2, 4., 6.3;
  558. // do the computation
  559. info = lm.minimize(x);
  560. // check return value
  561. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeErrorTooSmall);
  562. VERIFY_IS_EQUAL(lm.nfev(), 9);
  563. VERIFY_IS_EQUAL(lm.njev(), 8);
  564. // check norm^2
  565. // VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.428595533845e-25); // should be 1.4307867721E-25, but nist results are on 128-bit floats
  566. // check x
  567. VERIFY_IS_APPROX(x[0], 9.5100000027E-02);
  568. VERIFY_IS_APPROX(x[1], 1.0000000001E+00);
  569. VERIFY_IS_APPROX(x[2], 8.6070000013E-01);
  570. VERIFY_IS_APPROX(x[3], 3.0000000002E+00);
  571. VERIFY_IS_APPROX(x[4], 1.5575999998E+00);
  572. VERIFY_IS_APPROX(x[5], 5.0000000001E+00);
  573. }
  574. struct rat42_functor : DenseFunctor<double>
  575. {
  576. rat42_functor(void) : DenseFunctor<double>(3,9) {}
  577. static const double x[9];
  578. static const double y[9];
  579. int operator()(const VectorXd &b, VectorXd &fvec)
  580. {
  581. assert(b.size()==3);
  582. assert(fvec.size()==9);
  583. for(int i=0; i<9; i++) {
  584. fvec[i] = b[0] / (1.+exp(b[1]-b[2]*x[i])) - y[i];
  585. }
  586. return 0;
  587. }
  588. int df(const VectorXd &b, MatrixXd &fjac)
  589. {
  590. assert(b.size()==3);
  591. assert(fjac.rows()==9);
  592. assert(fjac.cols()==3);
  593. for(int i=0; i<9; i++) {
  594. double e = exp(b[1]-b[2]*x[i]);
  595. fjac(i,0) = 1./(1.+e);
  596. fjac(i,1) = -b[0]*e/(1.+e)/(1.+e);
  597. fjac(i,2) = +b[0]*e*x[i]/(1.+e)/(1.+e);
  598. }
  599. return 0;
  600. }
  601. };
  602. const double rat42_functor::x[9] = { 9.000E0, 14.000E0, 21.000E0, 28.000E0, 42.000E0, 57.000E0, 63.000E0, 70.000E0, 79.000E0 };
  603. const double rat42_functor::y[9] = { 8.930E0 ,10.800E0 ,18.590E0 ,22.330E0 ,39.350E0 ,56.110E0 ,61.730E0 ,64.620E0 ,67.080E0 };
  604. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky2.shtml
  605. void testNistRat42(void)
  606. {
  607. const int n=3;
  608. LevenbergMarquardtSpace::Status info;
  609. VectorXd x(n);
  610. /*
  611. * First try
  612. */
  613. x<< 100., 1., 0.1;
  614. // do the computation
  615. rat42_functor functor;
  616. LevenbergMarquardt<rat42_functor> lm(functor);
  617. info = lm.minimize(x);
  618. // check return value
  619. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeReductionTooSmall);
  620. VERIFY_IS_EQUAL(lm.nfev(), 10);
  621. VERIFY_IS_EQUAL(lm.njev(), 8);
  622. // check norm^2
  623. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.0565229338E+00);
  624. // check x
  625. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  626. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  627. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  628. /*
  629. * Second try
  630. */
  631. x<< 75., 2.5, 0.07;
  632. // do the computation
  633. info = lm.minimize(x);
  634. // check return value
  635. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeReductionTooSmall);
  636. VERIFY_IS_EQUAL(lm.nfev(), 6);
  637. VERIFY_IS_EQUAL(lm.njev(), 5);
  638. // check norm^2
  639. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.0565229338E+00);
  640. // check x
  641. VERIFY_IS_APPROX(x[0], 7.2462237576E+01);
  642. VERIFY_IS_APPROX(x[1], 2.6180768402E+00);
  643. VERIFY_IS_APPROX(x[2], 6.7359200066E-02);
  644. }
  645. struct MGH10_functor : DenseFunctor<double>
  646. {
  647. MGH10_functor(void) : DenseFunctor<double>(3,16) {}
  648. static const double x[16];
  649. static const double y[16];
  650. int operator()(const VectorXd &b, VectorXd &fvec)
  651. {
  652. assert(b.size()==3);
  653. assert(fvec.size()==16);
  654. for(int i=0; i<16; i++)
  655. fvec[i] = b[0] * exp(b[1]/(x[i]+b[2])) - y[i];
  656. return 0;
  657. }
  658. int df(const VectorXd &b, MatrixXd &fjac)
  659. {
  660. assert(b.size()==3);
  661. assert(fjac.rows()==16);
  662. assert(fjac.cols()==3);
  663. for(int i=0; i<16; i++) {
  664. double factor = 1./(x[i]+b[2]);
  665. double e = exp(b[1]*factor);
  666. fjac(i,0) = e;
  667. fjac(i,1) = b[0]*factor*e;
  668. fjac(i,2) = -b[1]*b[0]*factor*factor*e;
  669. }
  670. return 0;
  671. }
  672. };
  673. 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 };
  674. 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 };
  675. // http://www.itl.nist.gov/div898/strd/nls/data/mgh10.shtml
  676. void testNistMGH10(void)
  677. {
  678. const int n=3;
  679. LevenbergMarquardtSpace::Status info;
  680. VectorXd x(n);
  681. /*
  682. * First try
  683. */
  684. x<< 2., 400000., 25000.;
  685. // do the computation
  686. MGH10_functor functor;
  687. LevenbergMarquardt<MGH10_functor> lm(functor);
  688. info = lm.minimize(x);
  689. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeErrorTooSmall);
  690. // check norm^2
  691. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7945855171E+01);
  692. // check x
  693. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  694. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  695. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  696. // check return value
  697. //VERIFY_IS_EQUAL(info, 1);
  698. VERIFY_IS_EQUAL(lm.nfev(), 284 );
  699. VERIFY_IS_EQUAL(lm.njev(), 249 );
  700. /*
  701. * Second try
  702. */
  703. x<< 0.02, 4000., 250.;
  704. // do the computation
  705. info = lm.minimize(x);
  706. VERIFY_IS_EQUAL(info, LevenbergMarquardtSpace::RelativeReductionTooSmall);
  707. // check norm^2
  708. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7945855171E+01);
  709. // check x
  710. VERIFY_IS_APPROX(x[0], 5.6096364710E-03);
  711. VERIFY_IS_APPROX(x[1], 6.1813463463E+03);
  712. VERIFY_IS_APPROX(x[2], 3.4522363462E+02);
  713. // check return value
  714. //VERIFY_IS_EQUAL(info, 1);
  715. VERIFY_IS_EQUAL(lm.nfev(), 126);
  716. VERIFY_IS_EQUAL(lm.njev(), 116);
  717. }
  718. struct BoxBOD_functor : DenseFunctor<double>
  719. {
  720. BoxBOD_functor(void) : DenseFunctor<double>(2,6) {}
  721. static const double x[6];
  722. int operator()(const VectorXd &b, VectorXd &fvec)
  723. {
  724. static const double y[6] = { 109., 149., 149., 191., 213., 224. };
  725. assert(b.size()==2);
  726. assert(fvec.size()==6);
  727. for(int i=0; i<6; i++)
  728. fvec[i] = b[0]*(1.-exp(-b[1]*x[i])) - y[i];
  729. return 0;
  730. }
  731. int df(const VectorXd &b, MatrixXd &fjac)
  732. {
  733. assert(b.size()==2);
  734. assert(fjac.rows()==6);
  735. assert(fjac.cols()==2);
  736. for(int i=0; i<6; i++) {
  737. double e = exp(-b[1]*x[i]);
  738. fjac(i,0) = 1.-e;
  739. fjac(i,1) = b[0]*x[i]*e;
  740. }
  741. return 0;
  742. }
  743. };
  744. const double BoxBOD_functor::x[6] = { 1., 2., 3., 5., 7., 10. };
  745. // http://www.itl.nist.gov/div898/strd/nls/data/boxbod.shtml
  746. void testNistBoxBOD(void)
  747. {
  748. const int n=2;
  749. int info;
  750. VectorXd x(n);
  751. /*
  752. * First try
  753. */
  754. x<< 1., 1.;
  755. // do the computation
  756. BoxBOD_functor functor;
  757. LevenbergMarquardt<BoxBOD_functor> lm(functor);
  758. lm.setFtol(1.E6*NumTraits<double>::epsilon());
  759. lm.setXtol(1.E6*NumTraits<double>::epsilon());
  760. lm.setFactor(10);
  761. info = lm.minimize(x);
  762. // check norm^2
  763. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.1680088766E+03);
  764. // check x
  765. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  766. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  767. // check return value
  768. VERIFY_IS_EQUAL(info, 1);
  769. VERIFY(lm.nfev() < 31); // 31
  770. VERIFY(lm.njev() < 25); // 25
  771. /*
  772. * Second try
  773. */
  774. x<< 100., 0.75;
  775. // do the computation
  776. lm.resetParameters();
  777. lm.setFtol(NumTraits<double>::epsilon());
  778. lm.setXtol( NumTraits<double>::epsilon());
  779. info = lm.minimize(x);
  780. // check return value
  781. VERIFY_IS_EQUAL(info, 1);
  782. VERIFY_IS_EQUAL(lm.nfev(), 16 );
  783. VERIFY_IS_EQUAL(lm.njev(), 15 );
  784. // check norm^2
  785. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.1680088766E+03);
  786. // check x
  787. VERIFY_IS_APPROX(x[0], 2.1380940889E+02);
  788. VERIFY_IS_APPROX(x[1], 5.4723748542E-01);
  789. }
  790. struct MGH17_functor : DenseFunctor<double>
  791. {
  792. MGH17_functor(void) : DenseFunctor<double>(5,33) {}
  793. static const double x[33];
  794. static const double y[33];
  795. int operator()(const VectorXd &b, VectorXd &fvec)
  796. {
  797. assert(b.size()==5);
  798. assert(fvec.size()==33);
  799. for(int i=0; i<33; i++)
  800. fvec[i] = b[0] + b[1]*exp(-b[3]*x[i]) + b[2]*exp(-b[4]*x[i]) - y[i];
  801. return 0;
  802. }
  803. int df(const VectorXd &b, MatrixXd &fjac)
  804. {
  805. assert(b.size()==5);
  806. assert(fjac.rows()==33);
  807. assert(fjac.cols()==5);
  808. for(int i=0; i<33; i++) {
  809. fjac(i,0) = 1.;
  810. fjac(i,1) = exp(-b[3]*x[i]);
  811. fjac(i,2) = exp(-b[4]*x[i]);
  812. fjac(i,3) = -x[i]*b[1]*exp(-b[3]*x[i]);
  813. fjac(i,4) = -x[i]*b[2]*exp(-b[4]*x[i]);
  814. }
  815. return 0;
  816. }
  817. };
  818. 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 };
  819. 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 };
  820. // http://www.itl.nist.gov/div898/strd/nls/data/mgh17.shtml
  821. void testNistMGH17(void)
  822. {
  823. const int n=5;
  824. int info;
  825. VectorXd x(n);
  826. /*
  827. * First try
  828. */
  829. x<< 50., 150., -100., 1., 2.;
  830. // do the computation
  831. MGH17_functor functor;
  832. LevenbergMarquardt<MGH17_functor> lm(functor);
  833. lm.setFtol(NumTraits<double>::epsilon());
  834. lm.setXtol(NumTraits<double>::epsilon());
  835. lm.setMaxfev(1000);
  836. info = lm.minimize(x);
  837. // check norm^2
  838. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.4648946975E-05);
  839. // check x
  840. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  841. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  842. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  843. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  844. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  845. // check return value
  846. // VERIFY_IS_EQUAL(info, 2); //FIXME Use (lm.info() == Success)
  847. VERIFY(lm.nfev() < 700 ); // 602
  848. VERIFY(lm.njev() < 600 ); // 545
  849. /*
  850. * Second try
  851. */
  852. x<< 0.5 ,1.5 ,-1 ,0.01 ,0.02;
  853. // do the computation
  854. lm.resetParameters();
  855. info = lm.minimize(x);
  856. // check return value
  857. VERIFY_IS_EQUAL(info, 1);
  858. VERIFY_IS_EQUAL(lm.nfev(), 18);
  859. VERIFY_IS_EQUAL(lm.njev(), 15);
  860. // check norm^2
  861. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.4648946975E-05);
  862. // check x
  863. VERIFY_IS_APPROX(x[0], 3.7541005211E-01);
  864. VERIFY_IS_APPROX(x[1], 1.9358469127E+00);
  865. VERIFY_IS_APPROX(x[2], -1.4646871366E+00);
  866. VERIFY_IS_APPROX(x[3], 1.2867534640E-02);
  867. VERIFY_IS_APPROX(x[4], 2.2122699662E-02);
  868. }
  869. struct MGH09_functor : DenseFunctor<double>
  870. {
  871. MGH09_functor(void) : DenseFunctor<double>(4,11) {}
  872. static const double _x[11];
  873. static const double y[11];
  874. int operator()(const VectorXd &b, VectorXd &fvec)
  875. {
  876. assert(b.size()==4);
  877. assert(fvec.size()==11);
  878. for(int i=0; i<11; i++) {
  879. double x = _x[i], xx=x*x;
  880. fvec[i] = b[0]*(xx+x*b[1])/(xx+x*b[2]+b[3]) - y[i];
  881. }
  882. return 0;
  883. }
  884. int df(const VectorXd &b, MatrixXd &fjac)
  885. {
  886. assert(b.size()==4);
  887. assert(fjac.rows()==11);
  888. assert(fjac.cols()==4);
  889. for(int i=0; i<11; i++) {
  890. double x = _x[i], xx=x*x;
  891. double factor = 1./(xx+x*b[2]+b[3]);
  892. fjac(i,0) = (xx+x*b[1]) * factor;
  893. fjac(i,1) = b[0]*x* factor;
  894. fjac(i,2) = - b[0]*(xx+x*b[1]) * x * factor * factor;
  895. fjac(i,3) = - b[0]*(xx+x*b[1]) * factor * factor;
  896. }
  897. return 0;
  898. }
  899. };
  900. 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 };
  901. 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 };
  902. // http://www.itl.nist.gov/div898/strd/nls/data/mgh09.shtml
  903. void testNistMGH09(void)
  904. {
  905. const int n=4;
  906. int info;
  907. VectorXd x(n);
  908. /*
  909. * First try
  910. */
  911. x<< 25., 39, 41.5, 39.;
  912. // do the computation
  913. MGH09_functor functor;
  914. LevenbergMarquardt<MGH09_functor> lm(functor);
  915. lm.setMaxfev(1000);
  916. info = lm.minimize(x);
  917. // check norm^2
  918. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 3.0750560385E-04);
  919. // check x
  920. VERIFY_IS_APPROX(x[0], 0.1928077089); // should be 1.9280693458E-01
  921. VERIFY_IS_APPROX(x[1], 0.19126423573); // should be 1.9128232873E-01
  922. VERIFY_IS_APPROX(x[2], 0.12305309914); // should be 1.2305650693E-01
  923. VERIFY_IS_APPROX(x[3], 0.13605395375); // should be 1.3606233068E-01
  924. // check return value
  925. VERIFY_IS_EQUAL(info, 1);
  926. VERIFY(lm.nfev() < 510 ); // 490
  927. VERIFY(lm.njev() < 400 ); // 376
  928. /*
  929. * Second try
  930. */
  931. x<< 0.25, 0.39, 0.415, 0.39;
  932. // do the computation
  933. lm.resetParameters();
  934. info = lm.minimize(x);
  935. // check return value
  936. VERIFY_IS_EQUAL(info, 1);
  937. VERIFY_IS_EQUAL(lm.nfev(), 18);
  938. VERIFY_IS_EQUAL(lm.njev(), 16);
  939. // check norm^2
  940. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 3.0750560385E-04);
  941. // check x
  942. VERIFY_IS_APPROX(x[0], 0.19280781); // should be 1.9280693458E-01
  943. VERIFY_IS_APPROX(x[1], 0.19126265); // should be 1.9128232873E-01
  944. VERIFY_IS_APPROX(x[2], 0.12305280); // should be 1.2305650693E-01
  945. VERIFY_IS_APPROX(x[3], 0.13605322); // should be 1.3606233068E-01
  946. }
  947. struct Bennett5_functor : DenseFunctor<double>
  948. {
  949. Bennett5_functor(void) : DenseFunctor<double>(3,154) {}
  950. static const double x[154];
  951. static const double y[154];
  952. int operator()(const VectorXd &b, VectorXd &fvec)
  953. {
  954. assert(b.size()==3);
  955. assert(fvec.size()==154);
  956. for(int i=0; i<154; i++)
  957. fvec[i] = b[0]* pow(b[1]+x[i],-1./b[2]) - y[i];
  958. return 0;
  959. }
  960. int df(const VectorXd &b, MatrixXd &fjac)
  961. {
  962. assert(b.size()==3);
  963. assert(fjac.rows()==154);
  964. assert(fjac.cols()==3);
  965. for(int i=0; i<154; i++) {
  966. double e = pow(b[1]+x[i],-1./b[2]);
  967. fjac(i,0) = e;
  968. fjac(i,1) = - b[0]*e/b[2]/(b[1]+x[i]);
  969. fjac(i,2) = b[0]*e*log(b[1]+x[i])/b[2]/b[2];
  970. }
  971. return 0;
  972. }
  973. };
  974. 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,
  975. 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 };
  976. 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
  977. ,-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 ,
  978. -31.827299E0 ,-31.821600E0 ,-31.818701E0 ,-31.812901E0 ,-31.809999E0 ,-31.807100E0 ,-31.801300E0 ,-31.798401E0 ,-31.795500E0 ,-31.789700E0 ,-31.786800E0 };
  979. // http://www.itl.nist.gov/div898/strd/nls/data/bennett5.shtml
  980. void testNistBennett5(void)
  981. {
  982. const int n=3;
  983. int info;
  984. VectorXd x(n);
  985. /*
  986. * First try
  987. */
  988. x<< -2000., 50., 0.8;
  989. // do the computation
  990. Bennett5_functor functor;
  991. LevenbergMarquardt<Bennett5_functor> lm(functor);
  992. lm.setMaxfev(1000);
  993. info = lm.minimize(x);
  994. // check return value
  995. VERIFY_IS_EQUAL(info, 1);
  996. VERIFY_IS_EQUAL(lm.nfev(), 758);
  997. VERIFY_IS_EQUAL(lm.njev(), 744);
  998. // check norm^2
  999. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.2404744073E-04);
  1000. // check x
  1001. VERIFY_IS_APPROX(x[0], -2.5235058043E+03);
  1002. VERIFY_IS_APPROX(x[1], 4.6736564644E+01);
  1003. VERIFY_IS_APPROX(x[2], 9.3218483193E-01);
  1004. /*
  1005. * Second try
  1006. */
  1007. x<< -1500., 45., 0.85;
  1008. // do the computation
  1009. lm.resetParameters();
  1010. info = lm.minimize(x);
  1011. // check return value
  1012. VERIFY_IS_EQUAL(info, 1);
  1013. VERIFY_IS_EQUAL(lm.nfev(), 203);
  1014. VERIFY_IS_EQUAL(lm.njev(), 192);
  1015. // check norm^2
  1016. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.2404744073E-04);
  1017. // check x
  1018. VERIFY_IS_APPROX(x[0], -2523.3007865); // should be -2.5235058043E+03
  1019. VERIFY_IS_APPROX(x[1], 46.735705771); // should be 4.6736564644E+01);
  1020. VERIFY_IS_APPROX(x[2], 0.93219881891); // should be 9.3218483193E-01);
  1021. }
  1022. struct thurber_functor : DenseFunctor<double>
  1023. {
  1024. thurber_functor(void) : DenseFunctor<double>(7,37) {}
  1025. static const double _x[37];
  1026. static const double _y[37];
  1027. int operator()(const VectorXd &b, VectorXd &fvec)
  1028. {
  1029. // int called=0; printf("call hahn1_functor with iflag=%d, called=%d\n", iflag, called); if (iflag==1) called++;
  1030. assert(b.size()==7);
  1031. assert(fvec.size()==37);
  1032. for(int i=0; i<37; i++) {
  1033. double x=_x[i], xx=x*x, xxx=xx*x;
  1034. 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];
  1035. }
  1036. return 0;
  1037. }
  1038. int df(const VectorXd &b, MatrixXd &fjac)
  1039. {
  1040. assert(b.size()==7);
  1041. assert(fjac.rows()==37);
  1042. assert(fjac.cols()==7);
  1043. for(int i=0; i<37; i++) {
  1044. double x=_x[i], xx=x*x, xxx=xx*x;
  1045. double fact = 1./(1.+b[4]*x+b[5]*xx+b[6]*xxx);
  1046. fjac(i,0) = 1.*fact;
  1047. fjac(i,1) = x*fact;
  1048. fjac(i,2) = xx*fact;
  1049. fjac(i,3) = xxx*fact;
  1050. fact = - (b[0]+b[1]*x+b[2]*xx+b[3]*xxx) * fact * fact;
  1051. fjac(i,4) = x*fact;
  1052. fjac(i,5) = xx*fact;
  1053. fjac(i,6) = xxx*fact;
  1054. }
  1055. return 0;
  1056. }
  1057. };
  1058. 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 };
  1059. 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};
  1060. // http://www.itl.nist.gov/div898/strd/nls/data/thurber.shtml
  1061. void testNistThurber(void)
  1062. {
  1063. const int n=7;
  1064. int info;
  1065. VectorXd x(n);
  1066. /*
  1067. * First try
  1068. */
  1069. x<< 1000 ,1000 ,400 ,40 ,0.7,0.3,0.0 ;
  1070. // do the computation
  1071. thurber_functor functor;
  1072. LevenbergMarquardt<thurber_functor> lm(functor);
  1073. lm.setFtol(1.E4*NumTraits<double>::epsilon());
  1074. lm.setXtol(1.E4*NumTraits<double>::epsilon());
  1075. info = lm.minimize(x);
  1076. // check return value
  1077. VERIFY_IS_EQUAL(info, 1);
  1078. VERIFY_IS_EQUAL(lm.nfev(), 39);
  1079. VERIFY_IS_EQUAL(lm.njev(), 36);
  1080. // check norm^2
  1081. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6427082397E+03);
  1082. // check x
  1083. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1084. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1085. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1086. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1087. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1088. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1089. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1090. /*
  1091. * Second try
  1092. */
  1093. x<< 1300 ,1500 ,500 ,75 ,1 ,0.4 ,0.05 ;
  1094. // do the computation
  1095. lm.resetParameters();
  1096. lm.setFtol(1.E4*NumTraits<double>::epsilon());
  1097. lm.setXtol(1.E4*NumTraits<double>::epsilon());
  1098. info = lm.minimize(x);
  1099. // check return value
  1100. VERIFY_IS_EQUAL(info, 1);
  1101. VERIFY_IS_EQUAL(lm.nfev(), 29);
  1102. VERIFY_IS_EQUAL(lm.njev(), 28);
  1103. // check norm^2
  1104. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 5.6427082397E+03);
  1105. // check x
  1106. VERIFY_IS_APPROX(x[0], 1.2881396800E+03);
  1107. VERIFY_IS_APPROX(x[1], 1.4910792535E+03);
  1108. VERIFY_IS_APPROX(x[2], 5.8323836877E+02);
  1109. VERIFY_IS_APPROX(x[3], 7.5416644291E+01);
  1110. VERIFY_IS_APPROX(x[4], 9.6629502864E-01);
  1111. VERIFY_IS_APPROX(x[5], 3.9797285797E-01);
  1112. VERIFY_IS_APPROX(x[6], 4.9727297349E-02);
  1113. }
  1114. struct rat43_functor : DenseFunctor<double>
  1115. {
  1116. rat43_functor(void) : DenseFunctor<double>(4,15) {}
  1117. static const double x[15];
  1118. static const double y[15];
  1119. int operator()(const VectorXd &b, VectorXd &fvec)
  1120. {
  1121. assert(b.size()==4);
  1122. assert(fvec.size()==15);
  1123. for(int i=0; i<15; i++)
  1124. fvec[i] = b[0] * pow(1.+exp(b[1]-b[2]*x[i]),-1./b[3]) - y[i];
  1125. return 0;
  1126. }
  1127. int df(const VectorXd &b, MatrixXd &fjac)
  1128. {
  1129. assert(b.size()==4);
  1130. assert(fjac.rows()==15);
  1131. assert(fjac.cols()==4);
  1132. for(int i=0; i<15; i++) {
  1133. double e = exp(b[1]-b[2]*x[i]);
  1134. double power = -1./b[3];
  1135. fjac(i,0) = pow(1.+e, power);
  1136. fjac(i,1) = power*b[0]*e*pow(1.+e, power-1.);
  1137. fjac(i,2) = -power*b[0]*e*x[i]*pow(1.+e, power-1.);
  1138. fjac(i,3) = b[0]*power*power*log(1.+e)*pow(1.+e, power);
  1139. }
  1140. return 0;
  1141. }
  1142. };
  1143. const double rat43_functor::x[15] = { 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14., 15. };
  1144. 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 };
  1145. // http://www.itl.nist.gov/div898/strd/nls/data/ratkowsky3.shtml
  1146. void testNistRat43(void)
  1147. {
  1148. const int n=4;
  1149. int info;
  1150. VectorXd x(n);
  1151. /*
  1152. * First try
  1153. */
  1154. x<< 100., 10., 1., 1.;
  1155. // do the computation
  1156. rat43_functor functor;
  1157. LevenbergMarquardt<rat43_functor> lm(functor);
  1158. lm.setFtol(1.E6*NumTraits<double>::epsilon());
  1159. lm.setXtol(1.E6*NumTraits<double>::epsilon());
  1160. info = lm.minimize(x);
  1161. // check return value
  1162. VERIFY_IS_EQUAL(info, 1);
  1163. VERIFY_IS_EQUAL(lm.nfev(), 27);
  1164. VERIFY_IS_EQUAL(lm.njev(), 20);
  1165. // check norm^2
  1166. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7864049080E+03);
  1167. // check x
  1168. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1169. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1170. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1171. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1172. /*
  1173. * Second try
  1174. */
  1175. x<< 700., 5., 0.75, 1.3;
  1176. // do the computation
  1177. lm.resetParameters();
  1178. lm.setFtol(1.E5*NumTraits<double>::epsilon());
  1179. lm.setXtol(1.E5*NumTraits<double>::epsilon());
  1180. info = lm.minimize(x);
  1181. // check return value
  1182. VERIFY_IS_EQUAL(info, 1);
  1183. VERIFY_IS_EQUAL(lm.nfev(), 9);
  1184. VERIFY_IS_EQUAL(lm.njev(), 8);
  1185. // check norm^2
  1186. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 8.7864049080E+03);
  1187. // check x
  1188. VERIFY_IS_APPROX(x[0], 6.9964151270E+02);
  1189. VERIFY_IS_APPROX(x[1], 5.2771253025E+00);
  1190. VERIFY_IS_APPROX(x[2], 7.5962938329E-01);
  1191. VERIFY_IS_APPROX(x[3], 1.2792483859E+00);
  1192. }
  1193. struct eckerle4_functor : DenseFunctor<double>
  1194. {
  1195. eckerle4_functor(void) : DenseFunctor<double>(3,35) {}
  1196. static const double x[35];
  1197. static const double y[35];
  1198. int operator()(const VectorXd &b, VectorXd &fvec)
  1199. {
  1200. assert(b.size()==3);
  1201. assert(fvec.size()==35);
  1202. for(int i=0; i<35; i++)
  1203. fvec[i] = b[0]/b[1] * exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/(b[1]*b[1])) - y[i];
  1204. return 0;
  1205. }
  1206. int df(const VectorXd &b, MatrixXd &fjac)
  1207. {
  1208. assert(b.size()==3);
  1209. assert(fjac.rows()==35);
  1210. assert(fjac.cols()==3);
  1211. for(int i=0; i<35; i++) {
  1212. double b12 = b[1]*b[1];
  1213. double e = exp(-0.5*(x[i]-b[2])*(x[i]-b[2])/b12);
  1214. fjac(i,0) = e / b[1];
  1215. fjac(i,1) = ((x[i]-b[2])*(x[i]-b[2])/b12-1.) * b[0]*e/b12;
  1216. fjac(i,2) = (x[i]-b[2])*e*b[0]/b[1]/b12;
  1217. }
  1218. return 0;
  1219. }
  1220. };
  1221. 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};
  1222. 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 };
  1223. // http://www.itl.nist.gov/div898/strd/nls/data/eckerle4.shtml
  1224. void testNistEckerle4(void)
  1225. {
  1226. const int n=3;
  1227. int info;
  1228. VectorXd x(n);
  1229. /*
  1230. * First try
  1231. */
  1232. x<< 1., 10., 500.;
  1233. // do the computation
  1234. eckerle4_functor functor;
  1235. LevenbergMarquardt<eckerle4_functor> lm(functor);
  1236. info = lm.minimize(x);
  1237. // check return value
  1238. VERIFY_IS_EQUAL(info, 1);
  1239. VERIFY_IS_EQUAL(lm.nfev(), 18);
  1240. VERIFY_IS_EQUAL(lm.njev(), 15);
  1241. // check norm^2
  1242. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.4635887487E-03);
  1243. // check x
  1244. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1245. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1246. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1247. /*
  1248. * Second try
  1249. */
  1250. x<< 1.5, 5., 450.;
  1251. // do the computation
  1252. info = lm.minimize(x);
  1253. // check return value
  1254. VERIFY_IS_EQUAL(info, 1);
  1255. VERIFY_IS_EQUAL(lm.nfev(), 7);
  1256. VERIFY_IS_EQUAL(lm.njev(), 6);
  1257. // check norm^2
  1258. VERIFY_IS_APPROX(lm.fvec().squaredNorm(), 1.4635887487E-03);
  1259. // check x
  1260. VERIFY_IS_APPROX(x[0], 1.5543827178);
  1261. VERIFY_IS_APPROX(x[1], 4.0888321754);
  1262. VERIFY_IS_APPROX(x[2], 4.5154121844E+02);
  1263. }
  1264. void test_levenberg_marquardt()
  1265. {
  1266. // Tests using the examples provided by (c)minpack
  1267. CALL_SUBTEST(testLmder1());
  1268. CALL_SUBTEST(testLmder());
  1269. CALL_SUBTEST(testLmdif1());
  1270. // CALL_SUBTEST(testLmstr1());
  1271. // CALL_SUBTEST(testLmstr());
  1272. CALL_SUBTEST(testLmdif());
  1273. // NIST tests, level of difficulty = "Lower"
  1274. CALL_SUBTEST(testNistMisra1a());
  1275. CALL_SUBTEST(testNistChwirut2());
  1276. // NIST tests, level of difficulty = "Average"
  1277. CALL_SUBTEST(testNistHahn1());
  1278. CALL_SUBTEST(testNistMisra1d());
  1279. CALL_SUBTEST(testNistMGH17());
  1280. CALL_SUBTEST(testNistLanczos1());
  1281. // // NIST tests, level of difficulty = "Higher"
  1282. CALL_SUBTEST(testNistRat42());
  1283. CALL_SUBTEST(testNistMGH10());
  1284. CALL_SUBTEST(testNistBoxBOD());
  1285. // CALL_SUBTEST(testNistMGH09());
  1286. CALL_SUBTEST(testNistBennett5());
  1287. CALL_SUBTEST(testNistThurber());
  1288. CALL_SUBTEST(testNistRat43());
  1289. CALL_SUBTEST(testNistEckerle4());
  1290. }