Class xlifepp::ComplexSchur#
-
template<typename _MatrixType>
class ComplexSchur#
-
Inheritence diagram for xlifepp::ComplexSchur:
Collaboration diagram for xlifepp::ComplexSchur:
Performs a complex Schur decomposition of a real or complex square matrix.
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.
See also
class RealSchur, class EigenSolver, class ComplexEigenSolver
- 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.
Note
This code is inspired from Jampack
Public Types
-
typedef MatrixEigenDense<ComplexScalar> 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 MatrixType::type_t Scalar#
-
Scalar type for matrices of type
_MatrixType
.
Public Functions
-
inline ComplexSchur(const MatrixType &matrix, bool computeU = true)#
-
Constructor; computes Schur decomposition of given matrix.
This constructor calls compute() to compute the Schur decomposition.
- 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.
-
inline ComplexSchur(Index size)#
-
Default constructor.
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 wrongsize
, but it may impair performance.See also
compute() for an example.
- Parameters:
-
size – [in] Positive integer, size of the matrix whose Schur decomposition will be computed.
-
ComplexSchur &compute(const MatrixType &matrix, bool computeU = true)#
-
Computes Schur decomposition of given matrix.
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: \include test_EigenSolverDecomposition.cpp
- 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
-
inline ComputationInfo info() const#
-
Reports whether previous computation was successful.
- Returns:
-
_success
if computation was succesful,_noConvergence
otherwise.
-
inline const ComplexMatrixType &matrixT() const#
-
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: \include test_EigenSolverDecomposition.cpp
- Returns:
-
A const reference to the matrix T.
-
inline const ComplexMatrixType &matrixU() const#
-
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: \include test_EigenSolverDecomposition.cpp
- Returns:
-
A const reference to the matrix U.
-
void swapSchur(const std::vector<int> &order, bool computeU = true)#
-
Reorder the Schur decomposition according to permutation vector.
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 change upper triangular matrix T and unitary matrix U
Example: \include test_EigenSolverDecomposition.cpp
Public Static Attributes
-
static const int maxIterations_ = 30#
-
Maximum number of iterations.
Maximum number of iterations allowed for an eigenvalue to converge.