Eigen::HessenbergDecomposition< MatrixType_ > Class Template Reference

Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation. More...

#include <HessenbergDecomposition.h>

Public Types

enum  {
  Size = MatrixType::RowsAtCompileTime , SizeMinusOne = Size == Dynamic ? Dynamic : Size - 1 , Options = internal::traits<MatrixType>::Options , MaxSize = MatrixType::MaxRowsAtCompileTime ,
  MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : MaxSize - 1
}
 
typedef MatrixType_ MatrixType
 Synonym for the template parameter MatrixType_. More...
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type MatrixType. More...
 
typedef Eigen::Index Index
 
typedef Matrix< Scalar, SizeMinusOne, 1, Options &~RowMajor, MaxSizeMinusOne, 1 > CoeffVectorType
 Type for vector of Householder coefficients. More...
 
typedef HouseholderSequence< MatrixType, internal::remove_all_t< typename CoeffVectorType::ConjugateReturnType > > HouseholderSequenceType
 Return type of matrixQ() More...
 
typedef internal::HessenbergDecompositionMatrixHReturnType< MatrixTypeMatrixHReturnType
 

Public Member Functions

 HessenbergDecomposition (Index size=Size==Dynamic ? 2 :Size)
 Default constructor; the decomposition will be computed later. More...
 
template<typename InputType >
 HessenbergDecomposition (const EigenBase< InputType > &matrix)
 Constructor; computes Hessenberg decomposition of given matrix. More...
 
template<typename InputType >
HessenbergDecompositioncompute (const EigenBase< InputType > &matrix)
 Computes Hessenberg decomposition of given matrix. More...
 
const CoeffVectorTypehouseholderCoefficients () const
 Returns the Householder coefficients. More...
 
const MatrixTypepackedMatrix () const
 Returns the internal representation of the decomposition. More...
 
HouseholderSequenceType matrixQ () const
 Reconstructs the orthogonal matrix Q in the decomposition. More...
 
MatrixHReturnType matrixH () const
 Constructs the Hessenberg matrix H in the decomposition. More...
 

Protected Attributes

MatrixType m_matrix
 
CoeffVectorType m_hCoeffs
 
VectorType m_temp
 
bool m_isInitialized
 

Private Types

typedef Matrix< Scalar, 1, Size, int(Options)|int(RowMajor), 1, MaxSizeVectorType
 
typedef NumTraits< Scalar >::Real RealScalar
 

Static Private Member Functions

static void _compute (MatrixType &matA, CoeffVectorType &hCoeffs, VectorType &temp)
 

Detailed Description

template<typename MatrixType_>
class Eigen::HessenbergDecomposition< MatrixType_ >

Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation.

\eigenvalues_module

Template Parameters
MatrixType_the type of the matrix of which we are computing the Hessenberg decomposition

This class performs an Hessenberg decomposition of a matrix \( A \). In the real case, the Hessenberg decomposition consists of an orthogonal matrix \( Q \) and a Hessenberg matrix \( H \) such that \( A = Q H Q^T \). An orthogonal matrix is a matrix whose inverse equals its transpose ( \( Q^{-1} = Q^T \)). A Hessenberg matrix has zeros below the subdiagonal, so it is almost upper triangular. The Hessenberg decomposition of a complex matrix is \( A = Q H Q^* \) with \( Q \) unitary (that is, \( Q^{-1} = Q^* \)).

Call the function compute() to compute the Hessenberg decomposition of a given matrix. Alternatively, you can use the HessenbergDecomposition(const MatrixType&) constructor which computes the Hessenberg decomposition at construction time. Once the decomposition is computed, you can use the matrixH() and matrixQ() functions to construct the matrices H and Q in the decomposition.

The documentation for matrixH() contains an example of the typical use of this class.

See also
class ComplexSchur, class Tridiagonalization, QR Module

Member Typedef Documentation

◆ CoeffVectorType

template<typename MatrixType_ >
typedef Matrix<Scalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> Eigen::HessenbergDecomposition< MatrixType_ >::CoeffVectorType

