1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5
// Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
7
// Eigen is free software; you can redistribute it and/or
8
// modify it under the terms of the GNU Lesser General Public
9
// License as published by the Free Software Foundation; either
10
// version 3 of the License, or (at your option) any later version.
12
// Alternatively, you can redistribute it and/or
13
// modify it under the terms of the GNU General Public License as
14
// published by the Free Software Foundation; either version 2 of
15
// the License, or (at your option) any later version.
17
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20
// GNU General Public License for more details.
22
// You should have received a copy of the GNU Lesser General Public
23
// License and a copy of the GNU General Public License along with
24
// Eigen. If not, see <http://www.gnu.org/licenses/>.
26
#ifndef EIGEN_MATRIXBASEEIGENVALUES_H
27
#define EIGEN_MATRIXBASEEIGENVALUES_H
31
template<typename Derived, bool IsComplex>
32
struct eigenvalues_selector
34
// this is the implementation for the case IsComplex = true
35
static inline typename MatrixBase<Derived>::EigenvaluesReturnType const
36
run(const MatrixBase<Derived>& m)
38
typedef typename Derived::PlainObject PlainObject;
39
PlainObject m_eval(m);
40
return ComplexEigenSolver<PlainObject>(m_eval, false).eigenvalues();
44
template<typename Derived>
45
struct eigenvalues_selector<Derived, false>
47
static inline typename MatrixBase<Derived>::EigenvaluesReturnType const
48
run(const MatrixBase<Derived>& m)
50
typedef typename Derived::PlainObject PlainObject;
51
PlainObject m_eval(m);
52
return EigenSolver<PlainObject>(m_eval, false).eigenvalues();
56
} // end namespace internal
58
/** \brief Computes the eigenvalues of a matrix
59
* \returns Column vector containing the eigenvalues.
62
* This function computes the eigenvalues with the help of the EigenSolver
63
* class (for real matrices) or the ComplexEigenSolver class (for complex
66
* The eigenvalues are repeated according to their algebraic multiplicity,
67
* so there are as many eigenvalues as rows in the matrix.
69
* The SelfAdjointView class provides a better algorithm for selfadjoint
72
* Example: \include MatrixBase_eigenvalues.cpp
73
* Output: \verbinclude MatrixBase_eigenvalues.out
75
* \sa EigenSolver::eigenvalues(), ComplexEigenSolver::eigenvalues(),
76
* SelfAdjointView::eigenvalues()
78
template<typename Derived>
79
inline typename MatrixBase<Derived>::EigenvaluesReturnType
80
MatrixBase<Derived>::eigenvalues() const
82
typedef typename internal::traits<Derived>::Scalar Scalar;
83
return internal::eigenvalues_selector<Derived, NumTraits<Scalar>::IsComplex>::run(derived());
86
/** \brief Computes the eigenvalues of a matrix
87
* \returns Column vector containing the eigenvalues.
90
* This function computes the eigenvalues with the help of the
91
* SelfAdjointEigenSolver class. The eigenvalues are repeated according to
92
* their algebraic multiplicity, so there are as many eigenvalues as rows in
95
* Example: \include SelfAdjointView_eigenvalues.cpp
96
* Output: \verbinclude SelfAdjointView_eigenvalues.out
98
* \sa SelfAdjointEigenSolver::eigenvalues(), MatrixBase::eigenvalues()
100
template<typename MatrixType, unsigned int UpLo>
101
inline typename SelfAdjointView<MatrixType, UpLo>::EigenvaluesReturnType
102
SelfAdjointView<MatrixType, UpLo>::eigenvalues() const
104
typedef typename SelfAdjointView<MatrixType, UpLo>::PlainObject PlainObject;
105
PlainObject thisAsMatrix(*this);
106
return SelfAdjointEigenSolver<PlainObject>(thisAsMatrix, false).eigenvalues();
111
/** \brief Computes the L2 operator norm
112
* \returns Operator norm of the matrix.
114
* \eigenvalues_module
115
* This function computes the L2 operator norm of a matrix, which is also
116
* known as the spectral norm. The norm of a matrix \f$ A \f$ is defined to be
117
* \f[ \|A\|_2 = \max_x \frac{\|Ax\|_2}{\|x\|_2} \f]
118
* where the maximum is over all vectors and the norm on the right is the
119
* Euclidean vector norm. The norm equals the largest singular value, which is
120
* the square root of the largest eigenvalue of the positive semi-definite
121
* matrix \f$ A^*A \f$.
123
* The current implementation uses the eigenvalues of \f$ A^*A \f$, as computed
124
* by SelfAdjointView::eigenvalues(), to compute the operator norm of a
125
* matrix. The SelfAdjointView class provides a better algorithm for
126
* selfadjoint matrices.
128
* Example: \include MatrixBase_operatorNorm.cpp
129
* Output: \verbinclude MatrixBase_operatorNorm.out
131
* \sa SelfAdjointView::eigenvalues(), SelfAdjointView::operatorNorm()
133
template<typename Derived>
134
inline typename MatrixBase<Derived>::RealScalar
135
MatrixBase<Derived>::operatorNorm() const
137
typename Derived::PlainObject m_eval(derived());
138
// FIXME if it is really guaranteed that the eigenvalues are already sorted,
139
// then we don't need to compute a maxCoeff() here, comparing the 1st and last ones is enough.
140
return internal::sqrt((m_eval*m_eval.adjoint())
142
.template selfadjointView<Lower>()
148
/** \brief Computes the L2 operator norm
149
* \returns Operator norm of the matrix.
151
* \eigenvalues_module
152
* This function computes the L2 operator norm of a self-adjoint matrix. For a
153
* self-adjoint matrix, the operator norm is the largest eigenvalue.
155
* The current implementation uses the eigenvalues of the matrix, as computed
156
* by eigenvalues(), to compute the operator norm of the matrix.
158
* Example: \include SelfAdjointView_operatorNorm.cpp
159
* Output: \verbinclude SelfAdjointView_operatorNorm.out
161
* \sa eigenvalues(), MatrixBase::operatorNorm()
163
template<typename MatrixType, unsigned int UpLo>
164
inline typename SelfAdjointView<MatrixType, UpLo>::RealScalar
165
SelfAdjointView<MatrixType, UpLo>::operatorNorm() const
167
return eigenvalues().cwiseAbs().maxCoeff();