Eigen  5.0.1-dev+60122df6
 
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 <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_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 internal::make_complex_t<Scalar> 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 const Index n = m_eivalues.rows();
323 MatrixType matD = MatrixType::Zero(n, n);
324 Index i = 0;
325 for (; i < n - 1; ++i) {
326 RealScalar real = numext::real(m_eivalues.coeff(i));
327 RealScalar imag = numext::imag(m_eivalues.coeff(i));
328 matD.coeffRef(i, i) = real;
329 if (!internal::isMuchSmallerThan(imag, real, precision)) {
330 matD.coeffRef(i, i + 1) = imag;
331 matD.coeffRef(i + 1, i) = -imag;
332 matD.coeffRef(i + 1, i + 1) = real;
333 ++i;
334 }
335 }
336 if (i == n - 1) {
337 matD.coeffRef(i, i) = numext::real(m_eivalues.coeff(i));
338 }
339
340 return matD;
341}
342
343template <typename MatrixType>
345 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
346 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
347 const RealScalar precision = RealScalar(2) * NumTraits<RealScalar>::epsilon();
348 Index n = m_eivec.cols();
349 EigenvectorsType matV(n, n);
350 for (Index j = 0; j < n; ++j) {
351 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j)), precision) ||
352 j + 1 == n) {
353 // we have a real eigen value
354 matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
355 matV.col(j).normalize();
356 } else {
357 // we have a pair of complex eigen values
358 for (Index i = 0; i < n; ++i) {
359 matV.coeffRef(i, j) = ComplexScalar(m_eivec.coeff(i, j), m_eivec.coeff(i, j + 1));
360 matV.coeffRef(i, j + 1) = ComplexScalar(m_eivec.coeff(i, j), -m_eivec.coeff(i, j + 1));
361 }
362 matV.col(j).normalize();
363 matV.col(j + 1).normalize();
364 ++j;
365 }
366 }
367 return matV;
368}
369
370template <typename MatrixType>
371template <typename InputType>
373 bool computeEigenvectors) {
374 check_template_parameters();
375
376 using numext::isfinite;
377 using std::abs;
378 using std::sqrt;
379 eigen_assert(matrix.cols() == matrix.rows());
380
381 // Reduce to real Schur form.
382 m_realSchur.compute(matrix.derived(), computeEigenvectors);
383
384 m_info = m_realSchur.info();
385
386 if (m_info == Success) {
387 m_matT = m_realSchur.matrixT();
388 if (computeEigenvectors) m_eivec = m_realSchur.matrixU();
389
390 // Compute eigenvalues from matT
391 m_eivalues.resize(matrix.cols());
392 Index i = 0;
393 while (i < matrix.cols()) {
394 if (i == matrix.cols() - 1 || m_matT.coeff(i + 1, i) == Scalar(0)) {
395 m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
396 if (!(isfinite)(m_eivalues.coeffRef(i))) {
397 m_isInitialized = true;
398 m_eigenvectorsOk = false;
399 m_info = NumericalIssue;
400 return *this;
401 }
402 ++i;
403 } else {
404 Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i + 1, i + 1));
405 Scalar z;
406 // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
407 // without overflow
408 {
409 Scalar t0 = m_matT.coeff(i + 1, i);
410 Scalar t1 = m_matT.coeff(i, i + 1);
411 Scalar maxval = numext::maxi<Scalar>(abs(p), numext::maxi<Scalar>(abs(t0), abs(t1)));
412 t0 /= maxval;
413 t1 /= maxval;
414 Scalar p0 = p / maxval;
415 z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
416 }
417
418 m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i + 1, i + 1) + p, z);
419 m_eivalues.coeffRef(i + 1) = ComplexScalar(m_matT.coeff(i + 1, i + 1) + p, -z);
420 if (!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i + 1)))) {
421 m_isInitialized = true;
422 m_eigenvectorsOk = false;
423 m_info = NumericalIssue;
424 return *this;
425 }
426 i += 2;
427 }
428 }
429
430 // Compute eigenvectors.
431 if (computeEigenvectors) doComputeEigenvectors();
432 }
433
434 m_isInitialized = true;
435 m_eigenvectorsOk = computeEigenvectors;
436
437 return *this;
438}
439
440template <typename MatrixType>
441void EigenSolver<MatrixType>::doComputeEigenvectors() {
442 using std::abs;
443 const Index size = m_eivec.cols();
444 const Scalar eps = NumTraits<Scalar>::epsilon();
445
446 // inefficient! this is already computed in RealSchur
447 Scalar norm(0);
448 for (Index j = 0; j < size; ++j) {
449 norm += m_matT.row(j).segment((std::max)(j - 1, Index(0)), size - (std::max)(j - 1, Index(0))).cwiseAbs().sum();
450 }
451
452 // Backsubstitute to find vectors of upper triangular form
453 if (norm == Scalar(0)) {
454 return;
455 }
456
457 for (Index n = size - 1; n >= 0; n--) {
458 Scalar p = m_eivalues.coeff(n).real();
459 Scalar q = m_eivalues.coeff(n).imag();
460
461 // Scalar vector
462 if (q == Scalar(0)) {
463 Scalar lastr(0), lastw(0);
464 Index l = n;
465
466 m_matT.coeffRef(n, n) = Scalar(1);
467 for (Index i = n - 1; i >= 0; i--) {
468 Scalar w = m_matT.coeff(i, i) - p;
469 Scalar r = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n).segment(l, n - l + 1));
470
471 if (m_eivalues.coeff(i).imag() < Scalar(0)) {
472 lastw = w;
473 lastr = r;
474 } else {
475 l = i;
476 if (m_eivalues.coeff(i).imag() == Scalar(0)) {
477 if (w != Scalar(0))
478 m_matT.coeffRef(i, n) = -r / w;
479 else
480 m_matT.coeffRef(i, n) = -r / (eps * norm);
481 } else // Solve real equations
482 {
483 Scalar x = m_matT.coeff(i, i + 1);
484 Scalar y = m_matT.coeff(i + 1, i);
485 Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) +
486 m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
487 Scalar t = (x * lastr - lastw * r) / denom;
488 m_matT.coeffRef(i, n) = t;
489 if (abs(x) > abs(lastw))
490 m_matT.coeffRef(i + 1, n) = (-r - w * t) / x;
491 else
492 m_matT.coeffRef(i + 1, n) = (-lastr - y * t) / lastw;
493 }
494
495 // Overflow control
496 Scalar t = abs(m_matT.coeff(i, n));
497 if ((eps * t) * t > Scalar(1)) m_matT.col(n).tail(size - i) /= t;
498 }
499 }
500 } else if (q < Scalar(0) && n > 0) // Complex vector
501 {
502 Scalar lastra(0), lastsa(0), lastw(0);
503 Index l = n - 1;
504
505 // Last vector component imaginary so matrix is triangular
506 if (abs(m_matT.coeff(n, n - 1)) > abs(m_matT.coeff(n - 1, n))) {
507 m_matT.coeffRef(n - 1, n - 1) = q / m_matT.coeff(n, n - 1);
508 m_matT.coeffRef(n - 1, n) = -(m_matT.coeff(n, n) - p) / m_matT.coeff(n, n - 1);
509 } else {
510 ComplexScalar cc =
511 ComplexScalar(Scalar(0), -m_matT.coeff(n - 1, n)) / ComplexScalar(m_matT.coeff(n - 1, n - 1) - p, q);
512 m_matT.coeffRef(n - 1, n - 1) = numext::real(cc);
513 m_matT.coeffRef(n - 1, n) = numext::imag(cc);
514 }
515 m_matT.coeffRef(n, n - 1) = Scalar(0);
516 m_matT.coeffRef(n, n) = Scalar(1);
517 for (Index i = n - 2; i >= 0; i--) {
518 Scalar ra = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n - 1).segment(l, n - l + 1));
519 Scalar sa = m_matT.row(i).segment(l, n - l + 1).dot(m_matT.col(n).segment(l, n - l + 1));
520 Scalar w = m_matT.coeff(i, i) - p;
521
522 if (m_eivalues.coeff(i).imag() < Scalar(0)) {
523 lastw = w;
524 lastra = ra;
525 lastsa = sa;
526 } else {
527 l = i;
528 if (m_eivalues.coeff(i).imag() == RealScalar(0)) {
529 ComplexScalar cc = ComplexScalar(-ra, -sa) / ComplexScalar(w, q);
530 m_matT.coeffRef(i, n - 1) = numext::real(cc);
531 m_matT.coeffRef(i, n) = numext::imag(cc);
532 } else {
533 // Solve complex equations
534 Scalar x = m_matT.coeff(i, i + 1);
535 Scalar y = m_matT.coeff(i + 1, i);
536 Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) +
537 m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
538 Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
539 if ((vr == Scalar(0)) && (vi == Scalar(0)))
540 vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
541
542 ComplexScalar cc = ComplexScalar(x * lastra - lastw * ra + q * sa, x * lastsa - lastw * sa - q * ra) /
543 ComplexScalar(vr, vi);
544 m_matT.coeffRef(i, n - 1) = numext::real(cc);
545 m_matT.coeffRef(i, n) = numext::imag(cc);
546 if (abs(x) > (abs(lastw) + abs(q))) {
547 m_matT.coeffRef(i + 1, n - 1) = (-ra - w * m_matT.coeff(i, n - 1) + q * m_matT.coeff(i, n)) / x;
548 m_matT.coeffRef(i + 1, n) = (-sa - w * m_matT.coeff(i, n) - q * m_matT.coeff(i, n - 1)) / x;
549 } else {
550 cc = ComplexScalar(-lastra - y * m_matT.coeff(i, n - 1), -lastsa - y * m_matT.coeff(i, n)) /
551 ComplexScalar(lastw, q);
552 m_matT.coeffRef(i + 1, n - 1) = numext::real(cc);
553 m_matT.coeffRef(i + 1, n) = numext::imag(cc);
554 }
555 }
556
557 // Overflow control
558 Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i, n - 1)), abs(m_matT.coeff(i, n)));
559 if ((eps * t) * t > Scalar(1)) m_matT.block(i, n - 1, size - i, 2) /= t;
560 }
561 }
562
563 // We handled a pair of complex conjugate eigenvalues, so need to skip them both
564 n--;
565 } else {
566 eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
567 }
568 }
569
570 // Back transformation to get eigenvectors of original matrix
571 for (Index j = size - 1; j >= 0; j--) {
572 m_tmp.noalias() = m_eivec.leftCols(j + 1) * m_matT.col(j).segment(0, j + 1);
573 m_eivec.col(j) = m_tmp;
574 }
575}
576
577} // end namespace Eigen
578
579#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:344
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
internal::make_complex_t< Scalar > ComplexScalar
Complex scalar type for MatrixType.
Definition EigenSolver.h:92
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:191
ComputationInfo
Definition Constants.h:438
@ NumericalIssue
Definition Constants.h:442
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
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_real_op< typename Derived::Scalar >, const Derived > real(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:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_imag_op< typename Derived::Scalar >, const Derived > imag(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33
constexpr Index cols() const noexcept
Definition EigenBase.h:61
constexpr Derived & derived()
Definition EigenBase.h:49
constexpr Index rows() const noexcept
Definition EigenBase.h:59