aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/eigensolver_selfadjoint.cpp
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2016-05-19 13:07:33 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2016-05-19 13:07:33 +0200
commitdf9a5e13c67ae09ec1e942763f02d1a5500c091e (patch)
treed639c1651d4546b86c92ff2b0df04a7a361c43d3 /test/eigensolver_selfadjoint.cpp
parent6a2916df80037e49935281558ebe8a6a4d3f76cc (diff)
Fix SelfAdjointEigenSolver for some input expression types, and add new regression unit tests for sparse and selfadjointview inputs.
Diffstat (limited to 'test/eigensolver_selfadjoint.cpp')
-rw-r--r--test/eigensolver_selfadjoint.cpp29
1 files changed, 27 insertions, 2 deletions
diff --git a/test/eigensolver_selfadjoint.cpp b/test/eigensolver_selfadjoint.cpp
index f909761a1..cd0ae5c2f 100644
--- a/test/eigensolver_selfadjoint.cpp
+++ b/test/eigensolver_selfadjoint.cpp
@@ -12,6 +12,7 @@
#include "svd_fill.h"
#include <limits>
#include <Eigen/Eigenvalues>
+#include <Eigen/SparseCore>
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()
{
Matrix3d m;
@@ -173,6 +175,7 @@ void bug_854()
selfadjointeigensolver_essential_check(m);
}
+template<int>
void bug_1014()
{
Matrix3d m;
@@ -182,6 +185,26 @@ void bug_1014()
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()
{
int s = 0;
@@ -210,8 +233,10 @@ void test_eigensolver_selfadjoint()
CALL_SUBTEST_7( selfadjointeigensolver(Matrix<double,2,2>()) );
}
- CALL_SUBTEST_13( bug_854() );
- CALL_SUBTEST_13( bug_1014() );
+ CALL_SUBTEST_13( bug_854<0>() );
+ CALL_SUBTEST_13( bug_1014<0>() );
+ CALL_SUBTEST_13( bug_1204<0>() );
+ CALL_SUBTEST_13( bug_1225<0>() );
// Test problem size constructors
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);