Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
SelfAdjointEigenSolver.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
12#define EIGEN_SELFADJOINTEIGENSOLVER_H
13
14#include "./Tridiagonalization.h"
15
16// IWYU pragma: private
17#include "./InternalHeaderCheck.h"
18
19namespace Eigen {
20
21template <typename MatrixType_>
23
24namespace internal {
25template <typename SolverType, int Size, bool IsComplex>
26struct direct_selfadjoint_eigenvalues;
27
28template <typename MatrixType, typename DiagType, typename SubDiagType>
29EIGEN_DEVICE_FUNC ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
30 const Index maxIterations, bool computeEigenvectors,
31 MatrixType& eivec);
32} // namespace internal
33
81template <typename MatrixType_>
83 public:
84 typedef MatrixType_ MatrixType;
85 enum {
86 Size = MatrixType::RowsAtCompileTime,
87 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
88 Options = internal::traits<MatrixType>::Options,
89 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
90 };
91
93 typedef typename MatrixType::Scalar Scalar;
95
97
104 typedef typename NumTraits<Scalar>::Real RealScalar;
105
106 friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver, Size, NumTraits<Scalar>::IsComplex>;
107
113 typedef typename internal::plain_col_type<MatrixType, Scalar>::type VectorType;
114 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
115 typedef Tridiagonalization<MatrixType> TridiagonalizationType;
116 typedef typename TridiagonalizationType::SubDiagonalType SubDiagonalType;
117
128 EIGEN_DEVICE_FUNC SelfAdjointEigenSolver()
129 : m_eivec(),
130 m_workspace(),
131 m_eivalues(),
132 m_subdiag(),
133 m_hcoeffs(),
134 m_info(InvalidInput),
135 m_isInitialized(false),
136 m_eigenvectorsOk(false) {}
137
150 EIGEN_DEVICE_FUNC explicit SelfAdjointEigenSolver(Index size)
151 : m_eivec(size, size),
152 m_workspace(size),
153 m_eivalues(size),
154 m_subdiag(size > 1 ? size - 1 : 1),
155 m_hcoeffs(size > 1 ? size - 1 : 1),
156 m_isInitialized(false),
157 m_eigenvectorsOk(false) {}
158
174 template <typename InputType>
175 EIGEN_DEVICE_FUNC explicit SelfAdjointEigenSolver(const EigenBase<InputType>& matrix,
176 int options = ComputeEigenvectors)
177 : m_eivec(matrix.rows(), matrix.cols()),
178 m_workspace(matrix.cols()),
179 m_eivalues(matrix.cols()),
180 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
181 m_hcoeffs(matrix.cols() > 1 ? matrix.cols() - 1 : 1),
182 m_isInitialized(false),
183 m_eigenvectorsOk(false) {
184 compute(matrix.derived(), options);
185 }
186
217 template <typename InputType>
218 EIGEN_DEVICE_FUNC SelfAdjointEigenSolver& compute(const EigenBase<InputType>& matrix,
219 int options = ComputeEigenvectors);
220
239 EIGEN_DEVICE_FUNC SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
240
253 SelfAdjointEigenSolver& computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag,
254 int options = ComputeEigenvectors);
255
279 EIGEN_DEVICE_FUNC const EigenvectorsType& eigenvectors() const {
280 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
281 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
282 return m_eivec;
283 }
284
300 EIGEN_DEVICE_FUNC const RealVectorType& eigenvalues() const {
301 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
302 return m_eivalues;
303 }
304
322 EIGEN_DEVICE_FUNC MatrixType operatorSqrt() const {
323 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
324 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
325 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
326 }
327
338 EIGEN_DEVICE_FUNC MatrixType operatorExp() const {
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.array().exp().matrix().asDiagonal() * m_eivec.adjoint();
342 }
343
362 EIGEN_DEVICE_FUNC MatrixType operatorInverseSqrt() const {
363 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
364 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
365 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
366 }
367
372 EIGEN_DEVICE_FUNC ComputationInfo info() const {
373 eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
374 return m_info;
375 }
376
382 static const int m_maxIterations = 30;
383
384 protected:
385 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
386
387 EigenvectorsType m_eivec;
388 VectorType m_workspace;
389 RealVectorType m_eivalues;
390 typename TridiagonalizationType::SubDiagonalType m_subdiag;
391 typename TridiagonalizationType::CoeffVectorType m_hcoeffs;
392 ComputationInfo m_info;
393 bool m_isInitialized;
394 bool m_eigenvectorsOk;
395};
396
397namespace internal {
418template <int StorageOrder, typename RealScalar, typename Scalar, typename Index>
419EIGEN_DEVICE_FUNC static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end,
420 Scalar* matrixQ, Index n);
421} // namespace internal
422
423template <typename MatrixType>
424template <typename InputType>
425EIGEN_DEVICE_FUNC SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>::compute(
426 const EigenBase<InputType>& a_matrix, int options) {
427 const InputType& matrix(a_matrix.derived());
428
429 EIGEN_USING_STD(abs);
430 eigen_assert(matrix.cols() == matrix.rows());
431 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
432 "invalid option parameter");
433 bool computeEigenvectors = (options & ComputeEigenvectors) == ComputeEigenvectors;
434 Index n = matrix.cols();
435 m_eivalues.resize(n, 1);
436
437 if (n == 1) {
438 m_eivec = matrix;
439 m_eivalues.coeffRef(0, 0) = numext::real(m_eivec.coeff(0, 0));
440 if (computeEigenvectors) m_eivec.setOnes(n, n);
441 m_info = Success;
442 m_isInitialized = true;
443 m_eigenvectorsOk = computeEigenvectors;
444 return *this;
445 }
446
447 // declare some aliases
448 RealVectorType& diag = m_eivalues;
449 EigenvectorsType& mat = m_eivec;
450
451 // map the matrix coefficients to [-1:1] to avoid over- and underflow.
452 mat = matrix.template triangularView<Lower>();
453 RealScalar scale = mat.cwiseAbs().maxCoeff();
454 if (numext::is_exactly_zero(scale)) scale = RealScalar(1);
455 mat.template triangularView<Lower>() /= scale;
456 m_subdiag.resize(n - 1);
457 m_hcoeffs.resize(n - 1);
458 internal::tridiagonalization_inplace(mat, diag, m_subdiag, m_hcoeffs, m_workspace, computeEigenvectors);
459
460 m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
461
462 // scale back the eigen values
463 m_eivalues *= scale;
464
465 m_isInitialized = true;
466 m_eigenvectorsOk = computeEigenvectors;
467 return *this;
468}
469
470template <typename MatrixType>
472 const RealVectorType& diag, const SubDiagonalType& subdiag, int options) {
473 // TODO : Add an option to scale the values beforehand
474 bool computeEigenvectors = (options & ComputeEigenvectors) == ComputeEigenvectors;
475
476 m_eivalues = diag;
477 m_subdiag = subdiag;
478 if (computeEigenvectors) {
479 m_eivec.setIdentity(diag.size(), diag.size());
480 }
481 m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
482
483 m_isInitialized = true;
484 m_eigenvectorsOk = computeEigenvectors;
485 return *this;
486}
487
488namespace internal {
500template <typename MatrixType, typename DiagType, typename SubDiagType>
501EIGEN_DEVICE_FUNC ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag,
502 const Index maxIterations, bool computeEigenvectors,
503 MatrixType& eivec) {
504 ComputationInfo info;
505 typedef typename MatrixType::Scalar Scalar;
506
507 Index n = diag.size();
508 Index end = n - 1;
509 Index start = 0;
510 Index iter = 0; // total number of iterations
511
512 typedef typename DiagType::RealScalar RealScalar;
513 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
514 const RealScalar precision_inv = RealScalar(1) / NumTraits<RealScalar>::epsilon();
515 while (end > 0) {
516 for (Index i = start; i < end; ++i) {
517 if (numext::abs(subdiag[i]) < considerAsZero) {
518 subdiag[i] = RealScalar(0);
519 } else {
520 // abs(subdiag[i]) <= epsilon * sqrt(abs(diag[i]) + abs(diag[i+1]))
521 // Scaled to prevent underflows.
522 const RealScalar scaled_subdiag = precision_inv * subdiag[i];
523 if (scaled_subdiag * scaled_subdiag <= (numext::abs(diag[i]) + numext::abs(diag[i + 1]))) {
524 subdiag[i] = RealScalar(0);
525 }
526 }
527 }
528
529 // find the largest unreduced block at the end of the matrix.
530 while (end > 0 && numext::is_exactly_zero(subdiag[end - 1])) {
531 end--;
532 }
533 if (end <= 0) break;
534
535 // if we spent too many iterations, we give up
536 iter++;
537 if (iter > maxIterations * n) break;
538
539 start = end - 1;
540 while (start > 0 && !numext::is_exactly_zero(subdiag[start - 1])) start--;
541
542 internal::tridiagonal_qr_step<MatrixType::Flags & RowMajorBit ? RowMajor : ColMajor>(
543 diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
544 }
545 if (iter <= maxIterations * n)
546 info = Success;
547 else
548 info = NoConvergence;
549
550 // Sort eigenvalues and corresponding vectors.
551 // TODO make the sort optional ?
552 // TODO use a better sort algorithm !!
553 if (info == Success) {
554 for (Index i = 0; i < n - 1; ++i) {
555 Index k;
556 diag.segment(i, n - i).minCoeff(&k);
557 if (k > 0) {
558 numext::swap(diag[i], diag[k + i]);
559 if (computeEigenvectors) eivec.col(i).swap(eivec.col(k + i));
560 }
561 }
562 }
563 return info;
564}
565
566template <typename SolverType, int Size, bool IsComplex>
567struct direct_selfadjoint_eigenvalues {
568 EIGEN_DEVICE_FUNC static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options) {
569 eig.compute(A, options);
570 }
571};
572
573template <typename SolverType>
574struct direct_selfadjoint_eigenvalues<SolverType, 3, false> {
575 typedef typename SolverType::MatrixType MatrixType;
576 typedef typename SolverType::RealVectorType VectorType;
577 typedef typename SolverType::Scalar Scalar;
578 typedef typename SolverType::EigenvectorsType EigenvectorsType;
579
584 EIGEN_DEVICE_FUNC static inline void computeRoots(const MatrixType& m, VectorType& roots) {
585 EIGEN_USING_STD(sqrt)
586 EIGEN_USING_STD(atan2)
587 EIGEN_USING_STD(cos)
588 EIGEN_USING_STD(sin)
589 const Scalar s_inv3 = Scalar(1) / Scalar(3);
590 const Scalar s_sqrt3 = sqrt(Scalar(3));
591
592 // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
593 // eigenvalues are the roots to this equation, all guaranteed to be
594 // real-valued, because the matrix is symmetric.
595 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) -
596 m(1, 1) * m(2, 0) * m(2, 0) - m(2, 2) * m(1, 0) * m(1, 0);
597 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) -
598 m(2, 1) * m(2, 1);
599 Scalar c2 = m(0, 0) + m(1, 1) + m(2, 2);
600
601 // Construct the parameters used in classifying the roots of the equation
602 // and in solving the equation for the roots in closed form.
603 Scalar c2_over_3 = c2 * s_inv3;
604 Scalar a_over_3 = (c2 * c2_over_3 - c1) * s_inv3;
605 a_over_3 = numext::maxi(a_over_3, Scalar(0));
606
607 Scalar half_b = Scalar(0.5) * (c0 + c2_over_3 * (Scalar(2) * c2_over_3 * c2_over_3 - c1));
608
609 Scalar q = a_over_3 * a_over_3 * a_over_3 - half_b * half_b;
610 q = numext::maxi(q, Scalar(0));
611
612 // Compute the eigenvalues by solving for the roots of the polynomial.
613 Scalar rho = sqrt(a_over_3);
614 Scalar theta = atan2(sqrt(q), half_b) * s_inv3; // since sqrt(q) > 0, atan2 is in [0, pi] and theta is in [0, pi/3]
615 Scalar cos_theta = cos(theta);
616 Scalar sin_theta = sin(theta);
617 // roots are already sorted, since cos is monotonically decreasing on [0, pi]
618 roots(0) = c2_over_3 - rho * (cos_theta + s_sqrt3 * sin_theta); // == 2*rho*cos(theta+2pi/3)
619 roots(1) = c2_over_3 - rho * (cos_theta - s_sqrt3 * sin_theta); // == 2*rho*cos(theta+ pi/3)
620 roots(2) = c2_over_3 + Scalar(2) * rho * cos_theta;
621 }
622
623 EIGEN_DEVICE_FUNC static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res,
624 Ref<VectorType> representative) {
625 EIGEN_USING_STD(abs);
626 EIGEN_USING_STD(sqrt);
627 Index i0;
628 // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
629 mat.diagonal().cwiseAbs().maxCoeff(&i0);
630 // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
631 // so let's save it:
632 representative = mat.col(i0);
633 Scalar n0, n1;
634 VectorType c0, c1;
635 n0 = (c0 = representative.cross(mat.col((i0 + 1) % 3))).squaredNorm();
636 n1 = (c1 = representative.cross(mat.col((i0 + 2) % 3))).squaredNorm();
637 if (n0 > n1)
638 res = c0 / sqrt(n0);
639 else
640 res = c1 / sqrt(n1);
641
642 return true;
643 }
644
645 EIGEN_DEVICE_FUNC static inline void run(SolverType& solver, const MatrixType& mat, int options) {
646 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
647 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
648 "invalid option parameter");
649 bool computeEigenvectors = (options & ComputeEigenvectors) == ComputeEigenvectors;
650
651 EigenvectorsType& eivecs = solver.m_eivec;
652 VectorType& eivals = solver.m_eivalues;
653
654 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
655 Scalar shift = mat.trace() / Scalar(3);
656 // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for
657 // computing the eigenvectors later
658 MatrixType scaledMat = mat.template selfadjointView<Lower>();
659 scaledMat.diagonal().array() -= shift;
660 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
661 if (scale > 0) scaledMat /= scale; // TODO for scale==0 we could save the remaining operations
662
663 // compute the eigenvalues
664 computeRoots(scaledMat, eivals);
665
666 // compute the eigenvectors
667 if (computeEigenvectors) {
668 if ((eivals(2) - eivals(0)) <= Eigen::NumTraits<Scalar>::epsilon()) {
669 // All three eigenvalues are numerically the same
670 eivecs.setIdentity();
671 } else {
672 MatrixType tmp;
673 tmp = scaledMat;
674
675 // Compute the eigenvector of the most distinct eigenvalue
676 Scalar d0 = eivals(2) - eivals(1);
677 Scalar d1 = eivals(1) - eivals(0);
678 Index k(0), l(2);
679 if (d0 > d1) {
680 numext::swap(k, l);
681 d0 = d1;
682 }
683
684 // Compute the eigenvector of index k
685 {
686 tmp.diagonal().array() -= eivals(k);
687 // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
688 extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
689 }
690
691 // Compute eigenvector of index l
692 if (d0 <= 2 * Eigen::NumTraits<Scalar>::epsilon() * d1) {
693 // If d0 is too small, then the two other eigenvalues are numerically the same,
694 // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
695 eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l)) * eivecs.col(l);
696 eivecs.col(l).normalize();
697 } else {
698 tmp = scaledMat;
699 tmp.diagonal().array() -= eivals(l);
700
701 VectorType dummy;
702 extract_kernel(tmp, eivecs.col(l), dummy);
703 }
704
705 // Compute last eigenvector from the other two
706 eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
707 }
708 }
709
710 // Rescale back to the original size.
711 eivals *= scale;
712 eivals.array() += shift;
713
714 solver.m_info = Success;
715 solver.m_isInitialized = true;
716 solver.m_eigenvectorsOk = computeEigenvectors;
717 }
718};
719
720// 2x2 direct eigenvalues decomposition, code from Hauke Heibel
721template <typename SolverType>
722struct direct_selfadjoint_eigenvalues<SolverType, 2, false> {
723 typedef typename SolverType::MatrixType MatrixType;
724 typedef typename SolverType::RealVectorType VectorType;
725 typedef typename SolverType::Scalar Scalar;
726 typedef typename SolverType::EigenvectorsType EigenvectorsType;
727
728 EIGEN_DEVICE_FUNC static inline void computeRoots(const MatrixType& m, VectorType& roots) {
729 EIGEN_USING_STD(sqrt);
730 const Scalar t0 = Scalar(0.5) * sqrt(numext::abs2(m(0, 0) - m(1, 1)) + Scalar(4) * numext::abs2(m(1, 0)));
731 const Scalar t1 = Scalar(0.5) * (m(0, 0) + m(1, 1));
732 roots(0) = t1 - t0;
733 roots(1) = t1 + t0;
734 }
735
736 EIGEN_DEVICE_FUNC static inline void run(SolverType& solver, const MatrixType& mat, int options) {
737 EIGEN_USING_STD(sqrt);
738 EIGEN_USING_STD(abs);
739
740 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
741 eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
742 "invalid option parameter");
743 bool computeEigenvectors = (options & ComputeEigenvectors) == ComputeEigenvectors;
744
745 EigenvectorsType& eivecs = solver.m_eivec;
746 VectorType& eivals = solver.m_eivalues;
747
748 // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
749 Scalar shift = mat.trace() / Scalar(2);
750 MatrixType scaledMat = mat;
751 scaledMat.coeffRef(0, 1) = mat.coeff(1, 0);
752 scaledMat.diagonal().array() -= shift;
753 Scalar scale = scaledMat.cwiseAbs().maxCoeff();
754 if (scale > Scalar(0)) scaledMat /= scale;
755
756 // Compute the eigenvalues
757 computeRoots(scaledMat, eivals);
758
759 // compute the eigen vectors
760 if (computeEigenvectors) {
761 if ((eivals(1) - eivals(0)) <= abs(eivals(1)) * Eigen::NumTraits<Scalar>::epsilon()) {
762 eivecs.setIdentity();
763 } else {
764 scaledMat.diagonal().array() -= eivals(1);
765 Scalar a2 = numext::abs2(scaledMat(0, 0));
766 Scalar c2 = numext::abs2(scaledMat(1, 1));
767 Scalar b2 = numext::abs2(scaledMat(1, 0));
768 if (a2 > c2) {
769 eivecs.col(1) << -scaledMat(1, 0), scaledMat(0, 0);
770 eivecs.col(1) /= sqrt(a2 + b2);
771 } else {
772 eivecs.col(1) << -scaledMat(1, 1), scaledMat(1, 0);
773 eivecs.col(1) /= sqrt(c2 + b2);
774 }
775
776 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
777 }
778 }
779
780 // Rescale back to the original size.
781 eivals *= scale;
782 eivals.array() += shift;
783
784 solver.m_info = Success;
785 solver.m_isInitialized = true;
786 solver.m_eigenvectorsOk = computeEigenvectors;
787 }
788};
789
790} // namespace internal
791
792template <typename MatrixType>
794 const MatrixType& matrix, int options) {
795 internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver, Size, NumTraits<Scalar>::IsComplex>::run(
796 *this, matrix, options);
797 return *this;
798}
799
800namespace internal {
801
802// Francis implicit QR step.
803template <int StorageOrder, typename RealScalar, typename Scalar, typename Index>
804EIGEN_DEVICE_FUNC static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end,
805 Scalar* matrixQ, Index n) {
806 // Wilkinson Shift.
807 RealScalar td = (diag[end - 1] - diag[end]) * RealScalar(0.5);
808 RealScalar e = subdiag[end - 1];
809 // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
810 // underflow thus leading to inf/NaN values when using the following commented code:
811 // RealScalar e2 = numext::abs2(subdiag[end-1]);
812 // RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
813 // This explain the following, somewhat more complicated, version:
814 RealScalar mu = diag[end];
815 if (numext::is_exactly_zero(td)) {
816 mu -= numext::abs(e);
817 } else if (!numext::is_exactly_zero(e)) {
818 const RealScalar e2 = numext::abs2(e);
819 const RealScalar h = numext::hypot(td, e);
820 if (numext::is_exactly_zero(e2)) {
821 mu -= e / ((td + (td > RealScalar(0) ? h : -h)) / e);
822 } else {
823 mu -= e2 / (td + (td > RealScalar(0) ? h : -h));
824 }
825 }
826
827 RealScalar x = diag[start] - mu;
828 RealScalar z = subdiag[start];
829 // If z ever becomes zero, the Givens rotation will be the identity and
830 // z will stay zero for all future iterations.
831 for (Index k = start; k < end && !numext::is_exactly_zero(z); ++k) {
832 JacobiRotation<RealScalar> rot;
833 rot.makeGivens(x, z);
834
835 // do T = G' T G
836 RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
837 RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k + 1];
838
839 diag[k] =
840 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;
843
844 if (k > start) subdiag[k - 1] = rot.c() * subdiag[k - 1] - rot.s() * z;
845
846 // "Chasing the bulge" to return to triangular form.
847 x = subdiag[k];
848 if (k < end - 1) {
849 z = -rot.s() * subdiag[k + 1];
850 subdiag[k + 1] = rot.c() * subdiag[k + 1];
851 }
852
853 // apply the givens rotation to the unit matrix Q = Q * G
854 if (matrixQ) {
855 // FIXME if StorageOrder == RowMajor this operation is not very efficient
856 Map<Matrix<Scalar, Dynamic, Dynamic, StorageOrder> > q(matrixQ, n, n);
857 q.applyOnTheRight(k, k + 1, rot);
858 }
859 }
860}
861
862} // end namespace internal
863
864} // end namespace Eigen
865
866#endif // EIGEN_SELFADJOINTEIGENSOLVER_H
const std::enable_if_t< std::is_same< typename LhsDerived::Scalar, typename RhsDerived::Scalar >::value, Eigen::CwiseBinaryOp< Eigen::internal::scalar_atan2_op< typename LhsDerived::Scalar, typename RhsDerived::Scalar >, const LhsDerived, const RhsDerived > > atan2(const Eigen::ArrayBase< LhsDerived > &x, const Eigen::ArrayBase< RhsDerived > &exponents)
Definition GlobalFunctions.h:209
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition GeneralizedSelfAdjointEigenSolver.h:51
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
MatrixType operatorExp() const
Computes the matrix exponential the matrix.
Definition SelfAdjointEigenSolver.h:338
SelfAdjointEigenSolver & computeFromTridiagonal(const RealVectorType &diag, const SubDiagonalType &subdiag, int options=ComputeEigenvectors)
Computes the eigen decomposition from a tridiagonal symmetric matrix.
Definition SelfAdjointEigenSolver.h:471
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType_.
Definition SelfAdjointEigenSolver.h:104
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
Definition SelfAdjointEigenSolver.h:362
SelfAdjointEigenSolver()
Default constructor for fixed-size matrices.
Definition SelfAdjointEigenSolver.h:128
internal::plain_col_type< MatrixType, Scalar >::type VectorType
Type for vector of eigenvalues as returned by eigenvalues().
Definition SelfAdjointEigenSolver.h:113
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SelfAdjointEigenSolver.h:372
Eigen::Index Index
Definition SelfAdjointEigenSolver.h:94
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType_.
Definition SelfAdjointEigenSolver.h:93
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
Definition SelfAdjointEigenSolver.h:322
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition SelfAdjointEigenSolver.h:150
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.
Definition SelfAdjointEigenSolver.h:300
static const int m_maxIterations
Maximum number of iterations.
Definition SelfAdjointEigenSolver.h:382
const EigenvectorsType & eigenvectors() const
Returns the eigenvectors of given matrix.
Definition SelfAdjointEigenSolver.h:279
SelfAdjointEigenSolver(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Definition SelfAdjointEigenSolver.h:175
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a closed-form algorithm.
Definition SelfAdjointEigenSolver.h:793
constexpr const Scalar * data() const
Definition Transform.h:591
Tridiagonal decomposition of a selfadjoint matrix.
Definition Tridiagonalization.h:66
ComputationInfo
Definition Constants.h:438
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
@ ComputeEigenvectors
Definition Constants.h:401
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
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_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sin_op< typename Derived::Scalar >, const Derived > sin(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33
constexpr Derived & derived()
Definition EigenBase.h:49
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:232