Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
EigenSolver.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_EIGENSOLVER_H
12#define EIGEN_EIGENSOLVER_H
13
14#include "./RealSchur.h"
15
16// IWYU pragma: private
17#include "./InternalHeaderCheck.h"
18
19namespace Eigen {
20
67template <typename MatrixType_>
69 public:
71 typedef MatrixType_ MatrixType;
72
73 enum {
74 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
75 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
76 Options = internal::traits<MatrixType>::Options,
77 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
78 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
79 };
80
82 typedef typename MatrixType::Scalar Scalar;
83 typedef typename NumTraits<Scalar>::Real RealScalar;
85
92 typedef std::complex<RealScalar> ComplexScalar;
93
100
106 typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime,
107 MaxColsAtCompileTime>
109
118 : m_eivec(), m_eivalues(), m_isInitialized(false), m_eigenvectorsOk(false), m_realSchur(), m_matT(), m_tmp() {}
119
126 explicit EigenSolver(Index size)
127 : m_eivec(size, size),
128 m_eivalues(size),
129 m_isInitialized(false),
130 m_eigenvectorsOk(false),
131 m_realSchur(size),
132 m_matT(size, size),
133 m_tmp(size) {}
134
150 template <typename InputType>
151 explicit EigenSolver(const EigenBase<InputType>& matrix, bool computeEigenvectors = true)
152 : m_eivec(matrix.rows(), matrix.cols()),
153 m_eivalues(matrix.cols()),
154 m_isInitialized(false),
155 m_eigenvectorsOk(false),
156 m_realSchur(matrix.cols()),
157 m_matT(matrix.rows(), matrix.cols()),
158 m_tmp(matrix.cols()) {
159 compute(matrix.derived(), computeEigenvectors);
160 }
161
183
203 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
204 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
205 return m_eivec;
206 }
207
227
247 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
248 return m_eivalues;
249 }
250
278 template <typename InputType>
279 EigenSolver& compute(const EigenBase<InputType>& matrix, bool computeEigenvectors = true);
280
284 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
285 return m_info;
286 }
287
290 m_realSchur.setMaxIterations(maxIters);
291 return *this;
292 }
293
295 Index getMaxIterations() { return m_realSchur.getMaxIterations(); }
296
297 private:
298 void doComputeEigenvectors();
299
300 protected:
301 static void check_template_parameters() {
302 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
303 EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
304 }
305
306 MatrixType m_eivec;
307 EigenvalueType m_eivalues;
308 bool m_isInitialized;
309 bool m_eigenvectorsOk;
310 ComputationInfo m_info;
311 RealSchur<MatrixType> m_realSchur;
312 MatrixType m_matT;
313
314 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
315 ColumnVectorType m_tmp;
316};
317
318template <typename MatrixType>
320 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
321 const RealScalar precision = RealScalar(2) * NumTraits<RealScalar>::epsilon();
322 Index n = m_eivalues.rows();
323 MatrixType matD = MatrixType::Zero(n, n);
324 for (Index i = 0; i < n; ++i) {
325 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i)), precision))
326 matD.coeffRef(i, i) = numext::real(m_eivalues.coeff(i));
327 else {
328 matD.template block<2, 2>(i, i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
329 -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
330 ++i;
331 }
332 }
333 return matD;
334}
335
336template <typename MatrixType>
338 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
339 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
340 const RealScalar precision = RealScalar(2) * NumTraits<RealScalar>::epsilon();
341 Index n = m_eivec.cols();
342 EigenvectorsType matV(n, n);
343 for (Index j = 0; j < n; ++j) {
344 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j)), precision) ||
345 j + 1 == n) {
346 // we have a real eigen value
347 matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
348 matV.col(j).normalize();
349 } else {
350 // we have a pair of complex eigen values
351 for (Index i = 0; i < n; ++i) {
352 matV.coeffRef(i, j) = ComplexScalar(m_eivec.coeff(i, j), m_eivec.coeff(i, j + 1));
353 matV.coeffRef(i, j + 1) = ComplexScalar(m_eivec.coeff(i, j), -m_eivec.coeff(i, j + 1));
354 }
355 matV.col(j).normalize();
356 matV.col(j + 1).normalize();
357 ++j;
358 }
359 }
360 return matV;
361}
362
363template <typename MatrixType>
364template <typename InputType>
366 bool computeEigenvectors) {
367 check_template_parameters();
368
369 using numext::isfinite;
370 using std::abs;
371 using std::sqrt;
372 eigen_assert(matrix.cols() == matrix.rows());
373
374 // Reduce to real Schur form.
375 m_realSchur.compute(matrix.derived(), computeEigenvectors);
376
377 m_info = m_realSchur.info();
378
379 if (m_info == Success) {
380 m_matT = m_realSchur.matrixT();
381 if (computeEigenvectors) m_eivec = m_realSchur.matrixU();
382
383 // Compute eigenvalues from matT
384 m_eivalues.resize(matrix.cols());
385 Index i = 0;
386 while (i < matrix.cols()) {
387 if (i == matrix.cols() - 1 || m_matT.coeff(i + 1, i) == Scalar(0)) {
388 m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
389 if (!(isfinite)(m_eivalues.coeffRef(i))) {
390 m_isInitialized = true;
391 m_eigenvectorsOk = false;
392 m_info = NumericalIssue;
393 return *this;
394 }
395 ++i;
396 } else {
397 Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i + 1, i + 1));
398 Scalar z;
399 // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
400 // without overflow
401 {
402 Scalar t0 = m_matT.coeff(i + 1, i);
403 Scalar t1 = m_matT.coeff(i, i + 1);
404 Scalar maxval = numext::maxi<Scalar>(abs(p), numext::maxi<Scalar>(abs(t0), abs(t1)));
405 t0 /= maxval;
406 t1 /= maxval;
407 Scalar p0 = p / maxval;
408 z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
409 }
410
411 m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i + 1, i + 1) + p, z);
412 m_eivalues.coeffRef(i + 1) = ComplexScalar(m_matT.coeff(i + 1, i + 1) + p, -z);
413 if (!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i + 1)))) {
414 m_isInitialized = true;
415 m_eigenvectorsOk = false;
416 m_info = NumericalIssue;
417 return *this;
418 }
419 i += 2;
420 }
421 }
422
423 // Compute eigenvectors.
424 if (computeEigenvectors) doComputeEigenvectors();
425 }
426
427 m_isInitialized = true;
428 m_eigenvectorsOk = computeEigenvectors;
429
430 return *this;
431}
432
433template <typename MatrixType>
434void EigenSolver<MatrixType>::doComputeEigenvectors() {
435 using std::abs;
436 const Index size = m_eivec.cols();
437 const Scalar eps = NumTraits<Scalar>::epsilon();
438
439 // inefficient! this is already computed in RealSchur
440 Scalar norm(0);
441 for (Index j = 0; j < size; ++j) {
442 norm += m_matT.row(j).segment((std::max)(j - 1, Index(0)), size - (std::max)(j - 1, Index(0))).cwiseAbs().sum();
443 }
444
445 // Backsubstitute to find vectors of upper triangular form
446 if (norm == Scalar(0)) {
447 return;
448 }
449
450 for (Index n = size - 1; n >= 0; n--) {
451 Scalar p = m_eivalues.coeff(n).real();
452 Scalar q = m_eivalues.coeff(n).imag();
453
454 // Scalar vector
455 if (q == Scalar(0)) {
456 Scalar lastr(0), lastw(0);
457 Index l = n;
458
459 m_matT.coeffRef(n, n) = Scalar(1);
460 for (Index i = n - 1; i >= 0; i--) {
461 Scalar w = m_matT.coeff(i, i) - p;
462 Scalar r = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n).segment(l, n - l + 1));
463
464 if (m_eivalues.coeff(i).imag() < Scalar(0)) {
465 lastw = w;
466 lastr = r;
467 } else {
468 l = i;
469 if (m_eivalues.coeff(i).imag() == Scalar(0)) {
470 if (w != Scalar(0))
471 m_matT.coeffRef(i, n) = -r / w;
472 else
473 m_matT.coeffRef(i, n) = -r / (eps * norm);
474 } else // Solve real equations
475 {
476 Scalar x = m_matT.coeff(i, i + 1);
477 Scalar y = m_matT.coeff(i + 1, i);
478 Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) +
479 m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
480 Scalar t = (x * lastr - lastw * r) / denom;
481 m_matT.coeffRef(i, n) = t;
482 if (abs(x) > abs(lastw))
483 m_matT.coeffRef(i + 1, n) = (-r - w * t) / x;
484 else
485 m_matT.coeffRef(i + 1, n) = (-lastr - y * t) / lastw;
486 }
487
488 // Overflow control
489 Scalar t = abs(m_matT.coeff(i, n));
490 if ((eps * t) * t > Scalar(1)) m_matT.col(n).tail(size - i) /= t;
491 }
492 }
493 } else if (q < Scalar(0) && n > 0) // Complex vector
494 {
495 Scalar lastra(0), lastsa(0), lastw(0);
496 Index l = n - 1;
497
498 // Last vector component imaginary so matrix is triangular
499 if (abs(m_matT.coeff(n, n - 1)) > abs(m_matT.coeff(n - 1, n))) {
500 m_matT.coeffRef(n - 1, n - 1) = q / m_matT.coeff(n, n - 1);
501 m_matT.coeffRef(n - 1, n) = -(m_matT.coeff(n, n) - p) / m_matT.coeff(n, n - 1);
502 } else {
503 ComplexScalar cc =
504 ComplexScalar(Scalar(0), -m_matT.coeff(n - 1, n)) / ComplexScalar(m_matT.coeff(n - 1, n - 1) - p, q);
505 m_matT.coeffRef(n - 1, n - 1) = numext::real(cc);
506 m_matT.coeffRef(n - 1, n) = numext::imag(cc);
507 }
508 m_matT.coeffRef(n, n - 1) = Scalar(0);
509 m_matT.coeffRef(n, n) = Scalar(1);
510 for (Index i = n - 2; i >= 0; i--) {
511 Scalar ra = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n - 1).segment(l, n - l + 1));
512 Scalar sa = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n).segment(l, n - l + 1));
513 Scalar w = m_matT.coeff(i, i) - p;
514
515 if (m_eivalues.coeff(i).imag() < Scalar(0)) {
516 lastw = w;
517 lastra = ra;
518 lastsa = sa;
519 } else {
520 l = i;
521 if (m_eivalues.coeff(i).imag() == RealScalar(0)) {
522 ComplexScalar cc = ComplexScalar(-ra, -sa) / ComplexScalar(w, q);
523 m_matT.coeffRef(i, n - 1) = numext::real(cc);
524 m_matT.coeffRef(i, n) = numext::imag(cc);
525 } else {
526 // Solve complex equations
527 Scalar x = m_matT.coeff(i, i + 1);
528 Scalar y = m_matT.coeff(i + 1, i);
529 Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) +
530 m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
531 Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
532 if ((vr == Scalar(0)) && (vi == Scalar(0)))
533 vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
534
535 ComplexScalar cc = ComplexScalar(x * lastra - lastw * ra + q * sa, x * lastsa - lastw * sa - q * ra) /
536 ComplexScalar(vr, vi);
537 m_matT.coeffRef(i, n - 1) = numext::real(cc);
538 m_matT.coeffRef(i, n) = numext::imag(cc);
539 if (abs(x) > (abs(lastw) + abs(q))) {
540 m_matT.coeffRef(i + 1, n - 1) = (-ra - w * m_matT.coeff(i, n - 1) + q * m_matT.coeff(i, n)) / x;
541 m_matT.coeffRef(i + 1, n) = (-sa - w * m_matT.coeff(i, n) - q * m_matT.coeff(i, n - 1)) / x;
542 } else {
543 cc = ComplexScalar(-lastra - y * m_matT.coeff(i, n - 1), -lastsa - y * m_matT.coeff(i, n)) /
544 ComplexScalar(lastw, q);
545 m_matT.coeffRef(i + 1, n - 1) = numext::real(cc);
546 m_matT.coeffRef(i + 1, n) = numext::imag(cc);
547 }
548 }
549
550 // Overflow control
551 Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i, n - 1)), abs(m_matT.coeff(i, n)));
552 if ((eps * t) * t > Scalar(1)) m_matT.block(i, n - 1, size - i, 2) /= t;
553 }
554 }
555
556 // We handled a pair of complex conjugate eigenvalues, so need to skip them both
557 n--;
558 } else {
559 eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
560 }
561 }
562
563 // Back transformation to get eigenvectors of original matrix
564 for (Index j = size - 1; j >= 0; j--) {
565 m_tmp.noalias() = m_eivec.leftCols(j + 1) * m_matT.col(j).segment(0, j + 1);
566 m_eivec.col(j) = m_tmp;
567 }
568}
569
570} // end namespace Eigen
571
572#endif // EIGEN_EIGENSOLVER_H
Computes eigenvalues and eigenvectors of general matrices.
Definition EigenSolver.h:68
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType.
Definition EigenSolver.h:82
EigenSolver(const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
Constructor; computes eigendecomposition of given matrix.
Definition EigenSolver.h:151
EigenSolver(Index size)
Default constructor with memory preallocation.
Definition EigenSolver.h:126
MatrixType pseudoEigenvalueMatrix() const
Returns the block-diagonal matrix in the pseudo-eigendecomposition.
Definition EigenSolver.h:319
ComputationInfo info() const
Definition EigenSolver.h:283
EigenvectorsType eigenvectors() const
Returns the eigenvectors of given matrix.
Definition EigenSolver.h:337
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType.
Definition EigenSolver.h:92
EigenSolver()
Default constructor.
Definition EigenSolver.h:117
EigenSolver & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition EigenSolver.h:289
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Type for vector of eigenvalues as returned by eigenvalues().
Definition EigenSolver.h:99
Eigen::Index Index
Definition EigenSolver.h:84
const MatrixType & pseudoEigenvectors() const
Returns the pseudo-eigenvectors of given matrix.
Definition EigenSolver.h:202
EigenSolver & compute(const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.
const EigenvalueType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition EigenSolver.h:246
Index getMaxIterations()
Returns the maximum number of iterations.
Definition EigenSolver.h:295
Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
Type for matrix of eigenvectors as returned by eigenvectors().
Definition EigenSolver.h:108
MatrixType_ MatrixType
Synonym for the template parameter MatrixType_.
Definition EigenSolver.h:71
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
Index getMaxIterations()
Returns the maximum number of iterations.
Definition RealSchur.h:210
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition RealSchur.h:204
ComputationInfo
Definition Constants.h:438
@ NumericalIssue
Definition Constants.h:442
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
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
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523