xref: /aosp_15_r20/external/eigen/test/sparseLM.cpp (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
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 #include <iostream>
11*bf2c3715SXin Li #include <fstream>
12*bf2c3715SXin Li #include <iomanip>
13*bf2c3715SXin Li 
14*bf2c3715SXin Li #include "main.h"
15*bf2c3715SXin Li #include <Eigen/LevenbergMarquardt>
16*bf2c3715SXin Li 
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 sparseGaussianTest : SparseFunctor<Scalar, int>
22*bf2c3715SXin Li {
23*bf2c3715SXin Li   typedef Matrix<Scalar,Dynamic,1> VectorType;
24*bf2c3715SXin Li   typedef SparseFunctor<Scalar,int> Base;
25*bf2c3715SXin Li   typedef typename Base::JacobianType JacobianType;
sparseGaussianTestsparseGaussianTest26*bf2c3715SXin Li   sparseGaussianTest(int inputs, int values) : SparseFunctor<Scalar,int>(inputs,values)
27*bf2c3715SXin Li   { }
28*bf2c3715SXin Li 
modelsparseGaussianTest29*bf2c3715SXin Li   VectorType model(const VectorType& uv, VectorType& x)
30*bf2c3715SXin Li   {
31*bf2c3715SXin Li     VectorType y; //Change this to use expression template
32*bf2c3715SXin Li     int m = Base::values();
33*bf2c3715SXin Li     int n = Base::inputs();
34*bf2c3715SXin Li     eigen_assert(uv.size()%2 == 0);
35*bf2c3715SXin Li     eigen_assert(uv.size() == n);
36*bf2c3715SXin Li     eigen_assert(x.size() == m);
37*bf2c3715SXin Li     y.setZero(m);
38*bf2c3715SXin Li     int half = n/2;
39*bf2c3715SXin Li     VectorBlock<const VectorType> u(uv, 0, half);
40*bf2c3715SXin Li     VectorBlock<const VectorType> v(uv, half, half);
41*bf2c3715SXin Li     Scalar coeff;
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       {
46*bf2c3715SXin Li         coeff = (x(j)-i)/v(i);
47*bf2c3715SXin Li         coeff *= coeff;
48*bf2c3715SXin Li         if (coeff < 1. && coeff > 0.)
49*bf2c3715SXin Li           y(j) += u(i)*std::pow((1-coeff), 2);
50*bf2c3715SXin Li       }
51*bf2c3715SXin Li     }
52*bf2c3715SXin Li     return y;
53*bf2c3715SXin Li   }
initPointssparseGaussianTest54*bf2c3715SXin Li   void initPoints(VectorType& uv_ref, VectorType& x)
55*bf2c3715SXin Li   {
56*bf2c3715SXin Li     m_x = x;
57*bf2c3715SXin Li     m_y = this->model(uv_ref,x);
58*bf2c3715SXin Li   }
operator ()sparseGaussianTest59*bf2c3715SXin Li   int operator()(const VectorType& uv, VectorType& fvec)
60*bf2c3715SXin Li   {
61*bf2c3715SXin Li     int m = Base::values();
62*bf2c3715SXin Li     int n = Base::inputs();
63*bf2c3715SXin Li     eigen_assert(uv.size()%2 == 0);
64*bf2c3715SXin Li     eigen_assert(uv.size() == n);
65*bf2c3715SXin Li     int half = n/2;
66*bf2c3715SXin Li     VectorBlock<const VectorType> u(uv, 0, half);
67*bf2c3715SXin Li     VectorBlock<const VectorType> v(uv, half, half);
68*bf2c3715SXin Li     fvec = m_y;
69*bf2c3715SXin Li     Scalar coeff;
70*bf2c3715SXin Li     for (int j = 0; j < m; j++)
71*bf2c3715SXin Li     {
72*bf2c3715SXin Li       for (int i = 0; i < half; i++)
73*bf2c3715SXin Li       {
74*bf2c3715SXin Li         coeff = (m_x(j)-i)/v(i);
75*bf2c3715SXin Li         coeff *= coeff;
76*bf2c3715SXin Li         if (coeff < 1. && coeff > 0.)
77*bf2c3715SXin Li           fvec(j) -= u(i)*std::pow((1-coeff), 2);
78*bf2c3715SXin Li       }
79*bf2c3715SXin Li     }
80*bf2c3715SXin Li     return 0;
81*bf2c3715SXin Li   }
82*bf2c3715SXin Li 
dfsparseGaussianTest83*bf2c3715SXin Li   int df(const VectorType& uv, JacobianType& fjac)
84*bf2c3715SXin Li   {
85*bf2c3715SXin Li     int m = Base::values();
86*bf2c3715SXin Li     int n = Base::inputs();
87*bf2c3715SXin Li     eigen_assert(n == uv.size());
88*bf2c3715SXin Li     eigen_assert(fjac.rows() == m);
89*bf2c3715SXin Li     eigen_assert(fjac.cols() == n);
90*bf2c3715SXin Li     int half = n/2;
91*bf2c3715SXin Li     VectorBlock<const VectorType> u(uv, 0, half);
92*bf2c3715SXin Li     VectorBlock<const VectorType> v(uv, half, half);
93*bf2c3715SXin Li     Scalar coeff;
94*bf2c3715SXin Li 
95*bf2c3715SXin Li     //Derivatives with respect to u
96*bf2c3715SXin Li     for (int col = 0; col < half; col++)
97*bf2c3715SXin Li     {
98*bf2c3715SXin Li       for (int row = 0; row < m; row++)
99*bf2c3715SXin Li       {
100*bf2c3715SXin Li         coeff = (m_x(row)-col)/v(col);
101*bf2c3715SXin Li           coeff = coeff*coeff;
102*bf2c3715SXin Li         if(coeff < 1. && coeff > 0.)
103*bf2c3715SXin Li         {
104*bf2c3715SXin Li           fjac.coeffRef(row,col) = -(1-coeff)*(1-coeff);
105*bf2c3715SXin Li         }
106*bf2c3715SXin Li       }
107*bf2c3715SXin Li     }
108*bf2c3715SXin Li     //Derivatives with respect to v
109*bf2c3715SXin Li     for (int col = 0; col < half; col++)
110*bf2c3715SXin Li     {
111*bf2c3715SXin Li       for (int row = 0; row < m; row++)
112*bf2c3715SXin Li       {
113*bf2c3715SXin Li         coeff = (m_x(row)-col)/v(col);
114*bf2c3715SXin Li         coeff = coeff*coeff;
115*bf2c3715SXin Li         if(coeff < 1. && coeff > 0.)
116*bf2c3715SXin Li         {
117*bf2c3715SXin Li           fjac.coeffRef(row,col+half) = -4 * (u(col)/v(col))*coeff*(1-coeff);
118*bf2c3715SXin Li         }
119*bf2c3715SXin Li       }
120*bf2c3715SXin Li     }
121*bf2c3715SXin Li     return 0;
122*bf2c3715SXin Li   }
123*bf2c3715SXin Li 
124*bf2c3715SXin Li   VectorType m_x, m_y; //Data points
125*bf2c3715SXin Li };
126*bf2c3715SXin Li 
127*bf2c3715SXin Li 
128*bf2c3715SXin Li template<typename T>
test_sparseLM_T()129*bf2c3715SXin Li void test_sparseLM_T()
130*bf2c3715SXin Li {
131*bf2c3715SXin Li   typedef Matrix<T,Dynamic,1> VectorType;
132*bf2c3715SXin Li 
133*bf2c3715SXin Li   int inputs = 10;
134*bf2c3715SXin Li   int values = 2000;
135*bf2c3715SXin Li   sparseGaussianTest<T> sparse_gaussian(inputs, values);
136*bf2c3715SXin Li   VectorType uv(inputs),uv_ref(inputs);
137*bf2c3715SXin Li   VectorType x(values);
138*bf2c3715SXin Li   // Generate the reference solution
139*bf2c3715SXin Li   uv_ref << -2, 1, 4 ,8, 6, 1.8, 1.2, 1.1, 1.9 , 3;
140*bf2c3715SXin Li   //Generate the reference data points
141*bf2c3715SXin Li   x.setRandom();
142*bf2c3715SXin Li   x = 10*x;
143*bf2c3715SXin Li   x.array() += 10;
144*bf2c3715SXin Li   sparse_gaussian.initPoints(uv_ref, x);
145*bf2c3715SXin Li 
146*bf2c3715SXin Li 
147*bf2c3715SXin Li   // Generate the initial parameters
148*bf2c3715SXin Li   VectorBlock<VectorType> u(uv, 0, inputs/2);
149*bf2c3715SXin Li   VectorBlock<VectorType> v(uv, inputs/2, inputs/2);
150*bf2c3715SXin Li   v.setOnes();
151*bf2c3715SXin Li   //Generate u or Solve for u from v
152*bf2c3715SXin Li   u.setOnes();
153*bf2c3715SXin Li 
154*bf2c3715SXin Li   // Solve the optimization problem
155*bf2c3715SXin Li   LevenbergMarquardt<sparseGaussianTest<T> > lm(sparse_gaussian);
156*bf2c3715SXin Li   int info;
157*bf2c3715SXin Li //   info = lm.minimize(uv);
158*bf2c3715SXin Li 
159*bf2c3715SXin Li   VERIFY_IS_EQUAL(info,1);
160*bf2c3715SXin Li     // Do a step by step solution and save the residual
161*bf2c3715SXin Li   int maxiter = 200;
162*bf2c3715SXin Li   int iter = 0;
163*bf2c3715SXin Li   MatrixXd Err(values, maxiter);
164*bf2c3715SXin Li   MatrixXd Mod(values, maxiter);
165*bf2c3715SXin Li   LevenbergMarquardtSpace::Status status;
166*bf2c3715SXin Li   status = lm.minimizeInit(uv);
167*bf2c3715SXin Li   if (status==LevenbergMarquardtSpace::ImproperInputParameters)
168*bf2c3715SXin Li       return ;
169*bf2c3715SXin Li 
170*bf2c3715SXin Li }
EIGEN_DECLARE_TEST(sparseLM)171*bf2c3715SXin Li EIGEN_DECLARE_TEST(sparseLM)
172*bf2c3715SXin Li {
173*bf2c3715SXin Li   CALL_SUBTEST_1(test_sparseLM_T<double>());
174*bf2c3715SXin Li 
175*bf2c3715SXin Li   // CALL_SUBTEST_2(test_sparseLM_T<std::complex<double>());
176*bf2c3715SXin Li }
177