xref: /aosp_15_r20/external/eigen/Eigen/src/Eigenvalues/MatrixBaseEigenvalues.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2010 Jitse Niesen <[email protected]>
6*bf2c3715SXin Li //
7*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
8*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
9*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10*bf2c3715SXin Li 
11*bf2c3715SXin Li #ifndef EIGEN_MATRIXBASEEIGENVALUES_H
12*bf2c3715SXin Li #define EIGEN_MATRIXBASEEIGENVALUES_H
13*bf2c3715SXin Li 
14*bf2c3715SXin Li namespace Eigen {
15*bf2c3715SXin Li 
16*bf2c3715SXin Li namespace internal {
17*bf2c3715SXin Li 
18*bf2c3715SXin Li template<typename Derived, bool IsComplex>
19*bf2c3715SXin Li struct eigenvalues_selector
20*bf2c3715SXin Li {
21*bf2c3715SXin Li   // this is the implementation for the case IsComplex = true
22*bf2c3715SXin Li   static inline typename MatrixBase<Derived>::EigenvaluesReturnType const
runeigenvalues_selector23*bf2c3715SXin Li   run(const MatrixBase<Derived>& m)
24*bf2c3715SXin Li   {
25*bf2c3715SXin Li     typedef typename Derived::PlainObject PlainObject;
26*bf2c3715SXin Li     PlainObject m_eval(m);
27*bf2c3715SXin Li     return ComplexEigenSolver<PlainObject>(m_eval, false).eigenvalues();
28*bf2c3715SXin Li   }
29*bf2c3715SXin Li };
30*bf2c3715SXin Li 
31*bf2c3715SXin Li template<typename Derived>
32*bf2c3715SXin Li struct eigenvalues_selector<Derived, false>
33*bf2c3715SXin Li {
34*bf2c3715SXin Li   static inline typename MatrixBase<Derived>::EigenvaluesReturnType const
35*bf2c3715SXin Li   run(const MatrixBase<Derived>& m)
36*bf2c3715SXin Li   {
37*bf2c3715SXin Li     typedef typename Derived::PlainObject PlainObject;
38*bf2c3715SXin Li     PlainObject m_eval(m);
39*bf2c3715SXin Li     return EigenSolver<PlainObject>(m_eval, false).eigenvalues();
40*bf2c3715SXin Li   }
41*bf2c3715SXin Li };
42*bf2c3715SXin Li 
43*bf2c3715SXin Li } // end namespace internal
44*bf2c3715SXin Li 
45*bf2c3715SXin Li /** \brief Computes the eigenvalues of a matrix
46*bf2c3715SXin Li   * \returns Column vector containing the eigenvalues.
47*bf2c3715SXin Li   *
48*bf2c3715SXin Li   * \eigenvalues_module
49*bf2c3715SXin Li   * This function computes the eigenvalues with the help of the EigenSolver
50*bf2c3715SXin Li   * class (for real matrices) or the ComplexEigenSolver class (for complex
51*bf2c3715SXin Li   * matrices).
52*bf2c3715SXin Li   *
53*bf2c3715SXin Li   * The eigenvalues are repeated according to their algebraic multiplicity,
54*bf2c3715SXin Li   * so there are as many eigenvalues as rows in the matrix.
55*bf2c3715SXin Li   *
56*bf2c3715SXin Li   * The SelfAdjointView class provides a better algorithm for selfadjoint
57*bf2c3715SXin Li   * matrices.
58*bf2c3715SXin Li   *
59*bf2c3715SXin Li   * Example: \include MatrixBase_eigenvalues.cpp
60*bf2c3715SXin Li   * Output: \verbinclude MatrixBase_eigenvalues.out
61*bf2c3715SXin Li   *
62*bf2c3715SXin Li   * \sa EigenSolver::eigenvalues(), ComplexEigenSolver::eigenvalues(),
63*bf2c3715SXin Li   *     SelfAdjointView::eigenvalues()
64*bf2c3715SXin Li   */
65*bf2c3715SXin Li template<typename Derived>
66*bf2c3715SXin Li inline typename MatrixBase<Derived>::EigenvaluesReturnType
67*bf2c3715SXin Li MatrixBase<Derived>::eigenvalues() const
68*bf2c3715SXin Li {
69*bf2c3715SXin Li   return internal::eigenvalues_selector<Derived, NumTraits<Scalar>::IsComplex>::run(derived());
70*bf2c3715SXin Li }
71*bf2c3715SXin Li 
72*bf2c3715SXin Li /** \brief Computes the eigenvalues of a matrix
73*bf2c3715SXin Li   * \returns Column vector containing the eigenvalues.
74*bf2c3715SXin Li   *
75*bf2c3715SXin Li   * \eigenvalues_module
76*bf2c3715SXin Li   * This function computes the eigenvalues with the help of the
77*bf2c3715SXin Li   * SelfAdjointEigenSolver class.  The eigenvalues are repeated according to
78*bf2c3715SXin Li   * their algebraic multiplicity, so there are as many eigenvalues as rows in
79*bf2c3715SXin Li   * the matrix.
80*bf2c3715SXin Li   *
81*bf2c3715SXin Li   * Example: \include SelfAdjointView_eigenvalues.cpp
82*bf2c3715SXin Li   * Output: \verbinclude SelfAdjointView_eigenvalues.out
83*bf2c3715SXin Li   *
84*bf2c3715SXin Li   * \sa SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()
85*bf2c3715SXin Li   */
86*bf2c3715SXin Li template<typename MatrixType, unsigned int UpLo>
87*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename SelfAdjointView<MatrixType, UpLo>::EigenvaluesReturnType
88*bf2c3715SXin Li SelfAdjointView<MatrixType, UpLo>::eigenvalues() const
89*bf2c3715SXin Li {
90*bf2c3715SXin Li   PlainObject thisAsMatrix(*this);
91*bf2c3715SXin Li   return SelfAdjointEigenSolver<PlainObject>(thisAsMatrix, false).eigenvalues();
92*bf2c3715SXin Li }
93*bf2c3715SXin Li 
94*bf2c3715SXin Li 
95*bf2c3715SXin Li 
96*bf2c3715SXin Li /** \brief Computes the L2 operator norm
97*bf2c3715SXin Li   * \returns Operator norm of the matrix.
98*bf2c3715SXin Li   *
99*bf2c3715SXin Li   * \eigenvalues_module
100*bf2c3715SXin Li   * This function computes the L2 operator norm of a matrix, which is also
101*bf2c3715SXin Li   * known as the spectral norm. The norm of a matrix \f$ A \f$ is defined to be
102*bf2c3715SXin Li   * \f[ \|A\|_2 = \max_x \frac{\|Ax\|_2}{\|x\|_2} \f]
103*bf2c3715SXin Li   * where the maximum is over all vectors and the norm on the right is the
104*bf2c3715SXin Li   * Euclidean vector norm. The norm equals the largest singular value, which is
105*bf2c3715SXin Li   * the square root of the largest eigenvalue of the positive semi-definite
106*bf2c3715SXin Li   * matrix \f$ A^*A \f$.
107*bf2c3715SXin Li   *
108*bf2c3715SXin Li   * The current implementation uses the eigenvalues of \f$ A^*A \f$, as computed
109*bf2c3715SXin Li   * by SelfAdjointView::eigenvalues(), to compute the operator norm of a
110*bf2c3715SXin Li   * matrix.  The SelfAdjointView class provides a better algorithm for
111*bf2c3715SXin Li   * selfadjoint matrices.
112*bf2c3715SXin Li   *
113*bf2c3715SXin Li   * Example: \include MatrixBase_operatorNorm.cpp
114*bf2c3715SXin Li   * Output: \verbinclude MatrixBase_operatorNorm.out
115*bf2c3715SXin Li   *
116*bf2c3715SXin Li   * \sa SelfAdjointView::eigenvalues(), SelfAdjointView::operatorNorm()
117*bf2c3715SXin Li   */
118*bf2c3715SXin Li template<typename Derived>
119*bf2c3715SXin Li inline typename MatrixBase<Derived>::RealScalar
120*bf2c3715SXin Li MatrixBase<Derived>::operatorNorm() const
121*bf2c3715SXin Li {
122*bf2c3715SXin Li   using std::sqrt;
123*bf2c3715SXin Li   typename Derived::PlainObject m_eval(derived());
124*bf2c3715SXin Li   // FIXME if it is really guaranteed that the eigenvalues are already sorted,
125*bf2c3715SXin Li   // then we don't need to compute a maxCoeff() here, comparing the 1st and last ones is enough.
126*bf2c3715SXin Li   return sqrt((m_eval*m_eval.adjoint())
127*bf2c3715SXin Li                  .eval()
128*bf2c3715SXin Li 		 .template selfadjointView<Lower>()
129*bf2c3715SXin Li 		 .eigenvalues()
130*bf2c3715SXin Li 		 .maxCoeff()
131*bf2c3715SXin Li 		 );
132*bf2c3715SXin Li }
133*bf2c3715SXin Li 
134*bf2c3715SXin Li /** \brief Computes the L2 operator norm
135*bf2c3715SXin Li   * \returns Operator norm of the matrix.
136*bf2c3715SXin Li   *
137*bf2c3715SXin Li   * \eigenvalues_module
138*bf2c3715SXin Li   * This function computes the L2 operator norm of a self-adjoint matrix. For a
139*bf2c3715SXin Li   * self-adjoint matrix, the operator norm is the largest eigenvalue.
140*bf2c3715SXin Li   *
141*bf2c3715SXin Li   * The current implementation uses the eigenvalues of the matrix, as computed
142*bf2c3715SXin Li   * by eigenvalues(), to compute the operator norm of the matrix.
143*bf2c3715SXin Li   *
144*bf2c3715SXin Li   * Example: \include SelfAdjointView_operatorNorm.cpp
145*bf2c3715SXin Li   * Output: \verbinclude SelfAdjointView_operatorNorm.out
146*bf2c3715SXin Li   *
147*bf2c3715SXin Li   * \sa eigenvalues(), MatrixBase::operatorNorm()
148*bf2c3715SXin Li   */
149*bf2c3715SXin Li template<typename MatrixType, unsigned int UpLo>
150*bf2c3715SXin Li EIGEN_DEVICE_FUNC inline typename SelfAdjointView<MatrixType, UpLo>::RealScalar
151*bf2c3715SXin Li SelfAdjointView<MatrixType, UpLo>::operatorNorm() const
152*bf2c3715SXin Li {
153*bf2c3715SXin Li   return eigenvalues().cwiseAbs().maxCoeff();
154*bf2c3715SXin Li }
155*bf2c3715SXin Li 
156*bf2c3715SXin Li } // end namespace Eigen
157*bf2c3715SXin Li 
158*bf2c3715SXin Li #endif
159