11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12 #define EIGEN_SELFADJOINTEIGENSOLVER_H
14 #include "./Tridiagonalization.h"
18 template<
typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
22 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
72 typedef _MatrixType MatrixType;
74 Size = MatrixType::RowsAtCompileTime,
75 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
76 Options = MatrixType::Options,
77 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
81 typedef typename MatrixType::Scalar
Scalar;
82 typedef typename MatrixType::Index Index;
101 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type
RealVectorType;
118 m_isInitialized(false)
134 : m_eivec(size, size),
136 m_subdiag(size > 1 ? size - 1 : 1),
137 m_isInitialized(false)
156 : m_eivec(matrix.rows(), matrix.cols()),
157 m_eivalues(matrix.cols()),
158 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
159 m_isInitialized(false)
232 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
233 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
254 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
278 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
279 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
280 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
303 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
304 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
305 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
314 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
325 #ifdef EIGEN2_SUPPORT
327 : m_eivec(matrix.rows(), matrix.cols()),
328 m_eivalues(matrix.cols()),
329 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
330 m_isInitialized(false)
332 compute(matrix, computeEigenvectors);
336 : m_eivec(matA.cols(), matA.cols()),
337 m_eivalues(matA.cols()),
338 m_subdiag(matA.cols() > 1 ? matA.cols() - 1 : 1),
339 m_isInitialized(false)
344 void compute(
const MatrixType& matrix,
bool computeEigenvectors)
349 void compute(
const MatrixType& matA,
const MatrixType& matB,
bool computeEigenvectors =
true)
353 #endif // EIGEN2_SUPPORT
356 static void check_template_parameters()
358 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar);
361 EigenvectorsType m_eivec;
363 typename TridiagonalizationType::SubDiagonalType m_subdiag;
365 bool m_isInitialized;
366 bool m_eigenvectorsOk;
386 template<
typename RealScalar,
typename Scalar,
typename Index>
387 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
390 template<
typename MatrixType>
394 check_template_parameters();
397 eigen_assert(matrix.cols() == matrix.rows());
398 eigen_assert((options&~(EigVecMask|GenEigMask))==0
399 && (options&EigVecMask)!=EigVecMask
400 &&
"invalid option parameter");
402 Index n = matrix.cols();
403 m_eivalues.resize(n,1);
407 m_eivalues.coeffRef(0,0) = numext::real(matrix.coeff(0,0));
408 if(computeEigenvectors)
409 m_eivec.setOnes(n,n);
411 m_isInitialized =
true;
412 m_eigenvectorsOk = computeEigenvectors;
421 mat = matrix.template triangularView<Lower>();
424 mat.template triangularView<Lower>() /= scale;
426 internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
434 for (Index i = start; i<end; ++i)
435 if (internal::isMuchSmallerThan(abs(m_subdiag[i]),(abs(diag[i])+abs(diag[i+1]))))
439 while (end>0 && m_subdiag[end-1]==0)
448 if(iter > m_maxIterations * n)
break;
451 while (start>0 && m_subdiag[start-1]!=0)
454 internal::tridiagonal_qr_step(diag.data(), m_subdiag.data(), start, end, computeEigenvectors ? m_eivec.data() : (
Scalar*)0, n);
457 if (iter <= m_maxIterations * n)
467 for (Index i = 0; i < n-1; ++i)
470 m_eivalues.segment(i,n-i).minCoeff(&k);
473 std::swap(m_eivalues[i], m_eivalues[k+i]);
474 if(computeEigenvectors)
475 m_eivec.col(i).swap(m_eivec.col(k+i));
483 m_isInitialized =
true;
484 m_eigenvectorsOk = computeEigenvectors;
491 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues
493 static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options)
494 { eig.compute(A,options); }
497 template<
typename SolverType>
struct direct_selfadjoint_eigenvalues<SolverType,3,false>
499 typedef typename SolverType::MatrixType MatrixType;
500 typedef typename SolverType::RealVectorType VectorType;
501 typedef typename SolverType::Scalar Scalar;
502 typedef typename MatrixType::Index Index;
503 typedef typename SolverType::EigenvectorsType EigenvectorsType;
509 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
515 const Scalar s_inv3 = Scalar(1.0)/Scalar(3.0);
516 const Scalar s_sqrt3 = sqrt(Scalar(3.0));
521 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);
522 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);
523 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
527 Scalar c2_over_3 = c2*s_inv3;
528 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
529 if(a_over_3<Scalar(0))
530 a_over_3 = Scalar(0);
532 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
534 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
539 Scalar rho = sqrt(a_over_3);
540 Scalar theta = atan2(sqrt(q),half_b)*s_inv3;
541 Scalar cos_theta = cos(theta);
542 Scalar sin_theta = sin(theta);
544 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
545 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
546 roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
549 static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
554 mat.diagonal().cwiseAbs().maxCoeff(&i0);
557 representative = mat.col(i0);
560 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
561 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
562 if(n0>n1) res = c0/std::sqrt(n0);
563 else res = c1/std::sqrt(n1);
568 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
570 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
571 eigen_assert((options&~(EigVecMask|GenEigMask))==0
572 && (options&EigVecMask)!=EigVecMask
573 &&
"invalid option parameter");
576 EigenvectorsType& eivecs = solver.m_eivec;
577 VectorType& eivals = solver.m_eivalues;
580 Scalar shift = mat.trace() / Scalar(3);
582 MatrixType scaledMat = mat.template selfadjointView<Lower>();
583 scaledMat.diagonal().array() -= shift;
584 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
585 if(scale > 0) scaledMat /= scale;
588 computeRoots(scaledMat,eivals);
591 if(computeEigenvectors)
596 eivecs.setIdentity();
604 Scalar d0 = eivals(2) - eivals(1);
605 Scalar d1 = eivals(1) - eivals(0);
615 tmp.diagonal().array () -= eivals(k);
617 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
625 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
626 eivecs.col(l).normalize();
631 tmp.diagonal().array () -= eivals(l);
634 extract_kernel(tmp, eivecs.col(l), dummy);
638 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
644 eivals.array() += shift;
647 solver.m_isInitialized =
true;
648 solver.m_eigenvectorsOk = computeEigenvectors;
653 template<
typename SolverType>
struct direct_selfadjoint_eigenvalues<SolverType,2,false>
655 typedef typename SolverType::MatrixType MatrixType;
656 typedef typename SolverType::RealVectorType VectorType;
657 typedef typename SolverType::Scalar Scalar;
658 typedef typename SolverType::EigenvectorsType EigenvectorsType;
660 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
663 const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
664 const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
669 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
674 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
675 eigen_assert((options&~(EigVecMask|GenEigMask))==0
676 && (options&EigVecMask)!=EigVecMask
677 &&
"invalid option parameter");
680 EigenvectorsType& eivecs = solver.m_eivec;
681 VectorType& eivals = solver.m_eivalues;
684 Scalar scale = mat.cwiseAbs().maxCoeff();
685 scale = (std::max)(scale,Scalar(1));
686 MatrixType scaledMat = mat / scale;
689 computeRoots(scaledMat,eivals);
692 if(computeEigenvectors)
696 eivecs.setIdentity();
700 scaledMat.diagonal().array () -= eivals(1);
701 Scalar a2 = numext::abs2(scaledMat(0,0));
702 Scalar c2 = numext::abs2(scaledMat(1,1));
703 Scalar b2 = numext::abs2(scaledMat(1,0));
706 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
707 eivecs.col(1) /= sqrt(a2+b2);
711 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
712 eivecs.col(1) /= sqrt(c2+b2);
715 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
723 solver.m_isInitialized =
true;
724 solver.m_eigenvectorsOk = computeEigenvectors;
730 template<
typename MatrixType>
734 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*
this,matrix,options);
739 template<
typename RealScalar,
typename Scalar,
typename Index>
740 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
743 RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
744 RealScalar e = subdiag[end-1];
750 RealScalar mu = diag[end];
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));
761 RealScalar x = diag[start] - mu;
762 RealScalar z = subdiag[start];
763 for (Index k = start; k < end; ++k)
765 JacobiRotation<RealScalar> rot;
766 rot.makeGivens(x, z);
769 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
770 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
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;
778 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
784 z = -rot.s() * subdiag[k+1];
785 subdiag[k + 1] = rot.c() * subdiag[k+1];
791 Map<Matrix<Scalar,Dynamic,Dynamic,ColMajor> > q(matrixQ,n,n);
792 q.applyOnTheRight(k,k+1,rot);
801 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SelfAdjointEigenSolver.h:312
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition: SelfAdjointEigenSolver.h:81
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition: SelfAdjointEigenSolver.h:301
Definition: Constants.h:339
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: SelfAdjointEigenSolver.h:252
SelfAdjointEigenSolver(const MatrixType &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Definition: SelfAdjointEigenSolver.h:155
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: SelfAdjointEigenSolver.h:133
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: SelfAdjointEigenSolver.h:68
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
Tridiagonal decomposition of a selfadjoint matrix.
Definition: Tridiagonalization.h:61
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a direct algorithm.
Definition: SelfAdjointEigenSolver.h:732
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition: GeneralizedSelfAdjointEigenSolver.h:48
static const int m_maxIterations
Maximum number of iterations.
Definition: SelfAdjointEigenSolver.h:323
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition: SelfAdjointEigenSolver.h:276
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
Definition: SelfAdjointEigenSolver.h:92
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition: SelfAdjointEigenSolver.h:101
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: SelfAdjointEigenSolver.h:230
Definition: Constants.h:380
Definition: Constants.h:376
void resize(Index nbRows, Index nbCols)
Definition: PlainObjectBase.h:235
SelfAdjointEigenSolver & compute(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
Definition: SelfAdjointEigenSolver.h:392
ComputationInfo
Definition: Constants.h:374
SelfAdjointEigenSolver()
Default constructor for fixed-size matrices.
Definition: SelfAdjointEigenSolver.h:114
Definition: Constants.h:336