00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #ifndef EIGEN_EIGENSOLVER_H
00012 #define EIGEN_EIGENSOLVER_H
00013
00014 #include "./RealSchur.h"
00015
00016 namespace Eigen {
00017
00064 template<typename _MatrixType> class EigenSolver
00065 {
00066 public:
00067
00069 typedef _MatrixType MatrixType;
00070
00071 enum {
00072 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00073 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00074 Options = MatrixType::Options,
00075 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00076 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00077 };
00078
00080 typedef typename MatrixType::Scalar Scalar;
00081 typedef typename NumTraits<Scalar>::Real RealScalar;
00082 typedef typename MatrixType::Index Index;
00083
00090 typedef std::complex<RealScalar> ComplexScalar;
00091
00097 typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
00098
00104 typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
00105
00113 EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
00114
00121 EigenSolver(Index size)
00122 : m_eivec(size, size),
00123 m_eivalues(size),
00124 m_isInitialized(false),
00125 m_eigenvectorsOk(false),
00126 m_realSchur(size),
00127 m_matT(size, size),
00128 m_tmp(size)
00129 {}
00130
00146 EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true)
00147 : m_eivec(matrix.rows(), matrix.cols()),
00148 m_eivalues(matrix.cols()),
00149 m_isInitialized(false),
00150 m_eigenvectorsOk(false),
00151 m_realSchur(matrix.cols()),
00152 m_matT(matrix.rows(), matrix.cols()),
00153 m_tmp(matrix.cols())
00154 {
00155 compute(matrix, computeEigenvectors);
00156 }
00157
00178 EigenvectorsType eigenvectors() const;
00179
00198 const MatrixType& pseudoEigenvectors() const
00199 {
00200 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00201 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00202 return m_eivec;
00203 }
00204
00223 MatrixType pseudoEigenvalueMatrix() const;
00224
00243 const EigenvalueType& eigenvalues() const
00244 {
00245 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00246 return m_eivalues;
00247 }
00248
00276 EigenSolver& compute(const MatrixType& matrix, bool computeEigenvectors = true);
00277
00278 ComputationInfo info() const
00279 {
00280 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00281 return m_realSchur.info();
00282 }
00283
00285 EigenSolver& setMaxIterations(Index maxIters)
00286 {
00287 m_realSchur.setMaxIterations(maxIters);
00288 return *this;
00289 }
00290
00292 Index getMaxIterations()
00293 {
00294 return m_realSchur.getMaxIterations();
00295 }
00296
00297 private:
00298 void doComputeEigenvectors();
00299
00300 protected:
00301 MatrixType m_eivec;
00302 EigenvalueType m_eivalues;
00303 bool m_isInitialized;
00304 bool m_eigenvectorsOk;
00305 RealSchur<MatrixType> m_realSchur;
00306 MatrixType m_matT;
00307
00308 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
00309 ColumnVectorType m_tmp;
00310 };
00311
00312 template<typename MatrixType>
00313 MatrixType EigenSolver<MatrixType>::pseudoEigenvalueMatrix() const
00314 {
00315 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00316 Index n = m_eivalues.rows();
00317 MatrixType matD = MatrixType::Zero(n,n);
00318 for (Index i=0; i<n; ++i)
00319 {
00320 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i))))
00321 matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
00322 else
00323 {
00324 matD.template block<2,2>(i,i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
00325 -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
00326 ++i;
00327 }
00328 }
00329 return matD;
00330 }
00331
00332 template<typename MatrixType>
00333 typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eigenvectors() const
00334 {
00335 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00336 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00337 Index n = m_eivec.cols();
00338 EigenvectorsType matV(n,n);
00339 for (Index j=0; j<n; ++j)
00340 {
00341 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j))) || j+1==n)
00342 {
00343
00344 matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
00345 matV.col(j).normalize();
00346 }
00347 else
00348 {
00349
00350 for (Index i=0; i<n; ++i)
00351 {
00352 matV.coeffRef(i,j) = ComplexScalar(m_eivec.coeff(i,j), m_eivec.coeff(i,j+1));
00353 matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
00354 }
00355 matV.col(j).normalize();
00356 matV.col(j+1).normalize();
00357 ++j;
00358 }
00359 }
00360 return matV;
00361 }
00362
00363 template<typename MatrixType>
00364 EigenSolver<MatrixType>&
00365 EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvectors)
00366 {
00367 using std::sqrt;
00368 using std::abs;
00369 eigen_assert(matrix.cols() == matrix.rows());
00370
00371
00372 m_realSchur.compute(matrix, computeEigenvectors);
00373
00374 if (m_realSchur.info() == Success)
00375 {
00376 m_matT = m_realSchur.matrixT();
00377 if (computeEigenvectors)
00378 m_eivec = m_realSchur.matrixU();
00379
00380
00381 m_eivalues.resize(matrix.cols());
00382 Index i = 0;
00383 while (i < matrix.cols())
00384 {
00385 if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0))
00386 {
00387 m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
00388 ++i;
00389 }
00390 else
00391 {
00392 Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
00393 Scalar z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
00394 m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
00395 m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
00396 i += 2;
00397 }
00398 }
00399
00400
00401 if (computeEigenvectors)
00402 doComputeEigenvectors();
00403 }
00404
00405 m_isInitialized = true;
00406 m_eigenvectorsOk = computeEigenvectors;
00407
00408 return *this;
00409 }
00410
00411
00412 template<typename Scalar>
00413 std::complex<Scalar> cdiv(const Scalar& xr, const Scalar& xi, const Scalar& yr, const Scalar& yi)
00414 {
00415 using std::abs;
00416 Scalar r,d;
00417 if (abs(yr) > abs(yi))
00418 {
00419 r = yi/yr;
00420 d = yr + r*yi;
00421 return std::complex<Scalar>((xr + r*xi)/d, (xi - r*xr)/d);
00422 }
00423 else
00424 {
00425 r = yr/yi;
00426 d = yi + r*yr;
00427 return std::complex<Scalar>((r*xr + xi)/d, (r*xi - xr)/d);
00428 }
00429 }
00430
00431
00432 template<typename MatrixType>
00433 void EigenSolver<MatrixType>::doComputeEigenvectors()
00434 {
00435 using std::abs;
00436 const Index size = m_eivec.cols();
00437 const Scalar eps = NumTraits<Scalar>::epsilon();
00438
00439
00440 Scalar norm(0);
00441 for (Index j = 0; j < size; ++j)
00442 {
00443 norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
00444 }
00445
00446
00447 if (norm == 0.0)
00448 {
00449 return;
00450 }
00451
00452 for (Index n = size-1; n >= 0; n--)
00453 {
00454 Scalar p = m_eivalues.coeff(n).real();
00455 Scalar q = m_eivalues.coeff(n).imag();
00456
00457
00458 if (q == Scalar(0))
00459 {
00460 Scalar lastr(0), lastw(0);
00461 Index l = n;
00462
00463 m_matT.coeffRef(n,n) = 1.0;
00464 for (Index i = n-1; i >= 0; i--)
00465 {
00466 Scalar w = m_matT.coeff(i,i) - p;
00467 Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
00468
00469 if (m_eivalues.coeff(i).imag() < 0.0)
00470 {
00471 lastw = w;
00472 lastr = r;
00473 }
00474 else
00475 {
00476 l = i;
00477 if (m_eivalues.coeff(i).imag() == 0.0)
00478 {
00479 if (w != 0.0)
00480 m_matT.coeffRef(i,n) = -r / w;
00481 else
00482 m_matT.coeffRef(i,n) = -r / (eps * norm);
00483 }
00484 else
00485 {
00486 Scalar x = m_matT.coeff(i,i+1);
00487 Scalar y = m_matT.coeff(i+1,i);
00488 Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
00489 Scalar t = (x * lastr - lastw * r) / denom;
00490 m_matT.coeffRef(i,n) = t;
00491 if (abs(x) > abs(lastw))
00492 m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
00493 else
00494 m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
00495 }
00496
00497
00498 Scalar t = abs(m_matT.coeff(i,n));
00499 if ((eps * t) * t > Scalar(1))
00500 m_matT.col(n).tail(size-i) /= t;
00501 }
00502 }
00503 }
00504 else if (q < Scalar(0) && n > 0)
00505 {
00506 Scalar lastra(0), lastsa(0), lastw(0);
00507 Index l = n-1;
00508
00509
00510 if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
00511 {
00512 m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
00513 m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
00514 }
00515 else
00516 {
00517 std::complex<Scalar> cc = cdiv<Scalar>(0.0,-m_matT.coeff(n-1,n),m_matT.coeff(n-1,n-1)-p,q);
00518 m_matT.coeffRef(n-1,n-1) = numext::real(cc);
00519 m_matT.coeffRef(n-1,n) = numext::imag(cc);
00520 }
00521 m_matT.coeffRef(n,n-1) = 0.0;
00522 m_matT.coeffRef(n,n) = 1.0;
00523 for (Index i = n-2; i >= 0; i--)
00524 {
00525 Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
00526 Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
00527 Scalar w = m_matT.coeff(i,i) - p;
00528
00529 if (m_eivalues.coeff(i).imag() < 0.0)
00530 {
00531 lastw = w;
00532 lastra = ra;
00533 lastsa = sa;
00534 }
00535 else
00536 {
00537 l = i;
00538 if (m_eivalues.coeff(i).imag() == RealScalar(0))
00539 {
00540 std::complex<Scalar> cc = cdiv(-ra,-sa,w,q);
00541 m_matT.coeffRef(i,n-1) = numext::real(cc);
00542 m_matT.coeffRef(i,n) = numext::imag(cc);
00543 }
00544 else
00545 {
00546
00547 Scalar x = m_matT.coeff(i,i+1);
00548 Scalar y = m_matT.coeff(i+1,i);
00549 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;
00550 Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
00551 if ((vr == 0.0) && (vi == 0.0))
00552 vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
00553
00554 std::complex<Scalar> cc = cdiv(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra,vr,vi);
00555 m_matT.coeffRef(i,n-1) = numext::real(cc);
00556 m_matT.coeffRef(i,n) = numext::imag(cc);
00557 if (abs(x) > (abs(lastw) + abs(q)))
00558 {
00559 m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
00560 m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
00561 }
00562 else
00563 {
00564 cc = cdiv(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n),lastw,q);
00565 m_matT.coeffRef(i+1,n-1) = numext::real(cc);
00566 m_matT.coeffRef(i+1,n) = numext::imag(cc);
00567 }
00568 }
00569
00570
00571 using std::max;
00572 Scalar t = (max)(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
00573 if ((eps * t) * t > Scalar(1))
00574 m_matT.block(i, n-1, size-i, 2) /= t;
00575
00576 }
00577 }
00578
00579
00580 n--;
00581 }
00582 else
00583 {
00584 eigen_assert(0 && "Internal bug in EigenSolver");
00585 }
00586 }
00587
00588
00589 for (Index j = size-1; j >= 0; j--)
00590 {
00591 m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
00592 m_eivec.col(j) = m_tmp;
00593 }
00594 }
00595
00596 }
00597
00598 #endif // EIGEN_EIGENSOLVER_H