11 #ifndef EIGEN_REAL_SCHUR_H
12 #define EIGEN_REAL_SCHUR_H
65 typedef typename MatrixType::Scalar
Scalar;
103 template<
typename InputType>
105 :
m_matT(matrix.rows(),matrix.cols()),
106 m_matU(matrix.rows(),matrix.cols()),
169 template<
typename InputType>
189 template<
typename HessMatrixType,
typename OrthMatrixType>
247 template<
typename MatrixType>
248 template<
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);
285 template<
typename MatrixType>
286 template<
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;
361 template<
typename MatrixType>
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>
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));
454 s = s * s + shiftInfo.
coeff(2);
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)
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);
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>
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);
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);
536 v.makeHouseholder(ess, tau, 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 AbsReturnType abs() const
Definition: ArrayCwiseUnaryOps.h:52
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
Definition: ArrayCwiseUnaryOps.h:80
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition: ArrayCwiseUnaryOps.h:187
#define eigen_assert(x)
Definition: Macros.h:1037
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:35
EIGEN_DEVICE_FUNC void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:162
EIGEN_DEVICE_FUNC JacobiRotation adjoint() const
Definition: Jacobi.h:67
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:175
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
Definition: CwiseNullaryOp.h:361
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:152
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:55
Index m_maxIters
Definition: RealSchur.h:234
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: RealSchur.h:66
_MatrixType MatrixType
Definition: RealSchur.h:57
@ 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
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
ColumnVectorType m_workspaceVector
Definition: RealSchur.h:229
MatrixType m_matT
Definition: RealSchur.h:227
void initFrancisQRStep(Index il, Index iu, const Vector3s &shiftInfo, Index &im, Vector3s &firstHouseholderVector)
Definition: RealSchur.h:472
void computeShift(Index iu, Index iter, Scalar &exshift, Vector3s &shiftInfo)
Definition: RealSchur.h:430
MatrixType m_matU
Definition: RealSchur.h:228
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
HessenbergDecomposition< MatrixType > m_hess
Definition: RealSchur.h:230
void splitOffTwoRows(Index iu, bool computeU, const Scalar &exshift)
Definition: RealSchur.h:395
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
bool m_isInitialized
Definition: RealSchur.h:232
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
Eigen::Index Index
Definition: RealSchur.h:67
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
MatrixType::Scalar Scalar
Definition: RealSchur.h:65
Scalar computeNormOfT()
Definition: RealSchur.h:362
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Definition: RealSchur.h:69
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
Index findSmallSubdiagEntry(Index iu, const Scalar &considerAsZero)
Definition: RealSchur.h:376
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
bool m_matUisUptodate
Definition: RealSchur.h:233
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: RealSchur.h:70
ComputationInfo m_info
Definition: RealSchur.h:231
Matrix< Scalar, 3, 1 > Vector3s
Definition: RealSchur.h:236
void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s &firstHouseholderVector, Scalar *workspace)
Definition: RealSchur.h:497
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
@ Success
Definition: Constants.h:442
@ NoConvergence
Definition: Constants.h:446
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:571
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: LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const int Dynamic
Definition: Constants.h:22
Definition: EigenBase.h:30
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:46
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:63
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:60
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:233