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 namespace Eigen {
17 
54 template<typename _MatrixType> class RealSchur
55 {
56  public:
57  typedef _MatrixType MatrixType;
58  enum {
59  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
60  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
61  Options = MatrixType::Options,
62  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
63  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
64  };
65  typedef typename MatrixType::Scalar Scalar;
66  typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
67  typedef Eigen::Index Index;
68 
71 
84  : m_matT(size, size),
85  m_matU(size, size),
87  m_hess(size),
88  m_isInitialized(false),
89  m_matUisUptodate(false),
90  m_maxIters(-1)
91  { }
92 
103  template<typename InputType>
104  explicit RealSchur(const EigenBase<InputType>& matrix, bool computeU = true)
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),
111  m_maxIters(-1)
112  {
113  compute(matrix.derived(), computeU);
114  }
115 
127  const MatrixType& matrixU() const
128  {
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.");
131  return m_matU;
132  }
133 
144  const MatrixType& matrixT() const
145  {
146  eigen_assert(m_isInitialized && "RealSchur is not initialized.");
147  return m_matT;
148  }
149 
169  template<typename InputType>
170  RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
171 
189  template<typename HessMatrixType, typename OrthMatrixType>
190  RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU);
196  {
197  eigen_assert(m_isInitialized && "RealSchur is not initialized.");
198  return m_info;
199  }
200 
206  RealSchur& setMaxIterations(Index maxIters)
207  {
208  m_maxIters = maxIters;
209  return *this;
210  }
211 
214  {
215  return m_maxIters;
216  }
217 
223  static const int m_maxIterationsPerRow = 40;
224 
225  private:
226 
227  MatrixType m_matT;
228  MatrixType m_matU;
229  ColumnVectorType m_workspaceVector;
234  Index m_maxIters;
235 
237 
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);
244 };
245 
246 
247 template<typename MatrixType>
248 template<typename InputType>
250 {
251  const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
252 
253  eigen_assert(matrix.cols() == matrix.rows());
254  Index maxIters = m_maxIters;
255  if (maxIters == -1)
256  maxIters = m_maxIterationsPerRow * matrix.rows();
257 
258  Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
259  if(scale<considerAsZero)
260  {
261  m_matT.setZero(matrix.rows(),matrix.cols());
262  if(computeU)
263  m_matU.setIdentity(matrix.rows(),matrix.cols());
264  m_info = Success;
265  m_isInitialized = true;
266  m_matUisUptodate = computeU;
267  return *this;
268  }
269 
270  // Step 1. Reduce to Hessenberg form
271  m_hess.compute(matrix.derived()/scale);
272 
273  // Step 2. Reduce to real Schur form
275 
276  m_matT *= scale;
277 
278  return *this;
279 }
280 template<typename MatrixType>
281 template<typename HessMatrixType, typename OrthMatrixType>
282 RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU)
283 {
284  using std::abs;
285 
286  m_matT = matrixH;
287  if(computeU)
288  m_matU = matrixQ;
289 
290  Index maxIters = m_maxIters;
291  if (maxIters == -1)
292  maxIters = m_maxIterationsPerRow * matrixH.rows();
294  Scalar* workspace = &m_workspaceVector.coeffRef(0);
295 
296  // The matrix m_matT is divided in three parts.
297  // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
298  // Rows il,...,iu is the part we are working on (the active window).
299  // Rows iu+1,...,end are already brought in triangular form.
300  Index iu = m_matT.cols() - 1;
301  Index iter = 0; // iteration count for current eigenvalue
302  Index totalIter = 0; // iteration count for whole matrix
303  Scalar exshift(0); // sum of exceptional shifts
304  Scalar norm = computeNormOfT();
305 
306  if(norm!=0)
307  {
308  while (iu >= 0)
309  {
310  Index il = findSmallSubdiagEntry(iu);
311 
312  // Check for convergence
313  if (il == iu) // One root found
314  {
315  m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
316  if (iu > 0)
317  m_matT.coeffRef(iu, iu-1) = Scalar(0);
318  iu--;
319  iter = 0;
320  }
321  else if (il == iu-1) // Two roots found
322  {
323  splitOffTwoRows(iu, computeU, exshift);
324  iu -= 2;
325  iter = 0;
326  }
327  else // No convergence yet
328  {
329  // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
330  Vector3s firstHouseholderVector(0,0,0), shiftInfo;
331  computeShift(iu, iter, exshift, shiftInfo);
332  iter = iter + 1;
333  totalIter = totalIter + 1;
334  if (totalIter > maxIters) break;
335  Index im;
336  initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
337  performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
338  }
339  }
340  }
341  if(totalIter <= maxIters)
342  m_info = Success;
343  else
345 
346  m_isInitialized = true;
347  m_matUisUptodate = computeU;
348  return *this;
349 }
350 
352 template<typename MatrixType>
353 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
354 {
355  const Index size = m_matT.cols();
356  // FIXME to be efficient the following would requires a triangular reduxion code
357  // Scalar norm = m_matT.upper().cwiseAbs().sum()
358  // + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
359  Scalar norm(0);
360  for (Index j = 0; j < size; ++j)
361  norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
362  return norm;
363 }
364 
366 template<typename MatrixType>
368 {
369  using std::abs;
370  Index res = iu;
371  while (res > 0)
372  {
373  Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
374  if (abs(m_matT.coeff(res,res-1)) <= NumTraits<Scalar>::epsilon() * s)
375  break;
376  res--;
377  }
378  return res;
379 }
380 
382 template<typename MatrixType>
383 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
384 {
385  using std::sqrt;
386  using std::abs;
387  const Index size = m_matT.cols();
388 
389  // The eigenvalues of the 2x2 matrix [a b; c d] are
390  // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
391  Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
392  Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu); // q = tr^2 / 4 - det = discr/4
393  m_matT.coeffRef(iu,iu) += exshift;
394  m_matT.coeffRef(iu-1,iu-1) += exshift;
395 
396  if (q >= Scalar(0)) // Two real eigenvalues
397  {
398  Scalar z = sqrt(abs(q));
400  if (p >= Scalar(0))
401  rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
402  else
403  rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
404 
405  m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
406  m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
407  m_matT.coeffRef(iu, iu-1) = Scalar(0);
408  if (computeU)
409  m_matU.applyOnTheRight(iu-1, iu, rot);
410  }
411 
412  if (iu > 1)
413  m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
414 }
415 
417 template<typename MatrixType>
418 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
419 {
420  using std::sqrt;
421  using std::abs;
422  shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
423  shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
424  shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
425 
426  // Wilkinson's original ad hoc shift
427  if (iter == 10)
428  {
429  exshift += shiftInfo.coeff(0);
430  for (Index i = 0; i <= iu; ++i)
431  m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
432  Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
433  shiftInfo.coeffRef(0) = Scalar(0.75) * s;
434  shiftInfo.coeffRef(1) = Scalar(0.75) * s;
435  shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
436  }
437 
438  // MATLAB's new ad hoc shift
439  if (iter == 30)
440  {
441  Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
442  s = s * s + shiftInfo.coeff(2);
443  if (s > Scalar(0))
444  {
445  s = sqrt(s);
446  if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
447  s = -s;
448  s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
449  s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
450  exshift += s;
451  for (Index i = 0; i <= iu; ++i)
452  m_matT.coeffRef(i,i) -= s;
453  shiftInfo.setConstant(Scalar(0.964));
454  }
455  }
456 }
457 
459 template<typename MatrixType>
460 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
461 {
462  using std::abs;
463  Vector3s& v = firstHouseholderVector; // alias to save typing
464 
465  for (im = iu-2; im >= il; --im)
466  {
467  const Scalar Tmm = m_matT.coeff(im,im);
468  const Scalar r = shiftInfo.coeff(0) - Tmm;
469  const Scalar s = shiftInfo.coeff(1) - Tmm;
470  v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
471  v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
472  v.coeffRef(2) = m_matT.coeff(im+2,im+1);
473  if (im == il) {
474  break;
475  }
476  const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
477  const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
478  if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
479  break;
480  }
481 }
482 
484 template<typename MatrixType>
485 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
486 {
487  eigen_assert(im >= il);
488  eigen_assert(im <= iu-2);
489 
490  const Index size = m_matT.cols();
491 
492  for (Index k = im; k <= iu-2; ++k)
493  {
494  bool firstIteration = (k == im);
495 
496  Vector3s v;
497  if (firstIteration)
498  v = firstHouseholderVector;
499  else
500  v = m_matT.template block<3,1>(k,k-1);
501 
502  Scalar tau, beta;
504  v.makeHouseholder(ess, tau, beta);
505 
506  if (beta != Scalar(0)) // if v is not zero
507  {
508  if (firstIteration && k > il)
509  m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
510  else if (!firstIteration)
511  m_matT.coeffRef(k,k-1) = beta;
512 
513  // These Householder transformations form the O(n^3) part of the algorithm
514  m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
515  m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
516  if (computeU)
517  m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
518  }
519  }
520 
521  Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
522  Scalar tau, beta;
524  v.makeHouseholder(ess, tau, beta);
525 
526  if (beta != Scalar(0)) // if v is not zero
527  {
528  m_matT.coeffRef(iu-1, iu-2) = beta;
529  m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
530  m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
531  if (computeU)
532  m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
533  }
534 
535  // clean up pollution due to round-off errors
536  for (Index i = im+2; i <= iu; ++i)
537  {
538  m_matT.coeffRef(i,i-2) = Scalar(0);
539  if (i > im+2)
540  m_matT.coeffRef(i,i-3) = Scalar(0);
541  }
542 }
543 
544 } // end namespace Eigen
545 
546 #endif // EIGEN_REAL_SCHUR_H
bool m_matUisUptodate
Definition: RealSchur.h:233
HouseholderSequenceType matrixQ() const
Reconstructs the orthogonal matrix Q in the decomposition.
MatrixHReturnType matrixH() const
Constructs the Hessenberg matrix H in the decomposition.
Index findSmallSubdiagEntry(Index iu)
Definition: RealSchur.h:367
ComputationInfo m_info
Definition: RealSchur.h:231
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
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
XmlRpcServer s
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition: LDLT.h:16
static constexpr size_t size(Tuple< Args... > &)
Provides access to the number of elements in a tuple as a compile-time constant expression.
Rotation given by a cosine-sine pair.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
MatrixType m_matT
Definition: RealSchur.h:227
Scalar computeNormOfT()
Definition: RealSchur.h:353
void computeShift(Index iu, Index iter, Scalar &exshift, Vector3s &shiftInfo)
Definition: RealSchur.h:418
MatrixType m_matU
Definition: RealSchur.h:228
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
HessenbergDecomposition< MatrixType > m_hess
Definition: RealSchur.h:230
Eigen::Index Index
Definition: RealSchur.h:67
bool m_isInitialized
Definition: RealSchur.h:232
HessenbergDecomposition & compute(const EigenBase< InputType > &matrix)
Computes Hessenberg decomposition of given matrix.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
_MatrixType MatrixType
Definition: RealSchur.h:57
MatrixType::Scalar Scalar
Definition: RealSchur.h:65
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Definition: RealSchur.h:69
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
#define eigen_assert(x)
Definition: Macros.h:577
JacobiRotation adjoint() const
Definition: Jacobi.h:62
EIGEN_DEVICE_FUNC const Scalar & q
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s &firstHouseholderVector, Scalar *workspace)
Definition: RealSchur.h:485
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
EIGEN_DEVICE_FUNC Index cols() const
Definition: EigenBase.h:62
TFSIMD_FORCE_INLINE const tfScalar & z() const
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: RealSchur.h:70
Matrix< Scalar, 3, 1 > Vector3s
Definition: RealSchur.h:236
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
EIGEN_DEVICE_FUNC Index rows() const
Definition: EigenBase.h:59
void initFrancisQRStep(Index il, Index iu, const Vector3s &shiftInfo, Index &im, Vector3s &firstHouseholderVector)
Definition: RealSchur.h:460
void splitOffTwoRows(Index iu, bool computeU, const Scalar &exshift)
Definition: RealSchur.h:383
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: RealSchur.h:66
const int Dynamic
Definition: Constants.h:21
ColumnVectorType m_workspaceVector
Definition: RealSchur.h:229
RealSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes real Schur decomposition of given matrix.
Definition: RealSchur.h:104
ComputationInfo
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:45
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:08:43