Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
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 <[email protected]>
5// Copyright (C) 2010,2012 Jitse Niesen <[email protected]>
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
16// IWYU pragma: private
17#include "./InternalHeaderCheck.h"
18
19namespace Eigen {
20
57template <typename MatrixType_>
58class RealSchur {
59 public:
60 typedef MatrixType_ MatrixType;
61 enum {
62 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
63 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
64 Options = internal::traits<MatrixType>::Options,
65 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
66 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
67 };
68 typedef typename MatrixType::Scalar Scalar;
69 typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
71
74
86 explicit RealSchur(Index size = RowsAtCompileTime == Dynamic ? 1 : RowsAtCompileTime)
87 : m_matT(size, size),
88 m_matU(size, size),
89 m_workspaceVector(size),
90 m_hess(size),
91 m_isInitialized(false),
92 m_matUisUptodate(false),
93 m_maxIters(-1) {}
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 compute(matrix.derived(), computeU);
115 }
116
128 const MatrixType& matrixU() const {
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 eigen_assert(m_isInitialized && "RealSchur is not initialized.");
146 return m_matT;
147 }
148
168 template <typename InputType>
169 RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
170
188 template <typename HessMatrixType, typename OrthMatrixType>
189 RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU);
195 eigen_assert(m_isInitialized && "RealSchur is not initialized.");
196 return m_info;
197 }
198
205 m_maxIters = maxIters;
206 return *this;
207 }
208
210 Index getMaxIterations() { return m_maxIters; }
211
217 static const int m_maxIterationsPerRow = 40;
218
219 private:
220 MatrixType m_matT;
221 MatrixType m_matU;
222 ColumnVectorType m_workspaceVector;
224 ComputationInfo m_info;
225 bool m_isInitialized;
226 bool m_matUisUptodate;
227 Index m_maxIters;
228
230
231 Scalar computeNormOfT();
232 Index findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero);
233 void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
234 void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
235 void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
236 void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector,
237 Scalar* workspace);
238};
239
240template <typename MatrixType>
241template <typename InputType>
243 const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
244
245 eigen_assert(matrix.cols() == matrix.rows());
246 Index maxIters = m_maxIters;
247 if (maxIters == -1) maxIters = m_maxIterationsPerRow * matrix.rows();
248
249 Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
250 if (scale < considerAsZero) {
251 m_matT.setZero(matrix.rows(), matrix.cols());
252 if (computeU) m_matU.setIdentity(matrix.rows(), matrix.cols());
253 m_info = Success;
254 m_isInitialized = true;
255 m_matUisUptodate = computeU;
256 return *this;
257 }
258
259 // Step 1. Reduce to Hessenberg form
260 m_hess.compute(matrix.derived() / scale);
261
262 // Step 2. Reduce to real Schur form
263 // Note: we copy m_hess.matrixQ() into m_matU here and not in computeFromHessenberg
264 // to be able to pass our working-space buffer for the Householder to Dense evaluation.
265 m_workspaceVector.resize(matrix.cols());
266 if (computeU) m_hess.matrixQ().evalTo(m_matU, m_workspaceVector);
267 computeFromHessenberg(m_hess.matrixH(), m_matU, computeU);
268
269 m_matT *= scale;
270
271 return *this;
272}
273template <typename MatrixType>
274template <typename HessMatrixType, typename OrthMatrixType>
275RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH,
276 const OrthMatrixType& matrixQ, bool computeU) {
277 using std::abs;
278
279 m_matT = matrixH;
280 m_workspaceVector.resize(m_matT.cols());
281 if (computeU && !internal::is_same_dense(m_matU, matrixQ)) m_matU = matrixQ;
282
283 Index maxIters = m_maxIters;
284 if (maxIters == -1) maxIters = m_maxIterationsPerRow * matrixH.rows();
285 Scalar* workspace = &m_workspaceVector.coeffRef(0);
286
287 // The matrix m_matT is divided in three parts.
288 // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
289 // Rows il,...,iu is the part we are working on (the active window).
290 // Rows iu+1,...,end are already brought in triangular form.
291 Index iu = m_matT.cols() - 1;
292 Index iter = 0; // iteration count for current eigenvalue
293 Index totalIter = 0; // iteration count for whole matrix
294 Scalar exshift(0); // sum of exceptional shifts
295 Scalar norm = computeNormOfT();
296 // sub-diagonal entries smaller than considerAsZero will be treated as zero.
297 // We use eps^2 to enable more precision in small eigenvalues.
298 Scalar considerAsZero =
299 numext::maxi<Scalar>(norm * numext::abs2(NumTraits<Scalar>::epsilon()), (std::numeric_limits<Scalar>::min)());
300
301 if (!numext::is_exactly_zero(norm)) {
302 while (iu >= 0) {
303 Index il = findSmallSubdiagEntry(iu, considerAsZero);
304
305 // Check for convergence
306 if (il == iu) // One root found
307 {
308 m_matT.coeffRef(iu, iu) = m_matT.coeff(iu, iu) + exshift;
309 if (iu > 0) m_matT.coeffRef(iu, iu - 1) = Scalar(0);
310 iu--;
311 iter = 0;
312 } else if (il == iu - 1) // Two roots found
313 {
314 splitOffTwoRows(iu, computeU, exshift);
315 iu -= 2;
316 iter = 0;
317 } else // No convergence yet
318 {
319 // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1
320 // -Wall -DNDEBUG )
321 Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
322 computeShift(iu, iter, exshift, shiftInfo);
323 iter = iter + 1;
324 totalIter = totalIter + 1;
325 if (totalIter > maxIters) break;
326 Index im;
327 initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
328 performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
329 }
330 }
331 }
332 if (totalIter <= maxIters)
333 m_info = Success;
334 else
335 m_info = NoConvergence;
336
337 m_isInitialized = true;
338 m_matUisUptodate = computeU;
339 return *this;
340}
341
343template <typename MatrixType>
344inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT() {
345 const Index size = m_matT.cols();
346 // FIXME to be efficient the following would requires a triangular reduxion code
347 // Scalar norm = m_matT.upper().cwiseAbs().sum()
348 // + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
349 Scalar norm(0);
350 for (Index j = 0; j < size; ++j) norm += m_matT.col(j).segment(0, (std::min)(size, j + 2)).cwiseAbs().sum();
351 return norm;
352}
353
355template <typename MatrixType>
356inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero) {
357 using std::abs;
358 Index res = iu;
359 while (res > 0) {
360 Scalar s = abs(m_matT.coeff(res - 1, res - 1)) + abs(m_matT.coeff(res, res));
361
362 s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
363
364 if (abs(m_matT.coeff(res, res - 1)) <= s) break;
365 res--;
366 }
367 return res;
368}
369
371template <typename MatrixType>
372inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift) {
373 using std::abs;
374 using std::sqrt;
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) m_matU.applyOnTheRight(iu - 1, iu, rot);
397 }
398
399 if (iu > 1) m_matT.coeffRef(iu - 1, iu - 2) = Scalar(0);
400}
401
403template <typename MatrixType>
404inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo) {
405 using std::abs;
406 using std::sqrt;
407 shiftInfo.coeffRef(0) = m_matT.coeff(iu, iu);
408 shiftInfo.coeffRef(1) = m_matT.coeff(iu - 1, iu - 1);
409 shiftInfo.coeffRef(2) = m_matT.coeff(iu, iu - 1) * m_matT.coeff(iu - 1, iu);
410
411 // Alternate exceptional shifting strategy every 16 iterations.
412 if (iter % 16 == 0) {
413 // Wilkinson's original ad hoc shift
414 if (iter % 32 != 0) {
415 exshift += shiftInfo.coeff(0);
416 for (Index i = 0; i <= iu; ++i) m_matT.coeffRef(i, i) -= shiftInfo.coeff(0);
417 Scalar s = abs(m_matT.coeff(iu, iu - 1)) + abs(m_matT.coeff(iu - 1, iu - 2));
418 shiftInfo.coeffRef(0) = Scalar(0.75) * s;
419 shiftInfo.coeffRef(1) = Scalar(0.75) * s;
420 shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
421 } else {
422 // MATLAB's new ad hoc shift
423 Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
424 s = s * s + shiftInfo.coeff(2);
425 if (s > Scalar(0)) {
426 s = sqrt(s);
427 if (shiftInfo.coeff(1) < shiftInfo.coeff(0)) s = -s;
428 s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
429 s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
430 exshift += s;
431 for (Index i = 0; i <= iu; ++i) m_matT.coeffRef(i, i) -= s;
432 shiftInfo.setConstant(Scalar(0.964));
433 }
434 }
435 }
436}
437
439template <typename MatrixType>
440inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im,
441 Vector3s& firstHouseholderVector) {
442 using std::abs;
443 Vector3s& v = firstHouseholderVector; // alias to save typing
444
445 for (im = iu - 2; im >= il; --im) {
446 const Scalar Tmm = m_matT.coeff(im, im);
447 const Scalar r = shiftInfo.coeff(0) - Tmm;
448 const Scalar s = shiftInfo.coeff(1) - Tmm;
449 v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im + 1, im) + m_matT.coeff(im, im + 1);
450 v.coeffRef(1) = m_matT.coeff(im + 1, im + 1) - Tmm - r - s;
451 v.coeffRef(2) = m_matT.coeff(im + 2, im + 1);
452 if (im == il) {
453 break;
454 }
455 const Scalar lhs = m_matT.coeff(im, im - 1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
456 const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im - 1, im - 1)) + abs(Tmm) + abs(m_matT.coeff(im + 1, im + 1)));
457 if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs) break;
458 }
459}
460
462template <typename MatrixType>
463inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU,
464 const Vector3s& firstHouseholderVector, Scalar* workspace) {
465 eigen_assert(im >= il);
466 eigen_assert(im <= iu - 2);
467
468 const Index size = m_matT.cols();
469
470 for (Index k = im; k <= iu - 2; ++k) {
471 bool firstIteration = (k == im);
472
473 Vector3s v;
474 if (firstIteration)
475 v = firstHouseholderVector;
476 else
477 v = m_matT.template block<3, 1>(k, k - 1);
478
479 Scalar tau, beta;
480 Matrix<Scalar, 2, 1> ess;
481 v.makeHouseholder(ess, tau, beta);
482
483 if (!numext::is_exactly_zero(beta)) // if v is not zero
484 {
485 if (firstIteration && k > il)
486 m_matT.coeffRef(k, k - 1) = -m_matT.coeff(k, k - 1);
487 else if (!firstIteration)
488 m_matT.coeffRef(k, k - 1) = beta;
489
490 // These Householder transformations form the O(n^3) part of the algorithm
491 m_matT.block(k, k, 3, size - k).applyHouseholderOnTheLeft(ess, tau, workspace);
492 m_matT.block(0, k, (std::min)(iu, k + 3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
493 if (computeU) m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
494 }
495 }
496
497 Matrix<Scalar, 2, 1> v = m_matT.template block<2, 1>(iu - 1, iu - 2);
498 Scalar tau, beta;
499 Matrix<Scalar, 1, 1> ess;
500 v.makeHouseholder(ess, tau, beta);
501
502 if (!numext::is_exactly_zero(beta)) // if v is not zero
503 {
504 m_matT.coeffRef(iu - 1, iu - 2) = beta;
505 m_matT.block(iu - 1, iu - 1, 2, size - iu + 1).applyHouseholderOnTheLeft(ess, tau, workspace);
506 m_matT.block(0, iu - 1, iu + 1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
507 if (computeU) m_matU.block(0, iu - 1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
508 }
509
510 // clean up pollution due to round-off errors
511 for (Index i = im + 2; i <= iu; ++i) {
512 m_matT.coeffRef(i, i - 2) = Scalar(0);
513 if (i > im + 2) m_matT.coeffRef(i, i - 3) = Scalar(0);
514 }
515}
516
517} // end namespace Eigen
518
519#endif // EIGEN_REAL_SCHUR_H
Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation.
Definition HessenbergDecomposition.h:61
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Performs a real Schur decomposition of a square matrix.
Definition RealSchur.h:58
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition RealSchur.h:128
ComputationInfo info() const
Reports whether previous computation was successful.
Definition RealSchur.h:194
Index getMaxIterations()
Returns the maximum number of iterations.
Definition RealSchur.h:210
RealSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition RealSchur.h:86
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:217
RealSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition RealSchur.h:204
RealSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T.
Eigen::Index Index
Definition RealSchur.h:70
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition RealSchur.h:144
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25
Definition EigenBase.h:33
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition EigenBase.h:61
Derived & derived()
Definition EigenBase.h:49
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition EigenBase.h:59