Type for vector of Householder coefficients.

This is column vector with entries of type Scalar. The length of the vector is one less than the size of MatrixType, if it is a fixed-side type.

◆ HouseholderSequenceType

template<typename MatrixType_ >
typedef HouseholderSequence<MatrixType, internal::remove_all_t<typename CoeffVectorType::ConjugateReturnType> > Eigen::HessenbergDecomposition< MatrixType_ >::HouseholderSequenceType

Return type of matrixQ()

◆ Index

template<typename MatrixType_ >
typedef Eigen::Index Eigen::HessenbergDecomposition< MatrixType_ >::Index
Deprecated:
since Eigen 3.3

◆ MatrixHReturnType

◆ MatrixType

template<typename MatrixType_ >
typedef MatrixType_ Eigen::HessenbergDecomposition< MatrixType_ >::MatrixType

Synonym for the template parameter MatrixType_.

◆ RealScalar

template<typename MatrixType_ >
typedef NumTraits<Scalar>::Real Eigen::HessenbergDecomposition< MatrixType_ >::RealScalar
private

◆ Scalar

template<typename MatrixType_ >
typedef MatrixType::Scalar Eigen::HessenbergDecomposition< MatrixType_ >::Scalar

Scalar type for matrices of type MatrixType.

◆ VectorType

template<typename MatrixType_ >
typedef Matrix<Scalar, 1, Size, int(Options) | int(RowMajor), 1, MaxSize> Eigen::HessenbergDecomposition< MatrixType_ >::VectorType
private

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
Size 
SizeMinusOne 
Options 
MaxSize 
MaxSizeMinusOne 
66  {
67  Size = MatrixType::RowsAtCompileTime,
68  SizeMinusOne = Size == Dynamic ? Dynamic : Size - 1,
69  Options = internal::traits<MatrixType>::Options,
70  MaxSize = MatrixType::MaxRowsAtCompileTime,
72  };
@ SizeMinusOne
Definition: HessenbergDecomposition.h:68
@ MaxSizeMinusOne
Definition: HessenbergDecomposition.h:71
@ MaxSize
Definition: HessenbergDecomposition.h:70
@ Options
Definition: HessenbergDecomposition.h:69
@ Size
Definition: HessenbergDecomposition.h:67
const int Dynamic
Definition: Constants.h:25

Constructor & Destructor Documentation

◆ HessenbergDecomposition() [1/2]

template<typename MatrixType_ >
Eigen::HessenbergDecomposition< MatrixType_ >::HessenbergDecomposition ( Index  size = Size == Dynamic ? 2 : Size)
inlineexplicit

Default constructor; the decomposition will be computed later.

Parameters
[in]sizeThe size of the matrix whose Hessenberg decomposition will be computed.

The default constructor is useful in cases in which the user intends to perform decompositions via compute(). The size parameter is only used as a hint. It is not an error to give a wrong size, but it may impair performance.

See also
compute() for an example.
103  : Size)
104  : m_matrix(size, size), m_temp(size), m_isInitialized(false) {
105  if (size > 1) m_hCoeffs.resize(size - 1);
106  }
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
MatrixType m_matrix
Definition: HessenbergDecomposition.h:261
VectorType m_temp
Definition: HessenbergDecomposition.h:263
bool m_isInitialized
Definition: HessenbergDecomposition.h:264
CoeffVectorType m_hCoeffs
Definition: HessenbergDecomposition.h:262
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294

References Eigen::HessenbergDecomposition< MatrixType_ >::m_hCoeffs, Eigen::PlainObjectBase< Derived >::resize(), and size.

◆ HessenbergDecomposition() [2/2]

template<typename MatrixType_ >
template<typename InputType >
Eigen::HessenbergDecomposition< MatrixType_ >::HessenbergDecomposition ( const EigenBase< InputType > &  matrix)
inlineexplicit

Constructor; computes Hessenberg decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Hessenberg decomposition is to be computed.

