11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12 #define EIGEN_SELFADJOINTEIGENSOLVER_H
20 template<
typename MatrixType_>
21 class GeneralizedSelfAdjointEigenSolver;
24 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
26 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
84 Size = MatrixType::RowsAtCompileTime,
111 typedef typename internal::plain_col_type<MatrixType, Scalar>::type
VectorType;
112 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type
RealVectorType;
176 template<
typename InputType>
220 template<
typename InputType>
414 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
419 template<
typename MatrixType>
420 template<
typename InputType>
425 const InputType &matrix(a_matrix.
derived());
431 &&
"invalid option parameter");
440 if(computeEigenvectors)
453 mat = matrix.template triangularView<Lower>();
456 mat.template triangularView<Lower>() /= scale;
471 template<
typename MatrixType>
480 if (computeEigenvectors)
482 m_eivec.setIdentity(diag.size(), diag.size());
503 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
508 typedef typename MatrixType::Scalar
Scalar;
515 typedef typename DiagType::RealScalar
RealScalar;
526 const RealScalar scaled_subdiag = precision_inv * subdiag[
i];
543 if(iter > maxIterations *
n)
break;
549 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start,
end, computeEigenvectors ? eivec.data() : (
Scalar*)0,
n);
551 if (iter <= maxIterations *
n)
564 diag.segment(
i,
n-
i).minCoeff(&k);
568 if(computeEigenvectors)
569 eivec.col(
i).swap(eivec.col(k+
i));
576 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues
580 { eig.compute(
A,options); }
583 template<
typename SolverType>
struct direct_selfadjoint_eigenvalues<SolverType,3,false>
586 typedef typename SolverType::RealVectorType
VectorType;
587 typedef typename SolverType::Scalar
Scalar;
608 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);
609 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);
614 Scalar c2_over_3 = c2*s_inv3;
615 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
620 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
629 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
630 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
631 roots(2) = c2_over_3 +
Scalar(2)*rho*cos_theta;
635 static inline bool extract_kernel(
MatrixType&
mat, Ref<VectorType>
res, Ref<VectorType> representative)
641 mat.diagonal().cwiseAbs().maxCoeff(&i0);
644 representative =
mat.col(i0);
647 n0 = (c0 = representative.cross(
mat.col((i0+1)%3))).squaredNorm();
648 n1 = (c1 = representative.cross(
mat.col((i0+2)%3))).squaredNorm();
661 &&
"invalid option parameter");
670 MatrixType scaledMat =
mat.template selfadjointView<Lower>();
671 scaledMat.diagonal().array() -= shift;
672 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
673 if(scale > 0) scaledMat /= scale;
676 computeRoots(scaledMat,
eivals);
679 if(computeEigenvectors)
684 eivecs.setIdentity();
703 tmp.diagonal().array () -=
eivals(k);
705 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
713 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
714 eivecs.col(l).normalize();
719 tmp.diagonal().array () -=
eivals(l);
722 extract_kernel(tmp, eivecs.col(l), dummy);
726 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
735 solver.m_isInitialized =
true;
736 solver.m_eigenvectorsOk = computeEigenvectors;
741 template<
typename SolverType>
742 struct direct_selfadjoint_eigenvalues<SolverType,2,false>
745 typedef typename SolverType::RealVectorType
VectorType;
746 typedef typename SolverType::Scalar
Scalar;
768 &&
"invalid option parameter");
777 scaledMat.coeffRef(0,1) =
mat.coeff(1,0);
778 scaledMat.diagonal().array() -= shift;
779 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
784 computeRoots(scaledMat,
eivals);
787 if(computeEigenvectors)
791 eivecs.setIdentity();
795 scaledMat.diagonal().array () -=
eivals(1);
801 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
802 eivecs.col(1) /=
sqrt(a2+b2);
806 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
807 eivecs.col(1) /=
sqrt(c2+b2);
810 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
819 solver.m_isInitialized =
true;
820 solver.m_eigenvectorsOk = computeEigenvectors;
826 template<
typename MatrixType>
831 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*
this,matrix,options);
838 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
859 mu -= e2 / (td + (td>
RealScalar(0) ? h : -h));
873 RealScalar sdk = rot.
s() * diag[k] + rot.
c() * subdiag[k];
874 RealScalar dkp1 = rot.
s() * subdiag[k] + rot.
c() * diag[k+1];
876 diag[k] = rot.
c() * (rot.
c() * diag[k] - rot.
s() * subdiag[k]) - rot.
s() * (rot.
c() * subdiag[k] - rot.
s() * diag[k+1]);
877 diag[k+1] = rot.
s() * sdk + rot.
c() * dkp1;
878 subdiag[k] = rot.
c() * sdk - rot.
s() * dkp1;
881 subdiag[k - 1] = rot.
c() * subdiag[k-1] - rot.
s() * z;
887 z = -rot.
s() * subdiag[k+1];
888 subdiag[k + 1] = rot.
c() * subdiag[k+1];
896 q.applyOnTheRight(k,k+1,rot);
const CwiseBinaryOp< atan2< Scalar >, const Derived, const OtherDerived > atan2(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
BiCGSTAB< SparseMatrix< double > > solver
RealReturnType real() const
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define EIGEN_USING_STD(FUNC)
#define EIGEN_DEVICE_FUNC
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Matrix< float, 1, Dynamic > MatrixType
internal::traits< Derived >::Scalar Scalar
Rotation given by a cosine-sine pair.
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
A matrix or vector expression mapping an existing array of data.
constexpr const Scalar & coeff(Index rowId, Index colId) const
Derived & setOnes(Index size)
constexpr void resize(Index rows, Index cols)
Computes eigenvalues and eigenvectors of selfadjoint matrices.
TridiagonalizationType::CoeffVectorType m_hcoeffs
internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
TridiagonalizationType::SubDiagonalType m_subdiag
RealVectorType m_eivalues
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType_.
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
ComputationInfo info() const
Reports whether previous computation was successful.
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
static const int m_maxIterations
Maximum number of iterations.
internal::plain_col_type< MatrixType, Scalar >::type VectorType
Type for vector of eigenvalues as returned by eigenvalues().
SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Tridiagonal decomposition of a selfadjoint matrix.
static const lastp1_t end
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
static void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
bool is_exactly_zero(const X &x)
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
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_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.