SelfAdjointEigenSolver.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-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 
18 namespace Eigen {
19 
20 template<typename MatrixType_>
21 class GeneralizedSelfAdjointEigenSolver;
22 
23 namespace internal {
24 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
25 
26 template<typename MatrixType, typename DiagType, typename SubDiagType>
28 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
29 }
30 
78 template<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;
92  typedef Eigen::Index Index;
93 
95 
103 
104  friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
105 
111  typedef typename internal::plain_col_type<MatrixType, Scalar>::type VectorType;
112  typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
115 
128  : m_eivec(),
129  m_workspace(),
130  m_eivalues(),
131  m_subdiag(),
132  m_hcoeffs(),
133  m_info(InvalidInput),
134  m_isInitialized(false),
135  m_eigenvectorsOk(false)
136  { }
137 
152  : m_eivec(size, size),
153  m_workspace(size),
154  m_eivalues(size),
155  m_subdiag(size > 1 ? size - 1 : 1),
156  m_hcoeffs(size > 1 ? size - 1 : 1),
157  m_isInitialized(false),
158  m_eigenvectorsOk(false)
159  {}
160 
176  template<typename InputType>
178  explicit SelfAdjointEigenSolver(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors)
179  : m_eivec(matrix.rows(), matrix.cols()),
180  m_workspace(matrix.cols()),
181  m_eivalues(matrix.cols()),
182  m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
183  m_hcoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1),
184  m_isInitialized(false),
185  m_eigenvectorsOk(false)
186  {
187  compute(matrix.derived(), options);
188  }
189 
220  template<typename InputType>
223 
244 
258 
284  {
285  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
286  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
287  return m_eivec;
288  }
289 
307  {
308  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
309  return m_eivalues;
310  }
311 
331  {
332  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
333  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
334  return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
335  }
336 
356  {
357  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
358  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
359  return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
360  }
361 
368  {
369  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
370  return m_info;
371  }
372 
378  static const int m_maxIterations = 30;
379 
380  protected:
382 
387  typename TridiagonalizationType::CoeffVectorType m_hcoeffs;
391 };
392 
393 namespace internal {
414 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
416 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
417 }
418 
419 template<typename MatrixType>
420 template<typename InputType>
423 ::compute(const EigenBase<InputType>& a_matrix, int options)
424 {
425  const InputType &matrix(a_matrix.derived());
426 
428  eigen_assert(matrix.cols() == matrix.rows());
429  eigen_assert((options&~(EigVecMask|GenEigMask))==0
430  && (options&EigVecMask)!=EigVecMask
431  && "invalid option parameter");
432  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
433  Index n = matrix.cols();
434  m_eivalues.resize(n,1);
435 
436  if(n==1)
437  {
438  m_eivec = matrix;
439  m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
440  if(computeEigenvectors)
441  m_eivec.setOnes(n,n);
442  m_info = Success;
443  m_isInitialized = true;
444  m_eigenvectorsOk = computeEigenvectors;
445  return *this;
446  }
447 
448  // declare some aliases
449  RealVectorType& diag = m_eivalues;
451 
452  // map the matrix coefficients to [-1:1] to avoid over- and underflow.
453  mat = matrix.template triangularView<Lower>();
454  RealScalar scale = mat.cwiseAbs().maxCoeff();
455  if(numext::is_exactly_zero(scale)) scale = RealScalar(1);
456  mat.template triangularView<Lower>() /= scale;
457  m_subdiag.resize(n-1);
458  m_hcoeffs.resize(n-1);
460 
462 
463  // scale back the eigen values
464  m_eivalues *= scale;
465 
466  m_isInitialized = true;
467  m_eigenvectorsOk = computeEigenvectors;
468  return *this;
469 }
470 
471 template<typename MatrixType>
473 ::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
474 {
475  //TODO : Add an option to scale the values beforehand
476  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
477 
478  m_eivalues = diag;
479  m_subdiag = subdiag;
480  if (computeEigenvectors)
481  {
482  m_eivec.setIdentity(diag.size(), diag.size());
483  }
485 
486  m_isInitialized = true;
487  m_eigenvectorsOk = computeEigenvectors;
488  return *this;
489 }
490 
491 namespace internal {
503 template<typename MatrixType, typename DiagType, typename SubDiagType>
505 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
506 {
508  typedef typename MatrixType::Scalar Scalar;
509 
510  Index n = diag.size();
511  Index end = n-1;
512  Index start = 0;
513  Index iter = 0; // total number of iterations
514 
515  typedef typename DiagType::RealScalar RealScalar;
516  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
517  const RealScalar precision_inv = RealScalar(1)/NumTraits<RealScalar>::epsilon();
518  while (end>0)
519  {
520  for (Index i = start; i<end; ++i) {
521  if (numext::abs(subdiag[i]) < considerAsZero) {
522  subdiag[i] = RealScalar(0);
523  } else {
524  // abs(subdiag[i]) <= epsilon * sqrt(abs(diag[i]) + abs(diag[i+1]))
525  // Scaled to prevent underflows.
526  const RealScalar scaled_subdiag = precision_inv * subdiag[i];
527  if (scaled_subdiag * scaled_subdiag <= (numext::abs(diag[i])+numext::abs(diag[i+1]))) {
528  subdiag[i] = RealScalar(0);
529  }
530  }
531  }
532 
533  // find the largest unreduced block at the end of the matrix.
534  while (end>0 && numext::is_exactly_zero(subdiag[end - 1]))
535  {
536  end--;
537  }
538  if (end<=0)
539  break;
540 
541  // if we spent too many iterations, we give up
542  iter++;
543  if(iter > maxIterations * n) break;
544 
545  start = end - 1;
546  while (start>0 && !numext::is_exactly_zero(subdiag[start - 1]))
547  start--;
548 
549  internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
550  }
551  if (iter <= maxIterations * n)
552  info = Success;
553  else
555 
556  // Sort eigenvalues and corresponding vectors.
557  // TODO make the sort optional ?
558  // TODO use a better sort algorithm !!
559  if (info == Success)
560  {
561  for (Index i = 0; i < n-1; ++i)
562  {
563  Index k;
564  diag.segment(i,n-i).minCoeff(&k);
565  if (k > 0)
566  {
567  numext::swap(diag[i], diag[k+i]);
568  if(computeEigenvectors)
569  eivec.col(i).swap(eivec.col(k+i));
570  }
571  }
572  }
573  return info;
574 }
575 
576 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
577 {
579  static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
580  { eig.compute(A,options); }
581 };
582 
583 template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
584 {
585  typedef typename SolverType::MatrixType MatrixType;
586  typedef typename SolverType::RealVectorType VectorType;
587  typedef typename SolverType::Scalar Scalar;
588  typedef typename SolverType::EigenvectorsType EigenvectorsType;
589 
590 
596  static inline void computeRoots(const MatrixType& m, VectorType& roots)
597  {
602  const Scalar s_inv3 = Scalar(1)/Scalar(3);
603  const Scalar s_sqrt3 = sqrt(Scalar(3));
604 
605  // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
606  // eigenvalues are the roots to this equation, all guaranteed to be
607  // real-valued, because the matrix is symmetric.
608  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);
609  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);
610  Scalar c2 = m(0,0) + m(1,1) + m(2,2);
611 
612  // Construct the parameters used in classifying the roots of the equation
613  // and in solving the equation for the roots in closed form.
614  Scalar c2_over_3 = c2*s_inv3;
615  Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
616  a_over_3 = numext::maxi(a_over_3, Scalar(0));
617 
618  Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
619 
620  Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
621  q = numext::maxi(q, Scalar(0));
622 
623  // Compute the eigenvalues by solving for the roots of the polynomial.
624  Scalar rho = sqrt(a_over_3);
625  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]
626  Scalar cos_theta = cos(theta);
627  Scalar sin_theta = sin(theta);
628  // roots are already sorted, since cos is monotonically decreasing on [0, pi]
629  roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
630  roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
631  roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
632  }
633 
635  static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
636  {
639  Index i0;
640  // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
641  mat.diagonal().cwiseAbs().maxCoeff(&i0);
642  // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
643  // so let's save it:
644  representative = mat.col(i0);
645  Scalar n0, n1;
646  VectorType c0, c1;
647  n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
648  n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
649  if(n0>n1) res = c0/sqrt(n0);
650  else res = c1/sqrt(n1);
651 
652  return true;
653  }
654 
656  static inline void run(SolverType& solver, const MatrixType& mat, int options)
657  {
658  eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
659  eigen_assert((options&~(EigVecMask|GenEigMask))==0
660  && (options&EigVecMask)!=EigVecMask
661  && "invalid option parameter");
662  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
663 
664  EigenvectorsType& eivecs = solver.m_eivec;
665  VectorType& eivals = solver.m_eivalues;
666 
667  // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
668  Scalar shift = mat.trace() / Scalar(3);
669  // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
670  MatrixType scaledMat = mat.template selfadjointView<Lower>();
671  scaledMat.diagonal().array() -= shift;
672  Scalar scale = scaledMat.cwiseAbs().maxCoeff();
673  if(scale > 0) scaledMat /= scale; // TODO for scale==0 we could save the remaining operations
674 
675  // compute the eigenvalues
676  computeRoots(scaledMat,eivals);
677 
678  // compute the eigenvectors
679  if(computeEigenvectors)
680  {
682  {
683  // All three eigenvalues are numerically the same
684  eivecs.setIdentity();
685  }
686  else
687  {
688  MatrixType tmp;
689  tmp = scaledMat;
690 
691  // Compute the eigenvector of the most distinct eigenvalue
692  Scalar d0 = eivals(2) - eivals(1);
693  Scalar d1 = eivals(1) - eivals(0);
694  Index k(0), l(2);
695  if(d0 > d1)
696  {
697  numext::swap(k,l);
698  d0 = d1;
699  }
700 
701  // Compute the eigenvector of index k
702  {
703  tmp.diagonal().array () -= eivals(k);
704  // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
705  extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
706  }
707 
708  // Compute eigenvector of index l
709  if(d0<=2*Eigen::NumTraits<Scalar>::epsilon()*d1)
710  {
711  // If d0 is too small, then the two other eigenvalues are numerically the same,
712  // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
713  eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
714  eivecs.col(l).normalize();
715  }
716  else
717  {
718  tmp = scaledMat;
719  tmp.diagonal().array () -= eivals(l);
720 
721  VectorType dummy;
722  extract_kernel(tmp, eivecs.col(l), dummy);
723  }
724 
725  // Compute last eigenvector from the other two
726  eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
727  }
728  }
729 
730  // Rescale back to the original size.
731  eivals *= scale;
732  eivals.array() += shift;
733 
734  solver.m_info = Success;
735  solver.m_isInitialized = true;
736  solver.m_eigenvectorsOk = computeEigenvectors;
737  }
738 };
739 
740 // 2x2 direct eigenvalues decomposition, code from Hauke Heibel
741 template<typename SolverType>
742 struct direct_selfadjoint_eigenvalues<SolverType,2,false>
743 {
744  typedef typename SolverType::MatrixType MatrixType;
745  typedef typename SolverType::RealVectorType VectorType;
746  typedef typename SolverType::Scalar Scalar;
747  typedef typename SolverType::EigenvectorsType EigenvectorsType;
748 
750  static inline void computeRoots(const MatrixType& m, VectorType& roots)
751  {
753  const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
754  const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
755  roots(0) = t1 - t0;
756  roots(1) = t1 + t0;
757  }
758 
760  static inline void run(SolverType& solver, const MatrixType& mat, int options)
761  {
764 
765  eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
766  eigen_assert((options&~(EigVecMask|GenEigMask))==0
767  && (options&EigVecMask)!=EigVecMask
768  && "invalid option parameter");
769  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
770 
771  EigenvectorsType& eivecs = solver.m_eivec;
772  VectorType& eivals = solver.m_eivalues;
773 
774  // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
775  Scalar shift = mat.trace() / Scalar(2);
776  MatrixType scaledMat = mat;
777  scaledMat.coeffRef(0,1) = mat.coeff(1,0);
778  scaledMat.diagonal().array() -= shift;
779  Scalar scale = scaledMat.cwiseAbs().maxCoeff();
780  if(scale > Scalar(0))
781  scaledMat /= scale;
782 
783  // Compute the eigenvalues
784  computeRoots(scaledMat,eivals);
785 
786  // compute the eigen vectors
787  if(computeEigenvectors)
788  {
790  {
791  eivecs.setIdentity();
792  }
793  else
794  {
795  scaledMat.diagonal().array () -= eivals(1);
796  Scalar a2 = numext::abs2(scaledMat(0,0));
797  Scalar c2 = numext::abs2(scaledMat(1,1));
798  Scalar b2 = numext::abs2(scaledMat(1,0));
799  if(a2>c2)
800  {
801  eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
802  eivecs.col(1) /= sqrt(a2+b2);
803  }
804  else
805  {
806  eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
807  eivecs.col(1) /= sqrt(c2+b2);
808  }
809 
810  eivecs.col(0) << eivecs.col(1).unitOrthogonal();
811  }
812  }
813 
814  // Rescale back to the original size.
815  eivals *= scale;
816  eivals.array() += shift;
817 
818  solver.m_info = Success;
819  solver.m_isInitialized = true;
820  solver.m_eigenvectorsOk = computeEigenvectors;
821  }
822 };
823 
824 }
825 
826 template<typename MatrixType>
828 SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
829 ::computeDirect(const MatrixType& matrix, int options)
830 {
831  internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
832  return *this;
833 }
834 
835 namespace internal {
836 
837 // Francis implicit QR step.
838 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
840 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
841 {
842  // Wilkinson Shift.
843  RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
844  RealScalar e = subdiag[end-1];
845  // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
846  // underflow thus leading to inf/NaN values when using the following commented code:
847  // RealScalar e2 = numext::abs2(subdiag[end-1]);
848  // RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
849  // This explain the following, somewhat more complicated, version:
850  RealScalar mu = diag[end];
851  if(numext::is_exactly_zero(td)) {
852  mu -= numext::abs(e);
853  } else if (!numext::is_exactly_zero(e)) {
854  const RealScalar e2 = numext::abs2(e);
855  const RealScalar h = numext::hypot(td,e);
856  if(numext::is_exactly_zero(e2)) {
857  mu -= e / ((td + (td>RealScalar(0) ? h : -h)) / e);
858  } else {
859  mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
860  }
861  }
862 
863  RealScalar x = diag[start] - mu;
864  RealScalar z = subdiag[start];
865  // If z ever becomes zero, the Givens rotation will be the identity and
866  // z will stay zero for all future iterations.
867  for (Index k = start; k < end && !numext::is_exactly_zero(z); ++k)
868  {
870  rot.makeGivens(x, z);
871 
872  // do T = G' T G
873  RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
874  RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
875 
876  diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
877  diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
878  subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
879 
880  if (k > start)
881  subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
882 
883  // "Chasing the bulge" to return to triangular form.
884  x = subdiag[k];
885  if (k < end - 1)
886  {
887  z = -rot.s() * subdiag[k+1];
888  subdiag[k + 1] = rot.c() * subdiag[k+1];
889  }
890 
891  // apply the givens rotation to the unit matrix Q = Q * G
892  if (matrixQ)
893  {
894  // FIXME if StorageOrder == RowMajor this operation is not very efficient
896  q.applyOnTheRight(k,k+1,rot);
897  }
898  }
899 }
900 
901 } // end namespace internal
902 
903 } // end namespace Eigen
904 
905 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
Matrix3f m
const CwiseBinaryOp< atan2< Scalar >, const Derived, const OtherDerived > atan2(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
BiCGSTAB< SparseMatrix< double > > solver
int n
RealReturnType real() const
MatrixXcf A
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1080
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
VectorXcd eivals
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:81
Matrix< float, 1, Dynamic > MatrixType
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:37
Scalar & s()
Definition: Jacobi.h:51
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition: Jacobi.h:164
Scalar & c()
Definition: Jacobi.h:49
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:98
constexpr const Scalar & coeff(Index rowId, Index colId) const
Derived & setOnes(Index size)
constexpr void resize(Index rows, Index cols)
Computes eigenvalues and eigenvectors of selfadjoint matrices.
TridiagonalizationType::CoeffVectorType m_hcoeffs
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
TridiagonalizationType::SubDiagonalType m_subdiag
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType_.
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
ComputationInfo info() const
Reports whether previous computation was successful.
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
static const int m_maxIterations
Maximum number of iterations.
internal::plain_col_type< MatrixType, Scalar >::type VectorType
Type for vector of eigenvalues as returned by eigenvalues().
SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Tridiagonal decomposition of a selfadjoint matrix.
static const lastp1_t end
ComputationInfo
Definition: Constants.h:444
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
@ GenEigMask
Definition: Constants.h:420
@ EigVecMask
Definition: Constants.h:409
@ ComputeEigenvectors
Definition: Constants.h:407
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
static void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
void swap(T &a, T &b)
Definition: Meta.h:419
bool abs2(bool x)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
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_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Derived & derived()
Definition: EigenBase.h:48
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231