This constructor calls compute() to compute the Hessenberg decomposition.

See also
matrixH() for an example.
119  : m_matrix(matrix.derived()), m_temp(matrix.rows()), m_isInitialized(false) {
120  if (matrix.rows() < 2) {
121  m_isInitialized = true;
122  return;
123  }
124  m_hCoeffs.resize(matrix.rows() - 1, 1);
126  m_isInitialized = true;
127  }
static void _compute(MatrixType &matA, CoeffVectorType &hCoeffs, VectorType &temp)
Definition: HessenbergDecomposition.h:280
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References Eigen::HessenbergDecomposition< MatrixType_ >::_compute(), Eigen::HessenbergDecomposition< MatrixType_ >::m_hCoeffs, Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized, Eigen::HessenbergDecomposition< MatrixType_ >::m_matrix, Eigen::HessenbergDecomposition< MatrixType_ >::m_temp, matrix(), and Eigen::PlainObjectBase< Derived >::resize().

Member Function Documentation

◆ _compute()

template<typename MatrixType >
void Eigen::HessenbergDecomposition< MatrixType >::_compute ( MatrixType matA,
CoeffVectorType hCoeffs,
VectorType temp 
)
staticprivate

Performs a tridiagonal decomposition of matA in place.

Parameters
matAthe input selfadjoint matrix
hCoeffsreturned Householder coefficients

The result is written in the lower triangular part of matA.

Implemented from Golub's "%Matrix Computations", algorithm 8.3.1.

See also
packedMatrix()
280  {
281  eigen_assert(matA.rows() == matA.cols());
282  Index n = matA.rows();
283  temp.resize(n);
284  for (Index i = 0; i < n - 1; ++i) {
285  // let's consider the vector v = i-th column starting at position i+1
286  Index remainingSize = n - i - 1;
288  Scalar h;
289  matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);
290  matA.col(i).coeffRef(i + 1) = beta;
291  hCoeffs.coeffRef(i) = h;
292 
293  // Apply similarity transformation to remaining columns,
294  // i.e., compute A = H A H'
295 
296  // A = H A
297  matA.bottomRightCorner(remainingSize, remainingSize)
298  .applyHouseholderOnTheLeft(matA.col(i).tail(remainingSize - 1), h, &temp.coeffRef(0));
299 
300  // A = A H'
301  matA.rightCols(remainingSize)
302  .applyHouseholderOnTheRight(matA.col(i).tail(remainingSize - 1), numext::conj(h), &temp.coeffRef(0));
303  }
304 }
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:133
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define eigen_assert(x)
Definition: Macros.h:910
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
Eigen::Index Index
Definition: HessenbergDecomposition.h:76
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:217
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:192
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:191
Scalar beta
Definition: level2_cplx_impl.h:36
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > matA(size, size)

References beta, Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >::coeffRef(), Eigen::PlainObjectBase< Derived >::cols(), conj(), eigen_assert, i, matA(), n, Eigen::PlainObjectBase< Derived >::resize(), and Eigen::PlainObjectBase< Derived >::rows().

Referenced by Eigen::HessenbergDecomposition< MatrixType_ >::compute(), and Eigen::HessenbergDecomposition< MatrixType_ >::HessenbergDecomposition().

◆ compute()

template<typename MatrixType_ >
template<typename InputType >
HessenbergDecomposition& Eigen::HessenbergDecomposition< MatrixType_ >::compute ( const EigenBase< InputType > &  matrix)
inline

Computes Hessenberg decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Hessenberg decomposition is to be computed.
Returns
Reference to *this

The Hessenberg decomposition is computed by bringing the columns of the matrix successively in the required form using Householder reflections (see, e.g., Algorithm 7.4.2 in Golub & Van Loan, Matrix Computations). The cost is \( 10n^3/3 \) flops, where \( n \) denotes the size of the given matrix.

This method reuses of the allocated data in the HessenbergDecomposition object.

Example:

