11#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
14#include "./Tridiagonalization.h"
18template<
typename _MatrixType>
19class GeneralizedSelfAdjointEigenSolver;
22template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
23template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
24ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
const Index maxIterations,
bool computeEigenvectors, MatrixType& eivec);
74 typedef _MatrixType MatrixType;
76 Size = MatrixType::RowsAtCompileTime,
77 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
78 Options = MatrixType::Options,
79 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
83 typedef typename MatrixType::Scalar
Scalar;
122 m_isInitialized(false)
139 : m_eivec(size, size),
141 m_subdiag(size > 1 ? size - 1 : 1),
142 m_isInitialized(
false)
160 template<
typename InputType>
163 : m_eivec(matrix.rows(), matrix.cols()),
164 m_eivalues(matrix.cols()),
165 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
166 m_isInitialized(
false)
168 compute(matrix.derived(), options);
201 template<
typename InputType>
260 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
261 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
283 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
308 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
309 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
310 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
334 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
335 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
336 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
346 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
358 static void check_template_parameters()
360 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar);
363 EigenvectorsType m_eivec;
365 typename TridiagonalizationType::SubDiagonalType m_subdiag;
367 bool m_isInitialized;
368 bool m_eigenvectorsOk;
388template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
390static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
393template<
typename MatrixType>
394template<
typename InputType>
396SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
397::compute(
const EigenBase<InputType>& a_matrix,
int options)
399 check_template_parameters();
401 const InputType &matrix(a_matrix.derived());
404 eigen_assert(matrix.cols() == matrix.rows());
405 eigen_assert((options&~(EigVecMask|GenEigMask))==0
406 && (options&EigVecMask)!=EigVecMask
407 &&
"invalid option parameter");
409 Index n = matrix.cols();
410 m_eivalues.resize(n,1);
414 m_eivalues.coeffRef(0,0) = numext::real(matrix(0,0));
415 if(computeEigenvectors)
416 m_eivec.setOnes(n,n);
418 m_isInitialized =
true;
419 m_eigenvectorsOk = computeEigenvectors;
424 RealVectorType& diag = m_eivalues;
425 EigenvectorsType& mat = m_eivec;
428 mat = matrix.template triangularView<Lower>();
429 RealScalar scale = mat.cwiseAbs().maxCoeff();
430 if(scale==RealScalar(0)) scale = RealScalar(1);
431 mat.template triangularView<Lower>() /= scale;
432 m_subdiag.resize(n-1);
433 internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
435 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
440 m_isInitialized =
true;
441 m_eigenvectorsOk = computeEigenvectors;
445template<
typename MatrixType>
446SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
456 m_eivec.setIdentity(
diag.size(),
diag.size());
458 m_info = computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations,
computeEigenvectors, m_eivec);
460 m_isInitialized =
true;
476template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
482 typedef typename MatrixType::Scalar Scalar;
484 Index n =
diag.size();
489 typedef typename DiagType::RealScalar RealScalar;
490 const RealScalar
considerAsZero = (std::numeric_limits<RealScalar>::min)();
494 for (Index i = start; i<end; ++i)
499 while (end>0 &&
subdiag[end-1]==0)
508 if(iter > maxIterations * n)
break;
511 while (start>0 && subdiag[start-1]!=0)
514 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
516 if (iter <= maxIterations * n)
526 for (Index i = 0; i < n-1; ++i)
529 diag.segment(i,n-i).minCoeff(&k);
533 if(computeEigenvectors)
534 eivec.col(i).swap(eivec.col(k+i));
544 static inline void run(
SolverType&
eig,
const typename SolverType::MatrixType& A,
int options)
545 {
eig.compute(A,options); }
550 typedef typename SolverType::MatrixType MatrixType;
551 typedef typename SolverType::RealVectorType VectorType;
552 typedef typename SolverType::Scalar Scalar;
553 typedef typename SolverType::EigenvectorsType EigenvectorsType;
561 static inline void computeRoots(
const MatrixType& m, VectorType&
roots)
567 const Scalar
s_inv3 = Scalar(1.0)/Scalar(3.0);
568 const Scalar
s_sqrt3 = sqrt(Scalar(3.0));
573 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);
574 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);
575 Scalar
c2 = m(0,0) + m(1,1) + m(2,2);
586 q = numext::maxi(
q, Scalar(0));
605 mat.diagonal().cwiseAbs().maxCoeff(&
i0);
614 else res =
c1/std::sqrt(
n1);
622 eigen_assert(
mat.cols() == 3 &&
mat.cols() ==
mat.rows());
623 eigen_assert((options&~(EigVecMask|GenEigMask))==0
624 && (options&EigVecMask)!=EigVecMask
625 &&
"invalid option parameter");
632 Scalar shift =
mat.trace() / Scalar(3);
636 Scalar scale =
scaledMat.cwiseAbs().maxCoeff();
667 tmp.diagonal().array () -=
eivals(k);
678 eivecs.col(l).normalize();
683 tmp.diagonal().array () -=
eivals(l);
686 extract_kernel(tmp,
eivecs.col(l), dummy);
699 solver.m_isInitialized =
true;
705template<
typename SolverType>
708 typedef typename SolverType::MatrixType MatrixType;
709 typedef typename SolverType::RealVectorType VectorType;
710 typedef typename SolverType::Scalar Scalar;
711 typedef typename SolverType::EigenvectorsType EigenvectorsType;
714 static inline void computeRoots(
const MatrixType& m, VectorType&
roots)
717 const Scalar
t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
718 const Scalar
t1 = Scalar(0.5) * (m(0,0) + m(1,1));
729 eigen_assert(
mat.cols() == 2 &&
mat.cols() ==
mat.rows());
730 eigen_assert((options&~(EigVecMask|GenEigMask))==0
731 && (options&EigVecMask)!=EigVecMask
732 &&
"invalid option parameter");
739 Scalar scale =
mat.cwiseAbs().maxCoeff();
740 scale = numext::maxi(scale,Scalar(1));
778 solver.m_isInitialized =
true;
785template<
typename MatrixType>
795template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
797static void tridiagonal_qr_step(RealScalar*
diag, RealScalar*
subdiag, Index start, Index end, Scalar* matrixQ, Index n)
800 RealScalar
td = (
diag[end-1] -
diag[end])*RealScalar(0.5);
807 RealScalar
mu =
diag[end];
812 RealScalar
e2 = numext::abs2(
subdiag[end-1]);
813 RealScalar h = numext::hypot(
td,e);
814 if(
e2==0)
mu -= (e / (
td + (
td>0 ? 1 : -1))) * (e / h);
815 else mu -=
e2 / (
td + (
td>0 ? h : -h));
818 RealScalar x = diag[start] - mu;
819 RealScalar z = subdiag[start];
820 for (Index k = start; k < end; ++k)
822 JacobiRotation<RealScalar> rot;
823 rot.makeGivens(x, z);
826 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
827 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
829 diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
830 diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
831 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
835 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
841 z = -rot.s() * subdiag[k+1];
842 subdiag[k + 1] = rot.c() * subdiag[k+1];
849 Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
850 q.applyOnTheRight(k,k+1,rot);
\eigenvalues_module
Definition SelfAdjointEigenSolver.h:71
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition SelfAdjointEigenSolver.h:83
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition SelfAdjointEigenSolver.h:447
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Definition SelfAdjointEigenSolver.h:788
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
Definition SelfAdjointEigenSolver.h:344
EIGEN_DEVICE_FUNC MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition SelfAdjointEigenSolver.h:306
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
Definition SelfAdjointEigenSolver.h:94
EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition SelfAdjointEigenSolver.h:332
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition SelfAdjointEigenSolver.h:103
Eigen::Index Index
Definition SelfAdjointEigenSolver.h:84
static const int m_maxIterations
Maximum number of iterations.
Definition SelfAdjointEigenSolver.h:355
EIGEN_DEVICE_FUNC const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition SelfAdjointEigenSolver.h:258
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Definition SelfAdjointEigenSolver.h:162
EIGEN_DEVICE_FUNC const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition SelfAdjointEigenSolver.h:281
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition SelfAdjointEigenSolver.h:138
Pseudo expression representing a solving operation.
Definition Solve.h:63
\eigenvalues_module
Definition Tridiagonalization.h:64
ComputationInfo
Enum for reporting the status of a computation.
Definition Constants.h:430
@ Success
Computation was successful.
Definition Constants.h:432
@ NoConvergence
Iterative procedure did not converge.
Definition Constants.h:436
@ ComputeEigenvectors
Used in SelfAdjointEigenSolver and GeneralizedSelfAdjointEigenSolver to specify that both the eigenva...
Definition Constants.h:395
NLOHMANN_BASIC_JSON_TPL_DECLARATION void swap(nlohmann::NLOHMANN_BASIC_JSON_TPL &j1, nlohmann::NLOHMANN_BASIC_JSON_TPL &j2) noexcept(//NOLINT(readability-inconsistent-declaration-parameter-name, cert-dcl58-cpp) is_nothrow_move_constructible< nlohmann::NLOHMANN_BASIC_JSON_TPL >::value &&//NOLINT(misc-redundant-expression) is_nothrow_move_assignable< nlohmann::NLOHMANN_BASIC_JSON_TPL >::value)
exchanges the values of two JSON objects
Definition json.hpp:24418
Holds information about the various numeric (i.e.
Definition NumTraits.h:108
Definition SelfAdjointEigenSolver.h:542
Definition XprHelper.h:597