2008-05-27 17:16:27 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-05-27 17:16:27 +08:00
|
|
|
//
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2010-06-04 05:59:57 +08:00
|
|
|
// Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
|
2008-05-27 17:16:27 +08:00
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2008-05-27 17:16:27 +08:00
|
|
|
|
|
|
|
#include "main.h"
|
2015-05-07 21:54:07 +08:00
|
|
|
#include "svd_fill.h"
|
2010-06-04 05:59:57 +08:00
|
|
|
#include <limits>
|
2009-09-04 15:23:38 +08:00
|
|
|
#include <Eigen/Eigenvalues>
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2015-05-07 21:54:07 +08:00
|
|
|
|
2015-05-13 00:45:39 +08:00
|
|
|
template<typename MatrixType> void selfadjointeigensolver_essential_check(const MatrixType& m)
|
|
|
|
{
|
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
2015-06-08 22:16:42 +08:00
|
|
|
RealScalar eival_eps = (std::min)(test_precision<RealScalar>(), NumTraits<Scalar>::dummy_precision()*20000);
|
2015-05-13 00:45:39 +08:00
|
|
|
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymm(m);
|
|
|
|
VERIFY_IS_EQUAL(eiSymm.info(), Success);
|
2015-06-08 22:16:42 +08:00
|
|
|
VERIFY_IS_APPROX(m.template selfadjointView<Lower>() * eiSymm.eigenvectors(),
|
|
|
|
eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal());
|
2015-05-13 00:45:39 +08:00
|
|
|
VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues(), eiSymm.eigenvalues());
|
|
|
|
VERIFY_IS_UNITARY(eiSymm.eigenvectors());
|
|
|
|
|
2015-06-08 22:16:42 +08:00
|
|
|
if(m.cols()<=4)
|
|
|
|
{
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eiDirect;
|
|
|
|
eiDirect.computeDirect(m);
|
|
|
|
VERIFY_IS_EQUAL(eiDirect.info(), Success);
|
|
|
|
VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiDirect.eigenvalues());
|
|
|
|
if(! eiSymm.eigenvalues().isApprox(eiDirect.eigenvalues(), eival_eps) )
|
|
|
|
{
|
|
|
|
std::cerr << "reference eigenvalues: " << eiSymm.eigenvalues().transpose() << "\n"
|
|
|
|
<< "obtained eigenvalues: " << eiDirect.eigenvalues().transpose() << "\n"
|
|
|
|
<< "diff: " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).transpose() << "\n"
|
|
|
|
<< "error (eps): " << (eiSymm.eigenvalues()-eiDirect.eigenvalues()).norm() / eiSymm.eigenvalues().norm() << " (" << eival_eps << ")\n";
|
|
|
|
}
|
|
|
|
VERIFY(eiSymm.eigenvalues().isApprox(eiDirect.eigenvalues(), eival_eps));
|
|
|
|
VERIFY_IS_APPROX(m.template selfadjointView<Lower>() * eiDirect.eigenvectors(),
|
|
|
|
eiDirect.eigenvectors() * eiDirect.eigenvalues().asDiagonal());
|
|
|
|
VERIFY_IS_APPROX(m.template selfadjointView<Lower>().eigenvalues(), eiDirect.eigenvalues());
|
|
|
|
VERIFY_IS_UNITARY(eiDirect.eigenvectors());
|
|
|
|
}
|
2015-05-13 00:45:39 +08:00
|
|
|
}
|
2015-05-07 21:54:07 +08:00
|
|
|
|
2008-10-01 18:17:08 +08:00
|
|
|
template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
|
2008-05-27 17:16:27 +08:00
|
|
|
{
|
2010-06-20 23:37:56 +08:00
|
|
|
typedef typename MatrixType::Index Index;
|
2008-05-27 17:16:27 +08:00
|
|
|
/* this test covers the following files:
|
2008-06-04 02:04:36 +08:00
|
|
|
EigenSolver.h, SelfAdjointEigenSolver.h (and indirectly: Tridiagonalization.h)
|
2008-05-27 17:16:27 +08:00
|
|
|
*/
|
2010-06-20 23:37:56 +08:00
|
|
|
Index rows = m.rows();
|
|
|
|
Index cols = m.cols();
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2008-08-23 23:14:20 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2008-08-23 23:14:20 +08:00
|
|
|
RealScalar largerEps = 10*test_precision<RealScalar>();
|
|
|
|
|
2008-09-02 01:31:21 +08:00
|
|
|
MatrixType a = MatrixType::Random(rows,cols);
|
|
|
|
MatrixType a1 = MatrixType::Random(rows,cols);
|
2008-08-23 23:14:20 +08:00
|
|
|
MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
|
2014-08-21 16:49:09 +08:00
|
|
|
MatrixType symmC = symmA;
|
|
|
|
|
2015-05-07 21:54:07 +08:00
|
|
|
svd_fill_random(symmA,Symmetric);
|
|
|
|
|
2010-06-10 22:39:46 +08:00
|
|
|
symmA.template triangularView<StrictlyUpper>().setZero();
|
2014-08-21 16:49:09 +08:00
|
|
|
symmC.template triangularView<StrictlyUpper>().setZero();
|
2008-08-23 23:14:20 +08:00
|
|
|
|
2008-09-02 01:31:21 +08:00
|
|
|
MatrixType b = MatrixType::Random(rows,cols);
|
|
|
|
MatrixType b1 = MatrixType::Random(rows,cols);
|
2008-08-23 23:14:20 +08:00
|
|
|
MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
|
2010-06-10 22:39:46 +08:00
|
|
|
symmB.template triangularView<StrictlyUpper>().setZero();
|
2015-05-13 00:45:39 +08:00
|
|
|
|
|
|
|
CALL_SUBTEST( selfadjointeigensolver_essential_check(symmA) );
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2008-06-15 03:42:12 +08:00
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
|
2008-08-23 23:14:20 +08:00
|
|
|
// generalized eigen pb
|
2014-08-21 16:49:09 +08:00
|
|
|
GeneralizedSelfAdjointEigenSolver<MatrixType> eiSymmGen(symmC, symmB);
|
2008-08-23 23:14:20 +08:00
|
|
|
|
2010-05-31 04:49:35 +08:00
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymmNoEivecs(symmA, false);
|
2010-06-04 05:59:57 +08:00
|
|
|
VERIFY_IS_EQUAL(eiSymmNoEivecs.info(), Success);
|
2010-05-31 04:49:35 +08:00
|
|
|
VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmNoEivecs.eigenvalues());
|
2011-07-22 01:07:52 +08:00
|
|
|
|
2008-06-15 03:42:12 +08:00
|
|
|
// generalized eigen problem Ax = lBx
|
2014-08-21 16:49:09 +08:00
|
|
|
eiSymmGen.compute(symmC, symmB,Ax_lBx);
|
2010-06-04 05:59:57 +08:00
|
|
|
VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
|
2014-08-21 16:49:09 +08:00
|
|
|
VERIFY((symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors()).isApprox(
|
2010-06-10 22:39:46 +08:00
|
|
|
symmB.template selfadjointView<Lower>() * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2010-06-17 16:16:15 +08:00
|
|
|
// generalized eigen problem BAx = lx
|
2014-08-21 16:49:09 +08:00
|
|
|
eiSymmGen.compute(symmC, symmB,BAx_lx);
|
2010-06-17 16:16:15 +08:00
|
|
|
VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
|
2014-08-21 16:49:09 +08:00
|
|
|
VERIFY((symmB.template selfadjointView<Lower>() * (symmC.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
|
2010-06-17 16:16:15 +08:00
|
|
|
(eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
|
|
|
|
|
|
|
|
// generalized eigen problem ABx = lx
|
2014-08-21 16:49:09 +08:00
|
|
|
eiSymmGen.compute(symmC, symmB,ABx_lx);
|
2010-06-17 16:16:15 +08:00
|
|
|
VERIFY_IS_EQUAL(eiSymmGen.info(), Success);
|
2014-08-21 16:49:09 +08:00
|
|
|
VERIFY((symmC.template selfadjointView<Lower>() * (symmB.template selfadjointView<Lower>() * eiSymmGen.eigenvectors())).isApprox(
|
2010-06-17 16:16:15 +08:00
|
|
|
(eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal()), largerEps));
|
|
|
|
|
|
|
|
|
2014-08-21 16:49:09 +08:00
|
|
|
eiSymm.compute(symmC);
|
2008-12-19 22:15:32 +08:00
|
|
|
MatrixType sqrtSymmA = eiSymm.operatorSqrt();
|
2014-08-21 16:49:09 +08:00
|
|
|
VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), sqrtSymmA*sqrtSymmA);
|
|
|
|
VERIFY_IS_APPROX(sqrtSymmA, symmC.template selfadjointView<Lower>()*eiSymm.operatorInverseSqrt());
|
2010-05-25 00:43:50 +08:00
|
|
|
|
|
|
|
MatrixType id = MatrixType::Identity(rows, cols);
|
|
|
|
VERIFY_IS_APPROX(id.template selfadjointView<Lower>().operatorNorm(), RealScalar(1));
|
2010-05-31 04:49:35 +08:00
|
|
|
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymmUninitialized;
|
2010-06-04 05:59:57 +08:00
|
|
|
VERIFY_RAISES_ASSERT(eiSymmUninitialized.info());
|
2010-05-31 04:49:35 +08:00
|
|
|
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());
|
2010-06-04 05:59:57 +08:00
|
|
|
|
2010-11-26 22:31:47 +08:00
|
|
|
// test Tridiagonalization's methods
|
2014-08-21 16:49:09 +08:00
|
|
|
Tridiagonalization<MatrixType> tridiag(symmC);
|
2014-10-01 20:33:55 +08:00
|
|
|
VERIFY_IS_APPROX(tridiag.diagonal(), tridiag.matrixT().diagonal());
|
2014-09-24 20:37:13 +08:00
|
|
|
VERIFY_IS_APPROX(tridiag.subDiagonal(), tridiag.matrixT().template diagonal<-1>());
|
2015-10-07 21:36:12 +08:00
|
|
|
Matrix<RealScalar,Dynamic,Dynamic> T = tridiag.matrixT();
|
2014-09-24 20:37:13 +08:00
|
|
|
if(rows>1 && cols>1) {
|
|
|
|
// FIXME check that upper and lower part are 0:
|
|
|
|
//VERIFY(T.topRightCorner(rows-2, cols-2).template triangularView<Upper>().isZero());
|
|
|
|
}
|
2015-10-07 21:36:12 +08:00
|
|
|
VERIFY_IS_APPROX(tridiag.diagonal(), T.diagonal());
|
|
|
|
VERIFY_IS_APPROX(tridiag.subDiagonal(), T.template diagonal<1>());
|
2014-08-21 16:49:09 +08:00
|
|
|
VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT().eval() * MatrixType(tridiag.matrixQ()).adjoint());
|
2014-09-24 20:37:13 +08:00
|
|
|
VERIFY_IS_APPROX(MatrixType(symmC.template selfadjointView<Lower>()), tridiag.matrixQ() * tridiag.matrixT() * tridiag.matrixQ().adjoint());
|
2010-11-26 22:31:47 +08:00
|
|
|
|
2013-07-18 16:32:31 +08:00
|
|
|
// Test computation of eigenvalues from tridiagonal matrix
|
|
|
|
if(rows > 1)
|
|
|
|
{
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymmTridiag;
|
|
|
|
eiSymmTridiag.computeFromTridiagonal(tridiag.matrixT().diagonal(), tridiag.matrixT().diagonal(-1), ComputeEigenvectors);
|
|
|
|
VERIFY_IS_APPROX(eiSymm.eigenvalues(), eiSymmTridiag.eigenvalues());
|
|
|
|
VERIFY_IS_APPROX(tridiag.matrixT(), eiSymmTridiag.eigenvectors().real() * eiSymmTridiag.eigenvalues().asDiagonal() * eiSymmTridiag.eigenvectors().real().transpose());
|
|
|
|
}
|
|
|
|
|
2015-10-30 21:44:22 +08:00
|
|
|
if (rows > 1 && rows < 20)
|
2010-06-04 05:59:57 +08:00
|
|
|
{
|
|
|
|
// Test matrix with NaN
|
2014-08-21 16:49:09 +08:00
|
|
|
symmC(0,0) = std::numeric_limits<typename MatrixType::RealScalar>::quiet_NaN();
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eiSymmNaN(symmC);
|
2010-06-04 05:59:57 +08:00
|
|
|
VERIFY_IS_EQUAL(eiSymmNaN.info(), NoConvergence);
|
|
|
|
}
|
2015-10-26 23:00:25 +08:00
|
|
|
|
|
|
|
// regression test for bug 1098
|
|
|
|
{
|
|
|
|
SelfAdjointEigenSolver<MatrixType> eig(a.adjoint() * a);
|
|
|
|
eig.compute(a.adjoint() * a);
|
|
|
|
}
|
2008-10-01 18:17:08 +08:00
|
|
|
}
|
2008-05-27 17:16:27 +08:00
|
|
|
|
2015-05-13 00:45:39 +08:00
|
|
|
void bug_854()
|
|
|
|
{
|
|
|
|
Matrix3d m;
|
|
|
|
m << 850.961, 51.966, 0,
|
|
|
|
51.966, 254.841, 0,
|
|
|
|
0, 0, 0;
|
|
|
|
selfadjointeigensolver_essential_check(m);
|
|
|
|
}
|
|
|
|
|
|
|
|
void bug_1014()
|
|
|
|
{
|
|
|
|
Matrix3d m;
|
|
|
|
m << 0.11111111111111114658, 0, 0,
|
|
|
|
0, 0.11111111111111109107, 0,
|
|
|
|
0, 0, 0.11111111111111107719;
|
|
|
|
selfadjointeigensolver_essential_check(m);
|
|
|
|
}
|
|
|
|
|
2009-03-23 22:38:59 +08:00
|
|
|
void test_eigensolver_selfadjoint()
|
2008-05-27 17:16:27 +08:00
|
|
|
{
|
2013-06-24 01:11:32 +08:00
|
|
|
int s = 0;
|
2008-08-23 23:14:20 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2014-10-30 00:46:54 +08:00
|
|
|
// trivial test for 1x1 matrices:
|
|
|
|
CALL_SUBTEST_1( selfadjointeigensolver(Matrix<float, 1, 1>()));
|
|
|
|
CALL_SUBTEST_1( selfadjointeigensolver(Matrix<double, 1, 1>()));
|
2011-12-02 01:17:19 +08:00
|
|
|
// very important to test 3x3 and 2x2 matrices since we provide special paths for them
|
2014-10-30 00:46:54 +08:00
|
|
|
CALL_SUBTEST_12( selfadjointeigensolver(Matrix2f()) );
|
|
|
|
CALL_SUBTEST_12( selfadjointeigensolver(Matrix2d()) );
|
|
|
|
CALL_SUBTEST_13( selfadjointeigensolver(Matrix3f()) );
|
|
|
|
CALL_SUBTEST_13( selfadjointeigensolver(Matrix3d()) );
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_2( selfadjointeigensolver(Matrix4d()) );
|
2015-02-18 18:30:44 +08:00
|
|
|
|
2011-07-12 20:41:00 +08:00
|
|
|
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
|
|
|
|
CALL_SUBTEST_3( selfadjointeigensolver(MatrixXf(s,s)) );
|
|
|
|
CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(s,s)) );
|
|
|
|
CALL_SUBTEST_5( selfadjointeigensolver(MatrixXcd(s,s)) );
|
|
|
|
CALL_SUBTEST_9( selfadjointeigensolver(Matrix<std::complex<double>,Dynamic,Dynamic,RowMajor>(s,s)) );
|
2015-02-18 18:30:44 +08:00
|
|
|
TEST_SET_BUT_UNUSED_VARIABLE(s)
|
2008-10-01 18:17:08 +08:00
|
|
|
|
2009-03-23 22:38:59 +08:00
|
|
|
// some trivial but implementation-wise tricky cases
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(1,1)) );
|
|
|
|
CALL_SUBTEST_4( selfadjointeigensolver(MatrixXd(2,2)) );
|
|
|
|
CALL_SUBTEST_6( selfadjointeigensolver(Matrix<double,1,1>()) );
|
|
|
|
CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
|
2008-05-27 17:16:27 +08:00
|
|
|
}
|
2015-05-13 00:45:39 +08:00
|
|
|
|
|
|
|
CALL_SUBTEST_13( bug_854() );
|
|
|
|
CALL_SUBTEST_13( bug_1014() );
|
2010-04-21 23:15:57 +08:00
|
|
|
|
|
|
|
// Test problem size constructors
|
2011-07-12 20:41:00 +08:00
|
|
|
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
|
2013-06-24 01:11:32 +08:00
|
|
|
CALL_SUBTEST_8(SelfAdjointEigenSolver<MatrixXf> tmp1(s));
|
|
|
|
CALL_SUBTEST_8(Tridiagonalization<MatrixXf> tmp2(s));
|
2011-10-31 11:55:20 +08:00
|
|
|
|
2013-06-25 17:42:04 +08:00
|
|
|
TEST_SET_BUT_UNUSED_VARIABLE(s)
|
2008-05-27 17:16:27 +08:00
|
|
|
}
|
2008-10-03 21:22:54 +08:00
|
|
|
|