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.

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