Eigen  3.3.0
 
Loading...
Searching...
No Matches
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
16namespace Eigen {
17
18template<typename _MatrixType>
19class GeneralizedSelfAdjointEigenSolver;
20
21namespace internal {
22template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
23template<typename MatrixType, typename DiagType, typename SubDiagType>
24ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
25}
26
70template<typename _MatrixType> class SelfAdjointEigenSolver
71{
72 public:
73
74 typedef _MatrixType MatrixType;
75 enum {
76 Size = MatrixType::RowsAtCompileTime,
77 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
78 Options = MatrixType::Options,
79 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
80 };
81
83 typedef typename MatrixType::Scalar Scalar;
85
87
94 typedef typename NumTraits<Scalar>::Real RealScalar;
95
96 friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
97
103 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
104 typedef Tridiagonalization<MatrixType> TridiagonalizationType;
106
117 EIGEN_DEVICE_FUNC
119 : m_eivec(),
120 m_eivalues(),
121 m_subdiag(),
122 m_isInitialized(false)
123 { }
124
137 EIGEN_DEVICE_FUNC
139 : m_eivec(size, size),
140 m_eivalues(size),
141 m_subdiag(size > 1 ? size - 1 : 1),
142 m_isInitialized(false)
143 {}
144
160 template<typename InputType>
161 EIGEN_DEVICE_FUNC
163 : m_eivec(matrix.rows(), matrix.cols()),
164 m_eivalues(matrix.cols()),
165 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
166 m_isInitialized(false)
167 {
168 compute(matrix.derived(), options);
169 }
170
201 template<typename InputType>
202 EIGEN_DEVICE_FUNC
204
223 EIGEN_DEVICE_FUNC
224 SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
225
239
258 EIGEN_DEVICE_FUNC
260 {
261 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
262 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
263 return m_eivec;
264 }
265
281 EIGEN_DEVICE_FUNC
283 {
284 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
285 return m_eivalues;
286 }
287
305 EIGEN_DEVICE_FUNC
306 MatrixType operatorSqrt() const
307 {
308 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
309 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
310 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
311 }
312
330 EIGEN_DEVICE_FUNC
331 MatrixType operatorInverseSqrt() const
332 {
333 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
334 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
335 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
336 }
337
342 EIGEN_DEVICE_FUNC
344 {
345 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
346 return m_info;
347 }
348
354 static const int m_maxIterations = 30;
355
356 protected:
357 static void check_template_parameters()
358 {
359 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
360 }
361
362 EigenvectorsType m_eivec;
363 RealVectorType m_eivalues;
364 typename TridiagonalizationType::SubDiagonalType m_subdiag;
365 ComputationInfo m_info;
366 bool m_isInitialized;
367 bool m_eigenvectorsOk;
368};
369
370namespace internal {
391template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
392EIGEN_DEVICE_FUNC
393static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
394}
395
396template<typename MatrixType>
397template<typename InputType>
398EIGEN_DEVICE_FUNC
399SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
400::compute(const EigenBase<InputType>& a_matrix, int options)
401{
402 check_template_parameters();
403
404 const InputType &matrix(a_matrix.derived());
405
406 using std::abs;
407 eigen_assert(matrix.cols() == matrix.rows());
408 eigen_assert((options&~(EigVecMask|GenEigMask))==0
409 && (options&EigVecMask)!=EigVecMask
410 && "invalid option parameter");
411 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
412 Index n = matrix.cols();
413 m_eivalues.resize(n,1);
414
415 if(n==1)
416 {
417 m_eivalues.coeffRef(0,0) = numext::real(matrix.diagonal()[0]);
418 if(computeEigenvectors)
419 m_eivec.setOnes(n,n);
420 m_info = Success;
421 m_isInitialized = true;
422 m_eigenvectorsOk = computeEigenvectors;
423 return *this;
424 }
425
426 // declare some aliases
427 RealVectorType& diag = m_eivalues;
428 EigenvectorsType& mat = m_eivec;
429
430 // map the matrix coefficients to [-1:1] to avoid over- and underflow.
431 mat = matrix.template triangularView<Lower>();
432 RealScalar scale = mat.cwiseAbs().maxCoeff();
433 if(scale==RealScalar(0)) scale = RealScalar(1);
434 mat.template triangularView<Lower>() /= scale;
435 m_subdiag.resize(n-1);
436 internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
437
438 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
439
440 // scale back the eigen values
441 m_eivalues *= scale;
442
443 m_isInitialized = true;
444 m_eigenvectorsOk = computeEigenvectors;
445 return *this;
446}
447
448template<typename MatrixType>
449SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
450::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
451{
452 //TODO : Add an option to scale the values beforehand
453 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
454
455 m_eivalues = diag;
456 m_subdiag = subdiag;
457 if (computeEigenvectors)
458 {
459 m_eivec.setIdentity(diag.size(), diag.size());
460 }
461 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
462
463 m_isInitialized = true;
464 m_eigenvectorsOk = computeEigenvectors;
465 return *this;
466}
467
468namespace internal {
480template<typename MatrixType, typename DiagType, typename SubDiagType>
481ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
482{
483 using std::abs;
484
485 ComputationInfo info;
486 typedef typename MatrixType::Scalar Scalar;
487
488 Index n = diag.size();
489 Index end = n-1;
490 Index start = 0;
491 Index iter = 0; // total number of iterations
492
493 typedef typename DiagType::RealScalar RealScalar;
494 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
495 const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
496
497 while (end>0)
498 {
499 for (Index i = start; i<end; ++i)
500 if (internal::isMuchSmallerThan(abs(subdiag[i]),(abs(diag[i])+abs(diag[i+1])),precision) || abs(subdiag[i]) <= considerAsZero)
501 subdiag[i] = 0;
502
503 // find the largest unreduced block
504 while (end>0 && subdiag[end-1]==RealScalar(0))
505 {
506 end--;
507 }
508 if (end<=0)
509 break;
510
511 // if we spent too many iterations, we give up
512 iter++;
513 if(iter > maxIterations * n) break;
514
515 start = end - 1;
516 while (start>0 && subdiag[start-1]!=0)
517 start--;
518
519 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
520 }
521 if (iter <= maxIterations * n)
522 info = Success;
523 else
524 info = NoConvergence;
525
526 // Sort eigenvalues and corresponding vectors.
527 // TODO make the sort optional ?
528 // TODO use a better sort algorithm !!
529 if (info == Success)
530 {
531 for (Index i = 0; i < n-1; ++i)
532 {
533 Index k;
534 diag.segment(i,n-i).minCoeff(&k);
535 if (k > 0)
536 {
537 std::swap(diag[i], diag[k+i]);
538 if(computeEigenvectors)
539 eivec.col(i).swap(eivec.col(k+i));
540 }
541 }
542 }
543 return info;
544}
545
546template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
547{
548 EIGEN_DEVICE_FUNC
549 static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
550 { eig.compute(A,options); }
551};
552
553template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
554{
555 typedef typename SolverType::MatrixType MatrixType;
556 typedef typename SolverType::RealVectorType VectorType;
557 typedef typename SolverType::Scalar Scalar;
558 typedef typename SolverType::EigenvectorsType EigenvectorsType;
559
560
565 EIGEN_DEVICE_FUNC
566 static inline void computeRoots(const MatrixType& m, VectorType& roots)
567 {
568 EIGEN_USING_STD_MATH(sqrt)
569 EIGEN_USING_STD_MATH(atan2)
570 EIGEN_USING_STD_MATH(cos)
571 EIGEN_USING_STD_MATH(sin)
572 const Scalar s_inv3 = Scalar(1)/Scalar(3);
573 const Scalar s_sqrt3 = sqrt(Scalar(3));
574
575 // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
576 // eigenvalues are the roots to this equation, all guaranteed to be
577 // real-valued, because the matrix is symmetric.
578 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);
579 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);
580 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
581
582 // Construct the parameters used in classifying the roots of the equation
583 // and in solving the equation for the roots in closed form.
584 Scalar c2_over_3 = c2*s_inv3;
585 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
586 a_over_3 = numext::maxi(a_over_3, Scalar(0));
587
588 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
589
590 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
591 q = numext::maxi(q, Scalar(0));
592
593 // Compute the eigenvalues by solving for the roots of the polynomial.
594 Scalar rho = sqrt(a_over_3);
595 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]
596 Scalar cos_theta = cos(theta);
597 Scalar sin_theta = sin(theta);
598 // roots are already sorted, since cos is monotonically decreasing on [0, pi]
599 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
600 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
601 roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
602 }
603
604 EIGEN_DEVICE_FUNC
605 static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
606 {
607 using std::abs;
608 Index i0;
609 // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
610 mat.diagonal().cwiseAbs().maxCoeff(&i0);
611 // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
612 // so let's save it:
613 representative = mat.col(i0);
614 Scalar n0, n1;
615 VectorType c0, c1;
616 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
617 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
618 if(n0>n1) res = c0/std::sqrt(n0);
619 else res = c1/std::sqrt(n1);
620
621 return true;
622 }
623
624 EIGEN_DEVICE_FUNC
625 static inline void run(SolverType& solver, const MatrixType& mat, int options)
626 {
627 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
628 eigen_assert((options&~(EigVecMask|GenEigMask))==0
629 && (options&EigVecMask)!=EigVecMask
630 && "invalid option parameter");
631 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
632
633 EigenvectorsType& eivecs = solver.m_eivec;
634 VectorType& eivals = solver.m_eivalues;
635
636 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
637 Scalar shift = mat.trace() / Scalar(3);
638 // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
639 MatrixType scaledMat = mat.template selfadjointView<Lower>();
640 scaledMat.diagonal().array() -= shift;
641 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
642 if(scale > 0) scaledMat /= scale; // TODO for scale==0 we could save the remaining operations
643
644 // compute the eigenvalues
645 computeRoots(scaledMat,eivals);
646
647 // compute the eigenvectors
648 if(computeEigenvectors)
649 {
650 if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
651 {
652 // All three eigenvalues are numerically the same
653 eivecs.setIdentity();
654 }
655 else
656 {
657 MatrixType tmp;
658 tmp = scaledMat;
659
660 // Compute the eigenvector of the most distinct eigenvalue
661 Scalar d0 = eivals(2) - eivals(1);
662 Scalar d1 = eivals(1) - eivals(0);
663 Index k(0), l(2);
664 if(d0 > d1)
665 {
666 numext::swap(k,l);
667 d0 = d1;
668 }
669
670 // Compute the eigenvector of index k
671 {
672 tmp.diagonal().array () -= eivals(k);
673 // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
674 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
675 }
676
677 // Compute eigenvector of index l
679 {
680 // If d0 is too small, then the two other eigenvalues are numerically the same,
681 // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
682 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
683 eivecs.col(l).normalize();
684 }
685 else
686 {
687 tmp = scaledMat;
688 tmp.diagonal().array () -= eivals(l);
689
690 VectorType dummy;
691 extract_kernel(tmp, eivecs.col(l), dummy);
692 }
693
694 // Compute last eigenvector from the other two
695 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
696 }
697 }
698
699 // Rescale back to the original size.
700 eivals *= scale;
701 eivals.array() += shift;
702
703 solver.m_info = Success;
704 solver.m_isInitialized = true;
705 solver.m_eigenvectorsOk = computeEigenvectors;
706 }
707};
708
709// 2x2 direct eigenvalues decomposition, code from Hauke Heibel
710template<typename SolverType>
711struct direct_selfadjoint_eigenvalues<SolverType,2,false>
712{
713 typedef typename SolverType::MatrixType MatrixType;
714 typedef typename SolverType::RealVectorType VectorType;
715 typedef typename SolverType::Scalar Scalar;
716 typedef typename SolverType::EigenvectorsType EigenvectorsType;
717
718 EIGEN_DEVICE_FUNC
719 static inline void computeRoots(const MatrixType& m, VectorType& roots)
720 {
721 using std::sqrt;
722 const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
723 const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
724 roots(0) = t1 - t0;
725 roots(1) = t1 + t0;
726 }
727
728 EIGEN_DEVICE_FUNC
729 static inline void run(SolverType& solver, const MatrixType& mat, int options)
730 {
731 EIGEN_USING_STD_MATH(sqrt);
732 EIGEN_USING_STD_MATH(abs);
733
734 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
735 eigen_assert((options&~(EigVecMask|GenEigMask))==0
736 && (options&EigVecMask)!=EigVecMask
737 && "invalid option parameter");
738 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
739
740 EigenvectorsType& eivecs = solver.m_eivec;
741 VectorType& eivals = solver.m_eivalues;
742
743 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
744 Scalar shift = mat.trace() / Scalar(2);
745 MatrixType scaledMat = mat;
746 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
747 scaledMat.diagonal().array() -= shift;
748 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
749 if(scale > Scalar(0))
750 scaledMat /= scale;
751
752 // Compute the eigenvalues
753 computeRoots(scaledMat,eivals);
754
755 // compute the eigen vectors
756 if(computeEigenvectors)
757 {
758 if((eivals(1)-eivals(0))<=abs(eivals(1))*Eigen::NumTraits<Scalar>::epsilon())
759 {
760 eivecs.setIdentity();
761 }
762 else
763 {
764 scaledMat.diagonal().array () -= eivals(1);
765 Scalar a2 = numext::abs2(scaledMat(0,0));
766 Scalar c2 = numext::abs2(scaledMat(1,1));
767 Scalar b2 = numext::abs2(scaledMat(1,0));
768 if(a2>c2)
769 {
770 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
771 eivecs.col(1) /= sqrt(a2+b2);
772 }
773 else
774 {
775 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
776 eivecs.col(1) /= sqrt(c2+b2);
777 }
778
779 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
780 }
781 }
782
783 // Rescale back to the original size.
784 eivals *= scale;
785 eivals.array() += shift;
786
787 solver.m_info = Success;
788 solver.m_isInitialized = true;
789 solver.m_eigenvectorsOk = computeEigenvectors;
790 }
791};
792
793}
794
795template<typename MatrixType>
796EIGEN_DEVICE_FUNC
797SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
798::computeDirect(const MatrixType& matrix, int options)
799{
800 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
801 return *this;
802}
803
804namespace internal {
805template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
806EIGEN_DEVICE_FUNC
807static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
808{
809 using std::abs;
810 RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
811 RealScalar e = subdiag[end-1];
812 // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
813 // underflow thus leading to inf/NaN values when using the following commented code:
814// RealScalar e2 = numext::abs2(subdiag[end-1]);
815// RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
816 // This explain the following, somewhat more complicated, version:
817 RealScalar mu = diag[end];
818 if(td==RealScalar(0))
819 mu -= abs(e);
820 else
821 {
822 RealScalar e2 = numext::abs2(subdiag[end-1]);
823 RealScalar h = numext::hypot(td,e);
824 if(e2==RealScalar(0)) mu -= (e / (td + (td>RealScalar(0) ? RealScalar(1) : RealScalar(-1)))) * (e / h);
825 else mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
826 }
827
828 RealScalar x = diag[start] - mu;
829 RealScalar z = subdiag[start];
830 for (Index k = start; k < end; ++k)
831 {
832 JacobiRotation<RealScalar> rot;
833 rot.makeGivens(x, z);
834
835 // do T = G' T G
836 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
837 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
838
839 diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
840 diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
841 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
842
843
844 if (k > start)
845 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
846
847 x = subdiag[k];
848
849 if (k < end - 1)
850 {
851 z = -rot.s() * subdiag[k+1];
852 subdiag[k + 1] = rot.c() * subdiag[k+1];
853 }
854
855 // apply the givens rotation to the unit matrix Q = Q * G
856 if (matrixQ)
857 {
858 // FIXME if StorageOrder == RowMajor this operation is not very efficient
859 Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
860 q.applyOnTheRight(k,k+1,rot);
861 }
862 }
863}
864
865} // end namespace internal
866
867} // end namespace Eigen
868
869#endif // EIGEN_SELFADJOINTEIGENSOLVER_H
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: SelfAdjointEigenSolver.h:71
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition: SelfAdjointEigenSolver.h:83
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition: SelfAdjointEigenSolver.h:450
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SelfAdjointEigenSolver.h:343
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: SelfAdjointEigenSolver.h:282
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
Definition: SelfAdjointEigenSolver.h:94
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: SelfAdjointEigenSolver.h:259
SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Definition: SelfAdjointEigenSolver.h:162
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition: SelfAdjointEigenSolver.h:103
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: SelfAdjointEigenSolver.h:138
Eigen::Index Index
Definition: SelfAdjointEigenSolver.h:84
static const int m_maxIterations
Maximum number of iterations.
Definition: SelfAdjointEigenSolver.h:354
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition: SelfAdjointEigenSolver.h:331
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition: SelfAdjointEigenSolver.h:306
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Definition: SelfAdjointEigenSolver.h:798
Tridiagonal decomposition of a selfadjoint matrix.
Definition: Tridiagonalization.h:64
ComputationInfo
Definition: Constants.h:430
@ Success
Definition: Constants.h:432
@ NoConvergence
Definition: Constants.h:436
@ ComputeEigenvectors
Definition: Constants.h:395
Namespace containing all symbols from the Eigen library.
Definition: Core:287
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Definition: EigenBase.h:29
Derived & derived()
Definition: EigenBase.h:44
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151