11 #ifndef EIGEN_REAL_SCHUR_H 12 #define EIGEN_REAL_SCHUR_H 14 #include "./HessenbergDecomposition.h" 57 typedef _MatrixType MatrixType;
59 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
60 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
61 Options = MatrixType::Options,
62 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
63 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
65 typedef typename MatrixType::Scalar Scalar;
66 typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
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;
229 ColumnVectorType m_workspaceVector;
232 bool m_isInitialized;
233 bool m_matUisUptodate;
238 Scalar computeNormOfT();
239 Index findSmallSubdiagEntry(Index iu,
const Scalar& considerAsZero);
240 void splitOffTwoRows(Index iu,
bool computeU,
const Scalar& exshift);
241 void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
242 void initFrancisQRStep(Index il, Index iu,
const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
243 void performFrancisQRStep(Index il, Index im, Index iu,
bool computeU,
const Vector3s& firstHouseholderVector, Scalar* workspace);
247 template<
typename MatrixType>
248 template<
typename InputType>
251 const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
253 eigen_assert(matrix.
cols() == matrix.
rows());
254 Index maxIters = m_maxIters;
258 Scalar scale = matrix.
derived().cwiseAbs().maxCoeff();
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;
278 m_hess.
matrixQ().evalTo(m_matU, m_workspaceVector);
285 template<
typename MatrixType>
286 template<
typename HessMatrixType,
typename OrthMatrixType>
292 m_workspaceVector.
resize(m_matT.cols());
293 if(computeU && !internal::is_same_dense(m_matU,matrixQ))
296 Index maxIters = m_maxIters;
299 Scalar* workspace = &m_workspaceVector.
coeffRef(0);
305 Index iu = m_matT.cols() - 1;
309 Scalar norm = computeNormOfT();
313 (std::numeric_limits<Scalar>::min)() );
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;
361 template<
typename MatrixType>
364 const Index size = m_matT.cols();
369 for (
Index j = 0; j < size; ++j)
370 norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
375 template<
typename MatrixType>
382 Scalar s =
abs(m_matT.coeff(res-1,res-1)) +
abs(m_matT.coeff(res,res));
386 if (
abs(m_matT.coeff(res,res-1)) <= s)
394 template<
typename MatrixType>
399 const Index size = m_matT.cols();
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;
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);
429 template<
typename MatrixType>
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);
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;
471 template<
typename MatrixType>
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)));
496 template<
typename MatrixType>
499 eigen_assert(im >= il);
500 eigen_assert(im <= iu-2);
502 const Index size = m_matT.cols();
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);
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);
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);
558 #endif // EIGEN_REAL_SCHUR_H void makeHouseholder(EssentialPart &essential, Scalar &tau, RealScalar &beta) const
Definition: Householder.h:67
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:54
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
Namespace containing all symbols from the Eigen library.
Definition: Core:141
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:34
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
Derived & derived()
Definition: EigenBase.h:46
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:271
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
Definition: EigenBase.h:29
Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:175
Eigen::Index Index
Definition: RealSchur.h:67
HessenbergDecomposition & compute(const EigenBase< InputType > &matrix)
Computes Hessenberg decomposition of given matrix.
Definition: HessenbergDecomposition.h:152
MatrixHReturnType matrixH() const
Constructs the Hessenberg matrix H in the decomposition.
Definition: HessenbergDecomposition.h:262
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:63
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
Derived & setConstant(Index size, const Scalar &val)
Definition: CwiseNullaryOp.h:361
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
Definition: Constants.h:442
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:152
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
HouseholderSequenceType matrixQ() const
Reconstructs the orthogonal matrix Q in the decomposition.
Definition: HessenbergDecomposition.h:234
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:60
const int Dynamic
Definition: Constants.h:22
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:104
ComputationInfo
Definition: Constants.h:440
Definition: Constants.h:446
JacobiRotation adjoint() const
Definition: Jacobi.h:67
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:162