Eigen::Tridiagonalization< MatrixType_ > Class Template Reference

Tridiagonal decomposition of a selfadjoint matrix. More...

#include <Tridiagonalization.h>

Public Types

enum  {
  Size = MatrixType::RowsAtCompileTime , SizeMinusOne = Size == Dynamic ? Dynamic : (Size > 1 ? Size - 1 : 1) , Options = internal::traits<MatrixType>::Options , MaxSize = MatrixType::MaxRowsAtCompileTime ,
  MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
}
 
typedef MatrixType_ MatrixType
 Synonym for the template parameter MatrixType_. More...
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Eigen::Index Index
 
typedef Matrix< Scalar, SizeMinusOne, 1, Options &~RowMajor, MaxSizeMinusOne, 1 > CoeffVectorType
 
typedef internal::plain_col_type< MatrixType, RealScalar >::type DiagonalType
 
typedef Matrix< RealScalar, SizeMinusOne, 1, Options &~RowMajor, MaxSizeMinusOne, 1 > SubDiagonalType
 
typedef internal::remove_all_t< typename MatrixType::RealReturnType > MatrixTypeRealView
 
typedef internal::TridiagonalizationMatrixTReturnType< MatrixTypeRealViewMatrixTReturnType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, internal::add_const_on_value_type_t< typename Diagonal< const MatrixType >::RealReturnType >, const Diagonal< const MatrixType > > DiagonalReturnType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, internal::add_const_on_value_type_t< typename Diagonal< const MatrixType, -1 >::RealReturnType >, const Diagonal< const MatrixType, -1 > > SubDiagonalReturnType
 
typedef HouseholderSequence< MatrixType, internal::remove_all_t< typename CoeffVectorType::ConjugateReturnType > > HouseholderSequenceType
 Return type of matrixQ() More...
 

Public Member Functions

 Tridiagonalization (Index size=Size==Dynamic ? 2 :Size)
 Default constructor. More...
 
template<typename InputType >
 Tridiagonalization (const EigenBase< InputType > &matrix)
 Constructor; computes tridiagonal decomposition of given matrix. More...
 
template<typename InputType >
Tridiagonalizationcompute (const EigenBase< InputType > &matrix)
 Computes tridiagonal decomposition of given matrix. More...
 
CoeffVectorType householderCoefficients () const
 Returns the Householder coefficients. More...
 
const MatrixTypepackedMatrix () const
 Returns the internal representation of the decomposition. More...
 
HouseholderSequenceType matrixQ () const
 Returns the unitary matrix Q in the decomposition. More...
 
MatrixTReturnType matrixT () const
 Returns an expression of the tridiagonal matrix T in the decomposition. More...
 
DiagonalReturnType diagonal () const
 Returns the diagonal of the tridiagonal matrix T in the decomposition. More...
 
SubDiagonalReturnType subDiagonal () const
 Returns the subdiagonal of the tridiagonal matrix T in the decomposition. More...
 

Protected Attributes

MatrixType m_matrix
 
CoeffVectorType m_hCoeffs
 
bool m_isInitialized
 

Detailed Description

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

Tridiagonal decomposition of a selfadjoint matrix.

\eigenvalues_module

Template Parameters
MatrixType_the type of the matrix of which we are computing the tridiagonal decomposition; this is expected to be an instantiation of the Matrix class template.

This class performs a tridiagonal decomposition of a selfadjoint matrix \( A \) such that: \( A = Q T Q^* \) where \( Q \) is unitary and \( T \) a real symmetric tridiagonal matrix.

A tridiagonal matrix is a matrix which has nonzero elements only on the main diagonal and the first diagonal below and above it. The Hessenberg decomposition of a selfadjoint matrix is in fact a tridiagonal decomposition. This class is used in SelfAdjointEigenSolver to compute the eigenvalues and eigenvectors of a selfadjoint matrix.

Call the function compute() to compute the tridiagonal decomposition of a given matrix. Alternatively, you can use the Tridiagonalization(const MatrixType&) constructor which computes the tridiagonal Schur decomposition at construction time. Once the decomposition is computed, you can use the matrixQ() and matrixT() functions to retrieve the matrices Q and T in the decomposition.

