source: pacpussensors/trunk/Vislab/lib3dv/eigen/test/eigensolver_selfadjoint.cpp@ 138

Last change on this file since 138 was 136, checked in by ldecherf, 8 years ago

Doc

File size: 6.8 KB
Line 
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#include "main.h"
12#include <limits>
13#include <Eigen/Eigenvalues>
14
15template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
16{
17 typedef typename MatrixType::Index Index;
18 /* this test covers the following files:
19 EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
20 */
21 Index rows = m.rows();
22 Index cols = m.cols();
23
24 typedef typename MatrixType::Scalar Scalar;
25 typedef typename NumTraits<Scalar>::Real RealScalar;
26
27 RealScalar largerEps = 10*test_precision<RealScalar>();
28
29 MatrixType a = MatrixType::Random(rows,cols);
30 MatrixType a1 = MatrixType::Random(rows,cols);
31 MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
32 MatrixType symmC = symmA;
33
34 // randomly nullify some rows/columns
35 {
36 Index count = 1;//internal::random<Index>(-cols,cols);
37 for(Index k=0; k<count; ++k)
38 {
39 Index i = internal::random<Index>(0,cols-1);
40 symmA.row(i).setZero();
41 symmA.col(i).setZero();
42 }
43 }
44
45 symmA.template triangularView<StrictlyUpper>().setZero();
46 symmC.template triangularView<StrictlyUpper>().setZero();
47
48 MatrixType b = MatrixType::Random(rows,cols);
49 MatrixType b1 = MatrixType::Random(rows,cols);
50 MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
51 symmB.template triangularView<StrictlyUpper>().setZero();
52
53 SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
54 SelfAdjointEigenSolver<MatrixType> eiDirect;
55 eiDirect.computeDirect(symmA);
56 // generalized eigen pb
57 GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmC, symmB);
58
59 VERIFY_IS_EQUAL(eiSymm.info(), Success);
60 VERIFY((symmA.template selfadjointView<Lower>() * eiSymm.eigenvectors()).isApprox(
61 eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal(), largerEps));
62 VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
63
64 VERIFY_IS_EQUAL(eiDirect.info(), Success);
65 VERIFY((symmA.template selfadjointView<Lower>() * eiDirect.eigenvectors()).isApprox(
66 eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal(), largerEps));
67 VERIFY_IS_APPROX(symmA.template selfadjointView<Lower>().eigenvalues(), eiDirect.eigenvalues());
68
69 SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
70 VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
71 VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
72
73 // generalized eigen problem Ax = lBx
74 eiSymmGen.compute(symmC, symmB,Ax_lBx);
75 VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
76 VERIFY((symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
77 symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
78
79 // generalized eigen problem BAx = lx
80 eiSymmGen.compute(symmC, symmB,BAx_lx);
81 VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
82 VERIFY((symmB.template selfadjointView<Lower>() * (symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
83 (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
84
85 // generalized eigen problem ABx = lx
86 eiSymmGen.compute(symmC, symmB,ABx_lx);
87 VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
88 VERIFY((symmC.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
89 (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
90
91
92 eiSymm.compute(symmC);
93 MatrixType sqrtSymmA = eiSymm.operatorSqrt();
94 VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
95 VERIFY_IS_APPROX(sqrtSymmA, symmC.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
96
97 MatrixType id = MatrixType::Identity(rows, cols);
98 VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
99
100 SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
101 VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
102 VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues());
103 VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
104 VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
105 VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
106
107 eiSymmUninitialized.compute(symmA, false);
108 VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors());
109 VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt());
110 VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt());
111
112 // test Tridiagonalization's methods
113 Tridiagonalization<MatrixType> tridiag(symmC);
114 // FIXME tridiag.matrixQ().adjoint() does not work
115 VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
116
117 if (rows > 1)
118 {
119 // Test matrix with NaN
120 symmC(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
121 SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmC);
122 VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
123 }
124}
125
126void test_eigensolver_selfadjoint()
127{
128 int s = 0;
129 for(int i = 0; i < g_repeat; i++) {
130 // very important to test 3x3 and 2x2 matrices since we provide special paths for them
131 CALL_SUBTEST_1( selfadjointeigensolver(Matrix2f()) );
132 CALL_SUBTEST_1( selfadjointeigensolver(Matrix2d()) );
133 CALL_SUBTEST_1( selfadjointeigensolver(Matrix3f()) );
134 CALL_SUBTEST_1( selfadjointeigensolver(Matrix3d()) );
135 CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
136 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
137 CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
138 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
139 CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
140 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
141 CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
142
143 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
144 CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
145
146 // some trivial but implementation-wise tricky cases
147 CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
148 CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
149 CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
150 CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
151 }
152
153 // Test problem size constructors
154 s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
155 CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
156 CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
157
158 TEST_SET_BUT_UNUSED_VARIABLE(s)
159}
160
Note: See TracBrowser for help on using the repository browser.