template<typename _MatrixType>
Eigen::ComplexSchur class

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

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.

This is defined in the Eigenvalues module. #include <Eigen/Eigenvalues>

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.

Public types

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

Public static variables

static const int m_maxIterationsPerRow
Maximum number of iterations per row.

Constructors, destructors, conversion operators

ComplexSchur(Index size = RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime) explicit
Default constructor.
template<typename InputType>
ComplexSchur(const EigenBase<InputType>& matrix, bool computeU = true) explicit
Constructor; computes Schur decomposition of given matrix.

Public functions

template<typename InputType>
auto compute(const EigenBase<InputType>& matrix, bool computeU = true) -> ComplexSchur&
Computes Schur decomposition of given matrix.
template<typename HessMatrixType, typename OrthMatrixType>
auto computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU = true) -> ComplexSchur&
Compute Schur decomposition from a given Hessenberg matrix.
auto getMaxIterations() -> Index
Returns the maximum number of iterations.
auto info() const -> ComputationInfo
Reports whether previous computation was successful.
auto matrixT() const -> const ComplexMatrixType&
Returns the triangular matrix in the Schur decomposition.
auto matrixU() const -> const ComplexMatrixType&
Returns the unitary matrix in the Schur decomposition.
auto setMaxIterations(Index maxIters) -> ComplexSchur&
Sets the maximum number of iterations allowed.

Typedef documentation

template<typename _MatrixType>
typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> Eigen::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.

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.

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

Function documentation

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

Default constructor.

Parameters
size in 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.

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

Constructor; computes Schur decomposition of given matrix.

Parameters
matrix in Square matrix whose Schur decomposition is to be computed.
computeU in If true, both T and U are computed; if false, only T is computed.

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

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

Computes Schur decomposition of given matrix.

Parameters
matrix in Square matrix whose Schur decomposition is to be computed.
computeU in 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)

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

Compute Schur decomposition from a given Hessenberg matrix.

Parameters
matrixH in Matrix in Hessenberg form H
matrixQ in orthogonal matrix Q that transform a matrix A to H : A = Q H Q^T
computeU Computes 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())

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

Reports whether previous computation was successful.

Returns Success if computation was successful, NoConvergence otherwise.

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

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)

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

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)

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

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.

Variable documentation

template<typename _MatrixType>
static const int Eigen::ComplexSchur<_MatrixType>::m_maxIterationsPerRow

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.