xref: /aosp_15_r20/external/eigen/test/eigensolver_generic.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) 2008 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2010,2012 Jitse Niesen <[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 "main.h"
12*bf2c3715SXin Li #include <limits>
13*bf2c3715SXin Li #include <Eigen/Eigenvalues>
14*bf2c3715SXin Li 
15*bf2c3715SXin Li template<typename EigType,typename MatType>
check_eigensolver_for_given_mat(const EigType & eig,const MatType & a)16*bf2c3715SXin Li void check_eigensolver_for_given_mat(const EigType &eig, const MatType& a)
17*bf2c3715SXin Li {
18*bf2c3715SXin Li   typedef typename NumTraits<typename MatType::Scalar>::Real RealScalar;
19*bf2c3715SXin Li   typedef Matrix<RealScalar, MatType::RowsAtCompileTime, 1> RealVectorType;
20*bf2c3715SXin Li   typedef typename std::complex<RealScalar> Complex;
21*bf2c3715SXin Li   Index n = a.rows();
22*bf2c3715SXin Li   VERIFY_IS_EQUAL(eig.info(), Success);
23*bf2c3715SXin Li   VERIFY_IS_APPROX(a * eig.pseudoEigenvectors(), eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix());
24*bf2c3715SXin Li   VERIFY_IS_APPROX(a.template cast<Complex>() * eig.eigenvectors(),
25*bf2c3715SXin Li                    eig.eigenvectors() * eig.eigenvalues().asDiagonal());
26*bf2c3715SXin Li   VERIFY_IS_APPROX(eig.eigenvectors().colwise().norm(), RealVectorType::Ones(n).transpose());
27*bf2c3715SXin Li   VERIFY_IS_APPROX(a.eigenvalues(), eig.eigenvalues());
28*bf2c3715SXin Li }
29*bf2c3715SXin Li 
eigensolver(const MatrixType & m)30*bf2c3715SXin Li template<typename MatrixType> void eigensolver(const MatrixType& m)
31*bf2c3715SXin Li {
32*bf2c3715SXin Li   /* this test covers the following files:
33*bf2c3715SXin Li      EigenSolver.h
34*bf2c3715SXin Li   */
35*bf2c3715SXin Li   Index rows = m.rows();
36*bf2c3715SXin Li   Index cols = m.cols();
37*bf2c3715SXin Li 
38*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
39*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
40*bf2c3715SXin Li   typedef typename std::complex<RealScalar> Complex;
41*bf2c3715SXin Li 
42*bf2c3715SXin Li   MatrixType a = MatrixType::Random(rows,cols);
43*bf2c3715SXin Li   MatrixType a1 = MatrixType::Random(rows,cols);
44*bf2c3715SXin Li   MatrixType symmA =  a.adjoint() * a + a1.adjoint() * a1;
45*bf2c3715SXin Li 
46*bf2c3715SXin Li   EigenSolver<MatrixType> ei0(symmA);
47*bf2c3715SXin Li   VERIFY_IS_EQUAL(ei0.info(), Success);
48*bf2c3715SXin Li   VERIFY_IS_APPROX(symmA * ei0.pseudoEigenvectors(), ei0.pseudoEigenvectors() * ei0.pseudoEigenvalueMatrix());
49*bf2c3715SXin Li   VERIFY_IS_APPROX((symmA.template cast<Complex>()) * (ei0.pseudoEigenvectors().template cast<Complex>()),
50*bf2c3715SXin Li     (ei0.pseudoEigenvectors().template cast<Complex>()) * (ei0.eigenvalues().asDiagonal()));
51*bf2c3715SXin Li 
52*bf2c3715SXin Li   EigenSolver<MatrixType> ei1(a);
53*bf2c3715SXin Li   CALL_SUBTEST( check_eigensolver_for_given_mat(ei1,a) );
54*bf2c3715SXin Li 
55*bf2c3715SXin Li   EigenSolver<MatrixType> ei2;
56*bf2c3715SXin Li   ei2.setMaxIterations(RealSchur<MatrixType>::m_maxIterationsPerRow * rows).compute(a);
57*bf2c3715SXin Li   VERIFY_IS_EQUAL(ei2.info(), Success);
58*bf2c3715SXin Li   VERIFY_IS_EQUAL(ei2.eigenvectors(), ei1.eigenvectors());
59*bf2c3715SXin Li   VERIFY_IS_EQUAL(ei2.eigenvalues(), ei1.eigenvalues());
60*bf2c3715SXin Li   if (rows > 2) {
61*bf2c3715SXin Li     ei2.setMaxIterations(1).compute(a);
62*bf2c3715SXin Li     VERIFY_IS_EQUAL(ei2.info(), NoConvergence);
63*bf2c3715SXin Li     VERIFY_IS_EQUAL(ei2.getMaxIterations(), 1);
64*bf2c3715SXin Li   }
65*bf2c3715SXin Li 
66*bf2c3715SXin Li   EigenSolver<MatrixType> eiNoEivecs(a, false);
67*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiNoEivecs.info(), Success);
68*bf2c3715SXin Li   VERIFY_IS_APPROX(ei1.eigenvalues(), eiNoEivecs.eigenvalues());
69*bf2c3715SXin Li   VERIFY_IS_APPROX(ei1.pseudoEigenvalueMatrix(), eiNoEivecs.pseudoEigenvalueMatrix());
70*bf2c3715SXin Li 
71*bf2c3715SXin Li   MatrixType id = MatrixType::Identity(rows, cols);
72*bf2c3715SXin Li   VERIFY_IS_APPROX(id.operatorNorm(), RealScalar(1));
73*bf2c3715SXin Li 
74*bf2c3715SXin Li   if (rows > 2 && rows < 20)
75*bf2c3715SXin Li   {
76*bf2c3715SXin Li     // Test matrix with NaN
77*bf2c3715SXin Li     a(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
78*bf2c3715SXin Li     EigenSolver<MatrixType> eiNaN(a);
79*bf2c3715SXin Li     VERIFY_IS_NOT_EQUAL(eiNaN.info(), Success);
80*bf2c3715SXin Li   }
81*bf2c3715SXin Li 
82*bf2c3715SXin Li   // regression test for bug 1098
83*bf2c3715SXin Li   {
84*bf2c3715SXin Li     EigenSolver<MatrixType> eig(a.adjoint() * a);
85*bf2c3715SXin Li     eig.compute(a.adjoint() * a);
86*bf2c3715SXin Li   }
87*bf2c3715SXin Li 
88*bf2c3715SXin Li   // regression test for bug 478
89*bf2c3715SXin Li   {
90*bf2c3715SXin Li     a.setZero();
91*bf2c3715SXin Li     EigenSolver<MatrixType> ei3(a);
92*bf2c3715SXin Li     VERIFY_IS_EQUAL(ei3.info(), Success);
93*bf2c3715SXin Li     VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
94*bf2c3715SXin Li     VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
95*bf2c3715SXin Li   }
96*bf2c3715SXin Li }
97*bf2c3715SXin Li 
eigensolver_verify_assert(const MatrixType & m)98*bf2c3715SXin Li template<typename MatrixType> void eigensolver_verify_assert(const MatrixType& m)
99*bf2c3715SXin Li {
100*bf2c3715SXin Li   EigenSolver<MatrixType> eig;
101*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.eigenvectors());
102*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
103*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.pseudoEigenvalueMatrix());
104*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.eigenvalues());
105*bf2c3715SXin Li 
106*bf2c3715SXin Li   MatrixType a = MatrixType::Random(m.rows(),m.cols());
107*bf2c3715SXin Li   eig.compute(a, false);
108*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.eigenvectors());
109*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eig.pseudoEigenvectors());
110*bf2c3715SXin Li }
111*bf2c3715SXin Li 
112*bf2c3715SXin Li 
113*bf2c3715SXin Li template<typename CoeffType>
114*bf2c3715SXin Li Matrix<typename CoeffType::Scalar,Dynamic,Dynamic>
make_companion(const CoeffType & coeffs)115*bf2c3715SXin Li make_companion(const CoeffType& coeffs)
116*bf2c3715SXin Li {
117*bf2c3715SXin Li   Index n = coeffs.size()-1;
118*bf2c3715SXin Li   Matrix<typename CoeffType::Scalar,Dynamic,Dynamic> res(n,n);
119*bf2c3715SXin Li   res.setZero();
120*bf2c3715SXin Li 	res.row(0) = -coeffs.tail(n) / coeffs(0);
121*bf2c3715SXin Li 	res.diagonal(-1).setOnes();
122*bf2c3715SXin Li   return res;
123*bf2c3715SXin Li }
124*bf2c3715SXin Li 
125*bf2c3715SXin Li template<int>
eigensolver_generic_extra()126*bf2c3715SXin Li void eigensolver_generic_extra()
127*bf2c3715SXin Li {
128*bf2c3715SXin Li   {
129*bf2c3715SXin Li     // regression test for bug 793
130*bf2c3715SXin Li     MatrixXd a(3,3);
131*bf2c3715SXin Li     a << 0,  0,  1,
132*bf2c3715SXin Li         1,  1, 1,
133*bf2c3715SXin Li         1, 1e+200,  1;
134*bf2c3715SXin Li     Eigen::EigenSolver<MatrixXd> eig(a);
135*bf2c3715SXin Li     double scale = 1e-200; // scale to avoid overflow during the comparisons
136*bf2c3715SXin Li     VERIFY_IS_APPROX(a * eig.pseudoEigenvectors()*scale, eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()*scale);
137*bf2c3715SXin Li     VERIFY_IS_APPROX(a * eig.eigenvectors()*scale, eig.eigenvectors() * eig.eigenvalues().asDiagonal()*scale);
138*bf2c3715SXin Li   }
139*bf2c3715SXin Li   {
140*bf2c3715SXin Li     // check a case where all eigenvalues are null.
141*bf2c3715SXin Li     MatrixXd a(2,2);
142*bf2c3715SXin Li     a << 1,  1,
143*bf2c3715SXin Li         -1, -1;
144*bf2c3715SXin Li     Eigen::EigenSolver<MatrixXd> eig(a);
145*bf2c3715SXin Li     VERIFY_IS_APPROX(eig.pseudoEigenvectors().squaredNorm(), 2.);
146*bf2c3715SXin Li     VERIFY_IS_APPROX((a * eig.pseudoEigenvectors()).norm()+1., 1.);
147*bf2c3715SXin Li     VERIFY_IS_APPROX((eig.pseudoEigenvectors() * eig.pseudoEigenvalueMatrix()).norm()+1., 1.);
148*bf2c3715SXin Li     VERIFY_IS_APPROX((a * eig.eigenvectors()).norm()+1., 1.);
149*bf2c3715SXin Li     VERIFY_IS_APPROX((eig.eigenvectors() * eig.eigenvalues().asDiagonal()).norm()+1., 1.);
150*bf2c3715SXin Li   }
151*bf2c3715SXin Li 
152*bf2c3715SXin Li   // regression test for bug 933
153*bf2c3715SXin Li   {
154*bf2c3715SXin Li     {
155*bf2c3715SXin Li       VectorXd coeffs(5); coeffs << 1, -3, -175, -225, 2250;
156*bf2c3715SXin Li       MatrixXd C = make_companion(coeffs);
157*bf2c3715SXin Li       EigenSolver<MatrixXd> eig(C);
158*bf2c3715SXin Li       CALL_SUBTEST( check_eigensolver_for_given_mat(eig,C) );
159*bf2c3715SXin Li     }
160*bf2c3715SXin Li     {
161*bf2c3715SXin Li       // this test is tricky because it requires high accuracy in smallest eigenvalues
162*bf2c3715SXin Li       VectorXd coeffs(5); coeffs << 6.154671e-15, -1.003870e-10, -9.819570e-01, 3.995715e+03, 2.211511e+08;
163*bf2c3715SXin Li       MatrixXd C = make_companion(coeffs);
164*bf2c3715SXin Li       EigenSolver<MatrixXd> eig(C);
165*bf2c3715SXin Li       CALL_SUBTEST( check_eigensolver_for_given_mat(eig,C) );
166*bf2c3715SXin Li       Index n = C.rows();
167*bf2c3715SXin Li       for(Index i=0;i<n;++i)
168*bf2c3715SXin Li       {
169*bf2c3715SXin Li         typedef std::complex<double> Complex;
170*bf2c3715SXin Li         MatrixXcd ac = C.cast<Complex>();
171*bf2c3715SXin Li         ac.diagonal().array() -= eig.eigenvalues()(i);
172*bf2c3715SXin Li         VectorXd sv = ac.jacobiSvd().singularValues();
173*bf2c3715SXin Li         // comparing to sv(0) is not enough here to catch the "bug",
174*bf2c3715SXin Li         // the hard-coded 1.0 is important!
175*bf2c3715SXin Li         VERIFY_IS_MUCH_SMALLER_THAN(sv(n-1), 1.0);
176*bf2c3715SXin Li       }
177*bf2c3715SXin Li     }
178*bf2c3715SXin Li   }
179*bf2c3715SXin Li   // regression test for bug 1557
180*bf2c3715SXin Li   {
181*bf2c3715SXin Li     // this test is interesting because it contains zeros on the diagonal.
182*bf2c3715SXin Li     MatrixXd A_bug1557(3,3);
183*bf2c3715SXin Li     A_bug1557 << 0, 0, 0, 1, 0, 0.5887907064808635127, 0, 1, 0;
184*bf2c3715SXin Li     EigenSolver<MatrixXd> eig(A_bug1557);
185*bf2c3715SXin Li     CALL_SUBTEST( check_eigensolver_for_given_mat(eig,A_bug1557) );
186*bf2c3715SXin Li   }
187*bf2c3715SXin Li 
188*bf2c3715SXin Li   // regression test for bug 1174
189*bf2c3715SXin Li   {
190*bf2c3715SXin Li     Index n = 12;
191*bf2c3715SXin Li     MatrixXf A_bug1174(n,n);
192*bf2c3715SXin Li     A_bug1174 <<  262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
193*bf2c3715SXin Li                   262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
194*bf2c3715SXin Li                   262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
195*bf2c3715SXin Li                   262144, 0, 0, 262144, 786432, 0, 0, 0, 0, 0, 0, 786432,
196*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
197*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
198*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
199*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
200*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
201*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
202*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0,
203*bf2c3715SXin Li                   0, 262144, 262144, 0, 0, 262144, 262144, 262144, 262144, 262144, 262144, 0;
204*bf2c3715SXin Li     EigenSolver<MatrixXf> eig(A_bug1174);
205*bf2c3715SXin Li     CALL_SUBTEST( check_eigensolver_for_given_mat(eig,A_bug1174) );
206*bf2c3715SXin Li   }
207*bf2c3715SXin Li }
208*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(eigensolver_generic)209*bf2c3715SXin Li EIGEN_DECLARE_TEST(eigensolver_generic)
210*bf2c3715SXin Li {
211*bf2c3715SXin Li   int s = 0;
212*bf2c3715SXin Li   for(int i = 0; i < g_repeat; i++) {
213*bf2c3715SXin Li     CALL_SUBTEST_1( eigensolver(Matrix4f()) );
214*bf2c3715SXin Li     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
215*bf2c3715SXin Li     CALL_SUBTEST_2( eigensolver(MatrixXd(s,s)) );
216*bf2c3715SXin Li     TEST_SET_BUT_UNUSED_VARIABLE(s)
217*bf2c3715SXin Li 
218*bf2c3715SXin Li     // some trivial but implementation-wise tricky cases
219*bf2c3715SXin Li     CALL_SUBTEST_2( eigensolver(MatrixXd(1,1)) );
220*bf2c3715SXin Li     CALL_SUBTEST_2( eigensolver(MatrixXd(2,2)) );
221*bf2c3715SXin Li     CALL_SUBTEST_3( eigensolver(Matrix<double,1,1>()) );
222*bf2c3715SXin Li     CALL_SUBTEST_4( eigensolver(Matrix2d()) );
223*bf2c3715SXin Li   }
224*bf2c3715SXin Li 
225*bf2c3715SXin Li   CALL_SUBTEST_1( eigensolver_verify_assert(Matrix4f()) );
226*bf2c3715SXin Li   s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
227*bf2c3715SXin Li   CALL_SUBTEST_2( eigensolver_verify_assert(MatrixXd(s,s)) );
228*bf2c3715SXin Li   CALL_SUBTEST_3( eigensolver_verify_assert(Matrix<double,1,1>()) );
229*bf2c3715SXin Li   CALL_SUBTEST_4( eigensolver_verify_assert(Matrix2d()) );
230*bf2c3715SXin Li 
231*bf2c3715SXin Li   // Test problem size constructors
232*bf2c3715SXin Li   CALL_SUBTEST_5(EigenSolver<MatrixXf> tmp(s));
233*bf2c3715SXin Li 
234*bf2c3715SXin Li   // regression test for bug 410
235*bf2c3715SXin Li   CALL_SUBTEST_2(
236*bf2c3715SXin Li   {
237*bf2c3715SXin Li      MatrixXd A(1,1);
238*bf2c3715SXin Li      A(0,0) = std::sqrt(-1.); // is Not-a-Number
239*bf2c3715SXin Li      Eigen::EigenSolver<MatrixXd> solver(A);
240*bf2c3715SXin Li      VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
241*bf2c3715SXin Li   }
242*bf2c3715SXin Li   );
243*bf2c3715SXin Li 
244*bf2c3715SXin Li   CALL_SUBTEST_2( eigensolver_generic_extra<0>() );
245*bf2c3715SXin Li 
246*bf2c3715SXin Li   TEST_SET_BUT_UNUSED_VARIABLE(s)
247*bf2c3715SXin Li }
248