MatrixXcf A = MatrixXcf::Random(4, 4);
HessenbergDecomposition<MatrixXcf> hd(4);
hd.compute(A);
cout << "The matrix H in the decomposition of A is:" << endl << hd.matrixH() << endl;
hd.compute(2 * A); // re-use hd to compute and store decomposition of 2A
cout << "The matrix H in the decomposition of 2A is:" << endl << hd.matrixH() << endl;
HessenbergDecomposition< MatrixXcf > hd(4)
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47

Output:

 
147  {
148  m_matrix = matrix.derived();
149  if (matrix.rows() < 2) {
150  m_isInitialized = true;
151  return *this;
152  }
153  m_hCoeffs.resize(matrix.rows() - 1, 1);
155  m_isInitialized = true;
156  return *this;
157  }

References Eigen::HessenbergDecomposition< MatrixType_ >::_compute(), Eigen::HessenbergDecomposition< MatrixType_ >::m_hCoeffs, Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized, Eigen::HessenbergDecomposition< MatrixType_ >::m_matrix, Eigen::HessenbergDecomposition< MatrixType_ >::m_temp, matrix(), and Eigen::PlainObjectBase< Derived >::resize().

Referenced by ctms_decompositions(), hessenberg(), Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, IsComplex >::run(), and Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, false >::run().

◆ householderCoefficients()

template<typename MatrixType_ >
const CoeffVectorType& Eigen::HessenbergDecomposition< MatrixType_ >::householderCoefficients ( ) const
inline

Returns the Householder coefficients.

