EigenSolver.h
Go to the documentation of this file.
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;
82  typedef typename MatrixType::Index Index;
83 
90  typedef std::complex<RealScalar> ComplexScalar;
91 
98 
105 
114 
121  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  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 
279  {
280  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
281  return m_realSchur.info();
282  }
283 
285  EigenSolver& setMaxIterations(Index maxIters)
286  {
287  m_realSchur.setMaxIterations(maxIters);
288  return *this;
289  }
290 
293  {
294  return m_realSchur.getMaxIterations();
295  }
296 
297  private:
298  void doComputeEigenvectors();
299 
300  protected:
301  MatrixType m_eivec;
302  EigenvalueType m_eivalues;
306  MatrixType m_matT;
307 
309  ColumnVectorType m_tmp;
310 };
311 
312 template<typename MatrixType>
314 {
315  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
316  Index n = m_eivalues.rows();
317  MatrixType matD = MatrixType::Zero(n,n);
318  for (Index i=0; i<n; ++i)
319  {
321  matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
322  else
323  {
324  matD.template block<2,2>(i,i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
326  ++i;
327  }
328  }
329  return matD;
330 }
331 
332 template<typename MatrixType>
334 {
335  eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
336  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
337  Index n = m_eivec.cols();
338  EigenvectorsType matV(n,n);
339  for (Index j=0; j<n; ++j)
340  {
342  {
343  // we have a real eigen value
344  matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
345  matV.col(j).normalize();
346  }
347  else
348  {
349  // we have a pair of complex eigen values
350  for (Index i=0; i<n; ++i)
351  {
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 
363 template<typename MatrixType>
365 EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvectors)
366 {
367  using std::sqrt;
368  using std::abs;
369  eigen_assert(matrix.cols() == matrix.rows());
370 
371  // Reduce to real Schur form.
372  m_realSchur.compute(matrix, computeEigenvectors);
373 
374  if (m_realSchur.info() == Success)
375  {
377  if (computeEigenvectors)
379 
380  // Compute eigenvalues from matT
381  m_eivalues.resize(matrix.cols());
382  Index i = 0;
383  while (i < matrix.cols())
384  {
385  if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0))
386  {
387  m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
388  ++i;
389  }
390  else
391  {
392  Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
393  Scalar z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
394  m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
395  m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
396  i += 2;
397  }
398  }
399 
400  // Compute eigenvectors.
401  if (computeEigenvectors)
403  }
404 
405  m_isInitialized = true;
406  m_eigenvectorsOk = computeEigenvectors;
407 
408  return *this;
409 }
410 
411 // Complex scalar division.
412 template<typename Scalar>
413 std::complex<Scalar> cdiv(const Scalar& xr, const Scalar& xi, const Scalar& yr, const Scalar& yi)
414 {
415  using std::abs;
416  Scalar r,d;
417  if (abs(yr) > abs(yi))
418  {
419  r = yi/yr;
420  d = yr + r*yi;
421  return std::complex<Scalar>((xr + r*xi)/d, (xi - r*xr)/d);
422  }
423  else
424  {
425  r = yr/yi;
426  d = yi + r*yr;
427  return std::complex<Scalar>((r*xr + xi)/d, (r*xi - xr)/d);
428  }
429 }
430 
431 
432 template<typename MatrixType>
434 {
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  {
443  norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
444  }
445 
446  // Backsubstitute to find vectors of upper triangular form
447  if (norm == 0.0)
448  {
449  return;
450  }
451 
452  for (Index n = size-1; n >= 0; n--)
453  {
454  Scalar p = m_eivalues.coeff(n).real();
455  Scalar q = m_eivalues.coeff(n).imag();
456 
457  // Scalar vector
458  if (q == Scalar(0))
459  {
460  Scalar lastr(0), lastw(0);
461  Index l = n;
462 
463  m_matT.coeffRef(n,n) = 1.0;
464  for (Index i = n-1; i >= 0; i--)
465  {
466  Scalar w = m_matT.coeff(i,i) - p;
467  Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
468 
469  if (m_eivalues.coeff(i).imag() < 0.0)
470  {
471  lastw = w;
472  lastr = r;
473  }
474  else
475  {
476  l = i;
477  if (m_eivalues.coeff(i).imag() == 0.0)
478  {
479  if (w != 0.0)
480  m_matT.coeffRef(i,n) = -r / w;
481  else
482  m_matT.coeffRef(i,n) = -r / (eps * norm);
483  }
484  else // Solve real equations
485  {
486  Scalar x = m_matT.coeff(i,i+1);
487  Scalar y = m_matT.coeff(i+1,i);
488  Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
489  Scalar t = (x * lastr - lastw * r) / denom;
490  m_matT.coeffRef(i,n) = t;
491  if (abs(x) > abs(lastw))
492  m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
493  else
494  m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
495  }
496 
497  // Overflow control
498  Scalar t = abs(m_matT.coeff(i,n));
499  if ((eps * t) * t > Scalar(1))
500  m_matT.col(n).tail(size-i) /= t;
501  }
502  }
503  }
504  else if (q < Scalar(0) && n > 0) // Complex vector
505  {
506  Scalar lastra(0), lastsa(0), lastw(0);
507  Index l = n-1;
508 
509  // Last vector component imaginary so matrix is triangular
510  if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
511  {
512  m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
513  m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
514  }
515  else
516  {
517  std::complex<Scalar> cc = cdiv<Scalar>(0.0,-m_matT.coeff(n-1,n),m_matT.coeff(n-1,n-1)-p,q);
518  m_matT.coeffRef(n-1,n-1) = numext::real(cc);
519  m_matT.coeffRef(n-1,n) = numext::imag(cc);
520  }
521  m_matT.coeffRef(n,n-1) = 0.0;
522  m_matT.coeffRef(n,n) = 1.0;
523  for (Index i = n-2; i >= 0; i--)
524  {
525  Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
526  Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
527  Scalar w = m_matT.coeff(i,i) - p;
528 
529  if (m_eivalues.coeff(i).imag() < 0.0)
530  {
531  lastw = w;
532  lastra = ra;
533  lastsa = sa;
534  }
535  else
536  {
537  l = i;
538  if (m_eivalues.coeff(i).imag() == RealScalar(0))
539  {
540  std::complex<Scalar> cc = cdiv(-ra,-sa,w,q);
541  m_matT.coeffRef(i,n-1) = numext::real(cc);
542  m_matT.coeffRef(i,n) = numext::imag(cc);
543  }
544  else
545  {
546  // Solve complex equations
547  Scalar x = m_matT.coeff(i,i+1);
548  Scalar y = m_matT.coeff(i+1,i);
549  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;
550  Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
551  if ((vr == 0.0) && (vi == 0.0))
552  vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
553 
554  std::complex<Scalar> cc = cdiv(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra,vr,vi);
555  m_matT.coeffRef(i,n-1) = numext::real(cc);
556  m_matT.coeffRef(i,n) = numext::imag(cc);
557  if (abs(x) > (abs(lastw) + abs(q)))
558  {
559  m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
560  m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
561  }
562  else
563  {
564  cc = cdiv(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n),lastw,q);
565  m_matT.coeffRef(i+1,n-1) = numext::real(cc);
566  m_matT.coeffRef(i+1,n) = numext::imag(cc);
567  }
568  }
569 
570  // Overflow control
571  using std::max;
572  Scalar t = (max)(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
573  if ((eps * t) * t > Scalar(1))
574  m_matT.block(i, n-1, size-i, 2) /= t;
575 
576  }
577  }
578 
579  // We handled a pair of complex conjugate eigenvalues, so need to skip them both
580  n--;
581  }
582  else
583  {
584  eigen_assert(0 && "Internal bug in EigenSolver"); // this should not happen
585  }
586  }
587 
588  // Back transformation to get eigenvectors of original matrix
589  for (Index j = size-1; j >= 0; j--)
590  {
591  m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
592  m_eivec.col(j) = m_tmp;
593  }
594 }
595 
596 } // end namespace Eigen
597 
598 #endif // EIGEN_EIGENSOLVER_H
d
MatrixType m_matT
Definition: EigenSolver.h:306
ColumnVectorType m_tmp
Definition: EigenSolver.h:309
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:333
EigenSolver(Index size)
Default constructor with memory preallocation.
Definition: EigenSolver.h:121
std::complex< Scalar > cdiv(const Scalar &xr, const Scalar &xi, const Scalar &yr, const Scalar &yi)
Definition: EigenSolver.h:413
const MatrixType & matrixU() const
Returns the orthogonal matrix in the Schur decomposition.
Definition: RealSchur.h:126
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: EigenSolver.h:292
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealSchur.h:193
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, typename NumTraits< Scalar >::Real precision=NumTraits< Scalar >::dummy_precision())
ComputationInfo info() const
Definition: EigenSolver.h:278
const MatrixType & matrixT() const
Returns the quasi-triangular matrix in the Schur decomposition.
Definition: RealSchur.h:143
RealSchur< MatrixType > m_realSchur
Definition: EigenSolver.h:305
const ImagReturnType imag() const
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::Index Index
Definition: EigenSolver.h:82
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: EigenSolver.h:308
EIGEN_STRONG_INLINE Index rows() const
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
RealReturnType real() const
MatrixType pseudoEigenvalueMatrix() const
Returns the block-diagonal matrix in the pseudo-eigendecomposition.
Definition: EigenSolver.h:313
NumTraits< Scalar >::Real RealScalar
Definition: EigenSolver.h:81
EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
const MatrixType & pseudoEigenvectors() const
Returns the pseudo-eigenvectors of given matrix.
Definition: EigenSolver.h:198
EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
EIGEN_STRONG_INLINE void resize(Index nbRows, Index nbCols)
TFSIMD_FORCE_INLINE const tfScalar & x() const
EigenSolver & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: EigenSolver.h:285
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType.
Definition: EigenSolver.h:90
RealSchur & compute(const MatrixType &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
Definition: RealSchur.h:246
_MatrixType MatrixType
Synonym for the template parameter _MatrixType.
Definition: EigenSolver.h:69
TFSIMD_FORCE_INLINE const tfScalar & z() const
EigenSolver()
Default constructor.
Definition: EigenSolver.h:113
TFSIMD_FORCE_INLINE const tfScalar & w() const
RealSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: RealSchur.h:204
void doComputeEigenvectors()
Definition: EigenSolver.h:433
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType.
Definition: EigenSolver.h:80
const CwiseUnaryOp< internal::scalar_sqrt_op< Scalar >, const Derived > sqrt() const
Computes eigenvalues and eigenvectors of general matrices.
Definition: EigenSolver.h:64
#define eigen_assert(x)
ComputationInfo
Definition: Constants.h:374
EigenvalueType m_eivalues
Definition: EigenSolver.h:302
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
MatrixType m_eivec
Definition: EigenSolver.h:301
EigenSolver & compute(const MatrixType &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.
Definition: EigenSolver.h:365


tuw_aruco
Author(s): Lukas Pfeifhofer
autogenerated on Mon Jun 10 2019 15:40:48