11#ifndef EIGEN_EIGENSOLVER_H
12#define EIGEN_EIGENSOLVER_H
14#include "./RealSchur.h"
17#include "./InternalHeaderCheck.h"
67template <
typename MatrixType_>
74 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
75 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
76 Options = internal::traits<MatrixType>::Options,
77 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
78 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
82 typedef typename MatrixType::Scalar
Scalar;
83 typedef typename NumTraits<Scalar>::Real RealScalar;
106 typedef Matrix<
ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime,
107 MaxColsAtCompileTime>
118 : m_eivec(), m_eivalues(), m_isInitialized(false), m_eigenvectorsOk(false), m_realSchur(), m_matT(), m_tmp() {}
127 : m_eivec(size, size),
129 m_isInitialized(false),
130 m_eigenvectorsOk(false),
150 template <
typename InputType>
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()) {
203 eigen_assert(m_isInitialized &&
"EigenSolver is not initialized.");
204 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
247 eigen_assert(m_isInitialized &&
"EigenSolver is not initialized.");
278 template <
typename InputType>
284 eigen_assert(m_isInitialized &&
"EigenSolver is not initialized.");
298 void doComputeEigenvectors();
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);
308 bool m_isInitialized;
309 bool m_eigenvectorsOk;
310 ComputationInfo m_info;
311 RealSchur<MatrixType> m_realSchur;
314 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
315 ColumnVectorType m_tmp;
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();
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;
337 matD.coeffRef(i, i) = numext::real(m_eivalues.coeff(i));
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();
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) ||
354 matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
355 matV.col(j).normalize();
358 for (
Index i = 0; i < n; ++i) {
362 matV.col(j).normalize();
363 matV.col(j + 1).normalize();
370template <
typename MatrixType>
371template <
typename InputType>
373 bool computeEigenvectors) {
374 check_template_parameters();
376 using numext::isfinite;
379 eigen_assert(matrix.
cols() == matrix.
rows());
382 m_realSchur.compute(matrix.
derived(), computeEigenvectors);
384 m_info = m_realSchur.info();
387 m_matT = m_realSchur.matrixT();
388 if (computeEigenvectors) m_eivec = m_realSchur.matrixU();
391 m_eivalues.resize(matrix.
cols());
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;
404 Scalar p =
Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i + 1, i + 1));
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)));
415 z = maxval *
sqrt(
abs(p0 * p0 + t0 * t1));
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;
431 if (computeEigenvectors) doComputeEigenvectors();
434 m_isInitialized =
true;
435 m_eigenvectorsOk = computeEigenvectors;
440template <
typename MatrixType>
441void EigenSolver<MatrixType>::doComputeEigenvectors() {
443 const Index size = m_eivec.cols();
444 const Scalar eps = NumTraits<Scalar>::epsilon();
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();
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();
463 Scalar lastr(0), lastw(0);
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));
471 if (m_eivalues.coeff(i).imag() <
Scalar(0)) {
476 if (m_eivalues.coeff(i).imag() ==
Scalar(0)) {
478 m_matT.coeffRef(i, n) = -r / w;
480 m_matT.coeffRef(i, n) = -r / (eps * norm);
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;
490 m_matT.coeffRef(i + 1, n) = (-r - w * t) / x;
492 m_matT.coeffRef(i + 1, n) = (-lastr - y * t) / lastw;
497 if ((eps * t) * t >
Scalar(1)) m_matT.col(n).tail(size - i) /= t;
500 }
else if (q <
Scalar(0) && n > 0)
502 Scalar lastra(0), lastsa(0), lastw(0);
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);
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);
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;
522 if (m_eivalues.coeff(i).imag() <
Scalar(0)) {
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);
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;
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);
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;
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);
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;
566 eigen_assert(0 &&
"Internal bug in EigenSolver (INF or NaN has not been detected)");
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;
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