The documentation of Tridiagonalization(const MatrixType&) contains an example of the typical use of this class.

See also
class HessenbergDecomposition, class SelfAdjointEigenSolver

Member Typedef Documentation

◆ CoeffVectorType

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

◆ DiagonalReturnType

template<typename MatrixType_ >
typedef std::conditional_t<NumTraits<Scalar>::IsComplex, internal::add_const_on_value_type_t<typename Diagonal<const MatrixType>::RealReturnType>, const Diagonal<const MatrixType> > Eigen::Tridiagonalization< MatrixType_ >::DiagonalReturnType

◆ DiagonalType

template<typename MatrixType_ >
typedef internal::plain_col_type<MatrixType, RealScalar>::type Eigen::Tridiagonalization< MatrixType_ >::DiagonalType

◆ HouseholderSequenceType

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

Return type of matrixQ()

◆ Index

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

◆ MatrixTReturnType

◆ MatrixType

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

Synonym for the template parameter MatrixType_.

◆ MatrixTypeRealView

template<typename MatrixType_ >
typedef internal::remove_all_t<typename MatrixType::RealReturnType> Eigen::Tridiagonalization< MatrixType_ >::MatrixTypeRealView

◆ RealScalar

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

◆ Scalar

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

◆ SubDiagonalReturnType

template<typename MatrixType_ >
typedef std::conditional_t< NumTraits<Scalar>::IsComplex, internal::add_const_on_value_type_t<typename Diagonal<const MatrixType, -1>::RealReturnType>, const Diagonal<const MatrixType, -1> > Eigen::Tridiagonalization< MatrixType_ >::SubDiagonalReturnType

◆ SubDiagonalType

template<typename MatrixType_ >
typedef Matrix<RealScalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> Eigen::Tridiagonalization< MatrixType_ >::SubDiagonalType

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
Size 
SizeMinusOne 
Options 
MaxSize 
MaxSizeMinusOne 
75  {
76  Size = MatrixType::RowsAtCompileTime,
77  SizeMinusOne = Size == Dynamic ? Dynamic : (Size > 1 ? Size - 1 : 1),
78  Options = internal::traits<MatrixType>::Options,
79  MaxSize = MatrixType::MaxRowsAtCompileTime,
80  MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
81  };
@ SizeMinusOne
Definition: Tridiagonalization.h:77
@ Size
Definition: Tridiagonalization.h:76
@ Options
Definition: Tridiagonalization.h:78
@ MaxSizeMinusOne
Definition: Tridiagonalization.h:80
@ MaxSize
Definition: Tridiagonalization.h:79
const int Dynamic
Definition: Constants.h:25

Constructor & Destructor Documentation

◆ Tridiagonalization() [1/2]

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

Default constructor.

Parameters
[in]sizePositive integer, size of the matrix whose tridiagonal 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.
116  : Size)
117  : m_matrix(size, size), m_hCoeffs(size > 1 ? size - 1 : 1), m_isInitialized(false) {}
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
bool m_isInitialized
Definition: Tridiagonalization.h:291
MatrixType m_matrix
Definition: Tridiagonalization.h:289
CoeffVectorType m_hCoeffs
Definition: Tridiagonalization.h:290

◆ Tridiagonalization() [2/2]

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

Constructor; computes tridiagonal decomposition of given matrix.

Parameters
[in]matrixSelfadjoint matrix whose tridiagonal decomposition is to be computed.

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

Example:

MatrixXd X = MatrixXd::Random(5, 5);
MatrixXd A = X + X.transpose();
cout << "Here is a random symmetric 5x5 matrix:" << endl << A << endl << endl;
Tridiagonalization<MatrixXd> triOfA(A);
MatrixXd Q = triOfA.matrixQ();
cout << "The orthogonal matrix Q is:" << endl << Q << endl;
MatrixXd T = triOfA.matrixT();
cout << "The tridiagonal matrix T is:" << endl << T << endl << endl;
cout << "Q * T * Q^T = " << endl << Q * T * Q.transpose() << endl;
MatrixXf Q
Definition: HouseholderQR_householderQ.cpp:1
cout<< "Here is a random self-adjoint 4x4 matrix:"<< endl<< A<< endl<< endl;Tridiagonalization< MatrixXcd > triOfA(A)
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
#define X
Definition: icosphere.cpp:20

