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& maxDiagEntry);
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 eigen_assert(matrix.
cols() == matrix.
rows());
252 Index maxIters = m_maxIters;
256 Scalar scale = matrix.
derived().cwiseAbs().maxCoeff();
268 template<
typename MatrixType>
269 template<
typename HessMatrixType,
typename OrthMatrixType>
278 Index maxIters = m_maxIters;
281 m_workspaceVector.
resize(m_matT.cols());
282 Scalar* workspace = &m_workspaceVector.
coeffRef(0);
288 Index iu = m_matT.cols() - 1;
292 Scalar norm = computeNormOfT();
296 Scalar maxDiagEntry = m_matT.cwiseAbs().diagonal().maxCoeff();
300 Index il = findSmallSubdiagEntry(iu,maxDiagEntry);
305 m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
307 maxDiagEntry = numext::maxi<Scalar>(maxDiagEntry,
abs(m_matT.coeffRef(iu,iu)));
309 m_matT.coeffRef(iu, iu-1) = Scalar(0);
315 splitOffTwoRows(iu, computeU, exshift);
317 maxDiagEntry = numext::maxi<Scalar>(maxDiagEntry,numext::maxi(
abs(m_matT.coeff(iu,iu)),
abs(m_matT.coeff(iu-1,iu-1))));
324 Vector3s firstHouseholderVector(0,0,0), shiftInfo;
325 computeShift(iu, iter, exshift, shiftInfo);
327 totalIter = totalIter + 1;
328 if (totalIter > maxIters)
break;
330 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
331 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
333 maxDiagEntry = numext::maxi(maxDiagEntry,m_matT.cwiseAbs().diagonal().segment(im,iu-im).maxCoeff());
337 if(totalIter <= maxIters)
342 m_isInitialized =
true;
343 m_matUisUptodate = computeU;
348 template<
typename MatrixType>
351 const Index size = m_matT.cols();
356 for (
Index j = 0; j < size; ++j)
357 norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
362 template<
typename MatrixType>
377 template<
typename MatrixType>
382 const Index size = m_matT.cols();
386 Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
387 Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
388 m_matT.coeffRef(iu,iu) += exshift;
389 m_matT.coeffRef(iu-1,iu-1) += exshift;
396 rot.
makeGivens(p + z, m_matT.coeff(iu, iu-1));
398 rot.
makeGivens(p - z, m_matT.coeff(iu, iu-1));
400 m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.
adjoint());
401 m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
402 m_matT.coeffRef(iu, iu-1) = Scalar(0);
404 m_matU.applyOnTheRight(iu-1, iu, rot);
408 m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
412 template<
typename MatrixType>
417 shiftInfo.
coeffRef(0) = m_matT.coeff(iu,iu);
418 shiftInfo.
coeffRef(1) = m_matT.coeff(iu-1,iu-1);
419 shiftInfo.
coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
424 exshift += shiftInfo.
coeff(0);
425 for (
Index i = 0; i <= iu; ++i)
426 m_matT.coeffRef(i,i) -= shiftInfo.
coeff(0);
427 Scalar s =
abs(m_matT.coeff(iu,iu-1)) +
abs(m_matT.coeff(iu-1,iu-2));
428 shiftInfo.
coeffRef(0) = Scalar(0.75) * s;
429 shiftInfo.
coeffRef(1) = Scalar(0.75) * s;
430 shiftInfo.
coeffRef(2) = Scalar(-0.4375) * s * s;
436 Scalar s = (shiftInfo.
coeff(1) - shiftInfo.
coeff(0)) / Scalar(2.0);
437 s = s * s + shiftInfo.
coeff(2);
443 s = s + (shiftInfo.
coeff(1) - shiftInfo.
coeff(0)) / Scalar(2.0);
444 s = shiftInfo.
coeff(0) - shiftInfo.
coeff(2) / s;
446 for (
Index i = 0; i <= iu; ++i)
447 m_matT.coeffRef(i,i) -= s;
454 template<
typename MatrixType>
458 Vector3s& v = firstHouseholderVector;
460 for (im = iu-2; im >= il; --im)
462 const Scalar Tmm = m_matT.coeff(im,im);
463 const Scalar r = shiftInfo.
coeff(0) - Tmm;
464 const Scalar s = shiftInfo.
coeff(1) - Tmm;
465 v.
coeffRef(0) = (r * s - shiftInfo.
coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
466 v.
coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
467 v.
coeffRef(2) = m_matT.coeff(im+2,im+1);
471 const Scalar lhs = m_matT.coeff(im,im-1) * (
abs(v.
coeff(1)) +
abs(v.
coeff(2)));
472 const Scalar rhs = v.
coeff(0) * (
abs(m_matT.coeff(im-1,im-1)) +
abs(Tmm) +
abs(m_matT.coeff(im+1,im+1)));
479 template<
typename MatrixType>
482 eigen_assert(im >= il);
483 eigen_assert(im <= iu-2);
485 const Index size = m_matT.cols();
487 for (
Index k = im; k <= iu-2; ++k)
489 bool firstIteration = (k == im);
493 v = firstHouseholderVector;
495 v = m_matT.template block<3,1>(k,k-1);
499 v.makeHouseholder(ess, tau, beta);
501 if (beta != Scalar(0))
503 if (firstIteration && k > il)
504 m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
505 else if (!firstIteration)
506 m_matT.coeffRef(k,k-1) = beta;
509 m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
510 m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
512 m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
521 if (beta != Scalar(0))
523 m_matT.coeffRef(iu-1, iu-2) = beta;
524 m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
525 m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
527 m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
531 for (
Index i = im+2; i <= iu; ++i)
533 m_matT.coeffRef(i,i-2) = Scalar(0);
535 m_matT.coeffRef(i,i-3) = Scalar(0);
541 #endif // EIGEN_REAL_SCHUR_H HouseholderSequenceType matrixQ() const
Reconstructs the orthogonal matrix Q in the decomposition.
Definition: HessenbergDecomposition.h:234
MatrixHReturnType matrixH() const
Constructs the Hessenberg matrix H in the decomposition.
Definition: HessenbergDecomposition.h:262
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:54
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Definition: Jacobi.h:149
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.
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
Namespace containing all symbols from the Eigen library.
Definition: Core:271
Rotation given by a cosine-sine pair.
Definition: ForwardDeclarations.h:263
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:167
void makeHouseholder(EssentialPart &essential, Scalar &tau, RealScalar &beta) const
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
Derived & derived()
Definition: EigenBase.h:44
void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:273
Index rows() const
Definition: EigenBase.h:58
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
Definition: EigenBase.h:28
Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:177
Eigen::Index Index
Definition: RealSchur.h:67
HessenbergDecomposition & compute(const EigenBase< InputType > &matrix)
Computes Hessenberg decomposition of given matrix.
Definition: HessenbergDecomposition.h:152
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: XprHelper.h:35
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
JacobiRotation adjoint() const
Definition: Jacobi.h:62
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:351
Definition: Constants.h:432
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Index cols() const
Definition: EigenBase.h:61
RealSchur(Index size=RowsAtCompileTime==Dynamic?1:RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
const int Dynamic
Definition: Constants.h:21
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:104
const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:154
ComputationInfo
Definition: Constants.h:430
Definition: Constants.h:436