aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Eigenvalues/HessenbergDecomposition.h8
-rw-r--r--test/CMakeLists.txt1
-rw-r--r--test/hessenberg.cpp46
3 files changed, 53 insertions, 2 deletions
diff --git a/Eigen/src/Eigenvalues/HessenbergDecomposition.h b/Eigen/src/Eigenvalues/HessenbergDecomposition.h
index 9f7df49bc..69597e77f 100644
--- a/Eigen/src/Eigenvalues/HessenbergDecomposition.h
+++ b/Eigen/src/Eigenvalues/HessenbergDecomposition.h
@@ -71,9 +71,11 @@ template<typename _MatrixType> class HessenbergDecomposition
{}
HessenbergDecomposition(const MatrixType& matrix)
- : m_matrix(matrix),
- m_hCoeffs(matrix.cols()-1)
+ : m_matrix(matrix)
{
+ if(matrix.rows()<=2)
+ return;
+ m_hCoeffs.resize(matrix.rows()-1,1);
_compute(m_matrix, m_hCoeffs);
}
@@ -84,6 +86,8 @@ template<typename _MatrixType> class HessenbergDecomposition
void compute(const MatrixType& matrix)
{
m_matrix = matrix;
+ if(matrix.rows()<=2)
+ return;
m_hCoeffs.resize(matrix.rows()-1,1);
_compute(m_matrix, m_hCoeffs);
}
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index b8efbcf51..ebb8a59d2 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -127,6 +127,7 @@ ei_add_test(inverse)
ei_add_test(qr)
ei_add_test(qr_colpivoting)
ei_add_test(qr_fullpivoting)
+ei_add_test(hessenberg " " "${GSL_LIBRARIES}")
ei_add_test(eigensolver_selfadjoint " " "${GSL_LIBRARIES}")
ei_add_test(eigensolver_generic " " "${GSL_LIBRARIES}")
ei_add_test(eigensolver_complex)
diff --git a/test/hessenberg.cpp b/test/hessenberg.cpp
new file mode 100644
index 000000000..d917be357
--- /dev/null
+++ b/test/hessenberg.cpp
@@ -0,0 +1,46 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
+//
+// 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 <http://www.gnu.org/licenses/>.
+
+#include "main.h"
+#include <Eigen/Eigenvalues>
+
+template<typename Scalar,int Size> void hessenberg(int size = Size)
+{
+ typedef Matrix<Scalar,Size,Size> MatrixType;
+ MatrixType m = MatrixType::Random(size,size);
+ HessenbergDecomposition<MatrixType> hess(m);
+
+ VERIFY_IS_APPROX(m, hess.matrixQ() * hess.matrixH() * hess.matrixQ().adjoint());
+}
+
+void test_hessenberg()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ CALL_SUBTEST_1(( hessenberg<std::complex<double>,1>() ));
+ CALL_SUBTEST_2(( hessenberg<std::complex<double>,2>() ));
+ CALL_SUBTEST_3(( hessenberg<std::complex<float>,4>() ));
+ CALL_SUBTEST_4(( hessenberg<float,Dynamic>(ei_random<int>(1,320)) ));
+ CALL_SUBTEST_5(( hessenberg<std::complex<double>,Dynamic>(ei_random<int>(1,320)) ));
+ }
+}