Output:

 
131  : m_matrix(matrix.derived()), m_hCoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1), m_isInitialized(false) {
133  m_isInitialized = true;
134  }
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
EIGEN_DEVICE_FUNC void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
Definition: Tridiagonalization.h:332

References Eigen::Tridiagonalization< MatrixType_ >::m_hCoeffs, Eigen::Tridiagonalization< MatrixType_ >::m_isInitialized, Eigen::Tridiagonalization< MatrixType_ >::m_matrix, and Eigen::internal::tridiagonalization_inplace().

Member Function Documentation

◆ compute()

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

Computes tridiagonal decomposition of given matrix.

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

The tridiagonal decomposition is computed by bringing the columns of the matrix successively in the required form using Householder reflections. The cost is \( 4n^3/3 \) flops, where \( n \) denotes the size of the given matrix.

This method reuses of the allocated data in the Tridiagonalization object, if the size of the matrix does not change.

Example:

Tridiagonalization<MatrixXf> tri;
MatrixXf X = MatrixXf::Random(4, 4);
MatrixXf A = X + X.transpose();
tri.compute(A);
cout << "The matrix T in the tridiagonal decomposition of A is: " << endl;
cout << tri.matrixT() << endl;
tri.compute(2 * A); // re-use tri to compute eigenvalues of 2A
cout << "The matrix T in the tridiagonal decomposition of 2A is: " << endl;
cout << tri.matrixT() << endl;
Tridiagonalization< MatrixXf > tri
Definition: Tridiagonalization_compute.cpp:1

Output:

 
154  {
155  m_matrix = matrix.derived();
156  m_hCoeffs.resize(matrix.rows() - 1, 1);
158  m_isInitialized = true;
159  return *this;
160  }
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294

References Eigen::Tridiagonalization< MatrixType_ >::m_hCoeffs, Eigen::Tridiagonalization< MatrixType_ >::m_isInitialized, Eigen::Tridiagonalization< MatrixType_ >::m_matrix, matrix(), Eigen::PlainObjectBase< Derived >::resize(), and Eigen::internal::tridiagonalization_inplace().

Referenced by ctms_decompositions().

◆ diagonal()

template<typename MatrixType >
Tridiagonalization< MatrixType >::DiagonalReturnType Eigen::Tridiagonalization< MatrixType >::diagonal

Returns the diagonal of the tridiagonal matrix T in the decomposition.

