11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H 12 #define EIGEN_SELFADJOINTEIGENSOLVER_H 18 template<
typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
24 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
82 Size = MatrixType::RowsAtCompileTime,
83 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
84 Options = MatrixType::Options,
85 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
129 m_info(InvalidInput),
130 m_isInitialized(false),
131 m_eigenvectorsOk(false)
148 : m_eivec(size, size),
150 m_subdiag(size > 1 ? size - 1 : 1),
151 m_hcoeffs(size > 1 ? size - 1 : 1),
152 m_isInitialized(false),
153 m_eigenvectorsOk(false)
171 template<
typename InputType>
174 : m_eivec(matrix.
rows(), matrix.
cols()),
175 m_eivalues(matrix.
cols()),
176 m_subdiag(matrix.
rows() > 1 ? matrix.
rows() - 1 : 1),
177 m_hcoeffs(matrix.
cols() > 1 ? matrix.
cols() - 1 : 1),
178 m_isInitialized(false),
179 m_eigenvectorsOk(false)
214 template<
typename InputType>
279 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
280 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
302 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
326 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
327 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
328 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
351 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
352 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
353 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
363 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
372 static const int m_maxIterations = 30;
411 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
416 template<
typename MatrixType>
417 template<
typename InputType>
422 check_template_parameters();
430 &&
"invalid option parameter");
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;
456 m_subdiag.resize(n-1);
457 m_hcoeffs.resize(n-1);
465 m_isInitialized =
true;
466 m_eigenvectorsOk = computeEigenvectors;
470 template<
typename MatrixType>
479 if (computeEigenvectors)
481 m_eivec.setIdentity(diag.size(), diag.size());
485 m_isInitialized =
true;
486 m_eigenvectorsOk = computeEigenvectors;
502 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
525 const RealScalar scaled_subdiag = precision_inv * subdiag[
i];
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)
563 diag.segment(
i,n-
i).minCoeff(&k);
567 if(computeEigenvectors)
568 eivec.col(i).swap(eivec.col(k+i));
579 { eig.compute(A,options); }
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;
617 Scalar half_b =
Scalar(0.5)*(c0 + c2_over_3*(
Scalar(2)*c2_over_3*c2_over_3 - c1));
619 Scalar
q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
623 Scalar rho =
sqrt(a_over_3);
624 Scalar theta =
atan2(
sqrt(q),half_b)*s_inv3;
625 Scalar cos_theta =
cos(theta);
626 Scalar sin_theta =
sin(theta);
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);
657 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
660 &&
"invalid option parameter");
663 EigenvectorsType& eivecs = solver.m_eivec;
664 VectorType&
eivals = solver.m_eivalues;
667 Scalar shift = mat.trace() /
Scalar(3);
669 MatrixType scaledMat = mat.template selfadjointView<Lower>();
670 scaledMat.diagonal().array() -= shift;
671 Scalar
scale = scaledMat.cwiseAbs().maxCoeff();
672 if(scale > 0) scaledMat /=
scale;
678 if(computeEigenvectors)
683 eivecs.setIdentity();
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;
740 template<
typename SolverType>
753 const Scalar t1 =
Scalar(0.5) * (
m(0,0) +
m(1,1));
764 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
767 &&
"invalid option parameter");
770 EigenvectorsType& eivecs = solver.m_eivec;
771 VectorType&
eivals = solver.m_eivalues;
774 Scalar shift = mat.trace() /
Scalar(2);
775 MatrixType scaledMat =
mat;
776 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
777 scaledMat.diagonal().array() -= shift;
778 Scalar
scale = scaledMat.cwiseAbs().maxCoeff();
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;
825 template<
typename MatrixType>
837 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
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);
904 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
static EIGEN_DEVICE_FUNC void run(SolverType &solver, const MatrixType &mat, int options)
int EIGEN_BLAS_FUNC() rot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, RealScalar *ps)
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Jet< T, N > cos(const Jet< T, N > &f)
SolverType::RealVectorType VectorType
Matrix diag(const std::vector< Matrix > &Hs)
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
SolverType::EigenvectorsType EigenvectorsType
SolverType::RealVectorType VectorType
EIGEN_DEVICE_FUNC Scalar & c()
SolverType::MatrixType MatrixType
EIGEN_DEVICE_FUNC ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
A matrix or vector expression mapping an existing array of data.
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Jet< T, N > sin(const Jet< T, N > &f)
static EIGEN_DEVICE_FUNC bool extract_kernel(MatrixType &mat, Ref< VectorType > res, Ref< VectorType > representative)
SolverType::Scalar Scalar
Namespace containing all symbols from the Eigen library.
SolverType::EigenvectorsType EigenvectorsType
Rotation given by a cosine-sine pair.
iterator iter(handle obj)
static EIGEN_DEVICE_FUNC void computeRoots(const MatrixType &m, VectorType &roots)
BiCGSTAB< SparseMatrix< double > > solver
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Tridiagonal decomposition of a selfadjoint matrix.
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
static EIGEN_DEVICE_FUNC void check_template_parameters()
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
EIGEN_STRONG_INLINE void swap(T &a, T &b)
static const Line3 l(Rot3(), 1, 1)
SolverType::MatrixType MatrixType
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE internal::enable_if< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real >::type abs(const T &x)
static EIGEN_DEVICE_FUNC void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
EIGEN_DEVICE_FUNC void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
static EIGEN_DEVICE_FUNC void computeRoots(const MatrixType &m, VectorType &roots)
static EIGEN_DEVICE_FUNC void run(SolverType &eig, const typename SolverType::MatrixType &A, int options)
EIGEN_DEVICE_FUNC const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
idx_t idx_t idx_t idx_t idx_t idx_t idx_t real_t real_t idx_t * options
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
SelfAdjointEigenSolver< PlainMatrixType > eig(mat, computeVectors?ComputeEigenvectors:EigenvaluesOnly)
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Array< double, 1, 3 > e(1./3., 0.5, 2.)
EIGEN_DEVICE_FUNC const Scalar & q
TridiagonalizationType::CoeffVectorType m_hcoeffs
RealVectorType m_eivalues
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
NumTraits< Scalar >::Real RealScalar
A matrix or vector expression mapping an existing expression.
EIGEN_CONSTEXPR Index size(const T &x)
#define EIGEN_DEVICE_FUNC
static EIGEN_DEVICE_FUNC void run(SolverType &solver, const MatrixType &mat, int options)
EIGEN_DEVICE_FUNC Scalar & s()
A triangularView< Lower >().adjoint().solveInPlace(B)
#define EIGEN_USING_STD(FUNC)
static EIGEN_DEPRECATED const end_t end
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
const AutoDiffScalar< Matrix< typename internal::traits< typename internal::remove_all< DerTypeA >::type >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics scale
void computeRoots(const Matrix &m, Roots &roots)
EIGEN_DEVICE_FUNC void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
Jet< T, N > sqrt(const Jet< T, N > &f)
EIGEN_DEVICE_FUNC const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
SolverType::Scalar Scalar
Generic expression where a coefficient-wise unary operator is applied to an expression.
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
TridiagonalizationType::SubDiagonalType m_subdiag
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
EIGEN_DEVICE_FUNC MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
EIGEN_DEVICE_FUNC bool abs2(bool x)
EIGEN_DEVICE_FUNC Derived & derived()
Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType