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) 2012 Desire Nuentsa <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2012 Gael Guennebaud <[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 #include <iostream>
12*bf2c3715SXin Li #include <fstream>
13*bf2c3715SXin Li #include <iomanip>
14*bf2c3715SXin Li
15*bf2c3715SXin Li #include "main.h"
16*bf2c3715SXin Li #include <Eigen/LevenbergMarquardt>
17*bf2c3715SXin Li using namespace std;
18*bf2c3715SXin Li using namespace Eigen;
19*bf2c3715SXin Li
20*bf2c3715SXin Li template<typename Scalar>
21*bf2c3715SXin Li struct DenseLM : DenseFunctor<Scalar>
22*bf2c3715SXin Li {
23*bf2c3715SXin Li typedef DenseFunctor<Scalar> Base;
24*bf2c3715SXin Li typedef typename Base::JacobianType JacobianType;
25*bf2c3715SXin Li typedef Matrix<Scalar,Dynamic,1> VectorType;
26*bf2c3715SXin Li
DenseLMDenseLM27*bf2c3715SXin Li DenseLM(int n, int m) : DenseFunctor<Scalar>(n,m)
28*bf2c3715SXin Li { }
29*bf2c3715SXin Li
modelDenseLM30*bf2c3715SXin Li VectorType model(const VectorType& uv, VectorType& x)
31*bf2c3715SXin Li {
32*bf2c3715SXin Li VectorType y; // Should change to use expression template
33*bf2c3715SXin Li int m = Base::values();
34*bf2c3715SXin Li int n = Base::inputs();
35*bf2c3715SXin Li eigen_assert(uv.size()%2 == 0);
36*bf2c3715SXin Li eigen_assert(uv.size() == n);
37*bf2c3715SXin Li eigen_assert(x.size() == m);
38*bf2c3715SXin Li y.setZero(m);
39*bf2c3715SXin Li int half = n/2;
40*bf2c3715SXin Li VectorBlock<const VectorType> u(uv, 0, half);
41*bf2c3715SXin Li VectorBlock<const VectorType> v(uv, half, half);
42*bf2c3715SXin Li for (int j = 0; j < m; j++)
43*bf2c3715SXin Li {
44*bf2c3715SXin Li for (int i = 0; i < half; i++)
45*bf2c3715SXin Li y(j) += u(i)*std::exp(-(x(j)-i)*(x(j)-i)/(v(i)*v(i)));
46*bf2c3715SXin Li }
47*bf2c3715SXin Li return y;
48*bf2c3715SXin Li
49*bf2c3715SXin Li }
initPointsDenseLM50*bf2c3715SXin Li void initPoints(VectorType& uv_ref, VectorType& x)
51*bf2c3715SXin Li {
52*bf2c3715SXin Li m_x = x;
53*bf2c3715SXin Li m_y = this->model(uv_ref, x);
54*bf2c3715SXin Li }
55*bf2c3715SXin Li
operator ()DenseLM56*bf2c3715SXin Li int operator()(const VectorType& uv, VectorType& fvec)
57*bf2c3715SXin Li {
58*bf2c3715SXin Li
59*bf2c3715SXin Li int m = Base::values();
60*bf2c3715SXin Li int n = Base::inputs();
61*bf2c3715SXin Li eigen_assert(uv.size()%2 == 0);
62*bf2c3715SXin Li eigen_assert(uv.size() == n);
63*bf2c3715SXin Li eigen_assert(fvec.size() == m);
64*bf2c3715SXin Li int half = n/2;
65*bf2c3715SXin Li VectorBlock<const VectorType> u(uv, 0, half);
66*bf2c3715SXin Li VectorBlock<const VectorType> v(uv, half, half);
67*bf2c3715SXin Li for (int j = 0; j < m; j++)
68*bf2c3715SXin Li {
69*bf2c3715SXin Li fvec(j) = m_y(j);
70*bf2c3715SXin Li for (int i = 0; i < half; i++)
71*bf2c3715SXin Li {
72*bf2c3715SXin Li fvec(j) -= u(i) *std::exp(-(m_x(j)-i)*(m_x(j)-i)/(v(i)*v(i)));
73*bf2c3715SXin Li }
74*bf2c3715SXin Li }
75*bf2c3715SXin Li
76*bf2c3715SXin Li return 0;
77*bf2c3715SXin Li }
dfDenseLM78*bf2c3715SXin Li int df(const VectorType& uv, JacobianType& fjac)
79*bf2c3715SXin Li {
80*bf2c3715SXin Li int m = Base::values();
81*bf2c3715SXin Li int n = Base::inputs();
82*bf2c3715SXin Li eigen_assert(n == uv.size());
83*bf2c3715SXin Li eigen_assert(fjac.rows() == m);
84*bf2c3715SXin Li eigen_assert(fjac.cols() == n);
85*bf2c3715SXin Li int half = n/2;
86*bf2c3715SXin Li VectorBlock<const VectorType> u(uv, 0, half);
87*bf2c3715SXin Li VectorBlock<const VectorType> v(uv, half, half);
88*bf2c3715SXin Li for (int j = 0; j < m; j++)
89*bf2c3715SXin Li {
90*bf2c3715SXin Li for (int i = 0; i < half; i++)
91*bf2c3715SXin Li {
92*bf2c3715SXin Li fjac.coeffRef(j,i) = -std::exp(-(m_x(j)-i)*(m_x(j)-i)/(v(i)*v(i)));
93*bf2c3715SXin Li fjac.coeffRef(j,i+half) = -2.*u(i)*(m_x(j)-i)*(m_x(j)-i)/(std::pow(v(i),3)) * std::exp(-(m_x(j)-i)*(m_x(j)-i)/(v(i)*v(i)));
94*bf2c3715SXin Li }
95*bf2c3715SXin Li }
96*bf2c3715SXin Li return 0;
97*bf2c3715SXin Li }
98*bf2c3715SXin Li VectorType m_x, m_y; //Data Points
99*bf2c3715SXin Li };
100*bf2c3715SXin Li
101*bf2c3715SXin Li template<typename FunctorType, typename VectorType>
test_minimizeLM(FunctorType & functor,VectorType & uv)102*bf2c3715SXin Li int test_minimizeLM(FunctorType& functor, VectorType& uv)
103*bf2c3715SXin Li {
104*bf2c3715SXin Li LevenbergMarquardt<FunctorType> lm(functor);
105*bf2c3715SXin Li LevenbergMarquardtSpace::Status info;
106*bf2c3715SXin Li
107*bf2c3715SXin Li info = lm.minimize(uv);
108*bf2c3715SXin Li
109*bf2c3715SXin Li VERIFY_IS_EQUAL(info, 1);
110*bf2c3715SXin Li //FIXME Check other parameters
111*bf2c3715SXin Li return info;
112*bf2c3715SXin Li }
113*bf2c3715SXin Li
114*bf2c3715SXin Li template<typename FunctorType, typename VectorType>
test_lmder(FunctorType & functor,VectorType & uv)115*bf2c3715SXin Li int test_lmder(FunctorType& functor, VectorType& uv)
116*bf2c3715SXin Li {
117*bf2c3715SXin Li typedef typename VectorType::Scalar Scalar;
118*bf2c3715SXin Li LevenbergMarquardtSpace::Status info;
119*bf2c3715SXin Li LevenbergMarquardt<FunctorType> lm(functor);
120*bf2c3715SXin Li info = lm.lmder1(uv);
121*bf2c3715SXin Li
122*bf2c3715SXin Li VERIFY_IS_EQUAL(info, 1);
123*bf2c3715SXin Li //FIXME Check other parameters
124*bf2c3715SXin Li return info;
125*bf2c3715SXin Li }
126*bf2c3715SXin Li
127*bf2c3715SXin Li template<typename FunctorType, typename VectorType>
test_minimizeSteps(FunctorType & functor,VectorType & uv)128*bf2c3715SXin Li int test_minimizeSteps(FunctorType& functor, VectorType& uv)
129*bf2c3715SXin Li {
130*bf2c3715SXin Li LevenbergMarquardtSpace::Status info;
131*bf2c3715SXin Li LevenbergMarquardt<FunctorType> lm(functor);
132*bf2c3715SXin Li info = lm.minimizeInit(uv);
133*bf2c3715SXin Li if (info==LevenbergMarquardtSpace::ImproperInputParameters)
134*bf2c3715SXin Li return info;
135*bf2c3715SXin Li do
136*bf2c3715SXin Li {
137*bf2c3715SXin Li info = lm.minimizeOneStep(uv);
138*bf2c3715SXin Li } while (info==LevenbergMarquardtSpace::Running);
139*bf2c3715SXin Li
140*bf2c3715SXin Li VERIFY_IS_EQUAL(info, 1);
141*bf2c3715SXin Li //FIXME Check other parameters
142*bf2c3715SXin Li return info;
143*bf2c3715SXin Li }
144*bf2c3715SXin Li
145*bf2c3715SXin Li template<typename T>
test_denseLM_T()146*bf2c3715SXin Li void test_denseLM_T()
147*bf2c3715SXin Li {
148*bf2c3715SXin Li typedef Matrix<T,Dynamic,1> VectorType;
149*bf2c3715SXin Li
150*bf2c3715SXin Li int inputs = 10;
151*bf2c3715SXin Li int values = 1000;
152*bf2c3715SXin Li DenseLM<T> dense_gaussian(inputs, values);
153*bf2c3715SXin Li VectorType uv(inputs),uv_ref(inputs);
154*bf2c3715SXin Li VectorType x(values);
155*bf2c3715SXin Li
156*bf2c3715SXin Li // Generate the reference solution
157*bf2c3715SXin Li uv_ref << -2, 1, 4 ,8, 6, 1.8, 1.2, 1.1, 1.9 , 3;
158*bf2c3715SXin Li
159*bf2c3715SXin Li //Generate the reference data points
160*bf2c3715SXin Li x.setRandom();
161*bf2c3715SXin Li x = 10*x;
162*bf2c3715SXin Li x.array() += 10;
163*bf2c3715SXin Li dense_gaussian.initPoints(uv_ref, x);
164*bf2c3715SXin Li
165*bf2c3715SXin Li // Generate the initial parameters
166*bf2c3715SXin Li VectorBlock<VectorType> u(uv, 0, inputs/2);
167*bf2c3715SXin Li VectorBlock<VectorType> v(uv, inputs/2, inputs/2);
168*bf2c3715SXin Li
169*bf2c3715SXin Li // Solve the optimization problem
170*bf2c3715SXin Li
171*bf2c3715SXin Li //Solve in one go
172*bf2c3715SXin Li u.setOnes(); v.setOnes();
173*bf2c3715SXin Li test_minimizeLM(dense_gaussian, uv);
174*bf2c3715SXin Li
175*bf2c3715SXin Li //Solve until the machine precision
176*bf2c3715SXin Li u.setOnes(); v.setOnes();
177*bf2c3715SXin Li test_lmder(dense_gaussian, uv);
178*bf2c3715SXin Li
179*bf2c3715SXin Li // Solve step by step
180*bf2c3715SXin Li v.setOnes(); u.setOnes();
181*bf2c3715SXin Li test_minimizeSteps(dense_gaussian, uv);
182*bf2c3715SXin Li
183*bf2c3715SXin Li }
184*bf2c3715SXin Li
EIGEN_DECLARE_TEST(denseLM)185*bf2c3715SXin Li EIGEN_DECLARE_TEST(denseLM)
186*bf2c3715SXin Li {
187*bf2c3715SXin Li CALL_SUBTEST_2(test_denseLM_T<double>());
188*bf2c3715SXin Li
189*bf2c3715SXin Li // CALL_SUBTEST_2(test_sparseLM_T<std::complex<double>());
190*bf2c3715SXin Li }
191