11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H 12 #define EIGEN_SELFADJOINTEIGENSOLVER_H 18 template<
typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
23 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
76 Size = MatrixType::RowsAtCompileTime,
77 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
78 Options = MatrixType::Options,
79 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
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)
201 template<
typename InputType>
261 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
262 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
284 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();
333 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
334 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
335 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
345 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
354 static const int m_maxIterations = 30;
391 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
396 template<
typename MatrixType>
397 template<
typename InputType>
402 check_template_parameters();
410 &&
"invalid option parameter");
413 m_eivalues.resize(n,1);
418 m_eivalues.coeffRef(0,0) =
numext::real(m_eivec.coeff(0,0));
419 if(computeEigenvectors)
420 m_eivec.setOnes(n,n);
422 m_isInitialized =
true;
423 m_eigenvectorsOk = computeEigenvectors;
436 m_subdiag.resize(n-1);
444 m_isInitialized =
true;
445 m_eigenvectorsOk = computeEigenvectors;
449 template<
typename MatrixType>
458 if (computeEigenvectors)
460 m_eivec.setIdentity(diag.size(), diag.size());
464 m_isInitialized =
true;
465 m_eigenvectorsOk = computeEigenvectors;
481 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
505 while (end>0 && subdiag[end-1]==
RealScalar(0))
514 if(iter > maxIterations * n)
break;
517 while (start>0 && subdiag[start-1]!=0)
520 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start,
end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
522 if (iter <= maxIterations * n)
535 diag.segment(
i,n-
i).minCoeff(&k);
539 if(computeEigenvectors)
540 eivec.col(i).swap(eivec.col(k+i));
551 { eig.compute(A,options); }
569 EIGEN_USING_STD_MATH(
sqrt)
570 EIGEN_USING_STD_MATH(
atan2)
571 EIGEN_USING_STD_MATH(
cos)
572 EIGEN_USING_STD_MATH(
sin)
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) -
m(1,1)*
m(2,0)*
m(2,0) -
m(2,2)*
m(1,0)*
m(1,0);
580 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);
581 Scalar
c2 =
m(0,0) +
m(1,1) +
m(2,2);
585 Scalar c2_over_3 = c2*s_inv3;
586 Scalar a_over_3 = (c2*c2_over_3 -
c1)*s_inv3;
589 Scalar half_b =
Scalar(0.5)*(c0 + c2_over_3*(
Scalar(2)*c2_over_3*c2_over_3 -
c1));
591 Scalar
q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
595 Scalar rho =
sqrt(a_over_3);
597 Scalar cos_theta =
cos(theta);
598 Scalar sin_theta =
sin(theta);
600 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
601 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
602 roots(2) = c2_over_3 +
Scalar(2)*rho*cos_theta;
611 mat.diagonal().cwiseAbs().maxCoeff(&i0);
614 representative = mat.col(i0);
617 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
618 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
628 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
631 &&
"invalid option parameter");
634 EigenvectorsType& eivecs = solver.m_eivec;
635 VectorType&
eivals = solver.m_eivalues;
638 Scalar shift = mat.trace() /
Scalar(3);
640 MatrixType scaledMat = mat.template selfadjointView<Lower>();
641 scaledMat.diagonal().array() -= shift;
642 Scalar
scale = scaledMat.cwiseAbs().maxCoeff();
643 if(scale > 0) scaledMat /=
scale;
649 if(computeEigenvectors)
654 eivecs.setIdentity();
673 tmp.diagonal().array () -=
eivals(k);
675 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
683 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
684 eivecs.col(l).normalize();
689 tmp.diagonal().array () -=
eivals(l);
692 extract_kernel(tmp, eivecs.col(l), dummy);
696 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
702 eivals.array() += shift;
705 solver.m_isInitialized =
true;
706 solver.m_eigenvectorsOk = computeEigenvectors;
711 template<
typename SolverType>
724 const Scalar t1 =
Scalar(0.5) * (
m(0,0) +
m(1,1));
732 EIGEN_USING_STD_MATH(
sqrt);
733 EIGEN_USING_STD_MATH(
abs);
735 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
738 &&
"invalid option parameter");
741 EigenvectorsType& eivecs = solver.m_eivec;
742 VectorType&
eivals = solver.m_eivalues;
745 Scalar shift = mat.trace() /
Scalar(2);
746 MatrixType scaledMat =
mat;
747 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
748 scaledMat.diagonal().array() -= shift;
749 Scalar
scale = scaledMat.cwiseAbs().maxCoeff();
757 if(computeEigenvectors)
761 eivecs.setIdentity();
765 scaledMat.diagonal().array () -=
eivals(1);
771 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
772 eivecs.col(1) /=
sqrt(a2+b2);
776 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
777 eivecs.col(1) /=
sqrt(c2+b2);
780 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
786 eivals.array() += shift;
789 solver.m_isInitialized =
true;
790 solver.m_eigenvectorsOk = computeEigenvectors;
796 template<
typename MatrixType>
806 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
826 else mu -= e2 / (td + (td>
RealScalar(0) ? h : -
h));
837 RealScalar sdk = rot.
s() * diag[k] + rot.
c() * subdiag[k];
838 RealScalar dkp1 = rot.
s() * subdiag[k] + rot.
c() * diag[k+1];
840 diag[k] = rot.
c() * (rot.
c() * diag[k] - rot.
s() * subdiag[k]) - rot.
s() * (rot.
c() * subdiag[k] - rot.
s() * diag[k+1]);
841 diag[k+1] = rot.
s() * sdk + rot.
c() * dkp1;
842 subdiag[k] = rot.
c() * sdk - rot.
s() * dkp1;
846 subdiag[k - 1] = rot.
c() * subdiag[k-1] - rot.
s() *
z;
852 z = -rot.
s() * subdiag[k+1];
853 subdiag[k + 1] = rot.
c() * subdiag[k+1];
861 q.applyOnTheRight(k,k+1,rot);
870 #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.
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
EIGEN_DEVICE_FUNC MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
SolverType::RealVectorType VectorType
Matrix diag(const std::vector< Matrix > &Hs)
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
SolverType::EigenvectorsType EigenvectorsType
SolverType::RealVectorType VectorType
SolverType::MatrixType MatrixType
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.
static EIGEN_DEVICE_FUNC bool extract_kernel(MatrixType &mat, Ref< VectorType > res, Ref< VectorType > representative)
SolverType::Scalar Scalar
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
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)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
EIGEN_DEVICE_FUNC const CosReturnType cos() const
EIGEN_STRONG_INLINE void swap(T &a, T &b)
static const Line3 l(Rot3(), 1, 1)
SolverType::MatrixType MatrixType
static EIGEN_DEVICE_FUNC void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
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)
ComputationInfo computeFromTridiagonal_impl(DiagType &diag, SubDiagType &subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType &eivec)
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
RealVectorType m_eivalues
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
EIGEN_DEVICE_FUNC const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
A matrix or vector expression mapping an existing expression.
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
static EIGEN_DEVICE_FUNC void run(SolverType &solver, const MatrixType &mat, int options)
static void check_template_parameters()
const mpreal hypot(const mpreal &x, const mpreal &y, mp_rnd_t rnd_mode=mpreal::get_default_rnd())
A triangularView< Lower >().adjoint().solveInPlace(B)
EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
EIGEN_DEVICE_FUNC const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
const AutoDiffScalar< Matrix< typename internal::traits< typename internal::remove_all< DerTypeA >::type >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
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 const SinReturnType sin() const
SolverType::Scalar Scalar
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
void run(Expr &expr, Dev &dev)
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 Derived & derived()
void swap(mpfr::mpreal &x, mpfr::mpreal &y)
Matrix< Scalar, Size, Size, ColMajor, MaxColsAtCompileTime, MaxColsAtCompileTime > EigenvectorsType
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)