11 #ifndef EIGEN_REAL_SCHUR_H
12 #define EIGEN_REAL_SCHUR_H
14 #include "./HessenbergDecomposition.h"
16 #include "./InternalHeaderCheck.h"
59 typedef MatrixType_ MatrixType;
61 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
62 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
63 Options = MatrixType::Options,
64 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
65 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
67 typedef typename MatrixType::Scalar Scalar;
68 typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
88 m_workspaceVector(size),
90 m_isInitialized(false),
91 m_matUisUptodate(false),
105 template<
typename InputType>
107 : m_matT(matrix.rows(),matrix.cols()),
108 m_matU(matrix.rows(),matrix.cols()),
109 m_workspaceVector(matrix.rows()),
110 m_hess(matrix.rows()),
111 m_isInitialized(false),
112 m_matUisUptodate(false),
131 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
132 eigen_assert(m_matUisUptodate &&
"The matrix U has not been computed during the RealSchur decomposition.");
148 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
171 template<
typename InputType>
191 template<
typename HessMatrixType,
typename OrthMatrixType>
199 eigen_assert(m_isInitialized &&
"RealSchur is not initialized.");
210 m_maxIters = maxIters;
234 bool m_isInitialized;
235 bool m_matUisUptodate;
240 Scalar computeNormOfT();
241 Index findSmallSubdiagEntry(
Index iu,
const Scalar& considerAsZero);
242 void splitOffTwoRows(
Index iu,
bool computeU,
const Scalar& exshift);
245 void performFrancisQRStep(
Index il,
Index im,
Index iu,
bool computeU,
const Vector3s& firstHouseholderVector, Scalar* workspace);
249 template<
typename MatrixType>
250 template<
typename InputType>
253 const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
255 eigen_assert(matrix.
cols() == matrix.
rows());
256 Index maxIters = m_maxIters;
258 maxIters = m_maxIterationsPerRow * matrix.
rows();
260 Scalar scale = matrix.
derived().cwiseAbs().maxCoeff();
261 if(scale<considerAsZero)
263 m_matT.setZero(matrix.
rows(),matrix.
cols());
265 m_matU.setIdentity(matrix.
rows(),matrix.
cols());
267 m_isInitialized =
true;
268 m_matUisUptodate = computeU;
273 m_hess.compute(matrix.
derived()/scale);
278 m_workspaceVector.resize(matrix.
cols());
280 m_hess.matrixQ().evalTo(m_matU, m_workspaceVector);
281 computeFromHessenberg(m_hess.matrixH(), m_matU, computeU);
287 template<
typename MatrixType>
288 template<
typename HessMatrixType,
typename OrthMatrixType>
294 m_workspaceVector.resize(m_matT.cols());
295 if(computeU && !internal::is_same_dense(m_matU,matrixQ))
298 Index maxIters = m_maxIters;
300 maxIters = m_maxIterationsPerRow * matrixH.rows();
301 Scalar* workspace = &m_workspaceVector.coeffRef(0);
307 Index iu = m_matT.cols() - 1;
311 Scalar norm = computeNormOfT();
314 Scalar considerAsZero = numext::maxi<Scalar>( norm * numext::abs2(NumTraits<Scalar>::epsilon()),
315 (std::numeric_limits<Scalar>::min)() );
317 if(!numext::is_exactly_zero(norm))
321 Index il = findSmallSubdiagEntry(iu,considerAsZero);
326 m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
328 m_matT.coeffRef(iu, iu-1) = Scalar(0);
334 splitOffTwoRows(iu, computeU, exshift);
341 Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
342 computeShift(iu, iter, exshift, shiftInfo);
344 totalIter = totalIter + 1;
345 if (totalIter > maxIters)
break;
347 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
348 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
352 if(totalIter <= maxIters)
357 m_isInitialized =
true;
358 m_matUisUptodate = computeU;
363 template<
typename MatrixType>
364 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
366 const Index size = m_matT.cols();
371 for (Index j = 0; j < size; ++j)
372 norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
377 template<
typename MatrixType>
378 inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu,
const Scalar& considerAsZero)
384 Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
386 s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
388 if (abs(m_matT.coeff(res,res-1)) <= s)
396 template<
typename MatrixType>
397 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu,
bool computeU,
const Scalar& exshift)
401 const Index size = m_matT.cols();
405 Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
406 Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
407 m_matT.coeffRef(iu,iu) += exshift;
408 m_matT.coeffRef(iu-1,iu-1) += exshift;
412 Scalar z = sqrt(abs(q));
413 JacobiRotation<Scalar> rot;
415 rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
417 rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
419 m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
420 m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
421 m_matT.coeffRef(iu, iu-1) = Scalar(0);
423 m_matU.applyOnTheRight(iu-1, iu, rot);
427 m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
431 template<
typename MatrixType>
432 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
436 shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
437 shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
438 shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
443 exshift += shiftInfo.coeff(0);
444 for (Index i = 0; i <= iu; ++i)
445 m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
446 Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
447 shiftInfo.coeffRef(0) = Scalar(0.75) * s;
448 shiftInfo.coeffRef(1) = Scalar(0.75) * s;
449 shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
455 Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
456 s = s * s + shiftInfo.coeff(2);
460 if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
462 s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
463 s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
465 for (Index i = 0; i <= iu; ++i)
466 m_matT.coeffRef(i,i) -= s;
467 shiftInfo.setConstant(Scalar(0.964));
473 template<
typename MatrixType>
474 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu,
const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
477 Vector3s& v = firstHouseholderVector;
479 for (im = iu-2; im >= il; --im)
481 const Scalar Tmm = m_matT.coeff(im,im);
482 const Scalar r = shiftInfo.coeff(0) - Tmm;
483 const Scalar s = shiftInfo.coeff(1) - Tmm;
484 v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
485 v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
486 v.coeffRef(2) = m_matT.coeff(im+2,im+1);
490 const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
491 const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
492 if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
498 template<
typename MatrixType>
499 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu,
bool computeU,
const Vector3s& firstHouseholderVector, Scalar* workspace)
501 eigen_assert(im >= il);
502 eigen_assert(im <= iu-2);
504 const Index size = m_matT.cols();
506 for (Index k = im; k <= iu-2; ++k)
508 bool firstIteration = (k == im);
512 v = firstHouseholderVector;
514 v = m_matT.template block<3,1>(k,k-1);
517 Matrix<Scalar, 2, 1> ess;
518 v.makeHouseholder(ess, tau, beta);
520 if (!numext::is_exactly_zero(beta))
522 if (firstIteration && k > il)
523 m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
524 else if (!firstIteration)
525 m_matT.coeffRef(k,k-1) = beta;
528 m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
529 m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
531 m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
535 Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
537 Matrix<Scalar, 1, 1> ess;
538 v.makeHouseholder(ess, tau, beta);
540 if (!numext::is_exactly_zero(beta))
542 m_matT.coeffRef(iu-1, iu-2) = beta;
543 m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
544 m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
546 m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
550 for (Index i = im+2; i <= iu; ++i)
552 m_matT.coeffRef(i,i-2) = Scalar(0);
554 m_matT.coeffRef(i,i-3) = Scalar(0);
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:57
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:208
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:129
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:197
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:215
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:85
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:106
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:225
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:146
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
Eigen::Index Index
Definition: RealSchur.h:69
ComputationInfo
Definition: Constants.h:442
@ Success
Definition: Constants.h:444
@ NoConvergence
Definition: Constants.h:448
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const int Dynamic
Definition: Constants.h:24
Definition: EigenBase.h:32
Derived & derived()
Definition: EigenBase.h:48
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62