12#ifndef EIGEN_COMPLEX_SCHUR_H
13#define EIGEN_COMPLEX_SCHUR_H
15#include "./HessenbergDecomposition.h"
18#include "./InternalHeaderCheck.h"
23template <
typename MatrixType,
bool IsComplex>
24struct complex_schur_reduce_to_hessenberg;
55template <
typename MatrixType_>
58 typedef MatrixType_ MatrixType;
60 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
61 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
62 Options = internal::traits<MatrixType>::Options,
63 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
68 typedef typename MatrixType::Scalar
Scalar;
69 typedef typename NumTraits<Scalar>::Real RealScalar;
85 typedef Matrix<
ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime,
101 : m_matT(size, size),
104 m_isInitialized(false),
105 m_matUisUptodate(false),
117 template <
typename InputType>
119 : m_matT(matrix.rows(), matrix.cols()),
120 m_matU(matrix.rows(), matrix.cols()),
121 m_hess(matrix.rows()),
122 m_isInitialized(false),
123 m_matUisUptodate(false),
143 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
144 eigen_assert(m_matUisUptodate &&
"The matrix U has not been computed during the ComplexSchur decomposition.");
166 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
192 template <
typename InputType>
212 template <
typename HessMatrixType,
typename OrthMatrixType>
214 bool computeU =
true);
221 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
231 m_maxIters = maxIters;
249 bool m_isInitialized;
250 bool m_matUisUptodate;
254 bool subdiagonalEntryIsNeglegible(
Index i);
256 void reduceToTriangularForm(
bool computeU);
257 friend struct internal::complex_schur_reduce_to_hessenberg<MatrixType,
NumTraits<
Scalar>::IsComplex>;
263template <typename MatrixType>
264inline bool
ComplexSchur<MatrixType>::subdiagonalEntryIsNeglegible(Index i) {
265 RealScalar d = numext::norm1(m_matT.
coeff(i, i)) + numext::norm1(m_matT.
coeff(i + 1, i + 1));
266 RealScalar sd = numext::norm1(m_matT.
coeff(i + 1, i));
275template <
typename MatrixType>
278 if ((iter == 10 || iter == 20) && iu > 1) {
280 return abs(numext::real(m_matT.
coeff(iu, iu - 1))) +
abs(numext::real(m_matT.
coeff(iu - 1, iu - 2)));
285 Matrix<ComplexScalar, 2, 2> t = m_matT.template block<2, 2>(iu - 1, iu - 1);
286 RealScalar normt = t.cwiseAbs().sum();
296 RealScalar eival1_norm = numext::norm1(eival1);
297 RealScalar eival2_norm = numext::norm1(eival2);
300 if (eival1_norm > eival2_norm)
301 eival2 = det / eival1;
302 else if (!numext::is_exactly_zero(eival2_norm))
303 eival1 = det / eival2;
306 if (numext::norm1(eival1 - t.coeff(1, 1)) < numext::norm1(eival2 - t.coeff(1, 1)))
307 return normt * eival1;
309 return normt * eival2;
312template <
typename MatrixType>
313template <
typename InputType>
315 m_matUisUptodate =
false;
316 eigen_assert(matrix.cols() == matrix.rows());
318 if (matrix.cols() == 1) {
319 m_matT = matrix.derived().template cast<ComplexScalar>();
320 if (computeU) m_matU = ComplexMatrixType::Identity(1, 1);
322 m_isInitialized =
true;
323 m_matUisUptodate = computeU;
327 internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>::run(*
this, matrix.derived(),
333template <
typename MatrixType>
334template <
typename HessMatrixType,
typename OrthMatrixType>
336 const OrthMatrixType& matrixQ,
339 if (computeU) m_matU = matrixQ;
340 reduceToTriangularForm(computeU);
346template <
typename MatrixType,
bool IsComplex>
347struct complex_schur_reduce_to_hessenberg {
349 static void run(ComplexSchur<MatrixType>& _this,
const MatrixType& matrix,
bool computeU) {
350 _this.m_hess.compute(matrix);
351 _this.m_matT = _this.m_hess.matrixH();
352 if (computeU) _this.m_matU = _this.m_hess.matrixQ();
356template <
typename MatrixType>
357struct complex_schur_reduce_to_hessenberg<MatrixType, false> {
358 static void run(ComplexSchur<MatrixType>& _this,
const MatrixType& matrix,
bool computeU) {
359 typedef typename ComplexSchur<MatrixType>::ComplexScalar ComplexScalar;
362 _this.m_hess.compute(matrix);
363 _this.m_matT = _this.m_hess.matrixH().template cast<ComplexScalar>();
366 MatrixType Q = _this.m_hess.matrixQ();
367 _this.m_matU = Q.template cast<ComplexScalar>();
375template <
typename MatrixType>
376void ComplexSchur<MatrixType>::reduceToTriangularForm(
bool computeU) {
377 Index maxIters = m_maxIters;
384 Index iu = m_matT.cols() - 1;
392 if (!subdiagonalEntryIsNeglegible(iu - 1))
break;
403 if (totalIter > maxIters)
break;
407 while (il > 0 && !subdiagonalEntryIsNeglegible(il - 1)) {
416 JacobiRotation<ComplexScalar> rot;
417 rot.makeGivens(m_matT.
coeff(il, il) - shift, m_matT.
coeff(il + 1, il));
418 m_matT.rightCols(m_matT.cols() - il).applyOnTheLeft(il, il + 1, rot.adjoint());
419 m_matT.topRows((std::min)(il + 2, iu) + 1).applyOnTheRight(il, il + 1, rot);
420 if (computeU) m_matU.applyOnTheRight(il, il + 1, rot);
422 for (Index i = il + 1; i < iu; i++) {
425 m_matT.rightCols(m_matT.cols() - i).applyOnTheLeft(i, i + 1, rot.adjoint());
426 m_matT.topRows((std::min)(i + 2, iu) + 1).applyOnTheRight(i, i + 1, rot);
427 if (computeU) m_matU.applyOnTheRight(i, i + 1, rot);
431 if (totalIter <= maxIters)
436 m_isInitialized =
true;
437 m_matUisUptodate = computeU;
Performs a complex Schur decomposition of a real or complex square matrix.
Definition ComplexSchur.h:56
ComputationInfo info() const
Reports whether previous computation was successful.
Definition ComplexSchur.h:220
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
Definition ComplexSchur.h:68
ComplexSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition ComplexSchur.h:100
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition ComplexSchur.h:165
ComplexSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
Compute Schur decomposition from a given Hessenberg matrix.
Index getMaxIterations()
Returns the maximum number of iterations.
Definition ComplexSchur.h:236
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition ComplexSchur.h:243
Eigen::Index Index
Definition ComplexSchur.h:70
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType_.
Definition ComplexSchur.h:78
ComplexSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes Schur decomposition of given matrix.
Definition ComplexSchur.h:118
Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > ComplexMatrixType
Type for the matrices in the Schur decomposition.
Definition ComplexSchur.h:87
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
ComplexSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition ComplexSchur.h:230
const ComplexMatrixType & matrixU() const
Returns the unitary matrix in the Schur decomposition.
Definition ComplexSchur.h:142
Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation.
Definition HessenbergDecomposition.h:61
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:217
constexpr const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:198
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25
Definition EigenBase.h:33
Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523