levenberg_marquardt.cpp 54 KB

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