Eigen  3.3.0
 
Loading...
Searching...
No Matches
RealSchur.h
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
14#include "./HessenbergDecomposition.h"
15
16namespace Eigen {
17
54template<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;
68
71
83 explicit RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
84 : m_matT(size, size),
85 m_matU(size, size),
86 m_workspaceVector(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
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;
231 ComputationInfo m_info;
232 bool m_isInitialized;
233 bool m_matUisUptodate;
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
247template<typename MatrixType>
248template<typename InputType>
250{
251 eigen_assert(matrix.cols() == matrix.rows());
252 Index maxIters = m_maxIters;
253 if (maxIters == -1)
254 maxIters = m_maxIterationsPerRow * matrix.rows();
255
256 Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
257
258 // Step 1. Reduce to Hessenberg form
259 m_hess.compute(matrix.derived()/scale);
260
261 // Step 2. Reduce to real Schur form
262 computeFromHessenberg(m_hess.matrixH(), m_hess.matrixQ(), computeU);
263
264 m_matT *= scale;
265
266 return *this;
267}
268template<typename MatrixType>
269template<typename HessMatrixType, typename OrthMatrixType>
270RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU)
271{
272 using std::abs;
273
274 m_matT = matrixH;
275 if(computeU)
276 m_matU = matrixQ;
277
278 Index maxIters = m_maxIters;
279 if (maxIters == -1)
280 maxIters = m_maxIterationsPerRow * matrixH.rows();
281 m_workspaceVector.resize(m_matT.cols());
282 Scalar* workspace = &m_workspaceVector.coeffRef(0);
283
284 // The matrix m_matT is divided in three parts.
285 // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
286 // Rows il,...,iu is the part we are working on (the active window).
287 // Rows iu+1,...,end are already brought in triangular form.
288 Index iu = m_matT.cols() - 1;
289 Index iter = 0; // iteration count for current eigenvalue
290 Index totalIter = 0; // iteration count for whole matrix
291 Scalar exshift(0); // sum of exceptional shifts
292 Scalar norm = computeNormOfT();
293
294 if(norm!=0)
295 {
296 while (iu >= 0)
297 {
298 Index il = findSmallSubdiagEntry(iu);
299
300 // Check for convergence
301 if (il == iu) // One root found
302 {
303 m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
304 if (iu > 0)
305 m_matT.coeffRef(iu, iu-1) = Scalar(0);
306 iu--;
307 iter = 0;
308 }
309 else if (il == iu-1) // Two roots found
310 {
311 splitOffTwoRows(iu, computeU, exshift);
312 iu -= 2;
313 iter = 0;
314 }
315 else // No convergence yet
316 {
317 // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
318 Vector3s firstHouseholderVector(0,0,0), shiftInfo;
319 computeShift(iu, iter, exshift, shiftInfo);
320 iter = iter + 1;
321 totalIter = totalIter + 1;
322 if (totalIter > maxIters) break;
323 Index im;
324 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
325 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
326 }
327 }
328 }
329 if(totalIter <= maxIters)
330 m_info = Success;
331 else
332 m_info = NoConvergence;
333
334 m_isInitialized = true;
335 m_matUisUptodate = computeU;
336 return *this;
337}
338
340template<typename MatrixType>
341inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
342{
343 const Index size = m_matT.cols();
344 // FIXME to be efficient the following would requires a triangular reduxion code
345 // Scalar norm = m_matT.upper().cwiseAbs().sum()
346 // + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
347 Scalar norm(0);
348 for (Index j = 0; j < size; ++j)
349 norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
350 return norm;
351}
352
354template<typename MatrixType>
355inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu)
356{
357 using std::abs;
358 Index res = iu;
359 while (res > 0)
360 {
361 Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
362 if (abs(m_matT.coeff(res,res-1)) <= NumTraits<Scalar>::epsilon() * s)
363 break;
364 res--;
365 }
366 return res;
367}
368
370template<typename MatrixType>
371inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
372{
373 using std::sqrt;
374 using std::abs;
375 const Index size = m_matT.cols();
376
377 // The eigenvalues of the 2x2 matrix [a b; c d] are
378 // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
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); // q = tr^2 / 4 - det = discr/4
381 m_matT.coeffRef(iu,iu) += exshift;
382 m_matT.coeffRef(iu-1,iu-1) += exshift;
383
384 if (q >= Scalar(0)) // Two real eigenvalues
385 {
386 Scalar z = sqrt(abs(q));
387 JacobiRotation<Scalar> rot;
388 if (p >= Scalar(0))
389 rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
390 else
391 rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
392
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);
396 if (computeU)
397 m_matU.applyOnTheRight(iu-1, iu, rot);
398 }
399
400 if (iu > 1)
401 m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
402}
403
405template<typename MatrixType>
406inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
407{
408 using std::sqrt;
409 using std::abs;
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);
413
414 // Wilkinson's original ad hoc shift
415 if (iter == 10)
416 {
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;
424 }
425
426 // MATLAB's new ad hoc shift
427 if (iter == 30)
428 {
429 Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
430 s = s * s + shiftInfo.coeff(2);
431 if (s > Scalar(0))
432 {
433 s = sqrt(s);
434 if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
435 s = -s;
436 s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
437 s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
438 exshift += s;
439 for (Index i = 0; i <= iu; ++i)
440 m_matT.coeffRef(i,i) -= s;
441 shiftInfo.setConstant(Scalar(0.964));
442 }
443 }
444}
445
447template<typename MatrixType>
448inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
449{
450 using std::abs;
451 Vector3s& v = firstHouseholderVector; // alias to save typing
452
453 for (im = iu-2; im >= il; --im)
454 {
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);
461 if (im == il) {
462 break;
463 }
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)));
466 if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
467 break;
468 }
469}
470
472template<typename MatrixType>
473inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
474{
475 eigen_assert(im >= il);
476 eigen_assert(im <= iu-2);
477
478 const Index size = m_matT.cols();
479
480 for (Index k = im; k <= iu-2; ++k)
481 {
482 bool firstIteration = (k == im);
483
484 Vector3s v;
485 if (firstIteration)
486 v = firstHouseholderVector;
487 else
488 v = m_matT.template block<3,1>(k,k-1);
489
490 Scalar tau, beta;
491 Matrix<Scalar, 2, 1> ess;
492 v.makeHouseholder(ess, tau, beta);
493
494 if (beta != Scalar(0)) // if v is not zero
495 {
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;
500
501 // These Householder transformations form the O(n^3) part of the algorithm
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);
504 if (computeU)
505 m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
506 }
507 }
508
509 Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
510 Scalar tau, beta;
511 Matrix<Scalar, 1, 1> ess;
512 v.makeHouseholder(ess, tau, beta);
513
514 if (beta != Scalar(0)) // if v is not zero
515 {
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);
519 if (computeU)
520 m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
521 }
522
523 // clean up pollution due to round-off errors
524 for (Index i = im+2; i <= iu; ++i)
525 {
526 m_matT.coeffRef(i,i-2) = Scalar(0);
527 if (i > im+2)
528 m_matT.coeffRef(i,i-3) = Scalar(0);
529 }
530}
531
532} // end namespace Eigen
533
534#endif // EIGEN_REAL_SCHUR_H
Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation.
Definition: HessenbergDecomposition.h:58
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:55
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:195
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:127
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition: RealSchur.h:223
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealSchur.h:83
Eigen::Index Index
Definition: RealSchur.h:67
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:144
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:213
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:206
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
@ Success
Definition: Constants.h:432
@ NoConvergence
Definition: Constants.h:436
Namespace containing all symbols from the Eigen library.
Definition: Core:287
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const int Dynamic
Definition: Constants.h:21
Definition: EigenBase.h:29
Index cols() const
Definition: EigenBase.h:61
Derived & derived()
Definition: EigenBase.h:44
Index rows() const
Definition: EigenBase.h:58