Eigen  3.4.90 (git rev a4098ac676528a83cfb73d4d26ce1b42ec05f47c)
SelfAdjointEigenSolver.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2010 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_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
13
14#include "./Tridiagonalization.h"
15
16#include "./InternalHeaderCheck.h"
17
18namespace Eigen {
19
20template<typename MatrixType_>
21class GeneralizedSelfAdjointEigenSolver;
22
23namespace internal {
24template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
25
26template<typename MatrixType, typename DiagType, typename SubDiagType>
27EIGEN_DEVICE_FUNC
28ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
29}
30
78template<typename MatrixType_> class SelfAdjointEigenSolver
79{
80 public:
81
82 typedef MatrixType_ MatrixType;
83 enum {
84 Size = MatrixType::RowsAtCompileTime,
85 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
86 Options = MatrixType::Options,
87 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
88 };
89
91 typedef typename MatrixType::Scalar Scalar;
93
95
103
104 friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
105
111 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
112 typedef Tridiagonalization<MatrixType> TridiagonalizationType;
114
125 EIGEN_DEVICE_FUNC
127 : m_eivec(),
128 m_eivalues(),
129 m_subdiag(),
130 m_hcoeffs(),
131 m_info(InvalidInput),
132 m_isInitialized(false),
133 m_eigenvectorsOk(false)
134 { }
135
148 EIGEN_DEVICE_FUNC
150 : m_eivec(size, size),
151 m_eivalues(size),
152 m_subdiag(size > 1 ? size - 1 : 1),
153 m_hcoeffs(size > 1 ? size - 1 : 1),
154 m_isInitialized(false),
155 m_eigenvectorsOk(false)
156 {}
157
173 template<typename InputType>
174 EIGEN_DEVICE_FUNC
176 : m_eivec(matrix.rows(), matrix.cols()),
177 m_eivalues(matrix.cols()),
178 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
179 m_hcoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1),
180 m_isInitialized(false),
181 m_eigenvectorsOk(false)
182 {
183 compute(matrix.derived(), options);
184 }
185
216 template<typename InputType>
217 EIGEN_DEVICE_FUNC
219
238 EIGEN_DEVICE_FUNC
239 SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
240
254
278 EIGEN_DEVICE_FUNC
280 {
281 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
282 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
283 return m_eivec;
284 }
285
301 EIGEN_DEVICE_FUNC
303 {
304 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
305 return m_eivalues;
306 }
307
325 EIGEN_DEVICE_FUNC
326 MatrixType operatorSqrt() const
327 {
328 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
329 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
330 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
331 }
332
350 EIGEN_DEVICE_FUNC
351 MatrixType operatorInverseSqrt() const
352 {
353 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
354 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
355 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
356 }
357
362 EIGEN_DEVICE_FUNC
364 {
365 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
366 return m_info;
367 }
368
374 static const int m_maxIterations = 30;
375
376 protected:
377 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
378
379 EigenvectorsType m_eivec;
380 RealVectorType m_eivalues;
383 ComputationInfo m_info;
384 bool m_isInitialized;
385 bool m_eigenvectorsOk;
386};
387
388namespace internal {
409template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
410EIGEN_DEVICE_FUNC
411static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
412}
413
414template<typename MatrixType>
415template<typename InputType>
416EIGEN_DEVICE_FUNC
417SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
418::compute(const EigenBase<InputType>& a_matrix, int options)
419{
420 const InputType &matrix(a_matrix.derived());
421
422 EIGEN_USING_STD(abs);
423 eigen_assert(matrix.cols() == matrix.rows());
424 eigen_assert((options&~(EigVecMask|GenEigMask))==0
425 && (options&EigVecMask)!=EigVecMask
426 && "invalid option parameter");
427 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
428 Index n = matrix.cols();
429 m_eivalues.resize(n,1);
430
431 if(n==1)
432 {
433 m_eivec = matrix;
434 m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
435 if(computeEigenvectors)
436 m_eivec.setOnes(n,n);
437 m_info = Success;
438 m_isInitialized = true;
439 m_eigenvectorsOk = computeEigenvectors;
440 return *this;
441 }
442
443 // declare some aliases
444 RealVectorType& diag = m_eivalues;
445 EigenvectorsType& mat = m_eivec;
446
447 // map the matrix coefficients to [-1:1] to avoid over- and underflow.
448 mat = matrix.template triangularView<Lower>();
449 RealScalar scale = mat.cwiseAbs().maxCoeff();
450 if(scale==RealScalar(0)) scale = RealScalar(1);
451 mat.template triangularView<Lower>() /= scale;
452 m_subdiag.resize(n-1);
453 m_hcoeffs.resize(n-1);
454 internal::tridiagonalization_inplace(mat, diag, m_subdiag, m_hcoeffs, computeEigenvectors);
455
456 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
457
458 // scale back the eigen values
459 m_eivalues *= scale;
460
461 m_isInitialized = true;
462 m_eigenvectorsOk = computeEigenvectors;
463 return *this;
464}
465
466template<typename MatrixType>
467SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
468::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
469{
470 //TODO : Add an option to scale the values beforehand
471 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
472
473 m_eivalues = diag;
474 m_subdiag = subdiag;
475 if (computeEigenvectors)
476 {
477 m_eivec.setIdentity(diag.size(), diag.size());
478 }
479 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
480
481 m_isInitialized = true;
482 m_eigenvectorsOk = computeEigenvectors;
483 return *this;
484}
485
486namespace internal {
498template<typename MatrixType, typename DiagType, typename SubDiagType>
499EIGEN_DEVICE_FUNC
500ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
501{
502 ComputationInfo info;
503 typedef typename MatrixType::Scalar Scalar;
504
505 Index n = diag.size();
506 Index end = n-1;
507 Index start = 0;
508 Index iter = 0; // total number of iterations
509
510 typedef typename DiagType::RealScalar RealScalar;
511 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
512 const RealScalar precision_inv = RealScalar(1)/NumTraits<RealScalar>::epsilon();
513 while (end>0)
514 {
515 for (Index i = start; i<end; ++i) {
516 if (numext::abs(subdiag[i]) < considerAsZero) {
517 subdiag[i] = RealScalar(0);
518 } else {
519 // abs(subdiag[i]) <= epsilon * sqrt(abs(diag[i]) + abs(diag[i+1]))
520 // Scaled to prevent underflows.
521 const RealScalar scaled_subdiag = precision_inv * subdiag[i];
522 if (scaled_subdiag * scaled_subdiag <= (numext::abs(diag[i])+numext::abs(diag[i+1]))) {
523 subdiag[i] = RealScalar(0);
524 }
525 }
526 }
527
528 // find the largest unreduced block at the end of the matrix.
529 while (end>0 && subdiag[end-1]==RealScalar(0))
530 {
531 end--;
532 }
533 if (end<=0)
534 break;
535
536 // if we spent too many iterations, we give up
537 iter++;
538 if(iter > maxIterations * n) break;
539
540 start = end - 1;
541 while (start>0 && subdiag[start-1]!=0)
542 start--;
543
544 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
545 }
546 if (iter <= maxIterations * n)
547 info = Success;
548 else
549 info = NoConvergence;
550
551 // Sort eigenvalues and corresponding vectors.
552 // TODO make the sort optional ?
553 // TODO use a better sort algorithm !!
554 if (info == Success)
555 {
556 for (Index i = 0; i < n-1; ++i)
557 {
558 Index k;
559 diag.segment(i,n-i).minCoeff(&k);
560 if (k > 0)
561 {
562 numext::swap(diag[i], diag[k+i]);
563 if(computeEigenvectors)
564 eivec.col(i).swap(eivec.col(k+i));
565 }
566 }
567 }
568 return info;
569}
570
571template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
572{
573 EIGEN_DEVICE_FUNC
574 static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
575 { eig.compute(A,options); }
576};
577
578template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
579{
580 typedef typename SolverType::MatrixType MatrixType;
581 typedef typename SolverType::RealVectorType VectorType;
582 typedef typename SolverType::Scalar Scalar;
583 typedef typename SolverType::EigenvectorsType EigenvectorsType;
584
585
590 EIGEN_DEVICE_FUNC
591 static inline void computeRoots(const MatrixType& m, VectorType& roots)
592 {
593 EIGEN_USING_STD(sqrt)
594 EIGEN_USING_STD(atan2)
595 EIGEN_USING_STD(cos)
596 EIGEN_USING_STD(sin)
597 const Scalar s_inv3 = Scalar(1)/Scalar(3);
598 const Scalar s_sqrt3 = sqrt(Scalar(3));
599
600 // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
601 // eigenvalues are the roots to this equation, all guaranteed to be
602 // real-valued, because the matrix is symmetric.
603 Scalar c0 = m(0,0)*m(1,1)*m(2,2) + Scalar(2)*m(1,0)*m(2,0)*m(2,1) - m(0,0)*m(2,1)*m(2,1) - m(1,1)*m(2,0)*m(2,0) - m(2,2)*m(1,0)*m(1,0);
604 Scalar c1 = m(0,0)*m(1,1) - m(1,0)*m(1,0) + m(0,0)*m(2,2) - m(2,0)*m(2,0) + m(1,1)*m(2,2) - m(2,1)*m(2,1);
605 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
606
607 // Construct the parameters used in classifying the roots of the equation
608 // and in solving the equation for the roots in closed form.
609 Scalar c2_over_3 = c2*s_inv3;
610 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
611 a_over_3 = numext::maxi(a_over_3, Scalar(0));
612
613 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
614
615 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
616 q = numext::maxi(q, Scalar(0));
617
618 // Compute the eigenvalues by solving for the roots of the polynomial.
619 Scalar rho = sqrt(a_over_3);
620 Scalar theta = atan2(sqrt(q),half_b)*s_inv3; // since sqrt(q) > 0, atan2 is in [0, pi] and theta is in [0, pi/3]
621 Scalar cos_theta = cos(theta);
622 Scalar sin_theta = sin(theta);
623 // roots are already sorted, since cos is monotonically decreasing on [0, pi]
624 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
625 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
626 roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
627 }
628
629 EIGEN_DEVICE_FUNC
630 static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
631 {
632 EIGEN_USING_STD(abs);
633 EIGEN_USING_STD(sqrt);
634 Index i0;
635 // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
636 mat.diagonal().cwiseAbs().maxCoeff(&i0);
637 // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
638 // so let's save it:
639 representative = mat.col(i0);
640 Scalar n0, n1;
641 VectorType c0, c1;
642 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
643 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
644 if(n0>n1) res = c0/sqrt(n0);
645 else res = c1/sqrt(n1);
646
647 return true;
648 }
649
650 EIGEN_DEVICE_FUNC
651 static inline void run(SolverType& solver, const MatrixType& mat, int options)
652 {
653 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
654 eigen_assert((options&~(EigVecMask|GenEigMask))==0
655 && (options&EigVecMask)!=EigVecMask
656 && "invalid option parameter");
657 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
658
659 EigenvectorsType& eivecs = solver.m_eivec;
660 VectorType& eivals = solver.m_eivalues;
661
662 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
663 Scalar shift = mat.trace() / Scalar(3);
664 // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
665 MatrixType scaledMat = mat.template selfadjointView<Lower>();
666 scaledMat.diagonal().array() -= shift;
667 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
668 if(scale > 0) scaledMat /= scale; // TODO for scale==0 we could save the remaining operations
669
670 // compute the eigenvalues
671 computeRoots(scaledMat,eivals);
672
673 // compute the eigenvectors
674 if(computeEigenvectors)
675 {
676 if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
677 {
678 // All three eigenvalues are numerically the same
679 eivecs.setIdentity();
680 }
681 else
682 {
683 MatrixType tmp;
684 tmp = scaledMat;
685
686 // Compute the eigenvector of the most distinct eigenvalue
687 Scalar d0 = eivals(2) - eivals(1);
688 Scalar d1 = eivals(1) - eivals(0);
689 Index k(0), l(2);
690 if(d0 > d1)
691 {
692 numext::swap(k,l);
693 d0 = d1;
694 }
695
696 // Compute the eigenvector of index k
697 {
698 tmp.diagonal().array () -= eivals(k);
699 // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
700 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
701 }
702
703 // Compute eigenvector of index l
705 {
706 // If d0 is too small, then the two other eigenvalues are numerically the same,
707 // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
708 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
709 eivecs.col(l).normalize();
710 }
711 else
712 {
713 tmp = scaledMat;
714 tmp.diagonal().array () -= eivals(l);
715
716 VectorType dummy;
717 extract_kernel(tmp, eivecs.col(l), dummy);
718 }
719
720 // Compute last eigenvector from the other two
721 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
722 }
723 }
724
725 // Rescale back to the original size.
726 eivals *= scale;
727 eivals.array() += shift;
728
729 solver.m_info = Success;
730 solver.m_isInitialized = true;
731 solver.m_eigenvectorsOk = computeEigenvectors;
732 }
733};
734
735// 2x2 direct eigenvalues decomposition, code from Hauke Heibel
736template<typename SolverType>
737struct direct_selfadjoint_eigenvalues<SolverType,2,false>
738{
739 typedef typename SolverType::MatrixType MatrixType;
740 typedef typename SolverType::RealVectorType VectorType;
741 typedef typename SolverType::Scalar Scalar;
742 typedef typename SolverType::EigenvectorsType EigenvectorsType;
743
744 EIGEN_DEVICE_FUNC
745 static inline void computeRoots(const MatrixType& m, VectorType& roots)
746 {
747 EIGEN_USING_STD(sqrt);
748 const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
749 const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
750 roots(0) = t1 - t0;
751 roots(1) = t1 + t0;
752 }
753
754 EIGEN_DEVICE_FUNC
755 static inline void run(SolverType& solver, const MatrixType& mat, int options)
756 {
757 EIGEN_USING_STD(sqrt);
758 EIGEN_USING_STD(abs);
759
760 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
761 eigen_assert((options&~(EigVecMask|GenEigMask))==0
762 && (options&EigVecMask)!=EigVecMask
763 && "invalid option parameter");
764 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
765
766 EigenvectorsType& eivecs = solver.m_eivec;
767 VectorType& eivals = solver.m_eivalues;
768
769 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
770 Scalar shift = mat.trace() / Scalar(2);
771 MatrixType scaledMat = mat;
772 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
773 scaledMat.diagonal().array() -= shift;
774 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
775 if(scale > Scalar(0))
776 scaledMat /= scale;
777
778 // Compute the eigenvalues
779 computeRoots(scaledMat,eivals);
780
781 // compute the eigen vectors
782 if(computeEigenvectors)
783 {
784 if((eivals(1)-eivals(0))<=abs(eivals(1))*Eigen::NumTraits<Scalar>::epsilon())
785 {
786 eivecs.setIdentity();
787 }
788 else
789 {
790 scaledMat.diagonal().array () -= eivals(1);
791 Scalar a2 = numext::abs2(scaledMat(0,0));
792 Scalar c2 = numext::abs2(scaledMat(1,1));
793 Scalar b2 = numext::abs2(scaledMat(1,0));
794 if(a2>c2)
795 {
796 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
797 eivecs.col(1) /= sqrt(a2+b2);
798 }
799 else
800 {
801 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
802 eivecs.col(1) /= sqrt(c2+b2);
803 }
804
805 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
806 }
807 }
808
809 // Rescale back to the original size.
810 eivals *= scale;
811 eivals.array() += shift;
812
813 solver.m_info = Success;
814 solver.m_isInitialized = true;
815 solver.m_eigenvectorsOk = computeEigenvectors;
816 }
817};
818
819}
820
821template<typename MatrixType>
822EIGEN_DEVICE_FUNC
823SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
824::computeDirect(const MatrixType& matrix, int options)
825{
826 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
827 return *this;
828}
829
830namespace internal {
831
832// Francis implicit QR step.
833template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
834EIGEN_DEVICE_FUNC
835static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
836{
837 // Wilkinson Shift.
838 RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
839 RealScalar e = subdiag[end-1];
840 // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
841 // underflow thus leading to inf/NaN values when using the following commented code:
842 // RealScalar e2 = numext::abs2(subdiag[end-1]);
843 // RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
844 // This explain the following, somewhat more complicated, version:
845 RealScalar mu = diag[end];
846 if(td==RealScalar(0)) {
847 mu -= numext::abs(e);
848 } else if (e != RealScalar(0)) {
849 const RealScalar e2 = numext::abs2(e);
850 const RealScalar h = numext::hypot(td,e);
851 if(e2 == RealScalar(0)) {
852 mu -= e / ((td + (td>RealScalar(0) ? h : -h)) / e);
853 } else {
854 mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
855 }
856 }
857
858 RealScalar x = diag[start] - mu;
859 RealScalar z = subdiag[start];
860 // If z ever becomes zero, the Givens rotation will be the identity and
861 // z will stay zero for all future iterations.
862 for (Index k = start; k < end && z != RealScalar(0); ++k)
863 {
864 JacobiRotation<RealScalar> rot;
865 rot.makeGivens(x, z);
866
867 // do T = G' T G
868 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
869 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
870
871 diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
872 diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
873 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
874
875 if (k > start)
876 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
877
878 // "Chasing the bulge" to return to triangular form.
879 x = subdiag[k];
880 if (k < end - 1)
881 {
882 z = -rot.s() * subdiag[k+1];
883 subdiag[k + 1] = rot.c() * subdiag[k+1];
884 }
885
886 // apply the givens rotation to the unit matrix Q = Q * G
887 if (matrixQ)
888 {
889 // FIXME if StorageOrder == RowMajor this operation is not very efficient
890 Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
891 q.applyOnTheRight(k,k+1,rot);
892 }
893 }
894}
895
896} // end namespace internal
897
898} // end namespace Eigen
899
900#endif // EIGEN_SELFADJOINTEIGENSOLVER_H
Derived & setIdentity()
Definition: CwiseNullaryOp.h:875
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: SelfAdjointEigenSolver.h:79
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition: SelfAdjointEigenSolver.h:468
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType_.
Definition: SelfAdjointEigenSolver.h:102
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition: SelfAdjointEigenSolver.h:111
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition: SelfAdjointEigenSolver.h:351
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SelfAdjointEigenSolver.h:363
Eigen::Index Index
Definition: SelfAdjointEigenSolver.h:92
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
Definition: SelfAdjointEigenSolver.h:91
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition: SelfAdjointEigenSolver.h:326
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: SelfAdjointEigenSolver.h:149
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: SelfAdjointEigenSolver.h:302
static const int m_maxIterations
Maximum number of iterations.
Definition: SelfAdjointEigenSolver.h:374
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: SelfAdjointEigenSolver.h:279
SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Definition: SelfAdjointEigenSolver.h:175
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Definition: SelfAdjointEigenSolver.h:824
Tridiagonal decomposition of a selfadjoint matrix.
Definition: Tridiagonalization.h:67
static const lastp1_t end
Definition: IndexedViewHelper.h:183
ComputationInfo
Definition: Constants.h:442
@ Success
Definition: Constants.h:444
@ NoConvergence
Definition: Constants.h:448
@ ComputeEigenvectors
Definition: Constants.h:407
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x)
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_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:59
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Definition: EigenBase.h:32
Derived & derived()
Definition: EigenBase.h:48
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:235