RealSchur.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2010,2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_REAL_SCHUR_H
12 #define EIGEN_REAL_SCHUR_H
13 
15 
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
56 template<typename MatrixType_> class RealSchur
57 {
58  public:
59  typedef MatrixType_ MatrixType;
60  enum {
61  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
62  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
63  Options = MatrixType::Options,
64  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
65  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
66  };
67  typedef typename MatrixType::Scalar Scalar;
68  typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
69  typedef Eigen::Index Index;
70 
73 
86  : m_matT(size, size),
87  m_matU(size, size),
89  m_hess(size),
90  m_isInitialized(false),
91  m_matUisUptodate(false),
92  m_maxIters(-1)
93  { }
94 
105  template<typename InputType>
106  explicit RealSchur(const EigenBase<InputType>& matrix, bool computeU = true)
107  : m_matT(matrix.rows(),matrix.cols()),
108  m_matU(matrix.rows(),matrix.cols()),
109  m_workspaceVector(matrix.rows()),
110  m_hess(matrix.rows()),
111  m_isInitialized(false),
112  m_matUisUptodate(false),
113  m_maxIters(-1)
114  {
115  compute(matrix.derived(), computeU);
116  }
117 
129  const MatrixType& matrixU() const
130  {
131  eigen_assert(m_isInitialized && "RealSchur is not initialized.");
132  eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
133  return m_matU;
134  }
135 
146  const MatrixType& matrixT() const
147  {
148  eigen_assert(m_isInitialized && "RealSchur is not initialized.");
149  return m_matT;
150  }
151 
171  template<typename InputType>
172  RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
173 
191  template<typename HessMatrixType, typename OrthMatrixType>
192  RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU);
198  {
199  eigen_assert(m_isInitialized && "RealSchur is not initialized.");
200  return m_info;
201  }
202 
209  {
210  m_maxIters = maxIters;
211  return *this;
212  }
213 
216  {
217  return m_maxIters;
218  }
219 
225  static const int m_maxIterationsPerRow = 40;
226 
227  private:
228 
237 
239 
241  Index findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero);
242  void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
243  void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
244  void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
245  void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
246 };
247 
248 
249 template<typename MatrixType>
250 template<typename InputType>
252 {
253  const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
254 
255  eigen_assert(matrix.cols() == matrix.rows());
256  Index maxIters = m_maxIters;
257  if (maxIters == -1)
258  maxIters = m_maxIterationsPerRow * matrix.rows();
259 
260  Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
261  if(scale<considerAsZero)
262  {
263  m_matT.setZero(matrix.rows(),matrix.cols());
264  if(computeU)
265  m_matU.setIdentity(matrix.rows(),matrix.cols());
266  m_info = Success;
267  m_isInitialized = true;
268  m_matUisUptodate = computeU;
269  return *this;
270  }
271 
272  // Step 1. Reduce to Hessenberg form
273  m_hess.compute(matrix.derived()/scale);
274 
275  // Step 2. Reduce to real Schur form
276  // Note: we copy m_hess.matrixQ() into m_matU here and not in computeFromHessenberg
277  // to be able to pass our working-space buffer for the Householder to Dense evaluation.
278  m_workspaceVector.resize(matrix.cols());
279  if(computeU)
280  m_hess.matrixQ().evalTo(m_matU, m_workspaceVector);
281  computeFromHessenberg(m_hess.matrixH(), m_matU, computeU);
282 
283  m_matT *= scale;
284 
285  return *this;
286 }
287 template<typename MatrixType>
288 template<typename HessMatrixType, typename OrthMatrixType>
289 RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU)
290 {
291  using std::abs;
292 
293  m_matT = matrixH;
294  m_workspaceVector.resize(m_matT.cols());
295  if(computeU && !internal::is_same_dense(m_matU,matrixQ))
296  m_matU = matrixQ;
297 
298  Index maxIters = m_maxIters;
299  if (maxIters == -1)
300  maxIters = m_maxIterationsPerRow * matrixH.rows();
301  Scalar* workspace = &m_workspaceVector.coeffRef(0);
302 
303  // The matrix m_matT is divided in three parts.
304  // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
305  // Rows il,...,iu is the part we are working on (the active window).
306  // Rows iu+1,...,end are already brought in triangular form.
307  Index iu = m_matT.cols() - 1;
308  Index iter = 0; // iteration count for current eigenvalue
309  Index totalIter = 0; // iteration count for whole matrix
310  Scalar exshift(0); // sum of exceptional shifts
311  Scalar norm = computeNormOfT();
312  // sub-diagonal entries smaller than considerAsZero will be treated as zero.
313  // We use eps^2 to enable more precision in small eigenvalues.
314  Scalar considerAsZero = numext::maxi<Scalar>( norm * numext::abs2(NumTraits<Scalar>::epsilon()),
316 
317  if(!numext::is_exactly_zero(norm))
318  {
319  while (iu >= 0)
320  {
321  Index il = findSmallSubdiagEntry(iu,considerAsZero);
322 
323  // Check for convergence
324  if (il == iu) // One root found
325  {
326  m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
327  if (iu > 0)
328  m_matT.coeffRef(iu, iu-1) = Scalar(0);
329  iu--;
330  iter = 0;
331  }
332  else if (il == iu-1) // Two roots found
333  {
334  splitOffTwoRows(iu, computeU, exshift);
335  iu -= 2;
336  iter = 0;
337  }
338  else // No convergence yet
339  {
340  // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
341  Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
342  computeShift(iu, iter, exshift, shiftInfo);
343  iter = iter + 1;
344  totalIter = totalIter + 1;
345  if (totalIter > maxIters) break;
346  Index im;
347  initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
348  performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
349  }
350  }
351  }
352  if(totalIter <= maxIters)
353  m_info = Success;
354  else
355  m_info = NoConvergence;
356 
357  m_isInitialized = true;
358  m_matUisUptodate = computeU;
359  return *this;
360 }
361 
363 template<typename MatrixType>
364 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
365 {
366  const Index size = m_matT.cols();
367  // FIXME to be efficient the following would requires a triangular reduxion code
368  // Scalar norm = m_matT.upper().cwiseAbs().sum()
369  // + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
370  Scalar norm(0);
371  for (Index j = 0; j < size; ++j)
372  norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
373  return norm;
374 }
375 
377 template<typename MatrixType>
379 {
380  using std::abs;
381  Index res = iu;
382  while (res > 0)
383  {
384  Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
385 
386  s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
387 
388  if (abs(m_matT.coeff(res,res-1)) <= s)
389  break;
390  res--;
391  }
392  return res;
393 }
394 
396 template<typename MatrixType>
397 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
398 {
399  using std::sqrt;
400  using std::abs;
401  const Index size = m_matT.cols();
402 
403  // The eigenvalues of the 2x2 matrix [a b; c d] are
404  // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
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); // q = tr^2 / 4 - det = discr/4
407  m_matT.coeffRef(iu,iu) += exshift;
408  m_matT.coeffRef(iu-1,iu-1) += exshift;
409 
410  if (q >= Scalar(0)) // Two real eigenvalues
411  {
412  Scalar z = sqrt(abs(q));
414  if (p >= Scalar(0))
415  rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
416  else
417  rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
418 
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);
422  if (computeU)
423  m_matU.applyOnTheRight(iu-1, iu, rot);
424  }
425 
426  if (iu > 1)
427  m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
428 }
429 
431 template<typename MatrixType>
432 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
433 {
434  using std::sqrt;
435  using std::abs;
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);
439 
440  // Wilkinson's original ad hoc shift
441  if (iter == 10)
442  {
443  exshift += shiftInfo.coeff(0);
444  for (Index i = 0; i <= iu; ++i)
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));
447  shiftInfo.coeffRef(0) = Scalar(0.75) * s;
448  shiftInfo.coeffRef(1) = Scalar(0.75) * s;
449  shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
450  }
451 
452  // MATLAB's new ad hoc shift
453  if (iter == 30)
454  {
455  Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
456  s = s * s + shiftInfo.coeff(2);
457  if (s > Scalar(0))
458  {
459  s = sqrt(s);
460  if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
461  s = -s;
462  s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
463  s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
464  exshift += s;
465  for (Index i = 0; i <= iu; ++i)
466  m_matT.coeffRef(i,i) -= s;
467  shiftInfo.setConstant(Scalar(0.964));
468  }
469  }
470 }
471 
473 template<typename MatrixType>
474 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
475 {
476  using std::abs;
477  Vector3s& v = firstHouseholderVector; // alias to save typing
478 
479  for (im = iu-2; im >= il; --im)
480  {
481  const Scalar Tmm = m_matT.coeff(im,im);
482  const Scalar r = shiftInfo.coeff(0) - Tmm;
483  const Scalar s = shiftInfo.coeff(1) - Tmm;
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);
487  if (im == il) {
488  break;
489  }
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)));
492  if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
493  break;
494  }
495 }
496 
498 template<typename MatrixType>
499 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
500 {
501  eigen_assert(im >= il);
502  eigen_assert(im <= iu-2);
503 
504  const Index size = m_matT.cols();
505 
506  for (Index k = im; k <= iu-2; ++k)
507  {
508  bool firstIteration = (k == im);
509 
510  Vector3s v;
511  if (firstIteration)
512  v = firstHouseholderVector;
513  else
514  v = m_matT.template block<3,1>(k,k-1);
515 
516  Scalar tau, beta;
518  v.makeHouseholder(ess, tau, beta);
519 
520  if (!numext::is_exactly_zero(beta)) // if v is not zero
521  {
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;
526 
527  // These Householder transformations form the O(n^3) part of the algorithm
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);
530  if (computeU)
531  m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
532  }
533  }
534 
535  Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
536  Scalar tau, beta;
538  v.makeHouseholder(ess, tau, beta);
539 
540  if (!numext::is_exactly_zero(beta)) // if v is not zero
541  {
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);
545  if (computeU)
546  m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
547  }
548 
549  // clean up pollution due to round-off errors
550  for (Index i = im+2; i <= iu; ++i)
551  {
552  m_matT.coeffRef(i,i-2) = Scalar(0);
553  if (i > im+2)
554  m_matT.coeffRef(i,i-3) = Scalar(0);
555  }
556 }
557 
558 } // end namespace Eigen
559 
560 #endif // EIGEN_REAL_SCHUR_H
const Abs2ReturnType abs2() const
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
Array< int, Dynamic, 1 > v
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
float * p
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:37
JacobiRotation adjoint() const
Definition: Jacobi.h:69
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:164
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.
Definition: RealSchur.h:57
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:208
Matrix< Scalar, 3, 1 > Vector3s
Definition: RealSchur.h:238
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: RealSchur.h:68
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:129
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.
Definition: RealSchur.h:197
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Definition: RealSchur.h:71
void initFrancisQRStep(Index il, Index iu, const Vector3s &shiftInfo, Index &im, Vector3s &firstHouseholderVector)
Definition: RealSchur.h:474
void computeShift(Index iu, Index iter, Scalar &exshift, Vector3s &shiftInfo)
Definition: RealSchur.h:432
bool m_matUisUptodate
Definition: RealSchur.h:235
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:215
void splitOffTwoRows(Index iu, bool computeU, const Scalar &exshift)
Definition: RealSchur.h:397
MatrixType::Scalar Scalar
Definition: RealSchur.h:67
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:85
MatrixType m_matU
Definition: RealSchur.h:230
ColumnVectorType m_workspaceVector
Definition: RealSchur.h:231
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: RealSchur.h:72
MatrixType m_matT
Definition: RealSchur.h:229
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:106
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:225
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:146
Scalar computeNormOfT()
Definition: RealSchur.h:364
ComputationInfo m_info
Definition: RealSchur.h:233
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
HessenbergDecomposition< MatrixType > m_hess
Definition: RealSchur.h:232
Index findSmallSubdiagEntry(Index iu, const Scalar &considerAsZero)
Definition: RealSchur.h:378
bool m_isInitialized
Definition: RealSchur.h:234
MatrixType_ MatrixType
Definition: RealSchur.h:59
Eigen::Index Index
Definition: RealSchur.h:69
void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s &firstHouseholderVector, Scalar *workspace)
Definition: RealSchur.h:499
ComputationInfo
Definition: Constants.h:444
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
bool is_same_dense(const T1 &mat1, const T2 &mat2, std::enable_if_t< possibly_same_dense< T1, T2 >::value > *=0)
Definition: XprHelper.h:745
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24
Derived & derived()
Definition: EigenBase.h:48
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231
std::ptrdiff_t j