Public Types | Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | List of all members
Eigen::RealSchur< _MatrixType > Class Template Reference

Performs a real Schur decomposition of a square matrix. More...

#include <RealSchur.h>

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
 
typedef std::complex< typename NumTraits< Scalar >::RealComplexScalar
 
typedef Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
 
typedef Eigen::Index Index
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 

Public Member Functions

template<typename InputType >
RealSchur< MatrixType > & compute (const EigenBase< InputType > &matrix, bool computeU)
 
template<typename InputType >
RealSchurcompute (const EigenBase< InputType > &matrix, bool computeU=true)
 Computes Schur decomposition of given matrix. More...
 
template<typename HessMatrixType , typename OrthMatrixType >
RealSchurcomputeFromHessenberg (const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
 Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T. More...
 
template<typename HessMatrixType , typename OrthMatrixType >
RealSchur< MatrixType > & computeFromHessenberg (const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
 
Index getMaxIterations ()
 Returns the maximum number of iterations. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const MatrixTypematrixT () const
 Returns the quasi-triangular matrix in the Schur decomposition. More...
 
const MatrixTypematrixU () const
 Returns the orthogonal matrix in the Schur decomposition. More...
 
template<typename InputType >
 RealSchur (const EigenBase< InputType > &matrix, bool computeU=true)
 Constructor; computes real Schur decomposition of given matrix. More...
 
 RealSchur (Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
 Default constructor. More...
 
RealSchursetMaxIterations (Index maxIters)
 Sets the maximum number of iterations allowed. More...
 

Static Public Attributes

static const int m_maxIterationsPerRow = 40
 Maximum number of iterations per row. More...
 

Private Types

typedef Matrix< Scalar, 3, 1 > Vector3s
 

Private Member Functions

Scalar computeNormOfT ()
 
void computeShift (Index iu, Index iter, Scalar &exshift, Vector3s &shiftInfo)
 
Index findSmallSubdiagEntry (Index iu, const Scalar &considerAsZero)
 
void initFrancisQRStep (Index il, Index iu, const Vector3s &shiftInfo, Index &im, Vector3s &firstHouseholderVector)
 
void performFrancisQRStep (Index il, Index im, Index iu, bool computeU, const Vector3s &firstHouseholderVector, Scalar *workspace)
 
void splitOffTwoRows (Index iu, bool computeU, const Scalar &exshift)
 

Private Attributes

HessenbergDecomposition< MatrixTypem_hess
 
ComputationInfo m_info
 
bool m_isInitialized
 
MatrixType m_matT
 
MatrixType m_matU
 
bool m_matUisUptodate
 
Index m_maxIters
 
ColumnVectorType m_workspaceVector
 

Detailed Description

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

Performs a real Schur decomposition of a square matrix.

\eigenvalues_module

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

Given a real square matrix A, this class computes the real Schur decomposition: $ A = U T U^T $ where U is a real orthogonal matrix and T is a real quasi-triangular matrix. An orthogonal matrix is a matrix whose inverse is equal to its transpose, $ U^{-1} = U^T $. A quasi-triangular matrix is a block-triangular matrix whose diagonal consists of 1-by-1 blocks and 2-by-2 blocks with complex eigenvalues. The eigenvalues of the blocks on the diagonal of T are the same as the eigenvalues of the matrix A, and thus the real Schur decomposition is used in EigenSolver to compute the eigendecomposition of a matrix.

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

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

Note
The implementation is adapted from JAMA (public domain). Their code is based on EISPACK.
See also
class ComplexSchur, class EigenSolver, class ComplexEigenSolver

Definition at line 54 of file RealSchur.h.

Member Typedef Documentation

◆ ColumnVectorType

template<typename _MatrixType >
typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> Eigen::RealSchur< _MatrixType >::ColumnVectorType

Definition at line 70 of file RealSchur.h.

◆ ComplexScalar

template<typename _MatrixType >
typedef std::complex<typename NumTraits<Scalar>::Real> Eigen::RealSchur< _MatrixType >::ComplexScalar

Definition at line 66 of file RealSchur.h.

◆ EigenvalueType

template<typename _MatrixType >
typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> Eigen::RealSchur< _MatrixType >::EigenvalueType

Definition at line 69 of file RealSchur.h.

◆ Index

template<typename _MatrixType >
typedef Eigen::Index Eigen::RealSchur< _MatrixType >::Index
Deprecated:
since Eigen 3.3

Definition at line 67 of file RealSchur.h.

◆ MatrixType

template<typename _MatrixType >
typedef _MatrixType Eigen::RealSchur< _MatrixType >::MatrixType

Definition at line 57 of file RealSchur.h.

◆ Scalar

template<typename _MatrixType >
typedef MatrixType::Scalar Eigen::RealSchur< _MatrixType >::Scalar

Definition at line 65 of file RealSchur.h.

◆ Vector3s

template<typename _MatrixType >
typedef Matrix<Scalar,3,1> Eigen::RealSchur< _MatrixType >::Vector3s
private

Definition at line 236 of file RealSchur.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 58 of file RealSchur.h.

Constructor & Destructor Documentation

◆ RealSchur() [1/2]

template<typename _MatrixType >
Eigen::RealSchur< _MatrixType >::RealSchur ( Index  size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
inlineexplicit

Default constructor.

Parameters
[in]sizePositive integer, size of the matrix whose Schur 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.

Definition at line 83 of file RealSchur.h.

◆ RealSchur() [2/2]

template<typename _MatrixType >
template<typename InputType >
Eigen::RealSchur< _MatrixType >::RealSchur ( const EigenBase< InputType > &  matrix,
bool  computeU = true 
)
inlineexplicit

Constructor; computes real Schur decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Schur decomposition is to be computed.
[in]computeUIf true, both T and U are computed; if false, only T is computed.

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

Example:

MatrixXd A = MatrixXd::Random(6,6);
cout << "Here is a random 6x6 matrix, A:" << endl << A << endl << endl;
RealSchur<MatrixXd> schur(A);
cout << "The orthogonal matrix U is:" << endl << schur.matrixU() << endl;
cout << "The quasi-triangular matrix T is:" << endl << schur.matrixT() << endl << endl;
MatrixXd U = schur.matrixU();
MatrixXd T = schur.matrixT();
cout << "U * T * U^T = " << endl << U * T * U.transpose() << endl;

Output:

 

Definition at line 104 of file RealSchur.h.

Member Function Documentation

◆ compute() [1/2]

template<typename _MatrixType >
template<typename InputType >
RealSchur<MatrixType>& Eigen::RealSchur< _MatrixType >::compute ( const EigenBase< InputType > &  matrix,
bool  computeU 
)

Definition at line 249 of file RealSchur.h.

◆ compute() [2/2]

template<typename _MatrixType >
template<typename InputType >
RealSchur& Eigen::RealSchur< _MatrixType >::compute ( const EigenBase< InputType > &  matrix,
bool  computeU = true 
)

Computes Schur decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Schur decomposition is to be computed.
[in]computeUIf true, both T and U are computed; if false, only T is computed.
Returns
Reference to *this

The Schur decomposition is computed by first reducing the matrix to Hessenberg form using the class HessenbergDecomposition. The Hessenberg matrix is then reduced to triangular form by performing Francis QR iterations with implicit double shift. The cost of computing the Schur decomposition depends on the number of iterations; as a rough guide, it may be taken to be $25n^3$ flops if computeU is true and $10n^3$ flops if computeU is false.

Example:

MatrixXf A = MatrixXf::Random(4,4);
RealSchur<MatrixXf> schur(4);
schur.compute(A, /* computeU = */ false);
cout << "The matrix T in the decomposition of A is:" << endl << schur.matrixT() << endl;
schur.compute(A.inverse(), /* computeU = */ false);
cout << "The matrix T in the decomposition of A^(-1) is:" << endl << schur.matrixT() << endl;

Output:

See also
compute(const MatrixType&, bool, Index)

◆ computeFromHessenberg() [1/2]

template<typename _MatrixType >
template<typename HessMatrixType , typename OrthMatrixType >
RealSchur& Eigen::RealSchur< _MatrixType >::computeFromHessenberg ( const HessMatrixType &  matrixH,
const OrthMatrixType &  matrixQ,
bool  computeU 
)

Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.

Parameters
[in]matrixHMatrix in Hessenberg form H
[in]matrixQorthogonal matrix Q that transform a matrix A to H : A = Q H Q^T
computeUComputes the matriX U of the Schur vectors
Returns
Reference to *this

This routine assumes that the matrix is already reduced in Hessenberg form matrixH using either the class HessenbergDecomposition or another mean. It computes the upper quasi-triangular matrix T of the Schur decomposition of H When computeU is true, this routine computes the matrix U such that A = U T U^T = (QZ) T (QZ)^T = Q H Q^T where A is the initial matrix

NOTE Q is referenced if computeU is true; so, if the initial orthogonal matrix is not available, the user should give an identity matrix (Q.setIdentity())

See also
compute(const MatrixType&, bool)

◆ computeFromHessenberg() [2/2]

template<typename _MatrixType >
template<typename HessMatrixType , typename OrthMatrixType >
RealSchur<MatrixType>& Eigen::RealSchur< _MatrixType >::computeFromHessenberg ( const HessMatrixType &  matrixH,
const OrthMatrixType &  matrixQ,
bool  computeU 
)

Definition at line 287 of file RealSchur.h.

◆ computeNormOfT()

template<typename MatrixType >
MatrixType::Scalar Eigen::RealSchur< MatrixType >::computeNormOfT
inlineprivate

Definition at line 362 of file RealSchur.h.

◆ computeShift()

template<typename MatrixType >
void Eigen::RealSchur< MatrixType >::computeShift ( Index  iu,
Index  iter,
Scalar exshift,
Vector3s shiftInfo 
)
inlineprivate

Definition at line 430 of file RealSchur.h.

◆ findSmallSubdiagEntry()

template<typename MatrixType >
Index Eigen::RealSchur< MatrixType >::findSmallSubdiagEntry ( Index  iu,
const Scalar considerAsZero 
)
inlineprivate

Definition at line 376 of file RealSchur.h.

◆ getMaxIterations()

template<typename _MatrixType >
Index Eigen::RealSchur< _MatrixType >::getMaxIterations ( )
inline

Returns the maximum number of iterations.

Definition at line 213 of file RealSchur.h.

◆ info()

template<typename _MatrixType >
ComputationInfo Eigen::RealSchur< _MatrixType >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NoConvergence otherwise.

Definition at line 195 of file RealSchur.h.

◆ initFrancisQRStep()

template<typename MatrixType >
void Eigen::RealSchur< MatrixType >::initFrancisQRStep ( Index  il,
Index  iu,
const Vector3s shiftInfo,
Index im,
Vector3s firstHouseholderVector 
)
inlineprivate

Definition at line 472 of file RealSchur.h.

◆ matrixT()

template<typename _MatrixType >
const MatrixType& Eigen::RealSchur< _MatrixType >::matrixT ( ) const
inline

Returns the quasi-triangular matrix in the Schur decomposition.

Returns
A const reference to the matrix T.
Precondition
Either the constructor RealSchur(const MatrixType&, bool) or the member function compute(const MatrixType&, bool) has been called before to compute the Schur decomposition of a matrix.
See also
RealSchur(const MatrixType&, bool) for an example

Definition at line 144 of file RealSchur.h.

◆ matrixU()

template<typename _MatrixType >
const MatrixType& Eigen::RealSchur< _MatrixType >::matrixU ( ) const
inline

Returns the orthogonal matrix in the Schur decomposition.

Returns
A const reference to the matrix U.
Precondition
Either the constructor RealSchur(const MatrixType&, bool) or the member function compute(const MatrixType&, bool) has been called before to compute the Schur decomposition of a matrix, and computeU was set to true (the default value).
See also
RealSchur(const MatrixType&, bool) for an example

Definition at line 127 of file RealSchur.h.

◆ performFrancisQRStep()

template<typename MatrixType >
void Eigen::RealSchur< MatrixType >::performFrancisQRStep ( Index  il,
Index  im,
Index  iu,
bool  computeU,
const Vector3s firstHouseholderVector,
Scalar workspace 
)
inlineprivate

Definition at line 497 of file RealSchur.h.

◆ setMaxIterations()

template<typename _MatrixType >
RealSchur& Eigen::RealSchur< _MatrixType >::setMaxIterations ( Index  maxIters)
inline

Sets the maximum number of iterations allowed.

If not specified by the user, the maximum number of iterations is m_maxIterationsPerRow times the size of the matrix.

Definition at line 206 of file RealSchur.h.

◆ splitOffTwoRows()

template<typename MatrixType >
void Eigen::RealSchur< MatrixType >::splitOffTwoRows ( Index  iu,
bool  computeU,
const Scalar exshift 
)
inlineprivate

Definition at line 395 of file RealSchur.h.

Member Data Documentation

◆ m_hess

template<typename _MatrixType >
HessenbergDecomposition<MatrixType> Eigen::RealSchur< _MatrixType >::m_hess
private

Definition at line 230 of file RealSchur.h.

◆ m_info

template<typename _MatrixType >
ComputationInfo Eigen::RealSchur< _MatrixType >::m_info
private

Definition at line 231 of file RealSchur.h.

◆ m_isInitialized

template<typename _MatrixType >
bool Eigen::RealSchur< _MatrixType >::m_isInitialized
private

Definition at line 232 of file RealSchur.h.

◆ m_matT

template<typename _MatrixType >
MatrixType Eigen::RealSchur< _MatrixType >::m_matT
private

Definition at line 227 of file RealSchur.h.

◆ m_matU

template<typename _MatrixType >
MatrixType Eigen::RealSchur< _MatrixType >::m_matU
private

Definition at line 228 of file RealSchur.h.

◆ m_matUisUptodate

template<typename _MatrixType >
bool Eigen::RealSchur< _MatrixType >::m_matUisUptodate
private

Definition at line 233 of file RealSchur.h.

◆ m_maxIterationsPerRow

template<typename _MatrixType >
const int Eigen::RealSchur< _MatrixType >::m_maxIterationsPerRow = 40
static

Maximum number of iterations per row.

If not otherwise specified, the maximum number of iterations is this number times the size of the matrix. It is currently set to 40.

Definition at line 223 of file RealSchur.h.

◆ m_maxIters

template<typename _MatrixType >
Index Eigen::RealSchur< _MatrixType >::m_maxIters
private

Definition at line 234 of file RealSchur.h.

◆ m_workspaceVector

template<typename _MatrixType >
ColumnVectorType Eigen::RealSchur< _MatrixType >::m_workspaceVector
private

Definition at line 229 of file RealSchur.h.


The documentation for this class was generated from the following file:
A
Definition: test_numpy_dtypes.cpp:298
Eigen::Triplet
A small structure to hold a non zero as a triplet (i,j,value).
Definition: SparseUtil.h:162
schur
ComplexSchur< MatrixXcf > schur(4)
U
@ U
Definition: testDecisionTree.cpp:349


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:11:51