11#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
14#include "./Tridiagonalization.h"
16#include "./InternalHeaderCheck.h"
20template<
typename MatrixType_>
21class GeneralizedSelfAdjointEigenSolver;
24template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
26template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
28ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
const Index maxIterations,
bool computeEigenvectors, MatrixType& eivec);
82 typedef MatrixType_ MatrixType;
84 Size = MatrixType::RowsAtCompileTime,
85 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
86 Options = MatrixType::Options,
87 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
91 typedef typename MatrixType::Scalar
Scalar;
111 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type
RealVectorType;
131 m_info(InvalidInput),
132 m_isInitialized(false),
133 m_eigenvectorsOk(false)
150 : m_eivec(size, 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)
173 template<
typename InputType>
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)
216 template<
typename InputType>
281 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
282 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
304 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
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();
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();
365 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
377 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar)
384 bool m_isInitialized;
385 bool m_eigenvectorsOk;
409template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
411static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end, Scalar* matrixQ,
Index n);
414template<
typename MatrixType>
415template<
typename InputType>
418::compute(
const EigenBase<InputType>& a_matrix,
int options)
420 const InputType &matrix(a_matrix.derived());
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");
428 Index n = matrix.cols();
429 m_eivalues.resize(n,1);
434 m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
435 if(computeEigenvectors)
436 m_eivec.setOnes(n,n);
438 m_isInitialized =
true;
439 m_eigenvectorsOk = computeEigenvectors;
444 RealVectorType& diag = m_eivalues;
445 EigenvectorsType& mat = m_eivec;
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);
456 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
461 m_isInitialized =
true;
462 m_eigenvectorsOk = computeEigenvectors;
466template<
typename MatrixType>
475 if (computeEigenvectors)
479 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
481 m_isInitialized =
true;
482 m_eigenvectorsOk = computeEigenvectors;
498template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
500ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
const Index maxIterations,
bool computeEigenvectors, MatrixType& eivec)
503 typedef typename MatrixType::Scalar Scalar;
505 Index n = diag.size();
510 typedef typename DiagType::RealScalar RealScalar;
511 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
516 if (numext::abs(subdiag[i]) < considerAsZero) {
517 subdiag[i] = RealScalar(0);
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);
529 while (
end>0 && subdiag[
end-1]==RealScalar(0))
538 if(iter > maxIterations * n)
break;
541 while (start>0 && subdiag[start-1]!=0)
544 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start,
end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
546 if (iter <= maxIterations * n)
556 for (
Index i = 0; i < n-1; ++i)
559 diag.segment(i,n-i).minCoeff(&k);
562 numext::swap(diag[i], diag[k+i]);
563 if(computeEigenvectors)
564 eivec.col(i).swap(eivec.col(k+i));
571template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues
574 static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options)
575 { eig.compute(A,options); }
578template<
typename SolverType>
struct direct_selfadjoint_eigenvalues<SolverType,3,false>
580 typedef typename SolverType::MatrixType MatrixType;
581 typedef typename SolverType::RealVectorType VectorType;
582 typedef typename SolverType::Scalar Scalar;
583 typedef typename SolverType::EigenvectorsType EigenvectorsType;
591 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
593 EIGEN_USING_STD(
sqrt)
594 EIGEN_USING_STD(atan2)
597 const Scalar s_inv3 = Scalar(1)/Scalar(3);
598 const Scalar s_sqrt3 =
sqrt(Scalar(3));
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);
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));
613 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
615 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
616 q = numext::maxi(q, Scalar(0));
619 Scalar rho =
sqrt(a_over_3);
620 Scalar theta = atan2(
sqrt(q),half_b)*s_inv3;
621 Scalar cos_theta =
cos(theta);
622 Scalar sin_theta =
sin(theta);
624 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
625 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
626 roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
630 static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
632 EIGEN_USING_STD(
abs);
633 EIGEN_USING_STD(
sqrt);
636 mat.diagonal().cwiseAbs().maxCoeff(&i0);
639 representative = mat.col(i0);
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);
651 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
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");
659 EigenvectorsType& eivecs = solver.m_eivec;
660 VectorType& eivals = solver.m_eivalues;
663 Scalar shift = mat.trace() / Scalar(3);
665 MatrixType scaledMat = mat.template selfadjointView<Lower>();
666 scaledMat.diagonal().array() -= shift;
667 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
668 if(scale > 0) scaledMat /= scale;
671 computeRoots(scaledMat,eivals);
674 if(computeEigenvectors)
679 eivecs.setIdentity();
687 Scalar d0 = eivals(2) - eivals(1);
688 Scalar d1 = eivals(1) - eivals(0);
698 tmp.diagonal().array () -= eivals(k);
700 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
708 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
709 eivecs.col(l).normalize();
714 tmp.diagonal().array () -= eivals(l);
717 extract_kernel(tmp, eivecs.col(l), dummy);
721 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
727 eivals.array() += shift;
730 solver.m_isInitialized =
true;
731 solver.m_eigenvectorsOk = computeEigenvectors;
736template<
typename SolverType>
737struct direct_selfadjoint_eigenvalues<SolverType,2,false>
739 typedef typename SolverType::MatrixType MatrixType;
740 typedef typename SolverType::RealVectorType VectorType;
741 typedef typename SolverType::Scalar Scalar;
742 typedef typename SolverType::EigenvectorsType EigenvectorsType;
745 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
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));
755 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
757 EIGEN_USING_STD(
sqrt);
758 EIGEN_USING_STD(
abs);
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");
766 EigenvectorsType& eivecs = solver.m_eivec;
767 VectorType& eivals = solver.m_eivalues;
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))
779 computeRoots(scaledMat,eivals);
782 if(computeEigenvectors)
786 eivecs.setIdentity();
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));
796 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
797 eivecs.col(1) /=
sqrt(a2+b2);
801 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
802 eivecs.col(1) /=
sqrt(c2+b2);
805 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
811 eivals.array() += shift;
814 solver.m_isInitialized =
true;
815 solver.m_eigenvectorsOk = computeEigenvectors;
821template<
typename MatrixType>
826 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*
this,matrix,options);
833template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
835static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end, Scalar* matrixQ,
Index n)
838 RealScalar td = (diag[
end-1] - diag[
end])*RealScalar(0.5);
839 RealScalar e = subdiag[
end-1];
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);
854 mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
858 RealScalar x = diag[start] - mu;
859 RealScalar z = subdiag[start];
862 for (
Index k = start; k <
end && z != RealScalar(0); ++k)
864 JacobiRotation<RealScalar> rot;
865 rot.makeGivens(x, z);
868 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
869 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
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;
876 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
882 z = -rot.s() * subdiag[k+1];
883 subdiag[k + 1] = rot.c() * subdiag[k+1];
890 Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
891 q.applyOnTheRight(k,k+1,rot);
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