aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/hessenberg.cpp
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-12-23 12:18:54 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-12-23 12:18:54 +0100
commit791bac25f2feffd199f38da93deb54655eb07b13 (patch)
tree0780e04a17c498d3bf124bdc87b4a7bd0e428a6d /test/hessenberg.cpp
parentd65c8cb60acd34a0eb898194713ff45e604de0fe (diff)
fix #75, and add a basic unit test for Hessenberg
(it was indirectly tested by the eigenvalue decomposition)
Diffstat (limited to 'test/hessenberg.cpp')
-rw-r--r--test/hessenberg.cpp46
1 files changed, 46 insertions, 0 deletions
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)) ));
+ }
+}