Eigen  3.2.91
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 namespace Eigen {
17 
64 template<typename _MatrixType> class EigenSolver
65 {
66  public:
67 
69  typedef _MatrixType MatrixType;
70 
71  enum {
72  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
74  Options = MatrixType::Options,
75  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
76  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
77  };
78 
80  typedef typename MatrixType::Scalar Scalar;
81  typedef typename NumTraits<Scalar>::Real RealScalar;
82  typedef Eigen::Index Index;
83 
90  typedef std::complex<RealScalar> ComplexScalar;
91 
98 
105 
113  EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
114 
121  explicit EigenSolver(Index size)
122  : m_eivec(size, size),
123  m_eivalues(size),
124  m_isInitialized(false),
125  m_eigenvectorsOk(false),
126  m_realSchur(size),
127  m_matT(size, size),
128  m_tmp(size)
129  {}
130 
146  explicit EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true)
147  : m_eivec(matrix.rows(), matrix.cols()),
148  m_eivalues(matrix.cols()),
149  m_isInitialized(false),
150  m_eigenvectorsOk(false),
151  m_realSchur(matrix.cols()),
152  m_matT(matrix.rows(), matrix.cols()),
153  m_tmp(matrix.cols())
154  {
155  compute(matrix, computeEigenvectors);
156  }
157 
178  EigenvectorsType eigenvectors() const;
179 
198  const MatrixType& pseudoEigenvectors() const
199  {
200  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
201  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
202  return m_eivec;
203  }
204 
223  MatrixType pseudoEigenvalueMatrix() const;
224 
243  const EigenvalueType& eigenvalues() const
244  {
245  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
246  return m_eivalues;
247  }
248 
276  EigenSolver& compute(const MatrixType& matrix, bool computeEigenvectors = true);
277 
280  {
281  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
282  return m_info;
283  }
284 
286  EigenSolver& setMaxIterations(Index maxIters)
287  {
288  m_realSchur.setMaxIterations(maxIters);
289  return *this;
290  }
291 
294  {
295  return m_realSchur.getMaxIterations();
296  }
297 
298  private:
299  void doComputeEigenvectors();
300 
301  protected:
302 
303  static void check_template_parameters()
304  {
305  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
306  EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
307  }
308 
309  MatrixType m_eivec;
310  EigenvalueType m_eivalues;
311  bool m_isInitialized;
312  bool m_eigenvectorsOk;
313  ComputationInfo m_info;
314  RealSchur<MatrixType> m_realSchur;
315  MatrixType m_matT;
316 
317  typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
318  ColumnVectorType m_tmp;
319 };
320 
321 template<typename MatrixType>
323 {
324  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
325  Index n = m_eivalues.rows();
326  MatrixType matD = MatrixType::Zero(n,n);
327  for (Index i=0; i<n; ++i)
328  {
329  if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i))))
330  matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
331  else
332  {
333  matD.template block<2,2>(i,i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
334  -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
335  ++i;
336  }
337  }
338  return matD;
339 }
340 
341 template<typename MatrixType>
343 {
344  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
345  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
346  Index n = m_eivec.cols();
347  EigenvectorsType matV(n,n);
348  for (Index j=0; j<n; ++j)
349  {
350  if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j))) || j+1==n)
351  {
352  // we have a real eigen value
353  matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
354  matV.col(j).normalize();
355  }
356  else
357  {
358  // we have a pair of complex eigen values
359  for (Index i=0; i<n; ++i)
360  {
361  matV.coeffRef(i,j) = ComplexScalar(m_eivec.coeff(i,j), m_eivec.coeff(i,j+1));
362  matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
363  }
364  matV.col(j).normalize();
365  matV.col(j+1).normalize();
366  ++j;
367  }
368  }
369  return matV;
370 }
371 
372 template<typename MatrixType>
374 EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvectors)
375 {
376  check_template_parameters();
377 
378  using std::sqrt;
379  using std::abs;
380  using numext::isfinite;
381  eigen_assert(matrix.cols() == matrix.rows());
382 
383  // Reduce to real Schur form.
384  m_realSchur.compute(matrix, computeEigenvectors);
385 
386  m_info = m_realSchur.info();
387 
388  if (m_info == Success)
389  {
390  m_matT = m_realSchur.matrixT();
391  if (computeEigenvectors)
392  m_eivec = m_realSchur.matrixU();
393 
394  // Compute eigenvalues from matT
395  m_eivalues.resize(matrix.cols());
396  Index i = 0;
397  while (i < matrix.cols())
398  {
399  if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0))
400  {
401  m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
402  if(!(isfinite)(m_eivalues.coeffRef(i)))
403  {
404  m_isInitialized = true;
405  m_eigenvectorsOk = false;
406  m_info = NumericalIssue;
407  return *this;
408  }
409  ++i;
410  }
411  else
412  {
413  Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
414  Scalar z;
415  // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
416  // without overflow
417  {
418  Scalar t0 = m_matT.coeff(i+1, i);
419  Scalar t1 = m_matT.coeff(i, i+1);
420  Scalar maxval = numext::maxi<Scalar>(abs(p),numext::maxi<Scalar>(abs(t0),abs(t1)));
421  t0 /= maxval;
422  t1 /= maxval;
423  Scalar p0 = p/maxval;
424  z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
425  }
426 
427  m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
428  m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
429  if(!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i+1))))
430  {
431  m_isInitialized = true;
432  m_eigenvectorsOk = false;
433  m_info = NumericalIssue;
434  return *this;
435  }
436  i += 2;
437  }
438  }
439 
440  // Compute eigenvectors.
441  if (computeEigenvectors)
442  doComputeEigenvectors();
443  }
444 
445  m_isInitialized = true;
446  m_eigenvectorsOk = computeEigenvectors;
447 
448  return *this;
449 }
450 
451 // Complex scalar division.
452 template<typename Scalar>
453 std::complex<Scalar> cdiv(const Scalar& xr, const Scalar& xi, const Scalar& yr, const Scalar& yi)
454 {
455  using std::abs;
456  Scalar r,d;
457  if (abs(yr) > abs(yi))
458  {
459  r = yi/yr;
460  d = yr + r*yi;
461  return std::complex<Scalar>((xr + r*xi)/d, (xi - r*xr)/d);
462  }
463  else
464  {
465  r = yr/yi;
466  d = yi + r*yr;
467  return std::complex<Scalar>((r*xr + xi)/d, (r*xi - xr)/d);
468  }
469 }
470 
471 
472 template<typename MatrixType>
473 void EigenSolver<MatrixType>::doComputeEigenvectors()
474 {
475  using std::abs;
476  const Index size = m_eivec.cols();
477  const Scalar eps = NumTraits<Scalar>::epsilon();
478 
479  // inefficient! this is already computed in RealSchur
480  Scalar norm(0);
481  for (Index j = 0; j < size; ++j)
482  {
483  norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
484  }
485 
486  // Backsubstitute to find vectors of upper triangular form
487  if (norm == Scalar(0))
488  {
489  return;
490  }
491 
492  for (Index n = size-1; n >= 0; n--)
493  {
494  Scalar p = m_eivalues.coeff(n).real();
495  Scalar q = m_eivalues.coeff(n).imag();
496 
497  // Scalar vector
498  if (q == Scalar(0))
499  {
500  Scalar lastr(0), lastw(0);
501  Index l = n;
502 
503  m_matT.coeffRef(n,n) = 1.0;
504  for (Index i = n-1; i >= 0; i--)
505  {
506  Scalar w = m_matT.coeff(i,i) - p;
507  Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
508 
509  if (m_eivalues.coeff(i).imag() < Scalar(0))
510  {
511  lastw = w;
512  lastr = r;
513  }
514  else
515  {
516  l = i;
517  if (m_eivalues.coeff(i).imag() == Scalar(0))
518  {
519  if (w != Scalar(0))
520  m_matT.coeffRef(i,n) = -r / w;
521  else
522  m_matT.coeffRef(i,n) = -r / (eps * norm);
523  }
524  else // Solve real equations
525  {
526  Scalar x = m_matT.coeff(i,i+1);
527  Scalar y = m_matT.coeff(i+1,i);
528  Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
529  Scalar t = (x * lastr - lastw * r) / denom;
530  m_matT.coeffRef(i,n) = t;
531  if (abs(x) > abs(lastw))
532  m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
533  else
534  m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
535  }
536 
537  // Overflow control
538  Scalar t = abs(m_matT.coeff(i,n));
539  if ((eps * t) * t > Scalar(1))
540  m_matT.col(n).tail(size-i) /= t;
541  }
542  }
543  }
544  else if (q < Scalar(0) && n > 0) // Complex vector
545  {
546  Scalar lastra(0), lastsa(0), lastw(0);
547  Index l = n-1;
548 
549  // Last vector component imaginary so matrix is triangular
550  if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
551  {
552  m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
553  m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
554  }
555  else
556  {
557  std::complex<Scalar> cc = cdiv<Scalar>(Scalar(0),-m_matT.coeff(n-1,n),m_matT.coeff(n-1,n-1)-p,q);
558  m_matT.coeffRef(n-1,n-1) = numext::real(cc);
559  m_matT.coeffRef(n-1,n) = numext::imag(cc);
560  }
561  m_matT.coeffRef(n,n-1) = Scalar(0);
562  m_matT.coeffRef(n,n) = Scalar(1);
563  for (Index i = n-2; i >= 0; i--)
564  {
565  Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
566  Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
567  Scalar w = m_matT.coeff(i,i) - p;
568 
569  if (m_eivalues.coeff(i).imag() < Scalar(0))
570  {
571  lastw = w;
572  lastra = ra;
573  lastsa = sa;
574  }
575  else
576  {
577  l = i;
578  if (m_eivalues.coeff(i).imag() == RealScalar(0))
579  {
580  std::complex<Scalar> cc = cdiv(-ra,-sa,w,q);
581  m_matT.coeffRef(i,n-1) = numext::real(cc);
582  m_matT.coeffRef(i,n) = numext::imag(cc);
583  }
584  else
585  {
586  // Solve complex equations
587  Scalar x = m_matT.coeff(i,i+1);
588  Scalar y = m_matT.coeff(i+1,i);
589  Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
590  Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
591  if ((vr == Scalar(0)) && (vi == Scalar(0)))
592  vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
593 
594  std::complex<Scalar> cc = cdiv(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra,vr,vi);
595  m_matT.coeffRef(i,n-1) = numext::real(cc);
596  m_matT.coeffRef(i,n) = numext::imag(cc);
597  if (abs(x) > (abs(lastw) + abs(q)))
598  {
599  m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
600  m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
601  }
602  else
603  {
604  cc = cdiv(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n),lastw,q);
605  m_matT.coeffRef(i+1,n-1) = numext::real(cc);
606  m_matT.coeffRef(i+1,n) = numext::imag(cc);
607  }
608  }
609 
610  // Overflow control
611  Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
612  if ((eps * t) * t > Scalar(1))
613  m_matT.block(i, n-1, size-i, 2) /= t;
614 
615  }
616  }
617 
618  // We handled a pair of complex conjugate eigenvalues, so need to skip them both
619  n--;
620  }
621  else
622  {
623  eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
624  }
625  }
626 
627  // Back transformation to get eigenvectors of original matrix
628  for (Index j = size-1; j >= 0; j--)
629  {
630  m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
631  m_eivec.col(j) = m_tmp;
632  }
633 }
634 
635 } // end namespace Eigen
636 
637 #endif // EIGEN_EIGENSOLVER_H
Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
Type for matrix of eigenvectors as returned by eigenvectors().
Definition: EigenSolver.h:104
EigenvectorsType eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: EigenSolver.h:342
EigenSolver(Index size)
Default constructor with memory preallocation.
Definition: EigenSolver.h:121
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: EigenSolver.h:293
Eigen::Index Index
Definition: EigenSolver.h:82
ComputationInfo info() const
Definition: EigenSolver.h:279
EigenSolver(const MatrixType &matrix, bool computeEigenvectors=true)
Constructor; computes eigendecomposition of given matrix.
Definition: EigenSolver.h:146
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: RealSchur.h:211
MatrixType pseudoEigenvalueMatrix() const
Returns the block-diagonal matrix in the pseudo-eigendecomposition.
Definition: EigenSolver.h:322
const MatrixType & pseudoEigenvectors() const
Returns the pseudo-eigenvectors of given matrix.
Definition: EigenSolver.h:198
Definition: Constants.h:426
EigenSolver & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: EigenSolver.h:286
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType.
Definition: EigenSolver.h:90
Definition: Constants.h:424
_MatrixType MatrixType
Synonym for the template parameter _MatrixType.
Definition: EigenSolver.h:69
EigenSolver()
Default constructor.
Definition: EigenSolver.h:113
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:204
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType.
Definition: EigenSolver.h:80
Computes eigenvalues and eigenvectors of general matrices.
Definition: EigenSolver.h:64
ComputationInfo
Definition: Constants.h:422
const EigenvalueType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: EigenSolver.h:243
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Type for vector of eigenvalues as returned by eigenvalues().
Definition: EigenSolver.h:97
EigenSolver & compute(const MatrixType &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.
Definition: EigenSolver.h:374