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