Performs a real Schur decomposition of a square matrix. More...
#include <RealSchur.h>
Public Member Functions | |
| RealSchur & | compute (const MatrixType &matrix, bool computeU=true) |
| Computes Schur decomposition of given matrix. | |
| ComputationInfo | info () const |
| Reports whether previous computation was successful. | |
| const MatrixType & | matrixT () const |
| Returns the quasi-triangular matrix in the Schur decomposition. | |
| const MatrixType & | matrixU () const |
| Returns the orthogonal matrix in the Schur decomposition. | |
| RealSchur (const MatrixType &matrix, bool computeU=true) | |
| Constructor; computes real Schur decomposition of given matrix. | |
| RealSchur (Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime) | |
| Default constructor. | |
Static Public Attributes | |
| static const int | m_maxIterations |
| Maximum number of iterations. | |
Performs a real Schur decomposition of a square matrix.
This is defined in the Eigenvalues module.
| _MatrixType | the type of the matrix of which we are computing the real Schur decomposition; this is expected to be an instantiation of the Matrix class template. |
Given a real square matrix A, this class computes the real Schur decomposition: 

Call the function compute() to compute the real Schur decomposition of a given matrix. Alternatively, you can use the RealSchur(const MatrixType&, bool) constructor which computes the real Schur decomposition at construction time. Once the decomposition is computed, you can use the matrixU() and matrixT() functions to retrieve the matrices U and T in the decomposition.
The documentation of RealSchur(const MatrixType&, bool) contains an example of the typical use of this class.
|
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.
Referenced by compute().
|
inline |
Constructor; computes real 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.
Example:
Output:
Here is a random 6x6 matrix, A:
0.68 -0.33 -0.27 -0.717 -0.687 0.0259
-0.211 0.536 0.0268 0.214 -0.198 0.678
0.566 -0.444 0.904 -0.967 -0.74 0.225
0.597 0.108 0.832 -0.514 -0.782 -0.408
0.823 -0.0452 0.271 -0.726 0.998 0.275
-0.605 0.258 0.435 0.608 -0.563 0.0486
The orthogonal matrix U is:
0.348 -0.754 0.00435 -0.351 0.0145 0.432
-0.16 -0.266 -0.747 0.457 -0.366 0.0571
0.505 -0.157 0.0746 0.644 0.518 -0.177
0.703 0.324 -0.409 -0.349 -0.187 -0.275
0.296 0.372 0.24 0.324 -0.379 0.684
-0.126 0.305 -0.46 -0.161 0.647 0.485
The quasi-triangular matrix T is:
-0.2 -1.83 0.864 0.271 1.09 0.14
0.647 0.298 -0.0536 0.676 -0.288 0.023
0 0 0.967 -0.201 -0.429 0.847
0 0 0 0.353 0.602 0.694
0 0 0 0 0.572 -1.03
0 0 0 0 0.0184 0.664
U * T * U^T =
0.68 -0.33 -0.27 -0.717 -0.687 0.0259
-0.211 0.536 0.0268 0.214 -0.198 0.678
0.566 -0.444 0.904 -0.967 -0.74 0.225
0.597 0.108 0.832 -0.514 -0.782 -0.408
0.823 -0.0452 0.271 -0.726 0.998 0.275
-0.605 0.258 0.435 0.608 -0.563 0.0486
References compute().
| RealSchur< 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 Francis QR iterations with implicit double shift. The cost of computing the Schur decomposition depends on the number of iterations; as a rough guide, it may be taken to be 

Example:
Output:
The matrix T in the decomposition of A is:
0.523 -0.698 0.148 0.742
0.475 0.986 -0.793 0.721
0 0 -0.28 -0.77
0 0 0.0145 -0.367
The matrix T in the decomposition of A^(-1) is:
-3.06 -4.57 -6.05 5.39
0.168 -2.62 -3.33 3.86
0 0 0.434 0.56
0 0 -1.06 1.35
References m_maxIterations, Eigen::NoConvergence, RealSchur(), and Eigen::Success.
Referenced by RealSchur().
|
inline |
Reports whether previous computation was successful.
Success if computation was succesful, NoConvergence otherwise.
|
inline |
Returns the quasi-triangular matrix in the Schur decomposition.
|
inline |
Returns the orthogonal matrix in the Schur decomposition.
computeU was set to true (the default value).
|
static |
Maximum number of iterations.
Maximum number of iterations allowed for an eigenvalue to converge.
Referenced by compute().