xref: /aosp_15_r20/external/eigen/test/eigensolver_selfadjoint.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 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 "svd_fill.h"
13*bf2c3715SXin Li #include <limits>
14*bf2c3715SXin Li #include <Eigen/Eigenvalues>
15*bf2c3715SXin Li #include <Eigen/SparseCore>
16*bf2c3715SXin Li 
17*bf2c3715SXin Li 
selfadjointeigensolver_essential_check(const MatrixType & m)18*bf2c3715SXin Li template<typename MatrixType> void selfadjointeigensolver_essential_check(const MatrixType& m)
19*bf2c3715SXin Li {
20*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
21*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
22*bf2c3715SXin Li   RealScalar eival_eps = numext::mini<RealScalar>(test_precision<RealScalar>(),  NumTraits<Scalar>::dummy_precision()*20000);
23*bf2c3715SXin Li 
24*bf2c3715SXin Li   SelfAdjointEigenSolver<MatrixType> eiSymm(m);
25*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiSymm.info(), Success);
26*bf2c3715SXin Li 
27*bf2c3715SXin Li   RealScalar scaling = m.cwiseAbs().maxCoeff();
28*bf2c3715SXin Li 
29*bf2c3715SXin Li   if(scaling<(std::numeric_limits<RealScalar>::min)())
30*bf2c3715SXin Li   {
31*bf2c3715SXin Li     VERIFY(eiSymm.eigenvalues().cwiseAbs().maxCoeff() <= (std::numeric_limits<RealScalar>::min)());
32*bf2c3715SXin Li   }
33*bf2c3715SXin Li   else
34*bf2c3715SXin Li   {
35*bf2c3715SXin Li     VERIFY_IS_APPROX((m.template selfadjointView<Lower>() * eiSymm.eigenvectors())/scaling,
36*bf2c3715SXin Li                      (eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal())/scaling);
37*bf2c3715SXin Li   }
38*bf2c3715SXin Li   VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
39*bf2c3715SXin Li   VERIFY_IS_UNITARY(eiSymm.eigenvectors());
40*bf2c3715SXin Li 
41*bf2c3715SXin Li   if(m.cols()<=4)
42*bf2c3715SXin Li   {
43*bf2c3715SXin Li     SelfAdjointEigenSolver<MatrixType> eiDirect;
44*bf2c3715SXin Li     eiDirect.computeDirect(m);
45*bf2c3715SXin Li     VERIFY_IS_EQUAL(eiDirect.info(), Success);
46*bf2c3715SXin Li     if(! eiSymm.eigenvalues().isApprox(eiDirect.eigenvalues(), eival_eps) )
47*bf2c3715SXin Li     {
48*bf2c3715SXin Li       std::cerr << "reference eigenvalues: " << eiSymm.eigenvalues().transpose() << "\n"
49*bf2c3715SXin Li                 << "obtained eigenvalues:  " << eiDirect.eigenvalues().transpose() << "\n"
50*bf2c3715SXin Li                 << "diff:                  " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).transpose() << "\n"
51*bf2c3715SXin Li                 << "error (eps):           " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).norm() / eiSymm.eigenvalues().norm() << "  (" << eival_eps << ")\n";
52*bf2c3715SXin Li     }
53*bf2c3715SXin Li     if(scaling<(std::numeric_limits<RealScalar>::min)())
54*bf2c3715SXin Li     {
55*bf2c3715SXin Li       VERIFY(eiDirect.eigenvalues().cwiseAbs().maxCoeff() <= (std::numeric_limits<RealScalar>::min)());
56*bf2c3715SXin Li     }
57*bf2c3715SXin Li     else
58*bf2c3715SXin Li     {
59*bf2c3715SXin Li       VERIFY_IS_APPROX(eiSymm.eigenvalues()/scaling, eiDirect.eigenvalues()/scaling);
60*bf2c3715SXin Li       VERIFY_IS_APPROX((m.template selfadjointView<Lower>() * eiDirect.eigenvectors())/scaling,
61*bf2c3715SXin Li                        (eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal())/scaling);
62*bf2c3715SXin Li       VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues()/scaling, eiDirect.eigenvalues()/scaling);
63*bf2c3715SXin Li     }
64*bf2c3715SXin Li 
65*bf2c3715SXin Li     VERIFY_IS_UNITARY(eiDirect.eigenvectors());
66*bf2c3715SXin Li   }
67*bf2c3715SXin Li }
68*bf2c3715SXin Li 
selfadjointeigensolver(const MatrixType & m)69*bf2c3715SXin Li template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
70*bf2c3715SXin Li {
71*bf2c3715SXin Li   /* this test covers the following files:
72*bf2c3715SXin Li      EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
73*bf2c3715SXin Li   */
74*bf2c3715SXin Li   Index rows = m.rows();
75*bf2c3715SXin Li   Index cols = m.cols();
76*bf2c3715SXin Li 
77*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
78*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   RealScalar largerEps = 10*test_precision<RealScalar>();
81*bf2c3715SXin Li 
82*bf2c3715SXin Li   MatrixType a = MatrixType::Random(rows,cols);
83*bf2c3715SXin Li   MatrixType a1 = MatrixType::Random(rows,cols);
84*bf2c3715SXin Li   MatrixType symmA =  a.adjoint() * a + a1.adjoint() * a1;
85*bf2c3715SXin Li   MatrixType symmC = symmA;
86*bf2c3715SXin Li 
87*bf2c3715SXin Li   svd_fill_random(symmA,Symmetric);
88*bf2c3715SXin Li 
89*bf2c3715SXin Li   symmA.template triangularView<StrictlyUpper>().setZero();
90*bf2c3715SXin Li   symmC.template triangularView<StrictlyUpper>().setZero();
91*bf2c3715SXin Li 
92*bf2c3715SXin Li   MatrixType b = MatrixType::Random(rows,cols);
93*bf2c3715SXin Li   MatrixType b1 = MatrixType::Random(rows,cols);
94*bf2c3715SXin Li   MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
95*bf2c3715SXin Li   symmB.template triangularView<StrictlyUpper>().setZero();
96*bf2c3715SXin Li 
97*bf2c3715SXin Li   CALL_SUBTEST( selfadjointeigensolver_essential_check(symmA) );
98*bf2c3715SXin Li 
99*bf2c3715SXin Li   SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
100*bf2c3715SXin Li   // generalized eigen pb
101*bf2c3715SXin Li   GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmC, symmB);
102*bf2c3715SXin Li 
103*bf2c3715SXin Li   SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
104*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
105*bf2c3715SXin Li   VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
106*bf2c3715SXin Li 
107*bf2c3715SXin Li   // generalized eigen problem Ax = lBx
108*bf2c3715SXin Li   eiSymmGen.compute(symmC, symmB,Ax_lBx);
109*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
110*bf2c3715SXin Li   VERIFY((symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
111*bf2c3715SXin Li           symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
112*bf2c3715SXin Li 
113*bf2c3715SXin Li   // generalized eigen problem BAx = lx
114*bf2c3715SXin Li   eiSymmGen.compute(symmC, symmB,BAx_lx);
115*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
116*bf2c3715SXin Li   VERIFY((symmB.template selfadjointView<Lower>() * (symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
117*bf2c3715SXin Li          (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
118*bf2c3715SXin Li 
119*bf2c3715SXin Li   // generalized eigen problem ABx = lx
120*bf2c3715SXin Li   eiSymmGen.compute(symmC, symmB,ABx_lx);
121*bf2c3715SXin Li   VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
122*bf2c3715SXin Li   VERIFY((symmC.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
123*bf2c3715SXin Li          (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
124*bf2c3715SXin Li 
125*bf2c3715SXin Li 
126*bf2c3715SXin Li   eiSymm.compute(symmC);
127*bf2c3715SXin Li   MatrixType sqrtSymmA = eiSymm.operatorSqrt();
128*bf2c3715SXin Li   VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
129*bf2c3715SXin Li   VERIFY_IS_APPROX(sqrtSymmA, symmC.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
130*bf2c3715SXin Li 
131*bf2c3715SXin Li   MatrixType id = MatrixType::Identity(rows, cols);
132*bf2c3715SXin Li   VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
133*bf2c3715SXin Li 
134*bf2c3715SXin Li   SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
135*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
136*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues());
137*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
138*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
139*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
140*bf2c3715SXin Li 
141*bf2c3715SXin Li   eiSymmUninitialized.compute(symmA, false);
142*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
143*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
144*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
145*bf2c3715SXin Li 
146*bf2c3715SXin Li   // test Tridiagonalization's methods
147*bf2c3715SXin Li   Tridiagonalization<MatrixType> tridiag(symmC);
148*bf2c3715SXin Li   VERIFY_IS_APPROX(tridiag.diagonal(), tridiag.matrixT().diagonal());
149*bf2c3715SXin Li   VERIFY_IS_APPROX(tridiag.subDiagonal(), tridiag.matrixT().template diagonal<-1>());
150*bf2c3715SXin Li   Matrix<RealScalar,Dynamic,Dynamic> T = tridiag.matrixT();
151*bf2c3715SXin Li   if(rows>1 && cols>1) {
152*bf2c3715SXin Li     // FIXME check that upper and lower part are 0:
153*bf2c3715SXin Li     //VERIFY(T.topRightCorner(rows-2, cols-2).template triangularView<Upper>().isZero());
154*bf2c3715SXin Li   }
155*bf2c3715SXin Li   VERIFY_IS_APPROX(tridiag.diagonal(), T.diagonal());
156*bf2c3715SXin Li   VERIFY_IS_APPROX(tridiag.subDiagonal(), T.template diagonal<1>());
157*bf2c3715SXin Li   VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
158*bf2c3715SXin Li   VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT() * tridiag.matrixQ().adjoint());
159*bf2c3715SXin Li 
160*bf2c3715SXin Li   // Test computation of eigenvalues from tridiagonal matrix
161*bf2c3715SXin Li   if(rows > 1)
162*bf2c3715SXin Li   {
163*bf2c3715SXin Li     SelfAdjointEigenSolver<MatrixType> eiSymmTridiag;
164*bf2c3715SXin Li     eiSymmTridiag.computeFromTridiagonal(tridiag.matrixT().diagonal(), tridiag.matrixT().diagonal(-1), ComputeEigenvectors);
165*bf2c3715SXin Li     VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmTridiag.eigenvalues());
166*bf2c3715SXin Li     VERIFY_IS_APPROX(tridiag.matrixT(), eiSymmTridiag.eigenvectors().real() * eiSymmTridiag.eigenvalues().asDiagonal() * eiSymmTridiag.eigenvectors().real().transpose());
167*bf2c3715SXin Li   }
168*bf2c3715SXin Li 
169*bf2c3715SXin Li   if (rows > 1 && rows < 20)
170*bf2c3715SXin Li   {
171*bf2c3715SXin Li     // Test matrix with NaN
172*bf2c3715SXin Li     symmC(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
173*bf2c3715SXin Li     SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmC);
174*bf2c3715SXin Li     VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
175*bf2c3715SXin Li   }
176*bf2c3715SXin Li 
177*bf2c3715SXin Li   // regression test for bug 1098
178*bf2c3715SXin Li   {
179*bf2c3715SXin Li     SelfAdjointEigenSolver<MatrixType> eig(a.adjoint() * a);
180*bf2c3715SXin Li     eig.compute(a.adjoint() * a);
181*bf2c3715SXin Li   }
182*bf2c3715SXin Li 
183*bf2c3715SXin Li   // regression test for bug 478
184*bf2c3715SXin Li   {
185*bf2c3715SXin Li     a.setZero();
186*bf2c3715SXin Li     SelfAdjointEigenSolver<MatrixType> ei3(a);
187*bf2c3715SXin Li     VERIFY_IS_EQUAL(ei3.info(), Success);
188*bf2c3715SXin Li     VERIFY_IS_MUCH_SMALLER_THAN(ei3.eigenvalues().norm(),RealScalar(1));
189*bf2c3715SXin Li     VERIFY((ei3.eigenvectors().transpose()*ei3.eigenvectors().transpose()).eval().isIdentity());
190*bf2c3715SXin Li   }
191*bf2c3715SXin Li }
192*bf2c3715SXin Li 
193*bf2c3715SXin Li template<int>
bug_854()194*bf2c3715SXin Li void bug_854()
195*bf2c3715SXin Li {
196*bf2c3715SXin Li   Matrix3d m;
197*bf2c3715SXin Li   m << 850.961, 51.966, 0,
198*bf2c3715SXin Li        51.966, 254.841, 0,
199*bf2c3715SXin Li             0,       0, 0;
200*bf2c3715SXin Li   selfadjointeigensolver_essential_check(m);
201*bf2c3715SXin Li }
202*bf2c3715SXin Li 
203*bf2c3715SXin Li template<int>
bug_1014()204*bf2c3715SXin Li void bug_1014()
205*bf2c3715SXin Li {
206*bf2c3715SXin Li   Matrix3d m;
207*bf2c3715SXin Li   m <<        0.11111111111111114658, 0, 0,
208*bf2c3715SXin Li        0,     0.11111111111111109107, 0,
209*bf2c3715SXin Li        0, 0,  0.11111111111111107719;
210*bf2c3715SXin Li   selfadjointeigensolver_essential_check(m);
211*bf2c3715SXin Li }
212*bf2c3715SXin Li 
213*bf2c3715SXin Li template<int>
bug_1225()214*bf2c3715SXin Li void bug_1225()
215*bf2c3715SXin Li {
216*bf2c3715SXin Li   Matrix3d m1, m2;
217*bf2c3715SXin Li   m1.setRandom();
218*bf2c3715SXin Li   m1 = m1*m1.transpose();
219*bf2c3715SXin Li   m2 = m1.triangularView<Upper>();
220*bf2c3715SXin Li   SelfAdjointEigenSolver<Matrix3d> eig1(m1);
221*bf2c3715SXin Li   SelfAdjointEigenSolver<Matrix3d> eig2(m2.selfadjointView<Upper>());
222*bf2c3715SXin Li   VERIFY_IS_APPROX(eig1.eigenvalues(), eig2.eigenvalues());
223*bf2c3715SXin Li }
224*bf2c3715SXin Li 
225*bf2c3715SXin Li template<int>
bug_1204()226*bf2c3715SXin Li void bug_1204()
227*bf2c3715SXin Li {
228*bf2c3715SXin Li   SparseMatrix<double> A(2,2);
229*bf2c3715SXin Li   A.setIdentity();
230*bf2c3715SXin Li   SelfAdjointEigenSolver<Eigen::SparseMatrix<double> > eig(A);
231*bf2c3715SXin Li }
232*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(eigensolver_selfadjoint)233*bf2c3715SXin Li EIGEN_DECLARE_TEST(eigensolver_selfadjoint)
234*bf2c3715SXin Li {
235*bf2c3715SXin Li   int s = 0;
236*bf2c3715SXin Li   for(int i = 0; i < g_repeat; i++) {
237*bf2c3715SXin Li 
238*bf2c3715SXin Li     // trivial test for 1x1 matrices:
239*bf2c3715SXin Li     CALL_SUBTEST_1( selfadjointeigensolver(Matrix<float, 1, 1>()));
240*bf2c3715SXin Li     CALL_SUBTEST_1( selfadjointeigensolver(Matrix<double, 1, 1>()));
241*bf2c3715SXin Li     CALL_SUBTEST_1( selfadjointeigensolver(Matrix<std::complex<double>, 1, 1>()));
242*bf2c3715SXin Li 
243*bf2c3715SXin Li     // very important to test 3x3 and 2x2 matrices since we provide special paths for them
244*bf2c3715SXin Li     CALL_SUBTEST_12( selfadjointeigensolver(Matrix2f()) );
245*bf2c3715SXin Li     CALL_SUBTEST_12( selfadjointeigensolver(Matrix2d()) );
246*bf2c3715SXin Li     CALL_SUBTEST_12( selfadjointeigensolver(Matrix2cd()) );
247*bf2c3715SXin Li     CALL_SUBTEST_13( selfadjointeigensolver(Matrix3f()) );
248*bf2c3715SXin Li     CALL_SUBTEST_13( selfadjointeigensolver(Matrix3d()) );
249*bf2c3715SXin Li     CALL_SUBTEST_13( selfadjointeigensolver(Matrix3cd()) );
250*bf2c3715SXin Li     CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
251*bf2c3715SXin Li     CALL_SUBTEST_2( selfadjointeigensolver(Matrix4cd()) );
252*bf2c3715SXin Li 
253*bf2c3715SXin Li     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
254*bf2c3715SXin Li     CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
255*bf2c3715SXin Li     CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
256*bf2c3715SXin Li     CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
257*bf2c3715SXin Li     CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
258*bf2c3715SXin Li     TEST_SET_BUT_UNUSED_VARIABLE(s)
259*bf2c3715SXin Li 
260*bf2c3715SXin Li     // some trivial but implementation-wise tricky cases
261*bf2c3715SXin Li     CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
262*bf2c3715SXin Li     CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
263*bf2c3715SXin Li     CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(1,1)) );
264*bf2c3715SXin Li     CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(2,2)) );
265*bf2c3715SXin Li     CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
266*bf2c3715SXin Li     CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
267*bf2c3715SXin Li   }
268*bf2c3715SXin Li 
269*bf2c3715SXin Li   CALL_SUBTEST_13( bug_854<0>() );
270*bf2c3715SXin Li   CALL_SUBTEST_13( bug_1014<0>() );
271*bf2c3715SXin Li   CALL_SUBTEST_13( bug_1204<0>() );
272*bf2c3715SXin Li   CALL_SUBTEST_13( bug_1225<0>() );
273*bf2c3715SXin Li 
274*bf2c3715SXin Li   // Test problem size constructors
275*bf2c3715SXin Li   s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
276*bf2c3715SXin Li   CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
277*bf2c3715SXin Li   CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
278*bf2c3715SXin Li 
279*bf2c3715SXin Li   TEST_SET_BUT_UNUSED_VARIABLE(s)
280*bf2c3715SXin Li }
281*bf2c3715SXin Li 
282