11#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
14#include "./Tridiagonalization.h"
17#include "./InternalHeaderCheck.h"
21template <
typename MatrixType_>
25template <
typename SolverType,
int Size,
bool IsComplex>
26struct direct_selfadjoint_eigenvalues;
28template <
typename MatrixType,
typename DiagType,
typename SubDiagType>
29EIGEN_DEVICE_FUNC
ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
30 const Index maxIterations,
bool computeEigenvectors,
81template <
typename MatrixType_>
84 typedef MatrixType_ MatrixType;
86 Size = MatrixType::RowsAtCompileTime,
87 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
88 Options = internal::traits<MatrixType>::Options,
89 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
93 typedef typename MatrixType::Scalar
Scalar;
113 typedef typename internal::plain_col_type<MatrixType, Scalar>::type
VectorType;
114 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
116 typedef typename TridiagonalizationType::SubDiagonalType SubDiagonalType;
134 m_info(InvalidInput),
135 m_isInitialized(false),
136 m_eigenvectorsOk(false) {}
151 : m_eivec(size, size),
154 m_subdiag(size > 1 ? size - 1 : 1),
155 m_hcoeffs(size > 1 ? size - 1 : 1),
156 m_isInitialized(false),
157 m_eigenvectorsOk(false) {}
174 template <
typename InputType>
177 : m_eivec(matrix.rows(), matrix.cols()),
178 m_workspace(matrix.cols()),
179 m_eivalues(matrix.cols()),
180 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
181 m_hcoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1),
182 m_isInitialized(false),
183 m_eigenvectorsOk(false) {
217 template <
typename InputType>
280 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
281 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
301 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
323 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
324 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
325 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
339 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
340 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
341 return m_eivec * m_eivalues.array().exp().matrix().asDiagonal() * m_eivec.adjoint();
363 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
364 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
365 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
373 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
385 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar)
387 EigenvectorsType m_eivec;
389 RealVectorType m_eivalues;
390 typename TridiagonalizationType::SubDiagonalType m_subdiag;
391 typename TridiagonalizationType::CoeffVectorType m_hcoeffs;
393 bool m_isInitialized;
394 bool m_eigenvectorsOk;
418template <
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
419EIGEN_DEVICE_FUNC
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end,
423template <
typename MatrixType>
424template <
typename InputType>
426 const EigenBase<InputType>& a_matrix,
int options) {
427 const InputType& matrix(a_matrix.derived());
429 EIGEN_USING_STD(
abs);
430 eigen_assert(matrix.cols() == matrix.rows());
431 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
432 "invalid option parameter");
434 Index n = matrix.cols();
435 m_eivalues.resize(n, 1);
439 m_eivalues.coeffRef(0, 0) = numext::real(m_eivec.coeff(0, 0));
440 if (computeEigenvectors) m_eivec.setOnes(n, n);
442 m_isInitialized =
true;
443 m_eigenvectorsOk = computeEigenvectors;
448 RealVectorType& diag = m_eivalues;
449 EigenvectorsType& mat = m_eivec;
452 mat = matrix.template triangularView<Lower>();
453 RealScalar scale = mat.cwiseAbs().maxCoeff();
454 if (numext::is_exactly_zero(scale)) scale = RealScalar(1);
455 mat.template triangularView<Lower>() /= scale;
456 m_subdiag.resize(n - 1);
457 m_hcoeffs.resize(n - 1);
458 internal::tridiagonalization_inplace(mat, diag, m_subdiag, m_hcoeffs, m_workspace, computeEigenvectors);
460 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
465 m_isInitialized =
true;
466 m_eigenvectorsOk = computeEigenvectors;
470template <
typename MatrixType>
472 const RealVectorType& diag,
const SubDiagonalType& subdiag,
int options) {
478 if (computeEigenvectors) {
479 m_eivec.setIdentity(diag.size(), diag.size());
481 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag,
m_maxIterations, computeEigenvectors, m_eivec);
483 m_isInitialized =
true;
484 m_eigenvectorsOk = computeEigenvectors;
500template <
typename MatrixType,
typename DiagType,
typename SubDiagType>
501EIGEN_DEVICE_FUNC
ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
502 const Index maxIterations,
bool computeEigenvectors,
505 typedef typename MatrixType::Scalar
Scalar;
507 Index n = diag.size();
512 typedef typename DiagType::RealScalar RealScalar;
513 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
514 const RealScalar precision_inv = RealScalar(1) / NumTraits<RealScalar>::epsilon();
516 for (
Index i = start; i < end; ++i) {
517 if (numext::abs(subdiag[i]) < considerAsZero) {
518 subdiag[i] = RealScalar(0);
522 const RealScalar scaled_subdiag = precision_inv * subdiag[i];
523 if (scaled_subdiag * scaled_subdiag <= (numext::abs(diag[i]) + numext::abs(diag[i + 1]))) {
524 subdiag[i] = RealScalar(0);
530 while (end > 0 && numext::is_exactly_zero(subdiag[end - 1])) {
537 if (iter > maxIterations * n)
break;
540 while (start > 0 && !numext::is_exactly_zero(subdiag[start - 1])) start--;
542 internal::tridiagonal_qr_step<MatrixType::Flags & RowMajorBit ? RowMajor : ColMajor>(
543 diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.
data() : (Scalar*)0, n);
545 if (iter <= maxIterations * n)
554 for (
Index i = 0; i < n - 1; ++i) {
556 diag.segment(i, n - i).minCoeff(&k);
558 numext::swap(diag[i], diag[k + i]);
559 if (computeEigenvectors) eivec.col(i).swap(eivec.col(k + i));
566template <
typename SolverType,
int Size,
bool IsComplex>
567struct direct_selfadjoint_eigenvalues {
568 EIGEN_DEVICE_FUNC
static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options) {
569 eig.compute(A, options);
573template <
typename SolverType>
574struct direct_selfadjoint_eigenvalues<SolverType, 3, false> {
575 typedef typename SolverType::MatrixType MatrixType;
576 typedef typename SolverType::RealVectorType VectorType;
577 typedef typename SolverType::Scalar Scalar;
578 typedef typename SolverType::EigenvectorsType EigenvectorsType;
584 EIGEN_DEVICE_FUNC
static inline void computeRoots(
const MatrixType& m, VectorType& roots) {
585 EIGEN_USING_STD(
sqrt)
586 EIGEN_USING_STD(
atan2)
589 const Scalar s_inv3 = Scalar(1) / Scalar(3);
590 const Scalar s_sqrt3 =
sqrt(Scalar(3));
595 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) -
596 m(1, 1) * m(2, 0) * m(2, 0) - m(2, 2) * m(1, 0) * m(1, 0);
597 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) -
599 Scalar c2 = m(0, 0) + m(1, 1) + m(2, 2);
603 Scalar c2_over_3 = c2 * s_inv3;
604 Scalar a_over_3 = (c2 * c2_over_3 - c1) * s_inv3;
605 a_over_3 = numext::maxi(a_over_3, Scalar(0));
607 Scalar half_b = Scalar(0.5) * (c0 + c2_over_3 * (Scalar(2) * c2_over_3 * c2_over_3 - c1));
609 Scalar q = a_over_3 * a_over_3 * a_over_3 - half_b * half_b;
610 q = numext::maxi(q, Scalar(0));
613 Scalar rho =
sqrt(a_over_3);
614 Scalar theta =
atan2(
sqrt(q), half_b) * s_inv3;
615 Scalar cos_theta =
cos(theta);
616 Scalar sin_theta =
sin(theta);
618 roots(0) = c2_over_3 - rho * (cos_theta + s_sqrt3 * sin_theta);
619 roots(1) = c2_over_3 - rho * (cos_theta - s_sqrt3 * sin_theta);
620 roots(2) = c2_over_3 + Scalar(2) * rho * cos_theta;
623 EIGEN_DEVICE_FUNC
static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res,
624 Ref<VectorType> representative) {
625 EIGEN_USING_STD(
abs);
626 EIGEN_USING_STD(
sqrt);
629 mat.diagonal().cwiseAbs().maxCoeff(&i0);
632 representative = mat.col(i0);
635 n0 = (c0 = representative.cross(mat.col((i0 + 1) % 3))).squaredNorm();
636 n1 = (c1 = representative.cross(mat.col((i0 + 2) % 3))).squaredNorm();
645 EIGEN_DEVICE_FUNC
static inline void run(SolverType& solver,
const MatrixType& mat,
int options) {
646 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
647 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
648 "invalid option parameter");
651 EigenvectorsType& eivecs = solver.m_eivec;
652 VectorType& eivals = solver.m_eivalues;
655 Scalar shift = mat.trace() / Scalar(3);
658 MatrixType scaledMat = mat.template selfadjointView<Lower>();
659 scaledMat.diagonal().array() -= shift;
660 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
661 if (scale > 0) scaledMat /= scale;
664 computeRoots(scaledMat, eivals);
667 if (computeEigenvectors) {
668 if ((eivals(2) - eivals(0)) <= Eigen::NumTraits<Scalar>::epsilon()) {
670 eivecs.setIdentity();
676 Scalar d0 = eivals(2) - eivals(1);
677 Scalar d1 = eivals(1) - eivals(0);
686 tmp.diagonal().array() -= eivals(k);
688 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
692 if (d0 <= 2 * Eigen::NumTraits<Scalar>::epsilon() * d1) {
695 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l)) * eivecs.col(l);
696 eivecs.col(l).normalize();
699 tmp.diagonal().array() -= eivals(l);
702 extract_kernel(tmp, eivecs.col(l), dummy);
706 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
712 eivals.array() += shift;
715 solver.m_isInitialized =
true;
716 solver.m_eigenvectorsOk = computeEigenvectors;
721template <
typename SolverType>
722struct direct_selfadjoint_eigenvalues<SolverType, 2, false> {
723 typedef typename SolverType::MatrixType MatrixType;
724 typedef typename SolverType::RealVectorType VectorType;
725 typedef typename SolverType::Scalar Scalar;
726 typedef typename SolverType::EigenvectorsType EigenvectorsType;
728 EIGEN_DEVICE_FUNC
static inline void computeRoots(
const MatrixType& m, VectorType& roots) {
729 EIGEN_USING_STD(
sqrt);
730 const Scalar t0 = Scalar(0.5) *
sqrt(numext::abs2(m(0, 0) - m(1, 1)) + Scalar(4) * numext::abs2(m(1, 0)));
731 const Scalar t1 = Scalar(0.5) * (m(0, 0) + m(1, 1));
736 EIGEN_DEVICE_FUNC
static inline void run(SolverType& solver,
const MatrixType& mat,
int options) {
737 EIGEN_USING_STD(
sqrt);
738 EIGEN_USING_STD(
abs);
740 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
741 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
742 "invalid option parameter");
745 EigenvectorsType& eivecs = solver.m_eivec;
746 VectorType& eivals = solver.m_eivalues;
749 Scalar shift = mat.trace() / Scalar(2);
750 MatrixType scaledMat = mat;
751 scaledMat.coeffRef(0, 1) = mat.coeff(1, 0);
752 scaledMat.diagonal().array() -= shift;
753 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
754 if (scale > Scalar(0)) scaledMat /= scale;
757 computeRoots(scaledMat, eivals);
760 if (computeEigenvectors) {
761 if ((eivals(1) - eivals(0)) <=
abs(eivals(1)) * Eigen::NumTraits<Scalar>::epsilon()) {
762 eivecs.setIdentity();
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));
769 eivecs.col(1) << -scaledMat(1, 0), scaledMat(0, 0);
770 eivecs.col(1) /=
sqrt(a2 + b2);
772 eivecs.col(1) << -scaledMat(1, 1), scaledMat(1, 0);
773 eivecs.col(1) /=
sqrt(c2 + b2);
776 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
782 eivals.array() += shift;
785 solver.m_isInitialized =
true;
786 solver.m_eigenvectorsOk = computeEigenvectors;
792template <
typename MatrixType>
794 const MatrixType& matrix,
int options) {
795 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver, Size, NumTraits<Scalar>::IsComplex>::run(
796 *
this, matrix, options);
803template <
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
804EIGEN_DEVICE_FUNC
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end,
807 RealScalar td = (diag[end - 1] - diag[end]) * RealScalar(0.5);
808 RealScalar e = subdiag[end - 1];
814 RealScalar mu = diag[end];
815 if (numext::is_exactly_zero(td)) {
816 mu -= numext::abs(e);
817 }
else if (!numext::is_exactly_zero(e)) {
818 const RealScalar e2 = numext::abs2(e);
819 const RealScalar h = numext::hypot(td, e);
820 if (numext::is_exactly_zero(e2)) {
821 mu -= e / ((td + (td > RealScalar(0) ? h : -h)) / e);
823 mu -= e2 / (td + (td > RealScalar(0) ? h : -h));
827 RealScalar x = diag[start] - mu;
828 RealScalar z = subdiag[start];
831 for (
Index k = start; k < end && !numext::is_exactly_zero(z); ++k) {
832 JacobiRotation<RealScalar> rot;
833 rot.makeGivens(x, z);
836 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
837 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k + 1];
840 rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k + 1]);
841 diag[k + 1] = rot.s() * sdk + rot.c() * dkp1;
842 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
844 if (k > start) subdiag[k - 1] = rot.c() * subdiag[k - 1] - rot.s() * z;
849 z = -rot.s() * subdiag[k + 1];
850 subdiag[k + 1] = rot.c() * subdiag[k + 1];
856 Map<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > q(matrixQ, n, n);
857 q.applyOnTheRight(k, k + 1, rot);
const std::enable_if_t< std::is_same< typename LhsDerived::Scalar, typename RhsDerived::Scalar >::value, Eigen::CwiseBinaryOp< Eigen::internal::scalar_atan2_op< typename LhsDerived::Scalar, typename RhsDerived::Scalar >, const LhsDerived, const RhsDerived > > atan2(const Eigen::ArrayBase< LhsDerived > &x, const Eigen::ArrayBase< RhsDerived > &exponents)
Definition GlobalFunctions.h:209
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition GeneralizedSelfAdjointEigenSolver.h:51
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
MatrixType operatorExp() const
Computes the matrix exponential the matrix.
Definition SelfAdjointEigenSolver.h:338
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition SelfAdjointEigenSolver.h:471
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType_.
Definition SelfAdjointEigenSolver.h:104
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition SelfAdjointEigenSolver.h:362
SelfAdjointEigenSolver()
Default constructor for fixed-size matrices.
Definition SelfAdjointEigenSolver.h:128
internal::plain_col_type< MatrixType, Scalar >::type VectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition SelfAdjointEigenSolver.h:113
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SelfAdjointEigenSolver.h:372
Eigen::Index Index
Definition SelfAdjointEigenSolver.h:94
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
Definition SelfAdjointEigenSolver.h:93
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition SelfAdjointEigenSolver.h:322
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition SelfAdjointEigenSolver.h:150
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition SelfAdjointEigenSolver.h:300
static const int m_maxIterations
Maximum number of iterations.
Definition SelfAdjointEigenSolver.h:382
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:793
Tridiagonal decomposition of a selfadjoint matrix.
Definition Tridiagonalization.h:66
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
@ ComputeEigenvectors
Definition Constants.h:401
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:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33
constexpr Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:232