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 namespace Eigen {
17 
18 template<typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
20 
21 namespace internal {
22 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
23 }
24 
68 template<typename _MatrixType> class SelfAdjointEigenSolver
69 {
70  public:
71 
72  typedef _MatrixType MatrixType;
73  enum {
74  Size = MatrixType::RowsAtCompileTime,
75  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
76  Options = MatrixType::Options,
77  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
78  };
79 
81  typedef typename MatrixType::Scalar Scalar;
82  typedef typename MatrixType::Index Index;
83 
91 
93 
99  typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
101 
113  : m_eivec(),
114  m_eivalues(),
115  m_subdiag(),
116  m_isInitialized(false)
117  { }
118 
132  : m_eivec(size, size),
133  m_eivalues(size),
134  m_subdiag(size > 1 ? size - 1 : 1),
135  m_isInitialized(false)
136  {}
137 
153  SelfAdjointEigenSolver(const MatrixType& matrix, int options = ComputeEigenvectors)
154  : m_eivec(matrix.rows(), matrix.cols()),
155  m_eivalues(matrix.cols()),
156  m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
157  m_isInitialized(false)
158  {
159  compute(matrix, options);
160  }
161 
192  SelfAdjointEigenSolver& compute(const MatrixType& matrix, int options = ComputeEigenvectors);
193 
208  SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
209 
228  const MatrixType& eigenvectors() const
229  {
230  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
231  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
232  return m_eivec;
233  }
234 
251  {
252  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
253  return m_eivalues;
254  }
255 
274  MatrixType operatorSqrt() const
275  {
276  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
277  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
278  return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
279  }
280 
299  MatrixType operatorInverseSqrt() const
300  {
301  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
302  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
303  return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
304  }
305 
311  {
312  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
313  return m_info;
314  }
315 
321  static const int m_maxIterations = 30;
322 
323  #ifdef EIGEN2_SUPPORT
324  SelfAdjointEigenSolver(const MatrixType& matrix, bool computeEigenvectors)
325  : m_eivec(matrix.rows(), matrix.cols()),
326  m_eivalues(matrix.cols()),
327  m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
328  m_isInitialized(false)
329  {
330  compute(matrix, computeEigenvectors);
331  }
332 
333  SelfAdjointEigenSolver(const MatrixType& matA, const MatrixType& matB, bool computeEigenvectors = true)
334  : m_eivec(matA.cols(), matA.cols()),
335  m_eivalues(matA.cols()),
336  m_subdiag(matA.cols() > 1 ? matA.cols() - 1 : 1),
337  m_isInitialized(false)
338  {
339  static_cast<GeneralizedSelfAdjointEigenSolver<MatrixType>*>(this)->compute(matA, matB, computeEigenvectors ? ComputeEigenvectors : EigenvaluesOnly);
340  }
341 
342  void compute(const MatrixType& matrix, bool computeEigenvectors)
343  {
344  compute(matrix, computeEigenvectors ? ComputeEigenvectors : EigenvaluesOnly);
345  }
346 
347  void compute(const MatrixType& matA, const MatrixType& matB, bool computeEigenvectors = true)
348  {
349  compute(matA, matB, computeEigenvectors ? ComputeEigenvectors : EigenvaluesOnly);
350  }
351  #endif // EIGEN2_SUPPORT
352 
353  protected:
354  MatrixType m_eivec;
360 };
361 
378 namespace internal {
379 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
380 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
381 }
382 
383 template<typename MatrixType>
385 ::compute(const MatrixType& matrix, int options)
386 {
387  using std::abs;
388  eigen_assert(matrix.cols() == matrix.rows());
389  eigen_assert((options&~(EigVecMask|GenEigMask))==0
390  && (options&EigVecMask)!=EigVecMask
391  && "invalid option parameter");
392  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
393  Index n = matrix.cols();
394  m_eivalues.resize(n,1);
395 
396  if(n==1)
397  {
398  m_eivalues.coeffRef(0,0) = numext::real(matrix.coeff(0,0));
399  if(computeEigenvectors)
400  m_eivec.setOnes(n,n);
401  m_info = Success;
402  m_isInitialized = true;
403  m_eigenvectorsOk = computeEigenvectors;
404  return *this;
405  }
406 
407  // declare some aliases
408  RealVectorType& diag = m_eivalues;
409  MatrixType& mat = m_eivec;
410 
411  // map the matrix coefficients to [-1:1] to avoid over- and underflow.
412  mat = matrix.template triangularView<Lower>();
413  RealScalar scale = mat.cwiseAbs().maxCoeff();
414  if(scale==RealScalar(0)) scale = RealScalar(1);
415  mat.template triangularView<Lower>() /= scale;
416  m_subdiag.resize(n-1);
417  internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
418 
419  Index end = n-1;
420  Index start = 0;
421  Index iter = 0; // total number of iterations
422 
423  while (end>0)
424  {
425  for (Index i = start; i<end; ++i)
426  if (internal::isMuchSmallerThan(abs(m_subdiag[i]),(abs(diag[i])+abs(diag[i+1]))))
427  m_subdiag[i] = 0;
428 
429  // find the largest unreduced block
430  while (end>0 && m_subdiag[end-1]==0)
431  {
432  end--;
433  }
434  if (end<=0)
435  break;
436 
437  // if we spent too many iterations, we give up
438  iter++;
439  if(iter > m_maxIterations * n) break;
440 
441  start = end - 1;
442  while (start>0 && m_subdiag[start-1]!=0)
443  start--;
444 
445  internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), m_subdiag.data(), start, end, computeEigenvectors ? m_eivec.data() : (Scalar*)0, n);
446  }
447 
448  if (iter <= m_maxIterations * n)
449  m_info = Success;
450  else
451  m_info = NoConvergence;
452 
453  // Sort eigenvalues and corresponding vectors.
454  // TODO make the sort optional ?
455  // TODO use a better sort algorithm !!
456  if (m_info == Success)
457  {
458  for (Index i = 0; i < n-1; ++i)
459  {
460  Index k;
461  m_eivalues.segment(i,n-i).minCoeff(&k);
462  if (k > 0)
463  {
464  std::swap(m_eivalues[i], m_eivalues[k+i]);
465  if(computeEigenvectors)
466  m_eivec.col(i).swap(m_eivec.col(k+i));
467  }
468  }
469  }
470 
471  // scale back the eigen values
472  m_eivalues *= scale;
473 
474  m_isInitialized = true;
475  m_eigenvectorsOk = computeEigenvectors;
476  return *this;
477 }
478 
479 
480 namespace internal {
481 
482 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
483 {
484  static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
485  { eig.compute(A,options); }
486 };
487 
488 template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
489 {
490  typedef typename SolverType::MatrixType MatrixType;
491  typedef typename SolverType::RealVectorType VectorType;
492  typedef typename SolverType::Scalar Scalar;
493 
494  static inline void computeRoots(const MatrixType& m, VectorType& roots)
495  {
496  using std::sqrt;
497  using std::atan2;
498  using std::cos;
499  using std::sin;
500  const Scalar s_inv3 = Scalar(1.0)/Scalar(3.0);
501  const Scalar s_sqrt3 = sqrt(Scalar(3.0));
502 
503  // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
504  // eigenvalues are the roots to this equation, all guaranteed to be
505  // real-valued, because the matrix is symmetric.
506  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);
507  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);
508  Scalar c2 = m(0,0) + m(1,1) + m(2,2);
509 
510  // Construct the parameters used in classifying the roots of the equation
511  // and in solving the equation for the roots in closed form.
512  Scalar c2_over_3 = c2*s_inv3;
513  Scalar a_over_3 = (c1 - c2*c2_over_3)*s_inv3;
514  if (a_over_3 > Scalar(0))
515  a_over_3 = Scalar(0);
516 
517  Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
518 
519  Scalar q = half_b*half_b + a_over_3*a_over_3*a_over_3;
520  if (q > Scalar(0))
521  q = Scalar(0);
522 
523  // Compute the eigenvalues by solving for the roots of the polynomial.
524  Scalar rho = sqrt(-a_over_3);
525  Scalar theta = atan2(sqrt(-q),half_b)*s_inv3;
526  Scalar cos_theta = cos(theta);
527  Scalar sin_theta = sin(theta);
528  roots(0) = c2_over_3 + Scalar(2)*rho*cos_theta;
529  roots(1) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
530  roots(2) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
531 
532  // Sort in increasing order.
533  if (roots(0) >= roots(1))
534  std::swap(roots(0),roots(1));
535  if (roots(1) >= roots(2))
536  {
537  std::swap(roots(1),roots(2));
538  if (roots(0) >= roots(1))
539  std::swap(roots(0),roots(1));
540  }
541  }
542 
543  static inline void run(SolverType& solver, const MatrixType& mat, int options)
544  {
545  using std::sqrt;
546  eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
547  eigen_assert((options&~(EigVecMask|GenEigMask))==0
548  && (options&EigVecMask)!=EigVecMask
549  && "invalid option parameter");
550  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
551 
552  MatrixType& eivecs = solver.m_eivec;
553  VectorType& eivals = solver.m_eivalues;
554 
555  // map the matrix coefficients to [-1:1] to avoid over- and underflow.
556  Scalar scale = mat.cwiseAbs().maxCoeff();
557  MatrixType scaledMat = mat / scale;
558 
559  // compute the eigenvalues
560  computeRoots(scaledMat,eivals);
561 
562  // compute the eigen vectors
563  if(computeEigenvectors)
564  {
565  Scalar safeNorm2 = Eigen::NumTraits<Scalar>::epsilon();
566  safeNorm2 *= safeNorm2;
567  if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
568  {
569  eivecs.setIdentity();
570  }
571  else
572  {
573  scaledMat = scaledMat.template selfadjointView<Lower>();
574  MatrixType tmp;
575  tmp = scaledMat;
576 
577  Scalar d0 = eivals(2) - eivals(1);
578  Scalar d1 = eivals(1) - eivals(0);
579  int k = d0 > d1 ? 2 : 0;
580  d0 = d0 > d1 ? d1 : d0;
581 
582  tmp.diagonal().array () -= eivals(k);
583  VectorType cross;
584  Scalar n;
585  n = (cross = tmp.row(0).cross(tmp.row(1))).squaredNorm();
586 
587  if(n>safeNorm2)
588  eivecs.col(k) = cross / sqrt(n);
589  else
590  {
591  n = (cross = tmp.row(0).cross(tmp.row(2))).squaredNorm();
592 
593  if(n>safeNorm2)
594  eivecs.col(k) = cross / sqrt(n);
595  else
596  {
597  n = (cross = tmp.row(1).cross(tmp.row(2))).squaredNorm();
598 
599  if(n>safeNorm2)
600  eivecs.col(k) = cross / sqrt(n);
601  else
602  {
603  // the input matrix and/or the eigenvaues probably contains some inf/NaN,
604  // => exit
605  // scale back to the original size.
606  eivals *= scale;
607 
608  solver.m_info = NumericalIssue;
609  solver.m_isInitialized = true;
610  solver.m_eigenvectorsOk = computeEigenvectors;
611  return;
612  }
613  }
614  }
615 
616  tmp = scaledMat;
617  tmp.diagonal().array() -= eivals(1);
618 
620  eivecs.col(1) = eivecs.col(k).unitOrthogonal();
621  else
622  {
623  n = (cross = eivecs.col(k).cross(tmp.row(0).normalized())).squaredNorm();
624  if(n>safeNorm2)
625  eivecs.col(1) = cross / sqrt(n);
626  else
627  {
628  n = (cross = eivecs.col(k).cross(tmp.row(1))).squaredNorm();
629  if(n>safeNorm2)
630  eivecs.col(1) = cross / sqrt(n);
631  else
632  {
633  n = (cross = eivecs.col(k).cross(tmp.row(2))).squaredNorm();
634  if(n>safeNorm2)
635  eivecs.col(1) = cross / sqrt(n);
636  else
637  {
638  // we should never reach this point,
639  // if so the last two eigenvalues are likely to ve very closed to each other
640  eivecs.col(1) = eivecs.col(k).unitOrthogonal();
641  }
642  }
643  }
644 
645  // make sure that eivecs[1] is orthogonal to eivecs[2]
646  Scalar d = eivecs.col(1).dot(eivecs.col(k));
647  eivecs.col(1) = (eivecs.col(1) - d * eivecs.col(k)).normalized();
648  }
649 
650  eivecs.col(k==2 ? 0 : 2) = eivecs.col(k).cross(eivecs.col(1)).normalized();
651  }
652  }
653  // Rescale back to the original size.
654  eivals *= scale;
655 
656  solver.m_info = Success;
657  solver.m_isInitialized = true;
658  solver.m_eigenvectorsOk = computeEigenvectors;
659  }
660 };
661 
662 // 2x2 direct eigenvalues decomposition, code from Hauke Heibel
663 template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,2,false>
664 {
665  typedef typename SolverType::MatrixType MatrixType;
666  typedef typename SolverType::RealVectorType VectorType;
667  typedef typename SolverType::Scalar Scalar;
668 
669  static inline void computeRoots(const MatrixType& m, VectorType& roots)
670  {
671  using std::sqrt;
672  const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*m(1,0)*m(1,0));
673  const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
674  roots(0) = t1 - t0;
675  roots(1) = t1 + t0;
676  }
677 
678  static inline void run(SolverType& solver, const MatrixType& mat, int options)
679  {
680  using std::sqrt;
681  eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
682  eigen_assert((options&~(EigVecMask|GenEigMask))==0
683  && (options&EigVecMask)!=EigVecMask
684  && "invalid option parameter");
685  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
686 
687  MatrixType& eivecs = solver.m_eivec;
688  VectorType& eivals = solver.m_eivalues;
689 
690  // map the matrix coefficients to [-1:1] to avoid over- and underflow.
691  Scalar scale = mat.cwiseAbs().maxCoeff();
692  scale = (std::max)(scale,Scalar(1));
693  MatrixType scaledMat = mat / scale;
694 
695  // Compute the eigenvalues
696  computeRoots(scaledMat,eivals);
697 
698  // compute the eigen vectors
699  if(computeEigenvectors)
700  {
701  scaledMat.diagonal().array () -= eivals(1);
702  Scalar a2 = numext::abs2(scaledMat(0,0));
703  Scalar c2 = numext::abs2(scaledMat(1,1));
704  Scalar b2 = numext::abs2(scaledMat(1,0));
705  if(a2>c2)
706  {
707  eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
708  eivecs.col(1) /= sqrt(a2+b2);
709  }
710  else
711  {
712  eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
713  eivecs.col(1) /= sqrt(c2+b2);
714  }
715 
716  eivecs.col(0) << eivecs.col(1).unitOrthogonal();
717  }
718 
719  // Rescale back to the original size.
720  eivals *= scale;
721 
722  solver.m_info = Success;
723  solver.m_isInitialized = true;
724  solver.m_eigenvectorsOk = computeEigenvectors;
725  }
726 };
727 
728 }
729 
730 template<typename MatrixType>
732 ::computeDirect(const MatrixType& matrix, int options)
733 {
735  return *this;
736 }
737 
738 namespace internal {
739 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
740 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
741 {
742  using std::abs;
743  RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
744  RealScalar e = subdiag[end-1];
745  // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
746  // underflow thus leading to inf/NaN values when using the following commented code:
747 // RealScalar e2 = numext::abs2(subdiag[end-1]);
748 // RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
749  // This explain the following, somewhat more complicated, version:
750  RealScalar mu = diag[end];
751  if(td==0)
752  mu -= abs(e);
753  else
754  {
755  RealScalar e2 = numext::abs2(subdiag[end-1]);
756  RealScalar h = numext::hypot(td,e);
757  if(e2==0) mu -= (e / (td + (td>0 ? 1 : -1))) * (e / h);
758  else mu -= e2 / (td + (td>0 ? h : -h));
759  }
760 
761  RealScalar x = diag[start] - mu;
762  RealScalar z = subdiag[start];
763  for (Index k = start; k < end; ++k)
764  {
766  rot.makeGivens(x, z);
767 
768  // do T = G' T G
769  RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
770  RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
771 
772  diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
773  diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
774  subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
775 
776 
777  if (k > start)
778  subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
779 
780  x = subdiag[k];
781 
782  if (k < end - 1)
783  {
784  z = -rot.s() * subdiag[k+1];
785  subdiag[k + 1] = rot.c() * subdiag[k+1];
786  }
787 
788  // apply the givens rotation to the unit matrix Q = Q * G
789  if (matrixQ)
790  {
791  // FIXME if StorageOrder == RowMajor this operation is not very efficient
793  q.applyOnTheRight(k,k+1,rot);
794  }
795  }
796 }
797 
798 } // end namespace internal
799 
800 } // end namespace Eigen
801 
802 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
Scalar & c()
Definition: Jacobi.h:45
static void computeRoots(const MatrixType &m, VectorType &roots)
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
USING_NAMESPACE_ACADO IntermediateState sin(const Expression &arg)
IntermediateState sqrt(const Expression &arg)
static void run(SolverType &solver, const MatrixType &mat, int options)
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Definition: Jacobi.h:148
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:104
static void run(SolverType &solver, const MatrixType &mat, int options)
const AutoDiffScalar< Matrix< typename internal::traits< DerTypeA >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
Computes eigenvalues and eigenvectors of selfadjoint matrices.
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
SelfAdjointEigenSolver(const MatrixType &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Rotation given by a cosine-sine pair.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
Tridiagonal decomposition of a selfadjoint matrix.
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, typename NumTraits< Scalar >::Real precision=NumTraits< Scalar >::dummy_precision())
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs2_op< Scalar >, const Derived > abs2() const
IntermediateState cos(const Expression &arg)
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
RealReturnType real() const
ComputationInfo info() const
Reports whether previous computation was successful.
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Provides a generic way to set and pass user-specified options.
Definition: options.hpp:65
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
static void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
const MatrixType & eigenvectors() const
Returns the eigenvectors of given matrix.
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
SelfAdjointEigenSolver & compute(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
static void computeRoots(const MatrixType &m, VectorType &roots)
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a direct algorithm.
Scalar & s()
Definition: Jacobi.h:47
static void run(SolverType &eig, const typename SolverType::MatrixType &A, int options)
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
#define eigen_assert(x)
TridiagonalizationType::SubDiagonalType m_subdiag
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:35:04