diff options
author | Gael Guennebaud <g.gael@free.fr> | 2015-10-26 16:00:25 +0100 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2015-10-26 16:00:25 +0100 |
commit | f93654ae16b261e462ee00c5255072f8dd7d387b (patch) | |
tree | f9419de7849cd11f515a395e20081bcf684ac7de /test/eigensolver_generalized_real.cpp | |
parent | af2e25d482c4cca04701371785f3d08211ea4a2c (diff) |
bug #1098: fix regression introduced when generalizing some compute() methods in changeset 7031a851d45a8526474ac1ac972ad12a48e99f1a
.
Diffstat (limited to 'test/eigensolver_generalized_real.cpp')
-rw-r--r-- | test/eigensolver_generalized_real.cpp | 8 |
1 files changed, 8 insertions, 0 deletions
diff --git a/test/eigensolver_generalized_real.cpp b/test/eigensolver_generalized_real.cpp index 566a4bdc6..a46a2e50e 100644 --- a/test/eigensolver_generalized_real.cpp +++ b/test/eigensolver_generalized_real.cpp @@ -39,6 +39,14 @@ template<typename MatrixType> void generalized_eigensolver_real(const MatrixType VectorType realEigenvalues = eig.eigenvalues().real(); std::sort(realEigenvalues.data(), realEigenvalues.data()+realEigenvalues.size()); VERIFY_IS_APPROX(realEigenvalues, symmEig.eigenvalues()); + + // regression test for bug 1098 + { + GeneralizedSelfAdjointEigenSolver<MatrixType> eig1(a.adjoint() * a,b.adjoint() * b); + eig1.compute(a.adjoint() * a,b.adjoint() * b); + GeneralizedEigenSolver<MatrixType> eig2(a.adjoint() * a,b.adjoint() * b); + eig2.compute(a.adjoint() * a,b.adjoint() * b); + } } void test_eigensolver_generalized_real() |