11 #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H 12 #define EIGEN_SELFADJOINTEIGENSOLVER_H 18 template<
typename _MatrixType>
19 class GeneralizedSelfAdjointEigenSolver;
74 Size = MatrixType::RowsAtCompileTime,
75 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
77 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
81 typedef typename MatrixType::Scalar
Scalar;
82 typedef typename MatrixType::Index
Index;
116 m_isInitialized(false)
132 : m_eivec(size, size),
134 m_subdiag(size > 1 ? size - 1 : 1),
135 m_isInitialized(false)
154 : m_eivec(matrix.rows(), matrix.cols()),
155 m_eivalues(matrix.cols()),
156 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
157 m_isInitialized(false)
230 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
231 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
252 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
276 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
277 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
278 return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
301 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
302 eigen_assert(m_eigenvectorsOk &&
"The eigenvectors have not been computed together with the eigenvalues.");
303 return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
312 eigen_assert(m_isInitialized &&
"SelfAdjointEigenSolver is not initialized.");
321 static const int m_maxIterations = 30;
323 #ifdef EIGEN2_SUPPORT 325 : m_eivec(matrix.rows(), matrix.cols()),
326 m_eivalues(matrix.cols()),
327 m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
328 m_isInitialized(false)
330 compute(matrix, computeEigenvectors);
334 : m_eivec(matA.cols(), matA.cols()),
335 m_eivalues(matA.cols()),
336 m_subdiag(matA.cols() > 1 ? matA.cols() - 1 : 1),
337 m_isInitialized(
false)
342 void compute(
const MatrixType& matrix,
bool computeEigenvectors)
347 void compute(
const MatrixType& matA,
const MatrixType& matB,
bool computeEigenvectors =
true)
351 #endif // EIGEN2_SUPPORT 379 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
383 template<
typename MatrixType>
391 &&
"invalid option parameter");
393 Index n = matrix.cols();
394 m_eivalues.resize(n,1);
398 m_eivalues.coeffRef(0,0) =
numext::real(matrix.coeff(0,0));
399 if(computeEigenvectors)
400 m_eivec.setOnes(n,n);
402 m_isInitialized =
true;
403 m_eigenvectorsOk = computeEigenvectors;
412 mat = matrix.template triangularView<Lower>();
415 mat.template triangularView<Lower>() /= scale;
416 m_subdiag.resize(n-1);
425 for (
Index i = start; i<end; ++i)
430 while (end>0 && m_subdiag[end-1]==0)
439 if(iter > m_maxIterations * n)
break;
442 while (start>0 && m_subdiag[start-1]!=0)
445 internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), m_subdiag.data(), start, end, computeEigenvectors ? m_eivec.data() : (
Scalar*)0, n);
448 if (iter <= m_maxIterations * n)
458 for (
Index i = 0; i < n-1; ++i)
461 m_eivalues.segment(i,n-i).minCoeff(&k);
464 std::swap(m_eivalues[i], m_eivalues[k+i]);
465 if(computeEigenvectors)
466 m_eivec.col(i).swap(m_eivec.col(k+i));
474 m_isInitialized =
true;
475 m_eigenvectorsOk = computeEigenvectors;
482 template<
typename SolverType,
int Size,
bool IsComplex>
struct direct_selfadjoint_eigenvalues
484 static inline void run(SolverType& eig,
const typename SolverType::MatrixType&
A,
int options)
485 { eig.compute(A,options); }
492 typedef typename SolverType::Scalar
Scalar;
494 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
500 const Scalar s_inv3 = Scalar(1.0)/Scalar(3.0);
501 const Scalar s_sqrt3 =
sqrt(Scalar(3.0));
506 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);
507 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);
508 Scalar c2 = m(0,0) + m(1,1) + m(2,2);
512 Scalar c2_over_3 = c2*s_inv3;
513 Scalar a_over_3 = (c1 - c2*c2_over_3)*s_inv3;
514 if (a_over_3 > Scalar(0))
515 a_over_3 = Scalar(0);
517 Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
519 Scalar q = half_b*half_b + a_over_3*a_over_3*a_over_3;
525 Scalar theta =
atan2(
sqrt(-q),half_b)*s_inv3;
526 Scalar cos_theta =
cos(theta);
527 Scalar sin_theta =
sin(theta);
528 roots(0) = c2_over_3 + Scalar(2)*rho*cos_theta;
529 roots(1) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
530 roots(2) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
533 if (roots(0) >= roots(1))
534 std::swap(roots(0),roots(1));
535 if (roots(1) >= roots(2))
537 std::swap(roots(1),roots(2));
538 if (roots(0) >= roots(1))
539 std::swap(roots(0),roots(1));
543 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
546 eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
549 &&
"invalid option parameter");
552 MatrixType& eivecs = solver.m_eivec;
553 VectorType& eivals = solver.m_eivalues;
556 Scalar scale = mat.cwiseAbs().maxCoeff();
557 MatrixType scaledMat = mat / scale;
560 computeRoots(scaledMat,eivals);
563 if(computeEigenvectors)
566 safeNorm2 *= safeNorm2;
569 eivecs.setIdentity();
573 scaledMat = scaledMat.template selfadjointView<Lower>();
577 Scalar d0 = eivals(2) - eivals(1);
578 Scalar d1 = eivals(1) - eivals(0);
579 int k = d0 > d1 ? 2 : 0;
580 d0 = d0 > d1 ? d1 : d0;
582 tmp.diagonal().array () -= eivals(k);
585 n = (cross = tmp.row(0).cross(tmp.row(1))).squaredNorm();
588 eivecs.col(k) = cross /
sqrt(n);
591 n = (cross = tmp.row(0).cross(tmp.row(2))).squaredNorm();
594 eivecs.col(k) = cross /
sqrt(n);
597 n = (cross = tmp.row(1).cross(tmp.row(2))).squaredNorm();
600 eivecs.col(k) = cross /
sqrt(n);
609 solver.m_isInitialized =
true;
610 solver.m_eigenvectorsOk = computeEigenvectors;
617 tmp.diagonal().array() -= eivals(1);
620 eivecs.col(1) = eivecs.col(k).unitOrthogonal();
623 n = (cross = eivecs.col(k).cross(tmp.row(0).normalized())).squaredNorm();
625 eivecs.col(1) = cross /
sqrt(n);
628 n = (cross = eivecs.col(k).cross(tmp.row(1))).squaredNorm();
630 eivecs.col(1) = cross /
sqrt(n);
633 n = (cross = eivecs.col(k).cross(tmp.row(2))).squaredNorm();
635 eivecs.col(1) = cross /
sqrt(n);
640 eivecs.col(1) = eivecs.col(k).unitOrthogonal();
646 Scalar d = eivecs.col(1).dot(eivecs.col(k));
647 eivecs.col(1) = (eivecs.col(1) - d * eivecs.col(k)).normalized();
650 eivecs.col(k==2 ? 0 : 2) = eivecs.col(k).cross(eivecs.col(1)).normalized();
657 solver.m_isInitialized =
true;
658 solver.m_eigenvectorsOk = computeEigenvectors;
667 typedef typename SolverType::Scalar
Scalar;
669 static inline void computeRoots(
const MatrixType& m, VectorType& roots)
672 const Scalar t0 = Scalar(0.5) *
sqrt(
numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*m(1,0)*m(1,0));
673 const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
678 static inline void run(SolverType& solver,
const MatrixType& mat,
int options)
681 eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
684 &&
"invalid option parameter");
687 MatrixType& eivecs = solver.m_eivec;
688 VectorType& eivals = solver.m_eivalues;
691 Scalar scale = mat.cwiseAbs().maxCoeff();
692 scale = (std::max)(scale,Scalar(1));
693 MatrixType scaledMat = mat / scale;
696 computeRoots(scaledMat,eivals);
699 if(computeEigenvectors)
701 scaledMat.diagonal().array () -= eivals(1);
707 eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
708 eivecs.col(1) /=
sqrt(a2+b2);
712 eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
713 eivecs.col(1) /=
sqrt(c2+b2);
716 eivecs.col(0) << eivecs.col(1).unitOrthogonal();
723 solver.m_isInitialized =
true;
724 solver.m_eigenvectorsOk = computeEigenvectors;
730 template<
typename MatrixType>
739 template<
int StorageOrder,
typename RealScalar,
typename Scalar,
typename Index>
740 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
743 RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
744 RealScalar e = subdiag[end-1];
750 RealScalar mu = diag[end];
756 RealScalar h = numext::hypot(td,e);
757 if(e2==0) mu -= (e / (td + (td>0 ? 1 : -1))) * (e / h);
758 else mu -= e2 / (td + (td>0 ? h : -h));
761 RealScalar x = diag[start] - mu;
762 RealScalar z = subdiag[start];
763 for (Index k = start; k < end; ++k)
769 RealScalar sdk = rot.
s() * diag[k] + rot.
c() * subdiag[k];
770 RealScalar dkp1 = rot.
s() * subdiag[k] + rot.
c() * diag[k+1];
772 diag[k] = rot.
c() * (rot.
c() * diag[k] - rot.
s() * subdiag[k]) - rot.
s() * (rot.
c() * subdiag[k] - rot.
s() * diag[k+1]);
773 diag[k+1] = rot.
s() * sdk + rot.
c() * dkp1;
774 subdiag[k] = rot.
c() * sdk - rot.
s() * dkp1;
778 subdiag[k - 1] = rot.
c() * subdiag[k-1] - rot.
s() * z;
784 z = -rot.
s() * subdiag[k+1];
785 subdiag[k + 1] = rot.
c() * subdiag[k+1];
793 q.applyOnTheRight(k,k+1,rot);
802 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
static void computeRoots(const MatrixType &m, VectorType &roots)
MatrixType operatorSqrt() const
Computes the positive-definite square root of the matrix.
USING_NAMESPACE_ACADO IntermediateState sin(const Expression &arg)
IntermediateState sqrt(const Expression &arg)
SolverType::RealVectorType VectorType
static void run(SolverType &solver, const MatrixType &mat, int options)
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
SolverType::RealVectorType VectorType
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
SolverType::MatrixType MatrixType
A matrix or vector expression mapping an existing array of data.
static void run(SolverType &solver, const MatrixType &mat, int options)
const AutoDiffScalar< Matrix< typename internal::traits< DerTypeA >::Scalar, Dynamic, 1 > > atan2(const AutoDiffScalar< DerTypeA > &a, const AutoDiffScalar< DerTypeB > &b)
Computes eigenvalues and eigenvectors of selfadjoint matrices.
SolverType::Scalar Scalar
iterative scaling algorithm to equilibrate rows and column norms in matrices
SelfAdjointEigenSolver(const MatrixType &matrix, int options=ComputeEigenvectors)
Constructor; computes eigendecomposition of given matrix.
Rotation given by a cosine-sine pair.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Tridiagonal decomposition of a selfadjoint matrix.
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, typename NumTraits< Scalar >::Real precision=NumTraits< Scalar >::dummy_precision())
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs2_op< Scalar >, const Derived > abs2() const
IntermediateState cos(const Expression &arg)
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
RealReturnType real() const
SolverType::MatrixType MatrixType
ComputationInfo info() const
Reports whether previous computation was successful.
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Provides a generic way to set and pass user-specified options.
MatrixType operatorInverseSqrt() const
Computes the inverse square root of the matrix.
static void tridiagonal_qr_step(RealScalar *diag, RealScalar *subdiag, Index start, Index end, Scalar *matrixQ, Index n)
RealVectorType m_eivalues
const MatrixType & eigenvectors() const
Returns the eigenvectors of given matrix.
NumTraits< Scalar >::Real RealScalar
Real scalar type for _MatrixType.
SelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
void tridiagonalization_inplace(MatrixType &matA, CoeffVectorType &hCoeffs)
SelfAdjointEigenSolver & compute(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
static void computeRoots(const MatrixType &m, VectorType &roots)
SelfAdjointEigenSolver & computeDirect(const MatrixType &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix using a direct algorithm.
static void run(SolverType &eig, const typename SolverType::MatrixType &A, int options)
SolverType::Scalar Scalar
The matrix class, also used for vectors and row-vectors.
TridiagonalizationType::SubDiagonalType m_subdiag
const RealVectorType & eigenvalues() const
Returns the eigenvalues of given matrix.