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;
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();
347 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
348 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
349 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
357 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
369 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar)
373 RealVectorType m_eivalues;
377 bool m_isInitialized;
378 bool m_eigenvectorsOk;
402template <
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
403EIGEN_DEVICE_FUNC
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end,
404 Scalar* matrixQ,
Index n);
407template <
typename MatrixType>
408template <
typename InputType>
409EIGEN_DEVICE_FUNC SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>::compute(
410 const EigenBase<InputType>& a_matrix,
int options) {
411 const InputType& matrix(a_matrix.derived());
413 EIGEN_USING_STD(abs);
414 eigen_assert(matrix.cols() == matrix.rows());
415 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
416 "invalid option parameter");
417 bool computeEigenvectors = (options & ComputeEigenvectors) == ComputeEigenvectors;
418 Index n = matrix.cols();
419 m_eivalues.resize(n, 1);
423 m_eivalues.coeffRef(0, 0) = numext::real(m_eivec.coeff(0, 0));
424 if (computeEigenvectors) m_eivec.setOnes(n, n);
426 m_isInitialized =
true;
427 m_eigenvectorsOk = computeEigenvectors;
432 RealVectorType& diag = m_eivalues;
433 EigenvectorsType& mat = m_eivec;
436 mat = matrix.template triangularView<Lower>();
437 RealScalar scale = mat.cwiseAbs().maxCoeff();
438 if (numext::is_exactly_zero(scale)) scale = RealScalar(1);
439 mat.template triangularView<Lower>() /= scale;
440 m_subdiag.resize(n - 1);
441 m_hcoeffs.resize(n - 1);
442 internal::tridiagonalization_inplace(mat, diag, m_subdiag, m_hcoeffs, m_workspace, computeEigenvectors);
444 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
449 m_isInitialized =
true;
450 m_eigenvectorsOk = computeEigenvectors;
454template <
typename MatrixType>
456 const RealVectorType& diag,
const SubDiagonalType& subdiag,
int options) {
462 if (computeEigenvectors) {
465 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
467 m_isInitialized =
true;
468 m_eigenvectorsOk = computeEigenvectors;
484template <
typename MatrixType,
typename DiagType,
typename SubDiagType>
485EIGEN_DEVICE_FUNC
ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
486 const Index maxIterations,
bool computeEigenvectors,
489 typedef typename MatrixType::Scalar Scalar;
491 Index n = diag.size();
496 typedef typename DiagType::RealScalar RealScalar;
497 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
500 for (
Index i = start; i < end; ++i) {
501 if (numext::abs(subdiag[i]) < considerAsZero) {
502 subdiag[i] = RealScalar(0);
506 const RealScalar scaled_subdiag = precision_inv * subdiag[i];
507 if (scaled_subdiag * scaled_subdiag <= (numext::abs(diag[i]) + numext::abs(diag[i + 1]))) {
508 subdiag[i] = RealScalar(0);
514 while (end > 0 && numext::is_exactly_zero(subdiag[end - 1])) {
521 if (iter > maxIterations * n)
break;
524 while (start > 0 && !numext::is_exactly_zero(subdiag[start - 1])) start--;
526 internal::tridiagonal_qr_step<MatrixType::Flags & RowMajorBit ? RowMajor : ColMajor>(
527 diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
529 if (iter <= maxIterations * n)
538 for (Index i = 0; i < n - 1; ++i) {
540 diag.segment(i, n - i).minCoeff(&k);
542 numext::swap(diag[i], diag[k + i]);
543 if (computeEigenvectors) eivec.col(i).swap(eivec.col(k + i));
550template <
typename SolverType,
int Size,
bool IsComplex>
551struct direct_selfadjoint_eigenvalues {
552 EIGEN_DEVICE_FUNC
static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options) {
553 eig.compute(A, options);
557template <
typename SolverType>
558struct direct_selfadjoint_eigenvalues<SolverType, 3, false> {
559 typedef typename SolverType::MatrixType MatrixType;
560 typedef typename SolverType::RealVectorType VectorType;
561 typedef typename SolverType::Scalar Scalar;
562 typedef typename SolverType::EigenvectorsType EigenvectorsType;
568 EIGEN_DEVICE_FUNC
static inline void computeRoots(
const MatrixType& m, VectorType& roots) {
569 EIGEN_USING_STD(
sqrt)
570 EIGEN_USING_STD(atan2)
573 const Scalar s_inv3 = Scalar(1) / Scalar(3);
574 const Scalar s_sqrt3 =
sqrt(Scalar(3));
579 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) -
580 m(1, 1) * m(2, 0) * m(2, 0) - m(2, 2) * m(1, 0) * m(1, 0);
581 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) -
583 Scalar c2 = m(0, 0) + m(1, 1) + m(2, 2);
587 Scalar c2_over_3 = c2 * s_inv3;
588 Scalar a_over_3 = (c2 * c2_over_3 - c1) * s_inv3;
589 a_over_3 = numext::maxi(a_over_3, Scalar(0));
591 Scalar half_b = Scalar(0.5) * (c0 + c2_over_3 * (Scalar(2) * c2_over_3 * c2_over_3 - c1));
593 Scalar q = a_over_3 * a_over_3 * a_over_3 - half_b * half_b;
594 q = numext::maxi(q, Scalar(0));
597 Scalar rho =
sqrt(a_over_3);
598 Scalar theta = atan2(
sqrt(q), half_b) * s_inv3;
599 Scalar cos_theta =
cos(theta);
600 Scalar sin_theta =
sin(theta);
602 roots(0) = c2_over_3 - rho * (cos_theta + s_sqrt3 * sin_theta);
603 roots(1) = c2_over_3 - rho * (cos_theta - s_sqrt3 * sin_theta);
604 roots(2) = c2_over_3 + Scalar(2) * rho * cos_theta;
607 EIGEN_DEVICE_FUNC
static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res,
608 Ref<VectorType> representative) {
609 EIGEN_USING_STD(
abs);
610 EIGEN_USING_STD(
sqrt);
613 mat.diagonal().cwiseAbs().maxCoeff(&i0);
616 representative = mat.col(i0);
619 n0 = (c0 = representative.cross(mat.col((i0 + 1) % 3))).squaredNorm();
620 n1 = (c1 = representative.cross(mat.col((i0 + 2) % 3))).squaredNorm();
629 EIGEN_DEVICE_FUNC
static inline void run(SolverType& solver,
const MatrixType& mat,
int options) {
630 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
631 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
632 "invalid option parameter");
635 EigenvectorsType& eivecs = solver.m_eivec;
636 VectorType& eivals = solver.m_eivalues;
639 Scalar shift = mat.trace() / Scalar(3);
642 MatrixType scaledMat = mat.template selfadjointView<Lower>();
643 scaledMat.diagonal().array() -= shift;
644 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
645 if (scale > 0) scaledMat /= scale;
648 computeRoots(scaledMat, eivals);
651 if (computeEigenvectors) {
654 eivecs.setIdentity();
660 Scalar d0 = eivals(2) - eivals(1);
661 Scalar d1 = eivals(1) - eivals(0);
670 tmp.diagonal().array() -= eivals(k);
672 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
679 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l)) * eivecs.col(l);
680 eivecs.col(l).normalize();
683 tmp.diagonal().array() -= eivals(l);
686 extract_kernel(tmp, eivecs.col(l), dummy);
690 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
696 eivals.array() += shift;
699 solver.m_isInitialized =
true;
700 solver.m_eigenvectorsOk = computeEigenvectors;
705template <
typename SolverType>
706struct direct_selfadjoint_eigenvalues<SolverType, 2, false> {
707 typedef typename SolverType::MatrixType MatrixType;
708 typedef typename SolverType::RealVectorType VectorType;
709 typedef typename SolverType::Scalar Scalar;
710 typedef typename SolverType::EigenvectorsType EigenvectorsType;
712 EIGEN_DEVICE_FUNC
static inline void computeRoots(
const MatrixType& m, VectorType& roots) {
713 EIGEN_USING_STD(
sqrt);
714 const Scalar t0 = Scalar(0.5) *
sqrt(numext::abs2(m(0, 0) - m(1, 1)) + Scalar(4) * numext::abs2(m(1, 0)));
715 const Scalar t1 = Scalar(0.5) * (m(0, 0) + m(1, 1));
720 EIGEN_DEVICE_FUNC
static inline void run(SolverType& solver,
const MatrixType& mat,
int options) {
721 EIGEN_USING_STD(
sqrt);
722 EIGEN_USING_STD(
abs);
724 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
725 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
726 "invalid option parameter");
729 EigenvectorsType& eivecs = solver.m_eivec;
730 VectorType& eivals = solver.m_eivalues;
733 Scalar shift = mat.trace() / Scalar(2);
734 MatrixType scaledMat = mat;
735 scaledMat.coeffRef(0, 1) = mat.coeff(1, 0);
736 scaledMat.diagonal().array() -= shift;
737 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
738 if (scale > Scalar(0)) scaledMat /= scale;
741 computeRoots(scaledMat, eivals);
744 if (computeEigenvectors) {
746 eivecs.setIdentity();
748 scaledMat.diagonal().array() -= eivals(1);
749 Scalar a2 = numext::abs2(scaledMat(0, 0));
750 Scalar c2 = numext::abs2(scaledMat(1, 1));
751 Scalar b2 = numext::abs2(scaledMat(1, 0));
753 eivecs.col(1) << -scaledMat(1, 0), scaledMat(0, 0);
754 eivecs.col(1) /=
sqrt(a2 + b2);
756 eivecs.col(1) << -scaledMat(1, 1), scaledMat(1, 0);
757 eivecs.col(1) /=
sqrt(c2 + b2);
760 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
766 eivals.array() += shift;
769 solver.m_isInitialized =
true;
770 solver.m_eigenvectorsOk = computeEigenvectors;
776template <
typename MatrixType>
778 const MatrixType& matrix,
int options) {
779 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver, Size, NumTraits<Scalar>::IsComplex>::run(
780 *
this, matrix, options);
787template <
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
788EIGEN_DEVICE_FUNC
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end,
789 Scalar* matrixQ,
Index n) {
791 RealScalar td = (diag[end - 1] - diag[end]) * RealScalar(0.5);
792 RealScalar e = subdiag[end - 1];
798 RealScalar mu = diag[end];
799 if (numext::is_exactly_zero(td)) {
800 mu -= numext::abs(e);
801 }
else if (!numext::is_exactly_zero(e)) {
802 const RealScalar e2 = numext::abs2(e);
803 const RealScalar h = numext::hypot(td, e);
804 if (numext::is_exactly_zero(e2)) {
805 mu -= e / ((td + (td > RealScalar(0) ? h : -h)) / e);
807 mu -= e2 / (td + (td > RealScalar(0) ? h : -h));
811 RealScalar x = diag[start] - mu;
812 RealScalar z = subdiag[start];
815 for (Index k = start; k < end && !numext::is_exactly_zero(z); ++k) {
816 JacobiRotation<RealScalar> rot;
817 rot.makeGivens(x, z);
820 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
821 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k + 1];
824 rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k + 1]);
825 diag[k + 1] = rot.s() * sdk + rot.c() * dkp1;
826 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
828 if (k > start) subdiag[k - 1] = rot.c() * subdiag[k - 1] - rot.s() * z;
833 z = -rot.s() * subdiag[k + 1];
834 subdiag[k + 1] = rot.c() * subdiag[k + 1];
840 Map<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > q(matrixQ, n, n);
841 q.applyOnTheRight(k, k + 1, rot);
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition SelfAdjointEigenSolver.h:22
Derived & setIdentity()
Definition CwiseNullaryOp.h:839
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition SelfAdjointEigenSolver.h:82
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:455
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:346
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:356
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:366
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:777
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 Core:137
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:83
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33
Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523