Fix SelfAdjointEigenSolver for some input expression types, and add new regression unit tests for sparse and selfadjointview inputs.

This commit is contained in:
Gael Guennebaud 2016-05-19 13:07:33 +02:00
parent 6a2916df80
commit df9a5e13c6
2 changed files with 28 additions and 3 deletions

View File

@ -414,7 +414,7 @@ SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
if(n==1) if(n==1)
{ {
m_eivalues.coeffRef(0,0) = numext::real(matrix.coeff(0,0)); m_eivalues.coeffRef(0,0) = numext::real(matrix.diagonal()[0]);
if(computeEigenvectors) if(computeEigenvectors)
m_eivec.setOnes(n,n); m_eivec.setOnes(n,n);
m_info = Success; m_info = Success;

View File

@ -12,6 +12,7 @@
#include "svd_fill.h" #include "svd_fill.h"
#include <limits> #include <limits>
#include <Eigen/Eigenvalues> #include <Eigen/Eigenvalues>
#include <Eigen/SparseCore>
template<typename MatrixType> void selfadjointeigensolver_essential_check(const MatrixType& m) template<typename MatrixType> void selfadjointeigensolver_essential_check(const MatrixType& m)
@ -164,6 +165,7 @@ template<typename MatrixType> void selfadjointeigensolver(const MatrixType& m)
} }
} }
template<int>
void bug_854() void bug_854()
{ {
Matrix3d m; Matrix3d m;
@ -173,6 +175,7 @@ void bug_854()
selfadjointeigensolver_essential_check(m); selfadjointeigensolver_essential_check(m);
} }
template<int>
void bug_1014() void bug_1014()
{ {
Matrix3d m; Matrix3d m;
@ -182,6 +185,26 @@ void bug_1014()
selfadjointeigensolver_essential_check(m); selfadjointeigensolver_essential_check(m);
} }
template<int>
void bug_1225()
{
Matrix3d m1, m2;
m1.setRandom();
m1 = m1*m1.transpose();
m2 = m1.triangularView<Upper>();
SelfAdjointEigenSolver<Matrix3d> eig1(m1);
SelfAdjointEigenSolver<Matrix3d> eig2(m2.selfadjointView<Upper>());
VERIFY_IS_APPROX(eig1.eigenvalues(), eig2.eigenvalues());
}
template<int>
void bug_1204()
{
SparseMatrix<double> A(2,2);
A.setIdentity();
SelfAdjointEigenSolver<Eigen::SparseMatrix<double> > eig(A);
}
void test_eigensolver_selfadjoint() void test_eigensolver_selfadjoint()
{ {
int s = 0; int s = 0;
@ -210,8 +233,10 @@ void test_eigensolver_selfadjoint()
CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) ); CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
} }
CALL_SUBTEST_13( bug_854() ); CALL_SUBTEST_13( bug_854<0>() );
CALL_SUBTEST_13( bug_1014() ); CALL_SUBTEST_13( bug_1014<0>() );
CALL_SUBTEST_13( bug_1204<0>() );
CALL_SUBTEST_13( bug_1225<0>() );
// Test problem size constructors // Test problem size constructors
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4); s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);