Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
ComplexSchur.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009 Claire Maurice
5// Copyright (C) 2009 Gael Guennebaud <[email protected]>
6// Copyright (C) 2010,2012 Jitse Niesen <[email protected]>
7//
8// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11
12#ifndef EIGEN_COMPLEX_SCHUR_H
13#define EIGEN_COMPLEX_SCHUR_H
14
15#include "./HessenbergDecomposition.h"
16
17// IWYU pragma: private
18#include "./InternalHeaderCheck.h"
19
20namespace Eigen {
21
22namespace internal {
23template <typename MatrixType, bool IsComplex>
24struct complex_schur_reduce_to_hessenberg;
25}
26
55template <typename MatrixType_>
57 public:
58 typedef MatrixType_ MatrixType;
59 enum {
60 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
61 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
62 Options = internal::traits<MatrixType>::Options,
63 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
65 };
66
68 typedef typename MatrixType::Scalar Scalar;
69 typedef typename NumTraits<Scalar>::Real RealScalar;
71
78 typedef std::complex<RealScalar> ComplexScalar;
79
85 typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime,
86 MaxColsAtCompileTime>
88
100 explicit ComplexSchur(Index size = RowsAtCompileTime == Dynamic ? 1 : RowsAtCompileTime)
101 : m_matT(size, size),
102 m_matU(size, size),
103 m_hess(size),
104 m_isInitialized(false),
105 m_matUisUptodate(false),
106 m_maxIters(-1) {}
107
117 template <typename InputType>
118 explicit ComplexSchur(const EigenBase<InputType>& matrix, bool computeU = true)
119 : m_matT(matrix.rows(), matrix.cols()),
120 m_matU(matrix.rows(), matrix.cols()),
121 m_hess(matrix.rows()),
122 m_isInitialized(false),
123 m_matUisUptodate(false),
124 m_maxIters(-1) {
125 compute(matrix.derived(), computeU);
126 }
127
142 const ComplexMatrixType& matrixU() const {
143 eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
144 eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the ComplexSchur decomposition.");
145 return m_matU;
146 }
147
165 const ComplexMatrixType& matrixT() const {
166 eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
167 return m_matT;
168 }
169
192 template <typename InputType>
193 ComplexSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
194
212 template <typename HessMatrixType, typename OrthMatrixType>
213 ComplexSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,
214 bool computeU = true);
215
221 eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
222 return m_info;
223 }
224
231 m_maxIters = maxIters;
232 return *this;
233 }
234
236 Index getMaxIterations() { return m_maxIters; }
237
243 static const int m_maxIterationsPerRow = 30;
244
245 protected:
246 ComplexMatrixType m_matT, m_matU;
248 ComputationInfo m_info;
249 bool m_isInitialized;
250 bool m_matUisUptodate;
251 Index m_maxIters;
252
253 private:
254 bool subdiagonalEntryIsNeglegible(Index i);
255 ComplexScalar computeShift(Index iu, Index iter);
256 void reduceToTriangularForm(bool computeU);
257 friend struct internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>;
258};
259
263template <typename MatrixType>
264inline bool ComplexSchur<MatrixType>::subdiagonalEntryIsNeglegible(Index i) {
265 RealScalar d = numext::norm1(m_matT.coeff(i, i)) + numext::norm1(m_matT.coeff(i + 1, i + 1));
266 RealScalar sd = numext::norm1(m_matT.coeff(i + 1, i));
267 if (internal::isMuchSmallerThan(sd, d, NumTraits<RealScalar>::epsilon())) {
268 m_matT.coeffRef(i + 1, i) = ComplexScalar(0);
269 return true;
270 }
271 return false;
272}
273
275template <typename MatrixType>
276typename ComplexSchur<MatrixType>::ComplexScalar ComplexSchur<MatrixType>::computeShift(Index iu, Index iter) {
277 using std::abs;
278 if ((iter == 10 || iter == 20) && iu > 1) {
279 // exceptional shift, taken from http://www.netlib.org/eispack/comqr.f
280 return abs(numext::real(m_matT.coeff(iu, iu - 1))) + abs(numext::real(m_matT.coeff(iu - 1, iu - 2)));
281 }
282
283 // compute the shift as one of the eigenvalues of t, the 2x2
284 // diagonal block on the bottom of the active submatrix
285 Matrix<ComplexScalar, 2, 2> t = m_matT.template block<2, 2>(iu - 1, iu - 1);
286 RealScalar normt = t.cwiseAbs().sum();
287 t /= normt; // the normalization by sf is to avoid under/overflow
288
289 ComplexScalar b = t.coeff(0, 1) * t.coeff(1, 0);
290 ComplexScalar c = t.coeff(0, 0) - t.coeff(1, 1);
291 ComplexScalar disc = sqrt(c * c + RealScalar(4) * b);
292 ComplexScalar det = t.coeff(0, 0) * t.coeff(1, 1) - b;
293 ComplexScalar trace = t.coeff(0, 0) + t.coeff(1, 1);
294 ComplexScalar eival1 = (trace + disc) / RealScalar(2);
295 ComplexScalar eival2 = (trace - disc) / RealScalar(2);
296 RealScalar eival1_norm = numext::norm1(eival1);
297 RealScalar eival2_norm = numext::norm1(eival2);
298 // A division by zero can only occur if eival1==eival2==0.
299 // In this case, det==0, and all we have to do is checking that eival2_norm!=0
300 if (eival1_norm > eival2_norm)
301 eival2 = det / eival1;
302 else if (!numext::is_exactly_zero(eival2_norm))
303 eival1 = det / eival2;
304
305 // choose the eigenvalue closest to the bottom entry of the diagonal
306 if (numext::norm1(eival1 - t.coeff(1, 1)) < numext::norm1(eival2 - t.coeff(1, 1)))
307 return normt * eival1;
308 else
309 return normt * eival2;
310}
311
312template <typename MatrixType>
313template <typename InputType>
314ComplexSchur<MatrixType>& ComplexSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU) {
315 m_matUisUptodate = false;
316 eigen_assert(matrix.cols() == matrix.rows());
317
318 if (matrix.cols() == 1) {
319 m_matT = matrix.derived().template cast<ComplexScalar>();
320 if (computeU) m_matU = ComplexMatrixType::Identity(1, 1);
321 m_info = Success;
322 m_isInitialized = true;
323 m_matUisUptodate = computeU;
324 return *this;
325 }
326
327 internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>::run(*this, matrix.derived(),
328 computeU);
329 computeFromHessenberg(m_matT, m_matU, computeU);
330 return *this;
331}
332
333template <typename MatrixType>
334template <typename HessMatrixType, typename OrthMatrixType>
335ComplexSchur<MatrixType>& ComplexSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH,
336 const OrthMatrixType& matrixQ,
337 bool computeU) {
338 m_matT = matrixH;
339 if (computeU) m_matU = matrixQ;
340 reduceToTriangularForm(computeU);
341 return *this;
342}
343namespace internal {
344
345/* Reduce given matrix to Hessenberg form */
346template <typename MatrixType, bool IsComplex>
347struct complex_schur_reduce_to_hessenberg {
348 // this is the implementation for the case IsComplex = true
349 static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU) {
350 _this.m_hess.compute(matrix);
351 _this.m_matT = _this.m_hess.matrixH();
352 if (computeU) _this.m_matU = _this.m_hess.matrixQ();
353 }
354};
355
356template <typename MatrixType>
357struct complex_schur_reduce_to_hessenberg<MatrixType, false> {
358 static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU) {
359 typedef typename ComplexSchur<MatrixType>::ComplexScalar ComplexScalar;
360
361 // Note: m_hess is over RealScalar; m_matT and m_matU is over ComplexScalar
362 _this.m_hess.compute(matrix);
363 _this.m_matT = _this.m_hess.matrixH().template cast<ComplexScalar>();
364 if (computeU) {
365 // This may cause an allocation which seems to be avoidable
366 MatrixType Q = _this.m_hess.matrixQ();
367 _this.m_matU = Q.template cast<ComplexScalar>();
368 }
369 }
370};
371
372} // end namespace internal
373
374// Reduce the Hessenberg matrix m_matT to triangular form by QR iteration.
375template <typename MatrixType>
376void ComplexSchur<MatrixType>::reduceToTriangularForm(bool computeU) {
377 Index maxIters = m_maxIters;
378 if (maxIters == -1) maxIters = m_maxIterationsPerRow * m_matT.rows();
379
380 // The matrix m_matT is divided in three parts.
381 // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
382 // Rows il,...,iu is the part we are working on (the active submatrix).
383 // Rows iu+1,...,end are already brought in triangular form.
384 Index iu = m_matT.cols() - 1;
385 Index il;
386 Index iter = 0; // number of iterations we are working on the (iu,iu) element
387 Index totalIter = 0; // number of iterations for whole matrix
388
389 while (true) {
390 // find iu, the bottom row of the active submatrix
391 while (iu > 0) {
392 if (!subdiagonalEntryIsNeglegible(iu - 1)) break;
393 iter = 0;
394 --iu;
395 }
396
397 // if iu is zero then we are done; the whole matrix is triangularized
398 if (iu == 0) break;
399
400 // if we spent too many iterations, we give up
401 iter++;
402 totalIter++;
403 if (totalIter > maxIters) break;
404
405 // find il, the top row of the active submatrix
406 il = iu - 1;
407 while (il > 0 && !subdiagonalEntryIsNeglegible(il - 1)) {
408 --il;
409 }
410
411 /* perform the QR step using Givens rotations. The first rotation
412 creates a bulge; the (il+2,il) element becomes nonzero. This
413 bulge is chased down to the bottom of the active submatrix. */
414
415 ComplexScalar shift = computeShift(iu, iter);
416 JacobiRotation<ComplexScalar> rot;
417 rot.makeGivens(m_matT.coeff(il, il) - shift, m_matT.coeff(il + 1, il));
418 m_matT.rightCols(m_matT.cols() - il).applyOnTheLeft(il, il + 1, rot.adjoint());
419 m_matT.topRows((std::min)(il + 2, iu) + 1).applyOnTheRight(il, il + 1, rot);
420 if (computeU) m_matU.applyOnTheRight(il, il + 1, rot);
421
422 for (Index i = il + 1; i < iu; i++) {
423 rot.makeGivens(m_matT.coeffRef(i, i - 1), m_matT.coeffRef(i + 1, i - 1), &m_matT.coeffRef(i, i - 1));
424 m_matT.coeffRef(i + 1, i - 1) = ComplexScalar(0);
425 m_matT.rightCols(m_matT.cols() - i).applyOnTheLeft(i, i + 1, rot.adjoint());
426 m_matT.topRows((std::min)(i + 2, iu) + 1).applyOnTheRight(i, i + 1, rot);
427 if (computeU) m_matU.applyOnTheRight(i, i + 1, rot);
428 }
429 }
430
431 if (totalIter <= maxIters)
432 m_info = Success;
433 else
434 m_info = NoConvergence;
435
436 m_isInitialized = true;
437 m_matUisUptodate = computeU;
438}
439
440} // end namespace Eigen
441
442#endif // EIGEN_COMPLEX_SCHUR_H
Performs a complex Schur decomposition of a real or complex square matrix.
Definition ComplexSchur.h:56
ComputationInfo info() const
Reports whether previous computation was successful.
Definition ComplexSchur.h:220
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
Definition ComplexSchur.h:68
ComplexSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition ComplexSchur.h:100
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition ComplexSchur.h:165
ComplexSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
Compute Schur decomposition from a given Hessenberg matrix.
Index getMaxIterations()
Returns the maximum number of iterations.
Definition ComplexSchur.h:236
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition ComplexSchur.h:243
Eigen::Index Index
Definition ComplexSchur.h:70
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType_.
Definition ComplexSchur.h:78
ComplexSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes Schur decomposition of given matrix.
Definition ComplexSchur.h:118
Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > ComplexMatrixType
Type for the matrices in the Schur decomposition.
Definition ComplexSchur.h:87
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
ComplexSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition ComplexSchur.h:230
const ComplexMatrixType & matrixU() const
Returns the unitary matrix in the Schur decomposition.
Definition ComplexSchur.h:142
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
constexpr Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:217
constexpr const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:198
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
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
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
Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523