Returns
expression representing the diagonal of T
Precondition
Either the constructor Tridiagonalization(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the tridiagonal decomposition of a matrix.

Example:

MatrixXcd X = MatrixXcd::Random(4, 4);
MatrixXcd A = X + X.adjoint();
cout << "Here is a random self-adjoint 4x4 matrix:" << endl << A << endl << endl;
Tridiagonalization<MatrixXcd> triOfA(A);
MatrixXd T = triOfA.matrixT();
cout << "The tridiagonal matrix T is:" << endl << T << endl << endl;
cout << "We can also extract the diagonals of T directly ..." << endl;
VectorXd diag = triOfA.diagonal();
cout << "The diagonal is:" << endl << diag << endl;
VectorXd subdiag = triOfA.subDiagonal();
cout << "The subdiagonal is:" << endl << subdiag << endl;
const char const char const char * diag
Definition: level2_impl.h:86

Output:

See also
matrixT(), subDiagonal()
295  {
296  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
297  return m_matrix.diagonal().real();
298 }
#define eigen_assert(x)
Definition: Macros.h:910

References eigen_assert.

Referenced by selfadjointeigensolver().

◆ householderCoefficients()

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

Returns the Householder coefficients.

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

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

Example:

Matrix4d X = Matrix4d::Random(4, 4);
Matrix4d A = X + X.transpose();
cout << "Here is a random symmetric 4x4 matrix:" << endl << A << endl;
Tridiagonalization<Matrix4d> triOfA(A);
Vector3d hc = triOfA.householderCoefficients();
cout << "The vector of Householder coefficients is:" << endl << hc << endl;
Vector3d hc
Definition: Tridiagonalization_householderCoefficients.cpp:5

Output:

See also
packedMatrix(), Householder module
178  {
179  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
180  return m_hCoeffs;
181  }

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

◆ matrixQ()

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

Returns the unitary matrix Q in the decomposition.

Returns
object representing the matrix Q
Precondition
Either the constructor Tridiagonalization(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the tridiagonal 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
Tridiagonalization(const MatrixType&) for an example, matrixT(), class HouseholderSequence
234  {
235  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
236  return HouseholderSequenceType(m_matrix, m_hCoeffs.conjugate()).setLength(m_matrix.rows() - 1).setShift(1);
237  }
HouseholderSequence< MatrixType, internal::remove_all_t< typename CoeffVectorType::ConjugateReturnType > > HouseholderSequenceType
Return type of matrixQ()
Definition: Tridiagonalization.h:102

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

Referenced by selfadjointeigensolver().

◆ matrixT()

template<typename MatrixType_ >
MatrixTReturnType Eigen::Tridiagonalization< MatrixType_ >::matrixT ( ) const
inline

Returns an expression of the tridiagonal matrix T in the decomposition.

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

Currently, this function can be used to extract the matrix T from internal data and copy it to a dense matrix object. In most cases, it may be sufficient to directly use the packed matrix or the vector expressions returned by diagonal() and subDiagonal() instead of creating a new dense copy matrix with this function.

See also
Tridiagonalization(const MatrixType&) for an example, matrixQ(), packedMatrix(), diagonal(), subDiagonal()
256  {
257  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
258  return MatrixTReturnType(m_matrix.real());
259  }
internal::TridiagonalizationMatrixTReturnType< MatrixTypeRealView > MatrixTReturnType
Definition: Tridiagonalization.h:87

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

Referenced by selfadjointeigensolver().

◆ packedMatrix()

template<typename MatrixType_ >
const MatrixType& Eigen::Tridiagonalization< 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 Tridiagonalization(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the tridiagonal decomposition of a matrix.

The returned matrix contains the following information:

  • the strict upper triangular part is equal to the input matrix A.
  • the diagonal and lower sub-diagonal represent the real tridiagonal symmetric matrix T.
  • 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 X = Matrix4d::Random(4, 4);
Matrix4d A = X + X.transpose();
cout << "Here is a random symmetric 4x4 matrix:" << endl << A << endl;
Tridiagonalization<Matrix4d> triOfA(A);
Matrix4d pm = triOfA.packedMatrix();
cout << "The packed matrix M is:" << endl << pm << endl;
cout << "The diagonal and subdiagonal corresponds to the matrix T, which is:" << endl << triOfA.matrixT() << endl;
Matrix4d pm
Definition: HessenbergDecomposition_packedMatrix.cpp:4

Output:

See also
householderCoefficients()
214  {
215  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
216  return m_matrix;
217  }

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

◆ subDiagonal()

template<typename MatrixType >
Tridiagonalization< MatrixType >::SubDiagonalReturnType Eigen::Tridiagonalization< MatrixType >::subDiagonal

Returns the subdiagonal of the tridiagonal matrix T in the decomposition.

Returns
expression representing the subdiagonal of T
Precondition
Either the constructor Tridiagonalization(const MatrixType&) or the member function compute(const MatrixType&) has been called before to compute the tridiagonal decomposition of a matrix.
See also
diagonal() for an example, matrixT()
301  {
302  eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
303  return m_matrix.template diagonal<-1>().real();
304 }
DiagonalReturnType diagonal() const
Returns the diagonal of the tridiagonal matrix T in the decomposition.
Definition: Tridiagonalization.h:295
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:486

References diagonal(), eigen_assert, and Eigen::real().

Referenced by selfadjointeigensolver().

Member Data Documentation

◆ m_hCoeffs

◆ m_isInitialized

◆ m_matrix


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