Returns
a const reference to the vector of Householder coefficients
Precondition
Either the constructor HessenbergDecomposition(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the Hessenberg decomposition of a matrix.

The Householder coefficients allow the reconstruction of the matrix \( Q \) in the Hessenberg decomposition from the packed data.

See also
packedMatrix(), Householder module
172  {
173  eigen_assert(m_isInitialized && "HessenbergDecomposition is not initialized.");
174  return m_hCoeffs;
175  }

References eigen_assert, Eigen::HessenbergDecomposition< MatrixType_ >::m_hCoeffs, and Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized.

Referenced by hessenberg().

◆ matrixH()

template<typename MatrixType_ >
MatrixHReturnType Eigen::HessenbergDecomposition< MatrixType_ >::matrixH ( ) const
inline

Constructs the Hessenberg matrix H in the decomposition.

Returns
expression object representing the matrix H
Precondition
Either the constructor HessenbergDecomposition(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the Hessenberg decomposition of a matrix.

The object returned by this function constructs the Hessenberg matrix H when it is assigned to a matrix or otherwise evaluated. The matrix H is constructed from the packed matrix as returned by packedMatrix(): The upper part (including the subdiagonal) of the packed matrix contains the matrix H. It may sometimes be better to directly use the packed matrix instead of constructing the matrix H.

Example:

Matrix4f A = MatrixXf::Random(4, 4);
cout << "Here is a random 4x4 matrix:" << endl << A << endl;
HessenbergDecomposition<MatrixXf> hessOfA(A);
MatrixXf H = hessOfA.matrixH();
cout << "The Hessenberg matrix H is:" << endl << H << endl;
MatrixXf Q = hessOfA.matrixQ();
cout << "The orthogonal matrix Q is:" << endl << Q << endl;
cout << "Q H Q^T is:" << endl << Q * H * Q.transpose() << endl;
cout<< "Here is a random 4x4 matrix:"<< endl<< A<< endl;HessenbergDecomposition< MatrixXf > hessOfA(A)
MatrixXf H
Definition: HessenbergDecomposition_matrixH.cpp:4
MatrixXf Q
Definition: HouseholderQR_householderQ.cpp:1

Output:

See also
matrixQ(), packedMatrix()
250  {
251  eigen_assert(m_isInitialized && "HessenbergDecomposition is not initialized.");
252  return MatrixHReturnType(*this);
253  }
internal::HessenbergDecompositionMatrixHReturnType< MatrixType > MatrixHReturnType
Definition: HessenbergDecomposition.h:90

References eigen_assert, and Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized.

Referenced by hessenberg(), Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, IsComplex >::run(), and Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, false >::run().

◆ matrixQ()

template<typename MatrixType_ >
HouseholderSequenceType Eigen::HessenbergDecomposition< MatrixType_ >::matrixQ ( ) const
inline

Reconstructs the orthogonal matrix Q in the decomposition.

Returns
object representing the matrix Q
Precondition
Either the constructor HessenbergDecomposition(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the Hessenberg decomposition of a matrix.

This function returns a light-weight object of template class HouseholderSequence. You can either apply it directly to a matrix or you can convert it to a matrix of type MatrixType.

See also
matrixH() for an example, class HouseholderSequence
225  {
226  eigen_assert(m_isInitialized && "HessenbergDecomposition is not initialized.");
227  return HouseholderSequenceType(m_matrix, m_hCoeffs.conjugate()).setLength(m_matrix.rows() - 1).setShift(1);
228  }
HouseholderSequence< MatrixType, internal::remove_all_t< typename CoeffVectorType::ConjugateReturnType > > HouseholderSequenceType
Return type of matrixQ()
Definition: HessenbergDecomposition.h:88

References eigen_assert, Eigen::HessenbergDecomposition< MatrixType_ >::m_hCoeffs, Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized, and Eigen::HessenbergDecomposition< MatrixType_ >::m_matrix.

Referenced by hessenberg(), Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, IsComplex >::run(), and Eigen::internal::complex_schur_reduce_to_hessenberg< MatrixType, false >::run().

◆ packedMatrix()

template<typename MatrixType_ >
const MatrixType& Eigen::HessenbergDecomposition< MatrixType_ >::packedMatrix ( ) const
inline

Returns the internal representation of the decomposition.

Returns
a const reference to a matrix with the internal representation of the decomposition.
Precondition
Either the constructor HessenbergDecomposition(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the Hessenberg decomposition of a matrix.

The returned matrix contains the following information:

  • the upper part and lower sub-diagonal represent the Hessenberg matrix H
  • the rest of the lower part contains the Householder vectors that, combined with Householder coefficients returned by householderCoefficients(), allows to reconstruct the matrix Q as \( Q = H_{N-1} \ldots H_1 H_0 \). Here, the matrices \( H_i \) are the Householder transformations \( H_i = (I - h_i v_i v_i^T) \) where \( h_i \) is the \( i \)th Householder coefficient and \( v_i \) is the Householder vector defined by \( v_i = [ 0, \ldots, 0, 1, M(i+2,i), \ldots, M(N-1,i) ]^T \) with M the matrix returned by this function.

See LAPACK for further details on this packed storage.

Example:

Matrix4d A = Matrix4d::Random(4, 4);
cout << "Here is a random 4x4 matrix:" << endl << A << endl;
HessenbergDecomposition<Matrix4d> hessOfA(A);
Matrix4d pm = hessOfA.packedMatrix();
cout << "The packed matrix M is:" << endl << pm << endl;
cout << "The upper Hessenberg part corresponds to the matrix H, which is:" << endl << hessOfA.matrixH() << endl;
Vector3d hc = hessOfA.householderCoefficients();
cout << "The vector of Householder coefficients is:" << endl << hc << endl;
Matrix4d pm
Definition: HessenbergDecomposition_packedMatrix.cpp:4
Vector3d hc
Definition: Tridiagonalization_householderCoefficients.cpp:5

Output:

See also
householderCoefficients()
206  {
207  eigen_assert(m_isInitialized && "HessenbergDecomposition is not initialized.");
208  return m_matrix;
209  }

References eigen_assert, Eigen::HessenbergDecomposition< MatrixType_ >::m_isInitialized, and Eigen::HessenbergDecomposition< MatrixType_ >::m_matrix.

Referenced by hessenberg().

Member Data Documentation

◆ m_hCoeffs

◆ m_isInitialized

◆ m_matrix

◆ m_temp


The documentation for this class was generated from the following file: