Public Types | Public Member Functions | Protected Attributes
Tridiagonalization< _MatrixType > Class Template Reference

Tridiagonal decomposition of a selfadjoint matrix. More...

#include <Tridiagonalization.h>

List of all members.

Public Types

enum  {
  Size,
  SizeMinusOne,
  Options,
  MaxSize,
  MaxSizeMinusOne
}
typedef Matrix< Scalar,
SizeMinusOne, 1, Options
&~RowMajor, MaxSizeMinusOne, 1 > 
CoeffVectorType
typedef internal::conditional
< NumTraits< Scalar >
::IsComplex, typename
internal::add_const_on_value_type
< typename Diagonal< const
MatrixType >::RealReturnType >
::type, const Diagonal< const
MatrixType > >::type 
DiagonalReturnType
typedef
internal::plain_col_type
< MatrixType, RealScalar >
::type 
DiagonalType
typedef HouseholderSequence
< MatrixType, CoeffVectorType >
::ConjugateReturnType 
HouseholderSequenceType
 Return type of matrixQ()
typedef MatrixType::Index Index
typedef
internal::TridiagonalizationMatrixTReturnType
< MatrixTypeRealView
MatrixTReturnType
typedef _MatrixType MatrixType
 Synonym for the template parameter _MatrixType.
typedef internal::remove_all
< typename
MatrixType::RealReturnType >
::type 
MatrixTypeRealView
typedef NumTraits< Scalar >::Real RealScalar
typedef MatrixType::Scalar Scalar
typedef internal::conditional
< NumTraits< Scalar >
::IsComplex, typename
internal::add_const_on_value_type
< typename Diagonal< Block
< const MatrixType,
SizeMinusOne, SizeMinusOne >
>::RealReturnType >::type,
const Diagonal< Block< const
MatrixType, SizeMinusOne,
SizeMinusOne > > >::type 
SubDiagonalReturnType
typedef Matrix< RealScalar,
SizeMinusOne, 1, Options
&~RowMajor, MaxSizeMinusOne, 1 > 
SubDiagonalType

Public Member Functions

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

Protected Attributes

CoeffVectorType m_hCoeffs
bool m_isInitialized
MatrixType m_matrix

Detailed Description

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

Tridiagonal decomposition of a selfadjoint matrix.

This is defined in the Eigenvalues module.

#include <Eigen/Eigenvalues>
Template Parameters:
_MatrixTypethe 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

typedef internal::conditional<NumTraits<Scalar>::IsComplex, typename internal::add_const_on_value_type<typename Diagonal<const MatrixType>::RealReturnType>::type, const Diagonal<const MatrixType> >::type DiagonalReturnType
typedef internal::plain_col_type<MatrixType, RealScalar>::type DiagonalType

Return type of matrixQ()

typedef MatrixType::Index Index
typedef internal::TridiagonalizationMatrixTReturnType<MatrixTypeRealView> MatrixTReturnType
typedef _MatrixType MatrixType

Synonym for the template parameter _MatrixType.

typedef internal::remove_all<typename MatrixType::RealReturnType>::type MatrixTypeRealView
typedef NumTraits<Scalar>::Real RealScalar
typedef MatrixType::Scalar Scalar
typedef internal::conditional<NumTraits<Scalar>::IsComplex, typename internal::add_const_on_value_type<typename Diagonal< Block<const MatrixType,SizeMinusOne,SizeMinusOne> >::RealReturnType>::type, const Diagonal< Block<const MatrixType,SizeMinusOne,SizeMinusOne> > >::type SubDiagonalReturnType

Member Enumeration Documentation

anonymous enum
Enumerator:
Size 
SizeMinusOne 
Options 
MaxSize 
MaxSizeMinusOne 

Constructor & Destructor Documentation

Tridiagonalization ( Index  size = Size==Dynamic ? 2 : Size)
inline

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.
Tridiagonalization ( const MatrixType matrix)
inline

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;

Output:

Here is a random symmetric 5x5 matrix:
  1.36 -0.816  0.521   1.43 -0.144
-0.816 -0.659  0.794 -0.173 -0.406
 0.521  0.794 -0.541  0.461  0.179
  1.43 -0.173  0.461  -1.43  0.822
-0.144 -0.406  0.179  0.822  -1.37

The orthogonal matrix Q is:
       1        0        0        0        0
       0   -0.471    0.127   -0.671   -0.558
       0    0.301   -0.195    0.437   -0.825
       0    0.825   0.0459   -0.563 -0.00872
       0  -0.0832   -0.971   -0.202   0.0922
The tridiagonal matrix T is:
  1.36   1.73      0      0      0
  1.73   -1.2 -0.966      0      0
     0 -0.966  -1.28  0.214      0
     0      0  0.214  -1.69  0.345
     0      0      0  0.345  0.164

Q * T * Q^T = 
  1.36 -0.816  0.521   1.43 -0.144
-0.816 -0.659  0.794 -0.173 -0.406
 0.521  0.794 -0.541  0.461  0.179
  1.43 -0.173  0.461  -1.43  0.822
-0.144 -0.406  0.179  0.822  -1.37

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


Member Function Documentation

Tridiagonalization& compute ( const MatrixType 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 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;

Output:

The matrix T in the tridiagonal decomposition of A is: 
  1.36 -0.704      0      0
-0.704 0.0147   1.71      0
     0   1.71  0.856  0.641
     0      0  0.641 -0.506
The matrix T in the tridiagonal decomposition of 2A is: 
  2.72  -1.41      0      0
 -1.41 0.0294   3.43      0
     0   3.43   1.71   1.28
     0      0   1.28  -1.01

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

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 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;

Output:

Here is a random self-adjoint 4x4 matrix:
    (-0.422,0)  (0.705,-1.01) (-0.17,-0.552) (0.338,-0.357)
  (0.705,1.01)      (0.515,0) (0.241,-0.446)   (0.05,-1.64)
 (-0.17,0.552)  (0.241,0.446)      (-1.03,0)  (0.0449,1.72)
 (0.338,0.357)    (0.05,1.64) (0.0449,-1.72)       (1.36,0)

The tridiagonal matrix T is:
-0.422 -1.45     0     0
-1.45  1.01 -1.42     0
    0 -1.42   1.8  -1.2
    0     0  -1.2 -1.96

We can also extract the diagonals of T directly ...
The diagonal is:
-0.422
1.01
1.8
-1.96
The subdiagonal is:
-1.45
-1.42
-1.2
See also:
matrixT(), subDiagonal()

References eigen_assert.

CoeffVectorType householderCoefficients ( ) const
inline

Returns the Householder coefficients.

   \returns a const reference to the vector of Householder coefficients

   \pre 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: \include Tridiagonalization_householderCoefficients.cpp
   Output: \verbinclude Tridiagonalization_householderCoefficients.out

   \sa packedMatrix(), \ref Householder_Module "Householder module"

References eigen_assert, Tridiagonalization< _MatrixType >::m_hCoeffs, and Tridiagonalization< _MatrixType >::m_isInitialized.

HouseholderSequenceType 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

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

MatrixTReturnType 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()

References eigen_assert, Tridiagonalization< _MatrixType >::m_isInitialized, and Tridiagonalization< _MatrixType >::m_matrix.

const 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.

   \pre 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 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;

Output:

Here is a random symmetric 4x4 matrix:
   1.36   0.612   0.122   0.326
  0.612   -1.21  -0.222   0.563
  0.122  -0.222 -0.0904    1.16
  0.326   0.563    1.16    1.66
The packed matrix M is:
  1.36  0.612  0.122  0.326
-0.704 0.0147 -0.222  0.563
0.0925   1.71  0.856   1.16
 0.248  0.785  0.641 -0.506
The diagonal and subdiagonal corresponds to the matrix T, which is:
  1.36 -0.704      0      0
-0.704 0.0147   1.71      0
     0   1.71  0.856  0.641
     0      0  0.641 -0.506
See also:
householderCoefficients()

References eigen_assert, Tridiagonalization< _MatrixType >::m_isInitialized, and Tridiagonalization< _MatrixType >::m_matrix.

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()

References eigen_assert.


Member Data Documentation

CoeffVectorType m_hCoeffs
protected
bool m_isInitialized
protected
MatrixType m_matrix
protected

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