// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2008 Gael Guennebaud // Copyright (C) 2010 Jitse Niesen // // Eigen is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 3 of the License, or (at your option) any later version. // // Alternatively, you can redistribute it and/or // modify it under the terms of the GNU General Public License as // published by the Free Software Foundation; either version 2 of // the License, or (at your option) any later version. // // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the // GNU General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License and a copy of the GNU General Public License along with // Eigen. If not, see . #include "main.h" #include #include #ifdef HAS_GSL #include "gsl_helper.h" #endif template void selfadjointeigensolver(const MatrixType& m) { /* this test covers the following files: EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h) */ int rows = m.rows(); int cols = m.cols(); typedef typename MatrixType::Scalar Scalar; typedef typename NumTraits::Real RealScalar; typedef Matrix VectorType; typedef Matrix RealVectorType; typedef typename std::complex::Real> Complex; RealScalar largerEps = 10*test_precision(); MatrixType a = MatrixType::Random(rows,cols); MatrixType a1 = MatrixType::Random(rows,cols); MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1; symmA.template triangularView().setZero(); MatrixType b = MatrixType::Random(rows,cols); MatrixType b1 = MatrixType::Random(rows,cols); MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1; symmB.template triangularView().setZero(); SelfAdjointEigenSolver eiSymm(symmA); // generalized eigen pb GeneralizedSelfAdjointEigenSolver eiSymmGen(symmA, symmB); #ifdef HAS_GSL if (ei_is_same_type::ret) { // restore symmA and symmB. symmA = MatrixType(symmA.template selfadjointView()); symmB = MatrixType(symmB.template selfadjointView()); typedef GslTraits Gsl; typename Gsl::Matrix gEvec=0, gSymmA=0, gSymmB=0; typename GslTraits::Vector gEval=0; RealVectorType _eval; MatrixType _evec; convert(symmA, gSymmA); convert(symmB, gSymmB); convert(symmA, gEvec); gEval = GslTraits::createVector(rows); Gsl::eigen_symm(gSymmA, gEval, gEvec); convert(gEval, _eval); convert(gEvec, _evec); // test gsl itself ! VERIFY((symmA * _evec).isApprox(_evec * _eval.asDiagonal(), largerEps)); // compare with eigen VERIFY_IS_APPROX(_eval, eiSymm.eigenvalues()); VERIFY_IS_APPROX(_evec.cwiseAbs(), eiSymm.eigenvectors().cwiseAbs()); // generalized pb Gsl::eigen_symm_gen(gSymmA, gSymmB, gEval, gEvec); convert(gEval, _eval); convert(gEvec, _evec); // test GSL itself: VERIFY((symmA * _evec).isApprox(symmB * (_evec * _eval.asDiagonal()), largerEps)); // compare with eigen MatrixType normalized_eivec = eiSymmGen.eigenvectors()*eiSymmGen.eigenvectors().colwise().norm().asDiagonal().inverse(); VERIFY_IS_APPROX(_eval, eiSymmGen.eigenvalues()); VERIFY_IS_APPROX(_evec.cwiseAbs(), normalized_eivec.cwiseAbs()); Gsl::free(gSymmA); Gsl::free(gSymmB); GslTraits::free(gEval); Gsl::free(gEvec); } #endif VERIFY_IS_EQUAL(eiSymm.info(), Success); VERIFY((symmA.template selfadjointView() * eiSymm.eigenvectors()).isApprox( eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal(), largerEps)); VERIFY_IS_APPROX(symmA.template selfadjointView().eigenvalues(), eiSymm.eigenvalues()); SelfAdjointEigenSolver eiSymmNoEivecs(symmA, false); VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success); VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues()); // generalized eigen problem Ax = lBx VERIFY_IS_EQUAL(eiSymmGen.info(), Success); VERIFY((symmA.template selfadjointView() * eiSymmGen.eigenvectors()).isApprox( symmB.template selfadjointView() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps)); MatrixType sqrtSymmA = eiSymm.operatorSqrt(); VERIFY_IS_APPROX(MatrixType(symmA.template selfadjointView()), sqrtSymmA*sqrtSymmA); VERIFY_IS_APPROX(sqrtSymmA, symmA.template selfadjointView()*eiSymm.operatorInverseSqrt()); MatrixType id = MatrixType::Identity(rows, cols); VERIFY_IS_APPROX(id.template selfadjointView().operatorNorm(), RealScalar(1)); SelfAdjointEigenSolver eiSymmUninitialized; VERIFY_RAISES_ASSERT(eiSymmUninitialized.info()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvalues()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt()); eiSymmUninitialized.compute(symmA, false); VERIFY_RAISES_ASSERT(eiSymmUninitialized.eigenvectors()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorSqrt()); VERIFY_RAISES_ASSERT(eiSymmUninitialized.operatorInverseSqrt()); if (rows > 1) { // Test matrix with NaN symmA(0,0) = std::numeric_limits::quiet_NaN(); SelfAdjointEigenSolver eiSymmNaN(symmA); VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence); } } void test_eigensolver_selfadjoint() { for(int i = 0; i < g_repeat; i++) { // very important to test a 3x3 matrix since we provide a special path for it CALL_SUBTEST_1( selfadjointeigensolver(Matrix3f()) ); CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) ); CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(10,10)) ); CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(19,19)) ); CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(17,17)) ); // some trivial but implementation-wise tricky cases CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) ); CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) ); CALL_SUBTEST_6( selfadjointeigensolver(Matrix()) ); CALL_SUBTEST_7( selfadjointeigensolver(Matrix()) ); } // Test problem size constructors CALL_SUBTEST_8(SelfAdjointEigenSolver(10)); CALL_SUBTEST_8(Tridiagonalization(10)); }