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