Performs a complex Schur decomposition of a real or complex square matrix. More...
#include <ComplexSchur.h>
Public Types | |
enum | { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime } |
typedef _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
Scalar type for matrices of type _MatrixType . | |
typedef NumTraits< Scalar >::Real | RealScalar |
typedef MatrixType::Index | Index |
typedef std::complex< RealScalar > | ComplexScalar |
Complex scalar type for _MatrixType . | |
typedef Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > | ComplexMatrixType |
Type for the matrices in the Schur decomposition. | |
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. | |
const ComplexMatrixType & | matrixU () const |
Returns the unitary matrix in the Schur decomposition. | |
const ComplexMatrixType & | matrixT () const |
Returns the triangular matrix in the Schur decomposition. | |
ComplexSchur & | compute (const MatrixType &matrix, bool computeU=true) |
Computes Schur decomposition of given matrix. | |
ComputationInfo | info () const |
Reports whether previous computation was successful. | |
Static Public Attributes | |
static const int | m_maxIterations = 30 |
Maximum number of iterations. | |
Protected Attributes | |
ComplexMatrixType | m_matT |
ComplexMatrixType | m_matU |
HessenbergDecomposition < MatrixType > | m_hess |
ComputationInfo | m_info |
bool | m_isInitialized |
bool | m_matUisUptodate |
Friends | |
struct | internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex > |
Performs a complex Schur decomposition of a real or complex square matrix.
_MatrixType | the 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: 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.
typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> ComplexSchur< _MatrixType >::ComplexMatrixType |
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
.
typedef std::complex<RealScalar> ComplexSchur< _MatrixType >::ComplexScalar |
typedef MatrixType::Index ComplexSchur< _MatrixType >::Index |
typedef _MatrixType ComplexSchur< _MatrixType >::MatrixType |
typedef NumTraits<Scalar>::Real ComplexSchur< _MatrixType >::RealScalar |
typedef MatrixType::Scalar ComplexSchur< _MatrixType >::Scalar |
Scalar type for matrices of type _MatrixType
.
anonymous enum |
ComplexSchur< _MatrixType >::ComplexSchur | ( | Index | size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime |
) | [inline] |
Default constructor.
[in] | size | Positive 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.
ComplexSchur< _MatrixType >::ComplexSchur | ( | const MatrixType & | matrix, | |
bool | computeU = true | |||
) | [inline] |
Constructor; computes Schur decomposition of given matrix.
[in] | matrix | Square matrix whose Schur decomposition is to be computed. |
[in] | computeU | If true, both T and U are computed; if false, only T is computed. |
This constructor calls compute() to compute the Schur decomposition.
ComplexSchur& ComplexSchur< _MatrixType >::compute | ( | const MatrixType & | matrix, | |
bool | computeU = true | |||
) |
Computes Schur decomposition of given matrix.
[in] | matrix | Square matrix whose Schur decomposition is to be computed. |
[in] | computeU | If true, both T and U are computed; if false, only T is computed. |
*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 complex flops, or complex flops if computeU is false.
Example:
Output:
ComputationInfo ComplexSchur< _MatrixType >::info | ( | ) | const [inline] |
Reports whether previous computation was successful.
Success
if computation was succesful, NoConvergence
otherwise. const ComplexMatrixType& ComplexSchur< _MatrixType >::matrixT | ( | ) | const [inline] |
Returns the triangular matrix in the Schur decomposition.
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:
const ComplexMatrixType& ComplexSchur< _MatrixType >::matrixU | ( | ) | const [inline] |
Returns the unitary matrix in the Schur decomposition.
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:
friend struct internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex > [friend] |
HessenbergDecomposition<MatrixType> ComplexSchur< _MatrixType >::m_hess [protected] |
ComputationInfo ComplexSchur< _MatrixType >::m_info [protected] |
bool ComplexSchur< _MatrixType >::m_isInitialized [protected] |
ComplexMatrixType ComplexSchur< _MatrixType >::m_matT [protected] |
ComplexMatrixType ComplexSchur< _MatrixType >::m_matU [protected] |
bool ComplexSchur< _MatrixType >::m_matUisUptodate [protected] |
const int ComplexSchur< _MatrixType >::m_maxIterations = 30 [static] |
Maximum number of iterations.
Maximum number of iterations allowed for an eigenvalue to converge.