mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-11 11:19:02 +08:00
fix #75, and add a basic unit test for Hessenberg
(it was indirectly tested by the eigenvalue decomposition)
This commit is contained in:
parent
d65c8cb60a
commit
791bac25f2
@ -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);
|
||||
}
|
||||
|
@ -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)
|
||||
|
46
test/hessenberg.cpp
Normal file
46
test/hessenberg.cpp
Normal file
@ -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)) ));
|
||||
}
|
||||
}
|
Loading…
x
Reference in New Issue
Block a user