Public Types | Public Member Functions | Static Public Attributes | Protected Attributes | Private Member Functions | Friends
Eigen::ComplexSchur< _MatrixType > Class Template Reference

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

#include <ComplexSchur.h>

List of all members.

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
typedef Matrix< ComplexScalar,
RowsAtCompileTime,
ColsAtCompileTime, Options,
MaxRowsAtCompileTime,
MaxColsAtCompileTime
ComplexMatrixType
 Type for the matrices in the Schur decomposition.
typedef std::complex< RealScalarComplexScalar
 Complex scalar type for _MatrixType.
typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef NumTraits< Scalar >::Real RealScalar
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type _MatrixType.

Public Member Functions

 ComplexSchur (Index size=RowsAtCompileTime==Dynamic?1:RowsAtCompileTime)
 Default constructor.
 ComplexSchur (const MatrixType &matrix, bool computeU=true)
 Constructor; computes Schur decomposition of given matrix.
ComplexSchurcompute (const MatrixType &matrix, bool computeU=true)
 Computes Schur decomposition of given matrix.
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchurcomputeFromHessenberg (const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
 Compute Schur decomposition from a given Hessenberg matrix.
Index getMaxIterations ()
 Returns the maximum number of iterations.
ComputationInfo info () const
 Reports whether previous computation was successful.
const ComplexMatrixTypematrixT () const
 Returns the triangular matrix in the Schur decomposition.
const ComplexMatrixTypematrixU () const
 Returns the unitary matrix in the Schur decomposition.
ComplexSchursetMaxIterations (Index maxIters)
 Sets the maximum number of iterations allowed.

Static Public Attributes

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

Protected Attributes

HessenbergDecomposition
< MatrixType
m_hess
ComputationInfo m_info
bool m_isInitialized
ComplexMatrixType m_matT
ComplexMatrixType m_matU
bool m_matUisUptodate
Index m_maxIters

Private Member Functions

ComplexScalar computeShift (Index iu, Index iter)
void reduceToTriangularForm (bool computeU)
bool subdiagonalEntryIsNeglegible (Index i)

Friends

struct internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex >

Detailed Description

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

Performs a complex Schur decomposition of a real or complex square matrix.

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

Given a real or complex square matrix A, this class computes the Schur decomposition: $ A = U T U^*$ where U is a unitary complex matrix, and T is a complex upper triangular matrix. The diagonal of the matrix T corresponds to the eigenvalues of the matrix A.

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

Note:
This code is inspired from Jampack
See also:
class RealSchur, class EigenSolver, class ComplexEigenSolver

Definition at line 51 of file ComplexSchur.h.


Member Typedef Documentation

Type for the matrices in the Schur decomposition.

This is a square matrix with entries of type ComplexScalar. The size is the same as the size of _MatrixType.

Definition at line 81 of file ComplexSchur.h.

template<typename _MatrixType>
typedef std::complex<RealScalar> Eigen::ComplexSchur< _MatrixType >::ComplexScalar

Complex scalar type for _MatrixType.

This is std::complex<Scalar> if Scalar is real (e.g., float or double) and just Scalar if Scalar is complex.

Definition at line 74 of file ComplexSchur.h.

template<typename _MatrixType>
typedef MatrixType::Index Eigen::ComplexSchur< _MatrixType >::Index

Definition at line 66 of file ComplexSchur.h.

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

Definition at line 54 of file ComplexSchur.h.

template<typename _MatrixType>
typedef NumTraits<Scalar>::Real Eigen::ComplexSchur< _MatrixType >::RealScalar

Definition at line 65 of file ComplexSchur.h.

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

Scalar type for matrices of type _MatrixType.

Definition at line 64 of file ComplexSchur.h.


Member Enumeration Documentation

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

Definition at line 55 of file ComplexSchur.h.


Constructor & Destructor Documentation

template<typename _MatrixType>
Eigen::ComplexSchur< _MatrixType >::ComplexSchur ( Index  size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) [inline]

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 94 of file ComplexSchur.h.

template<typename _MatrixType>
Eigen::ComplexSchur< _MatrixType >::ComplexSchur ( const MatrixType matrix,
bool  computeU = true 
) [inline]

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

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

See also:
matrixT() and matrixU() for examples.

Definition at line 112 of file ComplexSchur.h.


Member Function Documentation

template<typename MatrixType >
ComplexSchur< MatrixType > & Eigen::ComplexSchur< MatrixType >::compute ( const MatrixType 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 QR iterations with a single shift. The cost of computing the Schur decomposition depends on the number of iterations; as a rough guide, it may be taken on the number of iterations; as a rough guide, it may be taken to be $25n^3$ complex flops, or $10n^3$ complex flops if computeU is false.

Example:

Output:

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

Definition at line 316 of file ComplexSchur.h.

template<typename MatrixType >
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchur< MatrixType > & Eigen::ComplexSchur< MatrixType >::computeFromHessenberg ( const HessMatrixType &  matrixH,
const OrthMatrixType &  matrixQ,
bool  computeU = true 
)

Compute Schur decomposition from a given Hessenberg matrix.

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)

Definition at line 338 of file ComplexSchur.h.

template<typename MatrixType >
ComplexSchur< MatrixType >::ComplexScalar Eigen::ComplexSchur< MatrixType >::computeShift ( Index  iu,
Index  iter 
) [private]

Compute the shift in the current QR iteration.

Definition at line 279 of file ComplexSchur.h.

template<typename _MatrixType>
Index Eigen::ComplexSchur< _MatrixType >::getMaxIterations ( ) [inline]

Returns the maximum number of iterations.

Definition at line 233 of file ComplexSchur.h.

template<typename _MatrixType>
ComputationInfo Eigen::ComplexSchur< _MatrixType >::info ( ) const [inline]

Reports whether previous computation was successful.

Returns:
Success if computation was succesful, NoConvergence otherwise.

Definition at line 215 of file ComplexSchur.h.

template<typename _MatrixType>
const ComplexMatrixType& Eigen::ComplexSchur< _MatrixType >::matrixT ( ) const [inline]

Returns the triangular matrix in the Schur decomposition.

Returns:
A const reference to the matrix T.

It is assumed that either the constructor ComplexSchur(const MatrixType& matrix, bool computeU) or the member function compute(const MatrixType& matrix, bool computeU) has been called before to compute the Schur decomposition of a matrix.

Note that this function returns a plain square matrix. If you want to reference only the upper triangular part, use:

 schur.matrixT().triangularView<Upper>() 

Example:

Output:

Definition at line 161 of file ComplexSchur.h.

template<typename _MatrixType>
const ComplexMatrixType& Eigen::ComplexSchur< _MatrixType >::matrixU ( ) const [inline]

Returns the unitary matrix in the Schur decomposition.

Returns:
A const reference to the matrix U.

It is assumed that either the constructor ComplexSchur(const MatrixType& matrix, bool computeU) or the member function compute(const MatrixType& matrix, bool computeU) has been called before to compute the Schur decomposition of a matrix, and that computeU was set to true (the default value).

Example:

Output:

Definition at line 137 of file ComplexSchur.h.

template<typename MatrixType >
void Eigen::ComplexSchur< MatrixType >::reduceToTriangularForm ( bool  computeU) [private]

Definition at line 384 of file ComplexSchur.h.

template<typename _MatrixType>
ComplexSchur& Eigen::ComplexSchur< _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 226 of file ComplexSchur.h.

template<typename MatrixType >
bool Eigen::ComplexSchur< MatrixType >::subdiagonalEntryIsNeglegible ( Index  i) [inline, private]

If m_matT(i+1,i) is neglegible in floating point arithmetic compared to m_matT(i,i) and m_matT(j,j), then set it to zero and return true, else return false.

Definition at line 264 of file ComplexSchur.h.


Friends And Related Function Documentation

template<typename _MatrixType>
friend struct internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex > [friend]

Definition at line 257 of file ComplexSchur.h.


Member Data Documentation

template<typename _MatrixType>
HessenbergDecomposition<MatrixType> Eigen::ComplexSchur< _MatrixType >::m_hess [protected]

Definition at line 247 of file ComplexSchur.h.

template<typename _MatrixType>
ComputationInfo Eigen::ComplexSchur< _MatrixType >::m_info [protected]

Definition at line 248 of file ComplexSchur.h.

template<typename _MatrixType>
bool Eigen::ComplexSchur< _MatrixType >::m_isInitialized [protected]

Definition at line 249 of file ComplexSchur.h.

template<typename _MatrixType>
ComplexMatrixType Eigen::ComplexSchur< _MatrixType >::m_matT [protected]

Definition at line 246 of file ComplexSchur.h.

template<typename _MatrixType>
ComplexMatrixType Eigen::ComplexSchur< _MatrixType >::m_matU [protected]

Definition at line 246 of file ComplexSchur.h.

template<typename _MatrixType>
bool Eigen::ComplexSchur< _MatrixType >::m_matUisUptodate [protected]

Definition at line 250 of file ComplexSchur.h.

template<typename _MatrixType>
const int Eigen::ComplexSchur< _MatrixType >::m_maxIterationsPerRow = 30 [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 30.

Definition at line 243 of file ComplexSchur.h.

template<typename _MatrixType>
Index Eigen::ComplexSchur< _MatrixType >::m_maxIters [protected]

Definition at line 251 of file ComplexSchur.h.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Sat Jun 8 2019 19:40:33