11#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
18template<
typename _MatrixType>
19class GeneralizedSelfAdjointEigenSolver;
22template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
24template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
82 Size = MatrixType::RowsAtCompileTime,
89 typedef typename MatrixType::Scalar
Scalar;
171 template<
typename InputType>
174 :
m_eivec(matrix.rows(), matrix.cols()),
176 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
177 m_hcoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1),
214 template<
typename InputType>
411template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
416template<
typename MatrixType>
417template<
typename InputType>
422 check_template_parameters();
424 const InputType &matrix(a_matrix.
derived());
430 &&
"invalid option parameter");
432 Index n = matrix.cols();
433 m_eivalues.resize(n,1);
438 m_eivalues.coeffRef(0,0) =
numext::real(m_eivec.coeff(0,0));
439 if(computeEigenvectors)
440 m_eivec.setOnes(n,n);
442 m_isInitialized =
true;
443 m_eigenvectorsOk = computeEigenvectors;
452 mat = matrix.template triangularView<Lower>();
455 mat.template triangularView<Lower>() /= scale;
457 m_hcoeffs.resize(n-1);
465 m_isInitialized =
true;
466 m_eigenvectorsOk = computeEigenvectors;
470template<
typename MatrixType>
479 if (computeEigenvectors)
481 m_eivec.setIdentity(diag.size(), diag.size());
485 m_isInitialized =
true;
486 m_eigenvectorsOk = computeEigenvectors;
502template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
507 typedef typename MatrixType::Scalar Scalar;
509 Index n = diag.size();
514 typedef typename DiagType::RealScalar RealScalar;
521 subdiag[i] = RealScalar(0);
525 const RealScalar scaled_subdiag = precision_inv * subdiag[i];
527 subdiag[i] = RealScalar(0);
533 while (
end>0 && subdiag[
end-1]==RealScalar(0))
542 if(iter > maxIterations * n)
break;
545 while (start>0 && subdiag[start-1]!=0)
548 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start,
end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
550 if (iter <= maxIterations * n)
560 for (
Index i = 0; i < n-1; ++i)
563 diag.segment(i,n-i).minCoeff(&k);
567 if(computeEigenvectors)
568 eivec.col(i).swap(eivec.col(k+i));
578 static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options)
579 { eig.compute(A,options); }
586 typedef typename SolverType::Scalar
Scalar;
607 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);
608 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);
609 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
613 Scalar c2_over_3 = c2*s_inv3;
614 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
619 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
628 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
629 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
630 roots(2) = c2_over_3 +
Scalar(2)*rho*cos_theta;
640 mat.diagonal().cwiseAbs().maxCoeff(&i0);
643 representative = mat.col(i0);
646 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
647 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
648 if(n0>n1) res = c0/
sqrt(n0);
649 else res = c1/
sqrt(n1);
655 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
657 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
660 &&
"invalid option parameter");
669 MatrixType scaledMat = mat.template selfadjointView<Lower>();
670 scaledMat.diagonal().array() -= shift;
671 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
672 if(scale > 0) scaledMat /= scale;
675 computeRoots(scaledMat,eivals);
678 if(computeEigenvectors)
683 eivecs.setIdentity();
691 Scalar d0 = eivals(2) - eivals(1);
692 Scalar d1 = eivals(1) - eivals(0);
702 tmp.diagonal().array () -= eivals(k);
704 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
712 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
713 eivecs.col(l).normalize();
718 tmp.diagonal().array () -= eivals(l);
721 extract_kernel(tmp, eivecs.col(l), dummy);
725 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
731 eivals.array() += shift;
734 solver.m_isInitialized =
true;
735 solver.m_eigenvectorsOk = computeEigenvectors;
740template<
typename SolverType>
745 typedef typename SolverType::Scalar
Scalar;
759 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
764 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
767 &&
"invalid option parameter");
776 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
777 scaledMat.diagonal().array() -= shift;
778 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
783 computeRoots(scaledMat,eivals);
786 if(computeEigenvectors)
790 eivecs.setIdentity();
794 scaledMat.diagonal().array () -= eivals(1);
800 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
801 eivecs.col(1) /=
sqrt(a2+b2);
805 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
806 eivecs.col(1) /=
sqrt(c2+b2);
809 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
815 eivals.array() += shift;
818 solver.m_isInitialized =
true;
819 solver.m_eigenvectorsOk = computeEigenvectors;
825template<
typename MatrixType>
837template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
842 RealScalar td = (diag[
end-1] - diag[
end])*RealScalar(0.5);
843 RealScalar
e = subdiag[
end-1];
849 RealScalar mu = diag[
end];
850 if(td==RealScalar(0)) {
852 }
else if (
e != RealScalar(0)) {
855 if(e2 == RealScalar(0)) {
856 mu -=
e / ((td + (td>RealScalar(0) ?
h : -
h)) /
e);
858 mu -= e2 / (td + (td>RealScalar(0) ?
h : -
h));
862 RealScalar x = diag[start] - mu;
863 RealScalar z = subdiag[start];
866 for (
Index k = start; k <
end && z != RealScalar(0); ++k)
872 RealScalar sdk = rot.
s() * diag[k] + rot.
c() * subdiag[k];
873 RealScalar dkp1 = rot.
s() * subdiag[k] + rot.
c() * diag[k+1];
875 diag[k] = rot.
c() * (rot.
c() * diag[k] - rot.
s() * subdiag[k]) - rot.
s() * (rot.
c() * subdiag[k] - rot.
s() * diag[k+1]);
876 diag[k+1] = rot.
s() * sdk + rot.
c() * dkp1;
877 subdiag[k] = rot.
c() * sdk - rot.
s() * dkp1;
880 subdiag[k - 1] = rot.
c() * subdiag[k-1] - rot.
s() * z;
886 z = -rot.
s() * subdiag[k+1];
887 subdiag[k + 1] = rot.
c() * subdiag[k+1];
895 q.applyOnTheRight(k,k+1,rot);
EIGEN_DEVICE_FUNC RealReturnType real() const
Definition: CommonCwiseUnaryOps.h:100
#define EIGEN_USING_STD(FUNC)
Definition: Macros.h:1195
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:986
#define eigen_assert(x)
Definition: Macros.h:1047
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:187
constexpr common_return_t< T1, T2 > atan2(const T1 y, const T2 x) noexcept
Compile-time two-argument arctangent function.
Definition: atan2.hpp:82
\jacobi_module
Definition: Jacobi.h:35
EIGEN_DEVICE_FUNC void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Makes *this as a Givens rotation G such that applying to the left of the vector yields: .
Definition: Jacobi.h:162
EIGEN_DEVICE_FUNC Scalar & s()
Definition: Jacobi.h:49
EIGEN_DEVICE_FUNC Scalar & c()
Definition: Jacobi.h:47
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:271
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:283
\eigenvalues_module
Definition: SelfAdjointEigenSolver.h:77
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition: SelfAdjointEigenSolver.h:89
Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
Definition: SelfAdjointEigenSolver.h:92
RealVectorType m_eivalues
Definition: SelfAdjointEigenSolver.h:382
bool m_eigenvectorsOk
Definition: SelfAdjointEigenSolver.h:387
TridiagonalizationType::CoeffVectorType m_hcoeffs
Definition: SelfAdjointEigenSolver.h:384
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition: SelfAdjointEigenSolver.h:472
@ Size
Definition: SelfAdjointEigenSolver.h:82
@ Options
Definition: SelfAdjointEigenSolver.h:84
@ ColsAtCompileTime
Definition: SelfAdjointEigenSolver.h:83
@ MaxColsAtCompileTime
Definition: SelfAdjointEigenSolver.h:85
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Definition: SelfAdjointEigenSolver.h:828
EigenvectorsType m_eivec
Definition: SelfAdjointEigenSolver.h:381
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SelfAdjointEigenSolver.h:361
EIGEN_DEVICE_FUNC MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition: SelfAdjointEigenSolver.h:324
TridiagonalizationType::SubDiagonalType m_subdiag
Definition: SelfAdjointEigenSolver.h:383
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
Definition: SelfAdjointEigenSolver.h:100
EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition: SelfAdjointEigenSolver.h:349
bool m_isInitialized
Definition: SelfAdjointEigenSolver.h:386
ComputationInfo m_info
Definition: SelfAdjointEigenSolver.h:385
Eigen::Index Index
Definition: SelfAdjointEigenSolver.h:90
static const int m_maxIterations
Maximum number of iterations.
Definition: SelfAdjointEigenSolver.h:372
EIGEN_DEVICE_FUNC const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition: SelfAdjointEigenSolver.h:277
_MatrixType MatrixType
Definition: SelfAdjointEigenSolver.h:80
static EIGEN_DEVICE_FUNC void check_template_parameters()
Definition: SelfAdjointEigenSolver.h:376
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:173
EIGEN_DEVICE_FUNC const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition: SelfAdjointEigenSolver.h:300
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: SelfAdjointEigenSolver.h:147
\eigenvalues_module
Definition: Tridiagonalization.h:65
type
Definition: core.h:575
UnitType abs(const UnitType x) noexcept
Compute absolute value.
Definition: math.h:721
dimensionless::scalar_t cos(const AngleUnit angle) noexcept
Compute cosine.
Definition: math.h:61
auto sqrt(const UnitType &value) noexcept -> unit_t< square_root< typename units::traits::unit_t_traits< UnitType >::unit_type >, typename units::traits::unit_t_traits< UnitType >::underlying_type, linear_scale >
computes the square root of value
Definition: math.h:483
dimensionless::scalar_t sin(const AngleUnit angle) noexcept
Compute sine.
Definition: math.h:81
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:440
@ Success
Computation was successful.
Definition: Constants.h:442
@ NoConvergence
Iterative procedure did not converge.
Definition: Constants.h:446
@ GenEigMask
Definition: Constants.h:418
@ EigVecMask
Definition: Constants.h:407
@ ComputeEigenvectors
Used in SelfAdjointEigenSolver and GeneralizedSelfAdjointEigenSolver to specify that both the eigenva...
Definition: Constants.h:405
constexpr common_return_t< T1, T2 > hypot(const T1 x, const T2 y) noexcept
Compile-time Pythagorean addition function.
Definition: hypot.hpp:84
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
static EIGEN_DEVICE_FUNC void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
Definition: SelfAdjointEigenSolver.h:839
EIGEN_DEVICE_FUNC void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
Definition: Tridiagonalization.h:349
EIGEN_CONSTEXPR Index size(const T &x)
Definition: Meta.h:479
EIGEN_DEVICE_FUNC ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
Definition: SelfAdjointEigenSolver.h:504
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:1091
EIGEN_STRONG_INLINE void swap(T &a, T &b)
Definition: Meta.h:766
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE internal::enable_if< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typenameNumTraits< T >::Real >::type abs(const T &x)
Definition: MathFunctions.h:1509
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1292
static EIGEN_DEPRECATED const end_t end
Definition: IndexedViewHelper.h:181
Namespace containing all symbols from the Eigen library.
Definition: Core:141
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Definition: Eigen_Colamd.h:50
static constexpr const unit_t< compound_unit< energy::joule, time::seconds > > h(6.626070040e-34)
Planck constant.
static constexpr const charge::coulomb_t e(1.6021766208e-19)
elementary charge.
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:30
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:46
Holds information about the various numeric (i.e.
Definition: NumTraits.h:233
static EIGEN_DEVICE_FUNC void computeRoots(const MatrixType &m, VectorType &roots)
Definition: SelfAdjointEigenSolver.h:749
SolverType::Scalar Scalar
Definition: SelfAdjointEigenSolver.h:745
SolverType::RealVectorType VectorType
Definition: SelfAdjointEigenSolver.h:744
SolverType::EigenvectorsType EigenvectorsType
Definition: SelfAdjointEigenSolver.h:746
SolverType::MatrixType MatrixType
Definition: SelfAdjointEigenSolver.h:743
static EIGEN_DEVICE_FUNC void run(SolverType &solver, const MatrixType &mat, int options)
Definition: SelfAdjointEigenSolver.h:759
SolverType::RealVectorType VectorType
Definition: SelfAdjointEigenSolver.h:585
SolverType::Scalar Scalar
Definition: SelfAdjointEigenSolver.h:586
SolverType::EigenvectorsType EigenvectorsType
Definition: SelfAdjointEigenSolver.h:587
SolverType::MatrixType MatrixType
Definition: SelfAdjointEigenSolver.h:584
static EIGEN_DEVICE_FUNC void run(SolverType &solver, const MatrixType &mat, int options)
Definition: SelfAdjointEigenSolver.h:655
static EIGEN_DEVICE_FUNC void computeRoots(const MatrixType &m, VectorType &roots)
Definition: SelfAdjointEigenSolver.h:595
static EIGEN_DEVICE_FUNC bool extract_kernel(MatrixType &mat, Ref< VectorType > res, Ref< VectorType > representative)
Definition: SelfAdjointEigenSolver.h:634
Definition: SelfAdjointEigenSolver.h:576
static EIGEN_DEVICE_FUNC void run(SolverType &eig, const typename SolverType::MatrixType &A, int options)
Definition: SelfAdjointEigenSolver.h:578
Definition: XprHelper.h:615