ComplexSchur< _MatrixType > Class Template Reference
[Eigenvalues module]

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

List of all members.

Public Types

typedef Matrix< ComplexScalar,
RowsAtCompileTime,
ColsAtCompileTime, Options,
MaxRowsAtCompileTime,
MaxColsAtCompileTime > 
ComplexMatrixType
 Type for the matrices in the Schur decomposition.
typedef std::complex< RealScalar > ComplexScalar
 Complex scalar type for _MatrixType.
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type _MatrixType.

Public Member Functions

 ComplexSchur (const MatrixType &matrix, bool computeU=true)
 Constructor; computes Schur decomposition of given matrix.
 ComplexSchur (Index size=RowsAtCompileTime==Dynamic?1:RowsAtCompileTime)
 Default constructor.
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.

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:
_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: $ 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

typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> 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> 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.


Constructor & Destructor Documentation

ComplexSchur ( Index  size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime  )  [inline]

Default constructor.

Parameters:
[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.

See also:
compute() for an example.
ComplexSchur ( const MatrixType &  matrix,
bool  computeU = true 
) [inline]

Constructor; computes Schur decomposition of given matrix.

Parameters:
[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.

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

Member Function Documentation

ComplexSchur< MatrixType > & compute ( const MatrixType &  matrix,
bool  computeU = true 
) [inline]

Computes Schur decomposition of given matrix.

Parameters:
[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.
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:

MatrixXcf A = MatrixXcf::Random(4,4);
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)
ComputationInfo info (  )  const [inline]

Reports whether previous computation was successful.

Returns:
Success if computation was succesful, NoConvergence otherwise.
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:

MatrixXcf A = MatrixXcf::Random(4,4);
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:

MatrixXcf A = MatrixXcf::Random(4,4);
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

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:
Generated on Sun Jul 3 00:55:26 2011 for Eigen by  doxygen 1.6.3