Public Types | Public Member Functions | Static Public Attributes | Protected Attributes | List of all members
ComplexSchur< _MatrixType > Class Template Reference

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

#include <ComplexSchur.h>

Public Types

enum  {
  RowsAtCompileTime,
  ColsAtCompileTime,
  Options,
  MaxRowsAtCompileTime,
  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.
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.

Static Public Attributes

static const int m_maxIterations
 Maximum number of iterations.

Protected Attributes

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

Detailed Description

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

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

This is defined in the Eigenvalues module.

#include <Eigen/Eigenvalues>
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

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.

typedef std::complex<RealScalar> 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.

typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef NumTraits<Scalar>::Real RealScalar
typedef MatrixType::Scalar Scalar

Scalar type for matrices of type _MatrixType.

Member Enumeration Documentation

anonymous enum
Enumerator:
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Constructor & Destructor Documentation

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

Member Function Documentation

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:

ComplexSchur<MatrixXcf> schur(4);
schur.compute(A);
cout << "The matrix T in the decomposition of A is:" << endl << schur.matrixT() << endl;
schur.compute(A.inverse());
cout << "The matrix T in the decomposition of A^(-1) is:" << endl << schur.matrixT() << endl;

Output:

The matrix T in the decomposition of A is:
 (-0.691,-1.63)  (0.763,-0.144) (-0.104,-0.836) (-0.462,-0.378)
          (0,0)   (-0.758,1.22)  (-0.65,-0.772)  (-0.244,0.113)
          (0,0)           (0,0)   (0.137,0.505) (0.0687,-0.404)
          (0,0)           (0,0)           (0,0)   (1.52,-0.402)
The matrix T in the decomposition of A^(-1) is:
    (0.501,-1.84)    (-1.01,-0.984)       (0.636,1.3)    (-0.676,0.352)
            (0,0)   (-0.369,-0.593)     (0.0733,0.18) (-0.0658,-0.0263)
            (0,0)             (0,0)    (-0.222,0.521)    (-0.191,0.121)
            (0,0)             (0,0)             (0,0)     (0.614,0.162)

References eigen_assert, Eigen::internal::IsComplex, ComplexSchur< _MatrixType >::m_info, ComplexSchur< _MatrixType >::m_isInitialized, ComplexSchur< _MatrixType >::m_matT, ComplexSchur< _MatrixType >::m_matU, ComplexSchur< _MatrixType >::m_matUisUptodate, and Eigen::Success.

Referenced by ComplexSchur< MatrixType >::ComplexSchur().

ComputationInfo info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NoConvergence otherwise.

Referenced by ComplexEigenSolver< _MatrixType >::info().

const ComplexMatrixType& 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:

cout << "Here is a random 4x4 matrix, A:" << endl << A << endl << endl;
ComplexSchur<MatrixXcf> schurOfA(A, false); // false means do not compute U
cout << "The triangular matrix T is:" << endl << schurOfA.matrixT() << endl;

Output:

Here is a random 4x4 matrix, A:
  (-0.211,0.68)  (0.108,-0.444)   (0.435,0.271) (-0.198,-0.687)
  (0.597,0.566) (0.258,-0.0452)  (0.214,-0.717)  (-0.782,-0.74)
 (-0.605,0.823)  (0.0268,-0.27) (-0.514,-0.967)  (-0.563,0.998)
  (0.536,-0.33)   (0.832,0.904)  (0.608,-0.726)  (0.678,0.0259)

The triangular matrix T is:
 (-0.691,-1.63)  (0.763,-0.144) (-0.104,-0.836) (-0.462,-0.378)
          (0,0)   (-0.758,1.22)  (-0.65,-0.772)  (-0.244,0.113)
          (0,0)           (0,0)   (0.137,0.505) (0.0687,-0.404)
          (0,0)           (0,0)           (0,0)   (1.52,-0.402)
const ComplexMatrixType& 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:

cout << "Here is a random 4x4 matrix, A:" << endl << A << endl << endl;
ComplexSchur<MatrixXcf> schurOfA(A);
cout << "The unitary matrix U is:" << endl << schurOfA.matrixU() << endl;

Output:

Here is a random 4x4 matrix, A:
  (-0.211,0.68)  (0.108,-0.444)   (0.435,0.271) (-0.198,-0.687)
  (0.597,0.566) (0.258,-0.0452)  (0.214,-0.717)  (-0.782,-0.74)
 (-0.605,0.823)  (0.0268,-0.27) (-0.514,-0.967)  (-0.563,0.998)
  (0.536,-0.33)   (0.832,0.904)  (0.608,-0.726)  (0.678,0.0259)

The unitary matrix U is:
 (-0.122,0.271)   (0.354,0.255)    (-0.7,0.321) (0.0909,-0.346)
   (0.247,0.23)  (0.435,-0.395)   (0.184,-0.38)  (0.492,-0.347)
(0.859,-0.0877)  (0.00469,0.21) (-0.256,0.0163)   (0.133,0.355)
 (-0.116,0.195) (-0.484,-0.432)  (-0.183,0.359)   (0.559,0.231)

Member Data Documentation

ComputationInfo m_info
protected
bool m_isInitialized
protected
ComplexMatrixType m_matT
protected
ComplexMatrixType m_matU
protected
bool m_matUisUptodate
protected
const int m_maxIterations
static

Maximum number of iterations.

Maximum number of iterations allowed for an eigenvalue to converge.


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