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);
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();
298 Index il = findSmallSubdiagEntry(iu);
303 m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
305 m_matT.coeffRef(iu, iu-1) = Scalar(0);
311 splitOffTwoRows(iu, computeU, exshift);
318 Vector3s firstHouseholderVector(0,0,0), shiftInfo;
319 computeShift(iu, iter, exshift, shiftInfo);
321 totalIter = totalIter + 1;
322 if (totalIter > maxIters)
break;
324 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
325 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
329 if(totalIter <= maxIters)
334 m_isInitialized =
true;
335 m_matUisUptodate = computeU;
340 template<
typename MatrixType>
343 const Index size = m_matT.cols();
348 for (
Index j = 0; j < size; ++j)
349 norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
354 template<
typename MatrixType>
361 Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
370 template<
typename MatrixType>
375 const Index size = m_matT.cols();
379 Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
380 Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
381 m_matT.coeffRef(iu,iu) += exshift;
382 m_matT.coeffRef(iu-1,iu-1) += exshift;
386 Scalar z = sqrt(abs(q));
389 rot.
makeGivens(p + z, m_matT.coeff(iu, iu-1));
391 rot.
makeGivens(p - z, m_matT.coeff(iu, iu-1));
393 m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.
adjoint());
394 m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
395 m_matT.coeffRef(iu, iu-1) = Scalar(0);
397 m_matU.applyOnTheRight(iu-1, iu, rot);
401 m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
405 template<
typename MatrixType>
410 shiftInfo.
coeffRef(0) = m_matT.coeff(iu,iu);
411 shiftInfo.
coeffRef(1) = m_matT.coeff(iu-1,iu-1);
412 shiftInfo.
coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
417 exshift += shiftInfo.
coeff(0);
418 for (
Index i = 0; i <= iu; ++i)
419 m_matT.coeffRef(i,i) -= shiftInfo.
coeff(0);
420 Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
421 shiftInfo.
coeffRef(0) = Scalar(0.75) * s;
422 shiftInfo.
coeffRef(1) = Scalar(0.75) * s;
423 shiftInfo.
coeffRef(2) = Scalar(-0.4375) * s * s;
429 Scalar s = (shiftInfo.
coeff(1) - shiftInfo.
coeff(0)) / Scalar(2.0);
430 s = s * s + shiftInfo.
coeff(2);
436 s = s + (shiftInfo.
coeff(1) - shiftInfo.
coeff(0)) / Scalar(2.0);
437 s = shiftInfo.
coeff(0) - shiftInfo.
coeff(2) / s;
439 for (
Index i = 0; i <= iu; ++i)
440 m_matT.coeffRef(i,i) -= s;
447 template<
typename MatrixType>
451 Vector3s& v = firstHouseholderVector;
453 for (im = iu-2; im >= il; --im)
455 const Scalar Tmm = m_matT.coeff(im,im);
456 const Scalar r = shiftInfo.
coeff(0) - Tmm;
457 const Scalar s = shiftInfo.
coeff(1) - Tmm;
458 v.
coeffRef(0) = (r * s - shiftInfo.
coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
459 v.
coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
460 v.
coeffRef(2) = m_matT.coeff(im+2,im+1);
464 const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.
coeff(1)) + abs(v.
coeff(2)));
465 const Scalar rhs = v.
coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
472 template<
typename MatrixType>
475 eigen_assert(im >= il);
476 eigen_assert(im <= iu-2);
478 const Index size = m_matT.cols();
480 for (
Index k = im; k <= iu-2; ++k)
482 bool firstIteration = (k == im);
486 v = firstHouseholderVector;
488 v = m_matT.template block<3,1>(k,k-1);
492 v.makeHouseholder(ess, tau, beta);
494 if (beta != Scalar(0))
496 if (firstIteration && k > il)
497 m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
498 else if (!firstIteration)
499 m_matT.coeffRef(k,k-1) = beta;
502 m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
503 m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
505 m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
512 v.makeHouseholder(ess, tau, beta);
514 if (beta != Scalar(0))
516 m_matT.coeffRef(iu-1, iu-2) = beta;
517 m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
518 m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
520 m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
524 for (
Index i = im+2; i <= iu; ++i)
526 m_matT.coeffRef(i,i-2) = Scalar(0);
528 m_matT.coeffRef(i,i-3) = Scalar(0);
534 #endif // EIGEN_REAL_SCHUR_H EIGEN_DEVICE_FUNC Index rows() const
Definition: EigenBase.h:58
Definition: RealSchur.h:54
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Makes *this as a Givens rotation G such that applying to the left of the vector yields: ...
Definition: Jacobi.h:149
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: bench_norm.cpp:85
Definition: ForwardDeclarations.h:263
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
This is an overloaded version of DenseCoeffsBase<Derived,WriteAccessors>::coeffRef(Index,Index) const provided to by-pass the creation of an evaluator of the expression, thus saving compilation efforts.
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_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:273
MatrixHReturnType matrixH() const
Constructs the Hessenberg matrix H in the decomposition.
Definition: HessenbergDecomposition.h:262
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
This is an overloaded version of DenseCoeffsBase<Derived,ReadOnlyAccessors>::coeff(Index,Index) const provided to by-pass the creation of an evaluator of the expression, thus saving compilation efforts.
Definition: PlainObjectBase.h:154
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
Computation was successful.
Definition: Constants.h:432
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
Resizes to the given size, and sets all coefficients in this expression to the given value val...
Definition: CwiseNullaryOp.h:341
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
EIGEN_DEVICE_FUNC Index cols() const
Definition: EigenBase.h:61
HouseholderSequenceType matrixQ() const
Reconstructs the orthogonal matrix Q in the decomposition.
Definition: HessenbergDecomposition.h:234
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
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
Enum for reporting the status of a computation.
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:44
Iterative procedure did not converge.
Definition: Constants.h:436
JacobiRotation adjoint() const
Returns the adjoint transformation.
Definition: Jacobi.h:62