11 #ifndef EIGEN_REAL_SCHUR_H
12 #define EIGEN_REAL_SCHUR_H
67 typedef typename MatrixType::Scalar
Scalar;
105 template<
typename InputType>
171 template<
typename InputType>
191 template<
typename HessMatrixType,
typename OrthMatrixType>
249 template<
typename MatrixType>
250 template<
typename InputType>
256 Index maxIters = m_maxIters;
258 maxIters = m_maxIterationsPerRow * matrix.
rows();
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());
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();
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>
372 norm += m_matT.col(
j).segment(0, (
std::min)(
size,
j+2)).cwiseAbs().sum();
377 template<
typename MatrixType>
396 template<
typename MatrixType>
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;
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>
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);
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));
456 s = s * s + shiftInfo.
coeff(2);
463 s = shiftInfo.
coeff(0) - shiftInfo.
coeff(2) / s;
466 m_matT.coeffRef(
i,
i) -= s;
473 template<
typename MatrixType>
479 for (im = iu-2; im >= il; --im)
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)));
498 template<
typename MatrixType>
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);
518 v.makeHouseholder(ess, tau, 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);
538 v.makeHouseholder(ess, tau, 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);
const Abs2ReturnType abs2() const
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
Array< int, Dynamic, 1 > v
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Rotation given by a cosine-sine pair.
JacobiRotation adjoint() const
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
constexpr const Scalar & coeff(Index rowId, Index colId) const
Derived & setConstant(Index size, const Scalar &val)
constexpr Scalar & coeffRef(Index rowId, Index colId)
Performs a real Schur decomposition of a square matrix.
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Matrix< Scalar, 3, 1 > Vector3s
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
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.
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
void initFrancisQRStep(Index il, Index iu, const Vector3s &shiftInfo, Index &im, Vector3s &firstHouseholderVector)
void computeShift(Index iu, Index iter, Scalar &exshift, Vector3s &shiftInfo)
Index getMaxIterations()
Returns the maximum number of iterations.
void splitOffTwoRows(Index iu, bool computeU, const Scalar &exshift)
MatrixType::Scalar Scalar
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
ColumnVectorType m_workspaceVector
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
HessenbergDecomposition< MatrixType > m_hess
Index findSmallSubdiagEntry(Index iu, const Scalar &considerAsZero)
void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s &firstHouseholderVector, Scalar *workspace)
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
bool is_same_dense(const T1 &mat1, const T2 &mat2, std::enable_if_t< possibly_same_dense< T1, T2 >::value > *=0)
bool is_exactly_zero(const X &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.