11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12 #define EIGEN_SELFADJOINTEIGENSOLVER_H
14 #include "./Tridiagonalization.h"
18 template<
typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
22 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues;
24 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
26 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
const Index maxIterations,
bool computeEigenvectors, MatrixType& eivec);
76 typedef _MatrixType MatrixType;
78 Size = MatrixType::RowsAtCompileTime,
79 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
80 Options = MatrixType::Options,
81 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
85 typedef typename MatrixType::Scalar
Scalar;
124 m_info(InvalidInput),
125 m_isInitialized(false),
126 m_eigenvectorsOk(false)
143 : m_eivec(size, size),
145 m_subdiag(size > 1 ? size - 1 : 1),
146 m_isInitialized(false),
147 m_eigenvectorsOk(false)
165 template<
typename InputType>
168 : m_eivec(matrix.rows(), matrix.cols()),
169 m_eivalues(matrix.cols()),
170 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
171 m_isInitialized(false),
172 m_eigenvectorsOk(false)
207 template<
typename InputType>
267 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
268 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
290 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
314 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
315 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
316 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
339 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
340 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
341 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
351 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
363 static EIGEN_DEVICE_FUNC
364 void check_template_parameters()
366 EIGEN_STATIC_ASSERT_NON_INTEGER(
Scalar);
369 EigenvectorsType m_eivec;
371 typename TridiagonalizationType::SubDiagonalType m_subdiag;
373 bool m_isInitialized;
374 bool m_eigenvectorsOk;
398 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
400 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end, Scalar* matrixQ,
Index n);
403 template<
typename MatrixType>
404 template<
typename InputType>
407 ::compute(
const EigenBase<InputType>& a_matrix,
int options)
409 check_template_parameters();
411 const InputType &matrix(a_matrix.derived());
413 EIGEN_USING_STD(abs);
414 eigen_assert(matrix.cols() == matrix.rows());
415 eigen_assert((options&~(EigVecMask|GenEigMask))==0
416 && (options&EigVecMask)!=EigVecMask
417 &&
"invalid option parameter");
419 Index n = matrix.cols();
420 m_eivalues.resize(n,1);
425 m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
426 if(computeEigenvectors)
427 m_eivec.setOnes(n,n);
429 m_isInitialized =
true;
430 m_eigenvectorsOk = computeEigenvectors;
435 RealVectorType& diag = m_eivalues;
436 EigenvectorsType& mat = m_eivec;
439 mat = matrix.template triangularView<Lower>();
440 RealScalar scale = mat.cwiseAbs().maxCoeff();
441 if(scale==RealScalar(0)) scale = RealScalar(1);
442 mat.template triangularView<Lower>() /= scale;
443 m_subdiag.resize(n-1);
444 internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
446 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
451 m_isInitialized =
true;
452 m_eigenvectorsOk = computeEigenvectors;
456 template<
typename MatrixType>
465 if (computeEigenvectors)
467 m_eivec.setIdentity(diag.size(), diag.size());
469 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
471 m_isInitialized =
true;
472 m_eigenvectorsOk = computeEigenvectors;
488 template<
typename MatrixType,
typename DiagType,
typename SubDiagType>
490 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
const Index maxIterations,
bool computeEigenvectors, MatrixType& eivec)
492 EIGEN_USING_STD(abs);
495 typedef typename MatrixType::Scalar Scalar;
497 Index n = diag.size();
502 typedef typename DiagType::RealScalar RealScalar;
503 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
508 for (
Index i = start; i<end; ++i)
509 if (internal::isMuchSmallerThan(abs(subdiag[i]),(abs(diag[i])+abs(diag[i+1])),precision) || abs(subdiag[i]) <= considerAsZero)
513 while (end>0 && subdiag[end-1]==RealScalar(0))
522 if(iter > maxIterations * n)
break;
525 while (start>0 && subdiag[start-1]!=0)
528 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
530 if (iter <= maxIterations * n)
540 for (
Index i = 0; i < n-1; ++i)
543 diag.segment(i,n-i).minCoeff(&k);
546 numext::swap(diag[i], diag[k+i]);
547 if(computeEigenvectors)
548 eivec.col(i).swap(eivec.col(k+i));
558 static inline void run(SolverType& eig,
const typename SolverType::MatrixType& A,
int options)
559 { eig.compute(A,options); }
564 typedef typename SolverType::MatrixType MatrixType;
565 typedef typename SolverType::RealVectorType VectorType;
566 typedef typename SolverType::Scalar Scalar;
567 typedef typename SolverType::EigenvectorsType EigenvectorsType;
575 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
577 EIGEN_USING_STD(sqrt)
578 EIGEN_USING_STD(atan2)
581 const Scalar s_inv3 = Scalar(1)/Scalar(3);
582 const Scalar s_sqrt3 = sqrt(Scalar(3));
587 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);
588 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);
589 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
593 Scalar c2_over_3 = c2*s_inv3;
594 Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
595 a_over_3 = numext::maxi(a_over_3, Scalar(0));
597 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
599 Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
600 q = numext::maxi(q, Scalar(0));
603 Scalar rho = sqrt(a_over_3);
604 Scalar theta = atan2(sqrt(q),half_b)*s_inv3;
605 Scalar cos_theta = cos(theta);
606 Scalar sin_theta = sin(theta);
608 roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
609 roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
610 roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
616 EIGEN_USING_STD(abs);
617 EIGEN_USING_STD(sqrt);
620 mat.diagonal().cwiseAbs().maxCoeff(&i0);
623 representative = mat.col(i0);
626 n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
627 n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
628 if(n0>n1) res = c0/sqrt(n0);
629 else res = c1/sqrt(n1);
635 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
637 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
638 eigen_assert((options&~(EigVecMask|GenEigMask))==0
639 && (options&EigVecMask)!=EigVecMask
640 &&
"invalid option parameter");
643 EigenvectorsType& eivecs = solver.m_eivec;
644 VectorType& eivals = solver.m_eivalues;
647 Scalar shift = mat.trace() / Scalar(3);
649 MatrixType scaledMat = mat.template selfadjointView<Lower>();
650 scaledMat.diagonal().array() -= shift;
651 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
652 if(scale > 0) scaledMat /= scale;
655 computeRoots(scaledMat,eivals);
658 if(computeEigenvectors)
663 eivecs.setIdentity();
671 Scalar d0 = eivals(2) - eivals(1);
672 Scalar d1 = eivals(1) - eivals(0);
682 tmp.diagonal().array () -= eivals(k);
684 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
692 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
693 eivecs.col(l).normalize();
698 tmp.diagonal().array () -= eivals(l);
701 extract_kernel(tmp, eivecs.col(l), dummy);
705 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
711 eivals.array() += shift;
714 solver.m_isInitialized =
true;
715 solver.m_eigenvectorsOk = computeEigenvectors;
720 template<
typename SolverType>
723 typedef typename SolverType::MatrixType MatrixType;
724 typedef typename SolverType::RealVectorType VectorType;
725 typedef typename SolverType::Scalar Scalar;
726 typedef typename SolverType::EigenvectorsType EigenvectorsType;
729 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
731 EIGEN_USING_STD(sqrt);
732 const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
733 const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
739 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
741 EIGEN_USING_STD(sqrt);
742 EIGEN_USING_STD(abs);
744 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
745 eigen_assert((options&~(EigVecMask|GenEigMask))==0
746 && (options&EigVecMask)!=EigVecMask
747 &&
"invalid option parameter");
750 EigenvectorsType& eivecs = solver.m_eivec;
751 VectorType& eivals = solver.m_eivalues;
754 Scalar shift = mat.trace() / Scalar(2);
755 MatrixType scaledMat = mat;
756 scaledMat.coeffRef(0,1) = mat.coeff(1,0);
757 scaledMat.diagonal().array() -= shift;
758 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
759 if(scale > Scalar(0))
763 computeRoots(scaledMat,eivals);
766 if(computeEigenvectors)
770 eivecs.setIdentity();
774 scaledMat.diagonal().array () -= eivals(1);
775 Scalar a2 = numext::abs2(scaledMat(0,0));
776 Scalar c2 = numext::abs2(scaledMat(1,1));
777 Scalar b2 = numext::abs2(scaledMat(1,0));
780 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
781 eivecs.col(1) /= sqrt(a2+b2);
785 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
786 eivecs.col(1) /= sqrt(c2+b2);
789 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
795 eivals.array() += shift;
798 solver.m_isInitialized =
true;
799 solver.m_eigenvectorsOk = computeEigenvectors;
805 template<
typename MatrixType>
815 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
817 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag,
Index start,
Index end, Scalar* matrixQ,
Index n)
819 EIGEN_USING_STD(abs);
820 RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
821 RealScalar e = subdiag[end-1];
827 RealScalar mu = diag[end];
828 if(td==RealScalar(0))
832 RealScalar e2 = numext::abs2(subdiag[end-1]);
833 RealScalar h = numext::hypot(td,e);
834 if(e2==RealScalar(0)) mu -= (e / (td + (td>RealScalar(0) ? RealScalar(1) : RealScalar(-1)))) * (e / h);
835 else mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
838 RealScalar x = diag[start] - mu;
839 RealScalar z = subdiag[start];
840 for (
Index k = start; k < end; ++k)
842 JacobiRotation<RealScalar> rot;
843 rot.makeGivens(x, z);
846 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
847 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
849 diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
850 diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
851 subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
855 subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
861 z = -rot.s() * subdiag[k+1];
862 subdiag[k + 1] = rot.c() * subdiag[k+1];
869 Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
870 q.applyOnTheRight(k,k+1,rot);
879 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H