11#ifndef EIGEN_REAL_SCHUR_H
12#define EIGEN_REAL_SCHUR_H
65 typedef typename MatrixType::Scalar
Scalar;
86 m_workspaceVector(
size),
88 m_isInitialized(false),
89 m_matUisUptodate(false),
103 template<
typename InputType>
105 : m_matT(matrix.rows(),matrix.cols()),
106 m_matU(matrix.rows(),matrix.cols()),
107 m_workspaceVector(matrix.rows()),
108 m_hess(matrix.rows()),
109 m_isInitialized(false),
110 m_matUisUptodate(false),
129 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
130 eigen_assert(m_matUisUptodate &&
"The matrix U has not been computed during the RealSchur decomposition.");
146 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
169 template<
typename InputType>
189 template<
typename HessMatrixType,
typename OrthMatrixType>
197 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
208 m_maxIters = maxIters;
232 bool m_isInitialized;
233 bool m_matUisUptodate;
240 void splitOffTwoRows(
Index iu,
bool computeU,
const Scalar& exshift);
247template<
typename MatrixType>
248template<
typename InputType>
254 Index maxIters = m_maxIters;
256 maxIters = m_maxIterationsPerRow * matrix.
rows();
259 if(scale<considerAsZero)
261 m_matT.setZero(matrix.
rows(),matrix.
cols());
263 m_matU.setIdentity(matrix.
rows(),matrix.
cols());
265 m_isInitialized =
true;
266 m_matUisUptodate = computeU;
271 m_hess.compute(matrix.
derived()/scale);
276 m_workspaceVector.resize(matrix.
cols());
278 m_hess.matrixQ().evalTo(m_matU, m_workspaceVector);
279 computeFromHessenberg(m_hess.matrixH(), m_matU, computeU);
285template<
typename MatrixType>
286template<
typename HessMatrixType,
typename OrthMatrixType>
292 m_workspaceVector.resize(m_matT.cols());
296 Index maxIters = m_maxIters;
298 maxIters = m_maxIterationsPerRow * matrixH.rows();
299 Scalar* workspace = &m_workspaceVector.coeffRef(0);
305 Index iu = m_matT.cols() - 1;
309 Scalar norm = computeNormOfT();
319 Index il = findSmallSubdiagEntry(iu,considerAsZero);
324 m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
326 m_matT.coeffRef(iu, iu-1) =
Scalar(0);
332 splitOffTwoRows(iu, computeU, exshift);
339 Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
340 computeShift(iu, iter, exshift, shiftInfo);
342 totalIter = totalIter + 1;
343 if (totalIter > maxIters)
break;
345 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
346 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
350 if(totalIter <= maxIters)
355 m_isInitialized =
true;
356 m_matUisUptodate = computeU;
361template<
typename MatrixType>
370 norm += m_matT.col(j).segment(0, (
std::min)(
size,j+2)).cwiseAbs().sum();
375template<
typename MatrixType>
376inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(
Index iu,
const Scalar& considerAsZero)
382 Scalar s =
abs(m_matT.coeff(res-1,res-1)) +
abs(m_matT.coeff(res,res));
384 s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
386 if (
abs(m_matT.coeff(res,res-1)) <= s)
394template<
typename MatrixType>
395inline void RealSchur<MatrixType>::splitOffTwoRows(
Index iu,
bool computeU,
const Scalar& exshift)
403 Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
404 Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
405 m_matT.coeffRef(iu,iu) += exshift;
406 m_matT.coeffRef(iu-1,iu-1) += exshift;
411 JacobiRotation<Scalar> rot;
413 rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
415 rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
417 m_matT.rightCols(
size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
418 m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
419 m_matT.coeffRef(iu, iu-1) = Scalar(0);
421 m_matU.applyOnTheRight(iu-1, iu, rot);
425 m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
429template<
typename MatrixType>
430inline void RealSchur<MatrixType>::computeShift(
Index iu,
Index iter, Scalar& exshift, Vector3s& shiftInfo)
434 shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
435 shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
436 shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
441 exshift += shiftInfo.coeff(0);
442 for (
Index i = 0; i <= iu; ++i)
443 m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
444 Scalar s =
abs(m_matT.coeff(iu,iu-1)) +
abs(m_matT.coeff(iu-1,iu-2));
445 shiftInfo.coeffRef(0) = Scalar(0.75) * s;
446 shiftInfo.coeffRef(1) = Scalar(0.75) * s;
447 shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
453 Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
454 s = s * s + shiftInfo.coeff(2);
458 if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
460 s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
461 s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
463 for (
Index i = 0; i <= iu; ++i)
464 m_matT.coeffRef(i,i) -= s;
465 shiftInfo.setConstant(Scalar(0.964));
471template<
typename MatrixType>
472inline void RealSchur<MatrixType>::initFrancisQRStep(
Index il,
Index iu,
const Vector3s& shiftInfo,
Index& im, Vector3s& firstHouseholderVector)
475 Vector3s& v = firstHouseholderVector;
477 for (im = iu-2; im >= il; --im)
479 const Scalar Tmm = m_matT.coeff(im,im);
480 const Scalar r = shiftInfo.coeff(0) - Tmm;
481 const Scalar s = shiftInfo.coeff(1) - Tmm;
482 v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
483 v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
484 v.coeffRef(2) = m_matT.coeff(im+2,im+1);
488 const Scalar lhs = m_matT.coeff(im,im-1) * (
abs(v.coeff(1)) +
abs(v.coeff(2)));
489 const Scalar rhs = v.coeff(0) * (
abs(m_matT.coeff(im-1,im-1)) +
abs(Tmm) +
abs(m_matT.coeff(im+1,im+1)));
490 if (
abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
496template<
typename MatrixType>
497inline void RealSchur<MatrixType>::performFrancisQRStep(
Index il,
Index im,
Index iu,
bool computeU,
const Vector3s& firstHouseholderVector, Scalar* workspace)
504 for (
Index k = im; k <= iu-2; ++k)
506 bool firstIteration = (k == im);
510 v = firstHouseholderVector;
512 v = m_matT.template block<3,1>(k,k-1);
515 Matrix<Scalar, 2, 1> ess;
516 v.makeHouseholder(ess, tau,
beta);
518 if (
beta != Scalar(0))
520 if (firstIteration && k > il)
521 m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
522 else if (!firstIteration)
523 m_matT.coeffRef(k,k-1) =
beta;
526 m_matT.block(k, k, 3,
size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
527 m_matT.block(0, k, (
std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
529 m_matU.block(0, k,
size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
533 Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
535 Matrix<Scalar, 1, 1> ess;
536 v.makeHouseholder(ess, tau,
beta);
538 if (
beta != Scalar(0))
540 m_matT.coeffRef(iu-1, iu-2) =
beta;
541 m_matT.block(iu-1, iu-1, 2,
size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
542 m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
544 m_matU.block(0, iu-1,
size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
548 for (
Index i = im+2; i <= iu; ++i)
550 m_matT.coeffRef(i,i-2) = Scalar(0);
552 m_matT.coeffRef(i,i-3) = Scalar(0);
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
Definition: ArrayCwiseUnaryOps.h:80
#define eigen_assert(x)
Definition: Macros.h:1047
constexpr common_return_t< T1, T2 > beta(const T1 a, const T2 b) noexcept
Compile-time beta function.
Definition: beta.hpp:36
\eigenvalues_module
Definition: RealSchur.h:55
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: RealSchur.h:66
_MatrixType MatrixType
Definition: RealSchur.h:57
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
Eigen::Index Index
Definition: RealSchur.h:67
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
MatrixType::Scalar Scalar
Definition: RealSchur.h:65
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Definition: RealSchur.h:69
@ MaxColsAtCompileTime
Definition: RealSchur.h:63
@ MaxRowsAtCompileTime
Definition: RealSchur.h:62
@ ColsAtCompileTime
Definition: RealSchur.h:60
@ Options
Definition: RealSchur.h:61
@ RowsAtCompileTime
Definition: RealSchur.h:59
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: RealSchur.h:70
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:104
UnitType abs(const UnitType x) noexcept
Compute absolute value.
Definition: math.h:721
auto sqrt(const UnitType &value) noexcept -> unit_t< square_root< typename units::traits::unit_t_traits< UnitType >::unit_type >, typename units::traits::unit_t_traits< UnitType >::underlying_type, linear_scale >
computes the square root of value
Definition: math.h:483
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:440
@ Success
Computation was successful.
Definition: Constants.h:442
@ NoConvergence
Iterative procedure did not converge.
Definition: Constants.h:446
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
EIGEN_DEVICE_FUNC bool is_same_dense(const T1 &mat1, const T2 &mat2, typename enable_if< possibly_same_dense< T1, T2 >::value >::type *=0)
Definition: XprHelper.h:695
EIGEN_CONSTEXPR Index size(const T &x)
Definition: Meta.h:479
Namespace containing all symbols from the Eigen library.
Definition: MatrixExponential.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time,...
Definition: Constants.h:22
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:30
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:63
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:46
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:60
Holds information about the various numeric (i.e.
Definition: NumTraits.h:233