SelfAdjointEigenSolver.h 32.9 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5 6
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
// Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
//
Don Gagne's avatar
Don Gagne committed
7 8 9
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
10 11 12 13 14 15

#ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
#define EIGEN_SELFADJOINTEIGENSOLVER_H

#include "./Tridiagonalization.h"

Don Gagne's avatar
Don Gagne committed
16 17
namespace Eigen { 

LM's avatar
LM committed
18 19 20
template<typename _MatrixType>
class GeneralizedSelfAdjointEigenSolver;

Don Gagne's avatar
Don Gagne committed
21 22
namespace internal {
template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
23 24
template<typename MatrixType, typename DiagType, typename SubDiagType>
ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
Don Gagne's avatar
Don Gagne committed
25 26
}

LM's avatar
LM committed
27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80
/** \eigenvalues_module \ingroup Eigenvalues_Module
  *
  *
  * \class SelfAdjointEigenSolver
  *
  * \brief Computes eigenvalues and eigenvectors of selfadjoint matrices
  *
  * \tparam _MatrixType the type of the matrix of which we are computing the
  * eigendecomposition; this is expected to be an instantiation of the Matrix
  * class template.
  *
  * A matrix \f$ A \f$ is selfadjoint if it equals its adjoint. For real
  * matrices, this means that the matrix is symmetric: it equals its
  * transpose. This class computes the eigenvalues and eigenvectors of a
  * selfadjoint matrix. These are the scalars \f$ \lambda \f$ and vectors
  * \f$ v \f$ such that \f$ Av = \lambda v \f$.  The eigenvalues of a
  * selfadjoint matrix are always real. If \f$ D \f$ is a diagonal matrix with
  * the eigenvalues on the diagonal, and \f$ V \f$ is a matrix with the
  * eigenvectors as its columns, then \f$ A = V D V^{-1} \f$ (for selfadjoint
  * matrices, the matrix \f$ V \f$ is always invertible). This is called the
  * eigendecomposition.
  *
  * The algorithm exploits the fact that the matrix is selfadjoint, making it
  * faster and more accurate than the general purpose eigenvalue algorithms
  * implemented in EigenSolver and ComplexEigenSolver.
  *
  * Only the \b lower \b triangular \b part of the input matrix is referenced.
  *
  * Call the function compute() to compute the eigenvalues and eigenvectors of
  * a given matrix. Alternatively, you can use the
  * SelfAdjointEigenSolver(const MatrixType&, int) constructor which computes
  * the eigenvalues and eigenvectors at construction time. Once the eigenvalue
  * and eigenvectors are computed, they can be retrieved with the eigenvalues()
  * and eigenvectors() functions.
  *
  * The documentation for SelfAdjointEigenSolver(const MatrixType&, int)
  * contains an example of the typical use of this class.
  *
  * To solve the \em generalized eigenvalue problem \f$ Av = \lambda Bv \f$ and
  * the likes, see the class GeneralizedSelfAdjointEigenSolver.
  *
  * \sa MatrixBase::eigenvalues(), class EigenSolver, class ComplexEigenSolver
  */
template<typename _MatrixType> class SelfAdjointEigenSolver
{
  public:

    typedef _MatrixType MatrixType;
    enum {
      Size = MatrixType::RowsAtCompileTime,
      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
      Options = MatrixType::Options,
      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
    };
Don Gagne's avatar
Don Gagne committed
81
    
LM's avatar
LM committed
82 83
    /** \brief Scalar type for matrices of type \p _MatrixType. */
    typedef typename MatrixType::Scalar Scalar;
84
    typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
85 86
    
    typedef Matrix<Scalar,Size,Size,ColMajor,MaxColsAtCompileTime,MaxColsAtCompileTime> EigenvectorsType;
LM's avatar
LM committed
87 88 89 90 91 92 93 94

    /** \brief Real scalar type for \p _MatrixType.
      *
      * This is just \c Scalar if #Scalar is real (e.g., \c float or
      * \c double), and the type of the real part of \c Scalar if #Scalar is
      * complex.
      */
    typedef typename NumTraits<Scalar>::Real RealScalar;
Don Gagne's avatar
Don Gagne committed
95 96
    
    friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
LM's avatar
LM committed
97 98 99 100 101 102 103 104

    /** \brief Type for vector of eigenvalues as returned by eigenvalues().
      *
      * This is a column vector with entries of type #RealScalar.
      * The length of the vector is the size of \p _MatrixType.
      */
    typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
    typedef Tridiagonalization<MatrixType> TridiagonalizationType;
105
    typedef typename TridiagonalizationType::SubDiagonalType SubDiagonalType;
LM's avatar
LM committed
106 107 108 109 110 111 112 113 114 115 116

    /** \brief Default constructor for fixed-size matrices.
      *
      * The default constructor is useful in cases in which the user intends to
      * perform decompositions via compute(). This constructor
      * can only be used if \p _MatrixType is a fixed-size matrix; use
      * SelfAdjointEigenSolver(Index) for dynamic-size matrices.
      *
      * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver.out
      */
117
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136
    SelfAdjointEigenSolver()
        : m_eivec(),
          m_eivalues(),
          m_subdiag(),
          m_isInitialized(false)
    { }

    /** \brief Constructor, pre-allocates memory for dynamic-size matrices.
      *
      * \param [in]  size  Positive integer, size of the matrix whose
      * eigenvalues and eigenvectors will be computed.
      *
      * This constructor is useful for dynamic-size matrices, when the user
      * intends to perform decompositions via compute(). The \p size
      * parameter is only used as a hint. It is not an error to give a wrong
      * \p size, but it may impair performance.
      *
      * \sa compute() for an example
      */
137 138
    EIGEN_DEVICE_FUNC
    explicit SelfAdjointEigenSolver(Index size)
LM's avatar
LM committed
139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159
        : m_eivec(size, size),
          m_eivalues(size),
          m_subdiag(size > 1 ? size - 1 : 1),
          m_isInitialized(false)
    {}

    /** \brief Constructor; computes eigendecomposition of given matrix.
      *
      * \param[in]  matrix  Selfadjoint matrix whose eigendecomposition is to
      *    be computed. Only the lower triangular part of the matrix is referenced.
      * \param[in]  options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
      *
      * This constructor calls compute(const MatrixType&, int) to compute the
      * eigenvalues of the matrix \p matrix. The eigenvectors are computed if
      * \p options equals #ComputeEigenvectors.
      *
      * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType.out
      *
      * \sa compute(const MatrixType&, int)
      */
160 161 162
    template<typename InputType>
    EIGEN_DEVICE_FUNC
    explicit SelfAdjointEigenSolver(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors)
LM's avatar
LM committed
163 164 165 166 167
      : m_eivec(matrix.rows(), matrix.cols()),
        m_eivalues(matrix.cols()),
        m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
        m_isInitialized(false)
    {
168
      compute(matrix.derived(), options);
LM's avatar
LM committed
169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200
    }

    /** \brief Computes eigendecomposition of given matrix.
      *
      * \param[in]  matrix  Selfadjoint matrix whose eigendecomposition is to
      *    be computed. Only the lower triangular part of the matrix is referenced.
      * \param[in]  options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
      * \returns    Reference to \c *this
      *
      * This function computes the eigenvalues of \p matrix.  The eigenvalues()
      * function can be used to retrieve them.  If \p options equals #ComputeEigenvectors,
      * then the eigenvectors are also computed and can be retrieved by
      * calling eigenvectors().
      *
      * This implementation uses a symmetric QR algorithm. The matrix is first
      * reduced to tridiagonal form using the Tridiagonalization class. The
      * tridiagonal matrix is then brought to diagonal form with implicit
      * symmetric QR steps with Wilkinson shift. Details can be found in
      * Section 8.3 of Golub \& Van Loan, <i>%Matrix Computations</i>.
      *
      * The cost of the computation is about \f$ 9n^3 \f$ if the eigenvectors
      * are required and \f$ 4n^3/3 \f$ if they are not required.
      *
      * This method reuses the memory in the SelfAdjointEigenSolver object that
      * was allocated when the object was constructed, if the size of the
      * matrix does not change.
      *
      * Example: \include SelfAdjointEigenSolver_compute_MatrixType.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_compute_MatrixType.out
      *
      * \sa SelfAdjointEigenSolver(const MatrixType&, int)
      */
201 202 203
    template<typename InputType>
    EIGEN_DEVICE_FUNC
    SelfAdjointEigenSolver& compute(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors);
Don Gagne's avatar
Don Gagne committed
204
    
205
    /** \brief Computes eigendecomposition of given matrix using a closed-form algorithm
Don Gagne's avatar
Don Gagne committed
206 207 208 209
      *
      * This is a variant of compute(const MatrixType&, int options) which
      * directly solves the underlying polynomial equation.
      * 
210
      * Currently only 2x2 and 3x3 matrices for which the sizes are known at compile time are supported (e.g., Matrix3d).
Don Gagne's avatar
Don Gagne committed
211
      * 
212
      * This method is usually significantly faster than the QR iterative algorithm
Don Gagne's avatar
Don Gagne committed
213 214 215
      * but it might also be less accurate. It is also worth noting that
      * for 3x3 matrices it involves trigonometric operations which are
      * not necessarily available for all scalar types.
216 217 218 219
      * 
      * For the 3x3 case, we observed the following worst case relative error regarding the eigenvalues:
      *   - double: 1e-8
      *   - float:  1e-3
Don Gagne's avatar
Don Gagne committed
220 221 222
      *
      * \sa compute(const MatrixType&, int options)
      */
223
    EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
224
    SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
LM's avatar
LM committed
225

226 227 228 229 230 231 232 233 234 235 236 237 238 239
    /**
      *\brief Computes the eigen decomposition from a tridiagonal symmetric matrix
      *
      * \param[in] diag The vector containing the diagonal of the matrix.
      * \param[in] subdiag The subdiagonal of the matrix.
      * \param[in] options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
      * \returns Reference to \c *this
      *
      * This function assumes that the matrix has been reduced to tridiagonal form.
      *
      * \sa compute(const MatrixType&, int) for more information
      */
    SelfAdjointEigenSolver& computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options=ComputeEigenvectors);

LM's avatar
LM committed
240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257
    /** \brief Returns the eigenvectors of given matrix.
      *
      * \returns  A const reference to the matrix whose columns are the eigenvectors.
      *
      * \pre The eigenvectors have been computed before.
      *
      * Column \f$ k \f$ of the returned matrix is an eigenvector corresponding
      * to eigenvalue number \f$ k \f$ as returned by eigenvalues().  The
      * eigenvectors are normalized to have (Euclidean) norm equal to one. If
      * this object was used to solve the eigenproblem for the selfadjoint
      * matrix \f$ A \f$, then the matrix returned by this function is the
      * matrix \f$ V \f$ in the eigendecomposition \f$ A = V D V^{-1} \f$.
      *
      * Example: \include SelfAdjointEigenSolver_eigenvectors.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_eigenvectors.out
      *
      * \sa eigenvalues()
      */
258
    EIGEN_DEVICE_FUNC
259
    const EigenvectorsType& eigenvectors() const
LM's avatar
LM committed
260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280
    {
      eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
      eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
      return m_eivec;
    }

    /** \brief Returns the eigenvalues of given matrix.
      *
      * \returns A const reference to the column vector containing the eigenvalues.
      *
      * \pre The eigenvalues have been computed before.
      *
      * The eigenvalues are repeated according to their algebraic multiplicity,
      * so there are as many eigenvalues as rows in the matrix. The eigenvalues
      * are sorted in increasing order.
      *
      * Example: \include SelfAdjointEigenSolver_eigenvalues.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_eigenvalues.out
      *
      * \sa eigenvectors(), MatrixBase::eigenvalues()
      */
281
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302
    const RealVectorType& eigenvalues() const
    {
      eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
      return m_eivalues;
    }

    /** \brief Computes the positive-definite square root of the matrix.
      *
      * \returns the positive-definite square root of the matrix
      *
      * \pre The eigenvalues and eigenvectors of a positive-definite matrix
      * have been computed before.
      *
      * The square root of a positive-definite matrix \f$ A \f$ is the
      * positive-definite matrix whose square equals \f$ A \f$. This function
      * uses the eigendecomposition \f$ A = V D V^{-1} \f$ to compute the
      * square root as \f$ A^{1/2} = V D^{1/2} V^{-1} \f$.
      *
      * Example: \include SelfAdjointEigenSolver_operatorSqrt.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_operatorSqrt.out
      *
303
      * \sa operatorInverseSqrt(), <a href="unsupported/group__MatrixFunctions__Module.html">MatrixFunctions Module</a>
LM's avatar
LM committed
304
      */
305
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327
    MatrixType operatorSqrt() const
    {
      eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
      eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
      return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
    }

    /** \brief Computes the inverse square root of the matrix.
      *
      * \returns the inverse positive-definite square root of the matrix
      *
      * \pre The eigenvalues and eigenvectors of a positive-definite matrix
      * have been computed before.
      *
      * This function uses the eigendecomposition \f$ A = V D V^{-1} \f$ to
      * compute the inverse square root as \f$ V D^{-1/2} V^{-1} \f$. This is
      * cheaper than first computing the square root with operatorSqrt() and
      * then its inverse with MatrixBase::inverse().
      *
      * Example: \include SelfAdjointEigenSolver_operatorInverseSqrt.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_operatorInverseSqrt.out
      *
328
      * \sa operatorSqrt(), MatrixBase::inverse(), <a href="unsupported/group__MatrixFunctions__Module.html">MatrixFunctions Module</a>
LM's avatar
LM committed
329
      */
330
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
331 332 333 334 335 336 337 338 339 340 341
    MatrixType operatorInverseSqrt() const
    {
      eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
      eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
      return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
    }

    /** \brief Reports whether previous computation was successful.
      *
      * \returns \c Success if computation was succesful, \c NoConvergence otherwise.
      */
342
    EIGEN_DEVICE_FUNC
LM's avatar
LM committed
343 344 345 346 347 348 349 350
    ComputationInfo info() const
    {
      eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
      return m_info;
    }

    /** \brief Maximum number of iterations.
      *
Don Gagne's avatar
Don Gagne committed
351 352
      * The algorithm terminates if it does not converge within m_maxIterations * n iterations, where n
      * denotes the size of the matrix. This value is currently set to 30 (copied from LAPACK).
LM's avatar
LM committed
353 354 355 356
      */
    static const int m_maxIterations = 30;

  protected:
357 358 359 360 361 362
    static void check_template_parameters()
    {
      EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
    }
    
    EigenvectorsType m_eivec;
LM's avatar
LM committed
363 364 365 366 367 368 369
    RealVectorType m_eivalues;
    typename TridiagonalizationType::SubDiagonalType m_subdiag;
    ComputationInfo m_info;
    bool m_isInitialized;
    bool m_eigenvectorsOk;
};

370
namespace internal {
LM's avatar
LM committed
371 372 373 374 375 376 377
/** \internal
  *
  * \eigenvalues_module \ingroup Eigenvalues_Module
  *
  * Performs a QR step on a tridiagonal symmetric matrix represented as a
  * pair of two vectors \a diag and \a subdiag.
  *
378 379 380 381 382 383
  * \param diag the diagonal part of the input selfadjoint tridiagonal matrix
  * \param subdiag the sub-diagonal part of the input selfadjoint tridiagonal matrix
  * \param start starting index of the submatrix to work on
  * \param end last+1 index of the submatrix to work on
  * \param matrixQ pointer to the column-major matrix holding the eigenvectors, can be 0
  * \param n size of the input matrix
LM's avatar
LM committed
384 385 386 387 388 389 390
  *
  * For compilation efficiency reasons, this procedure does not use eigen expression
  * for its arguments.
  *
  * Implemented from Golub's "Matrix Computations", algorithm 8.3.2:
  * "implicit symmetric QR step with Wilkinson shift"
  */
391 392
template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
EIGEN_DEVICE_FUNC
LM's avatar
LM committed
393 394 395 396
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
}

template<typename MatrixType>
397 398
template<typename InputType>
EIGEN_DEVICE_FUNC
LM's avatar
LM committed
399
SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
400
::compute(const EigenBase<InputType>& a_matrix, int options)
LM's avatar
LM committed
401
{
402 403
  check_template_parameters();
  
404 405
  const InputType &matrix(a_matrix.derived());
  
Don Gagne's avatar
Don Gagne committed
406
  using std::abs;
LM's avatar
LM committed
407 408 409 410 411 412 413 414 415 416
  eigen_assert(matrix.cols() == matrix.rows());
  eigen_assert((options&~(EigVecMask|GenEigMask))==0
          && (options&EigVecMask)!=EigVecMask
          && "invalid option parameter");
  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
  Index n = matrix.cols();
  m_eivalues.resize(n,1);

  if(n==1)
  {
417 418
    m_eivec = matrix;
    m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
LM's avatar
LM committed
419
    if(computeEigenvectors)
Don Gagne's avatar
Don Gagne committed
420
      m_eivec.setOnes(n,n);
LM's avatar
LM committed
421 422 423 424 425 426 427 428
    m_info = Success;
    m_isInitialized = true;
    m_eigenvectorsOk = computeEigenvectors;
    return *this;
  }

  // declare some aliases
  RealVectorType& diag = m_eivalues;
429
  EigenvectorsType& mat = m_eivec;
LM's avatar
LM committed
430 431

  // map the matrix coefficients to [-1:1] to avoid over- and underflow.
Don Gagne's avatar
Don Gagne committed
432 433 434 435
  mat = matrix.template triangularView<Lower>();
  RealScalar scale = mat.cwiseAbs().maxCoeff();
  if(scale==RealScalar(0)) scale = RealScalar(1);
  mat.template triangularView<Lower>() /= scale;
LM's avatar
LM committed
436 437
  m_subdiag.resize(n-1);
  internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
438 439

  m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
LM's avatar
LM committed
440
  
441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489
  // scale back the eigen values
  m_eivalues *= scale;

  m_isInitialized = true;
  m_eigenvectorsOk = computeEigenvectors;
  return *this;
}

template<typename MatrixType>
SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
{
  //TODO : Add an option to scale the values beforehand
  bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;

  m_eivalues = diag;
  m_subdiag = subdiag;
  if (computeEigenvectors)
  {
    m_eivec.setIdentity(diag.size(), diag.size());
  }
  m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);

  m_isInitialized = true;
  m_eigenvectorsOk = computeEigenvectors;
  return *this;
}

namespace internal {
/**
  * \internal
  * \brief Compute the eigendecomposition from a tridiagonal matrix
  *
  * \param[in,out] diag : On input, the diagonal of the matrix, on output the eigenvalues
  * \param[in,out] subdiag : The subdiagonal part of the matrix (entries are modified during the decomposition)
  * \param[in] maxIterations : the maximum number of iterations
  * \param[in] computeEigenvectors : whether the eigenvectors have to be computed or not
  * \param[out] eivec : The matrix to store the eigenvectors if computeEigenvectors==true. Must be allocated on input.
  * \returns \c Success or \c NoConvergence
  */
template<typename MatrixType, typename DiagType, typename SubDiagType>
ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
{
  using std::abs;

  ComputationInfo info;
  typedef typename MatrixType::Scalar Scalar;

  Index n = diag.size();
LM's avatar
LM committed
490 491
  Index end = n-1;
  Index start = 0;
Don Gagne's avatar
Don Gagne committed
492
  Index iter = 0; // total number of iterations
493 494 495 496 497
  
  typedef typename DiagType::RealScalar RealScalar;
  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
  const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
  
LM's avatar
LM committed
498 499 500
  while (end>0)
  {
    for (Index i = start; i<end; ++i)
501 502
      if (internal::isMuchSmallerThan(abs(subdiag[i]),(abs(diag[i])+abs(diag[i+1])),precision) || abs(subdiag[i]) <= considerAsZero)
        subdiag[i] = 0;
LM's avatar
LM committed
503 504

    // find the largest unreduced block
505
    while (end>0 && subdiag[end-1]==RealScalar(0))
LM's avatar
LM committed
506 507 508 509 510 511
    {
      end--;
    }
    if (end<=0)
      break;

Don Gagne's avatar
Don Gagne committed
512
    // if we spent too many iterations, we give up
LM's avatar
LM committed
513
    iter++;
514
    if(iter > maxIterations * n) break;
LM's avatar
LM committed
515 516

    start = end - 1;
517
    while (start>0 && subdiag[start-1]!=0)
LM's avatar
LM committed
518 519
      start--;

520
    internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
LM's avatar
LM committed
521
  }
522 523
  if (iter <= maxIterations * n)
    info = Success;
LM's avatar
LM committed
524
  else
525
    info = NoConvergence;
LM's avatar
LM committed
526 527 528 529

  // Sort eigenvalues and corresponding vectors.
  // TODO make the sort optional ?
  // TODO use a better sort algorithm !!
530
  if (info == Success)
LM's avatar
LM committed
531 532 533 534
  {
    for (Index i = 0; i < n-1; ++i)
    {
      Index k;
535
      diag.segment(i,n-i).minCoeff(&k);
LM's avatar
LM committed
536 537
      if (k > 0)
      {
538
        std::swap(diag[i], diag[k+i]);
LM's avatar
LM committed
539
        if(computeEigenvectors)
540
          eivec.col(i).swap(eivec.col(k+i));
LM's avatar
LM committed
541 542 543
      }
    }
  }
544
  return info;
LM's avatar
LM committed
545
}
Don Gagne's avatar
Don Gagne committed
546 547 548
  
template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
{
549
  EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
550 551 552 553 554 555 556 557 558
  static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
  { eig.compute(A,options); }
};

template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
{
  typedef typename SolverType::MatrixType MatrixType;
  typedef typename SolverType::RealVectorType VectorType;
  typedef typename SolverType::Scalar Scalar;
559
  typedef typename SolverType::EigenvectorsType EigenvectorsType;
Don Gagne's avatar
Don Gagne committed
560
  
561

562 563 564 565
  /** \internal
   * Computes the roots of the characteristic polynomial of \a m.
   * For numerical stability m.trace() should be near zero and to avoid over- or underflow m should be normalized.
   */
566
  EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
567 568
  static inline void computeRoots(const MatrixType& m, VectorType& roots)
  {
569 570 571 572 573 574
    EIGEN_USING_STD_MATH(sqrt)
    EIGEN_USING_STD_MATH(atan2)
    EIGEN_USING_STD_MATH(cos)
    EIGEN_USING_STD_MATH(sin)
    const Scalar s_inv3 = Scalar(1)/Scalar(3);
    const Scalar s_sqrt3 = sqrt(Scalar(3));
Don Gagne's avatar
Don Gagne committed
575 576 577 578 579 580 581 582 583 584 585

    // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0.  The
    // eigenvalues are the roots to this equation, all guaranteed to be
    // real-valued, because the matrix is symmetric.
    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);
    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);
    Scalar c2 = m(0,0) + m(1,1) + m(2,2);

    // Construct the parameters used in classifying the roots of the equation
    // and in solving the equation for the roots in closed form.
    Scalar c2_over_3 = c2*s_inv3;
586
    Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
587
    a_over_3 = numext::maxi(a_over_3, Scalar(0));
Don Gagne's avatar
Don Gagne committed
588 589 590

    Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));

591
    Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
592
    q = numext::maxi(q, Scalar(0));
Don Gagne's avatar
Don Gagne committed
593 594

    // Compute the eigenvalues by solving for the roots of the polynomial.
595 596
    Scalar rho = sqrt(a_over_3);
    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]
Don Gagne's avatar
Don Gagne committed
597 598
    Scalar cos_theta = cos(theta);
    Scalar sin_theta = sin(theta);
599 600 601 602
    // roots are already sorted, since cos is monotonically decreasing on [0, pi]
    roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
    roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
    roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
Don Gagne's avatar
Don Gagne committed
603
  }
604

605
  EIGEN_DEVICE_FUNC
606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624
  static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
  {
    using std::abs;
    Index i0;
    // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
    mat.diagonal().cwiseAbs().maxCoeff(&i0);
    // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
    // so let's save it:
    representative = mat.col(i0);
    Scalar n0, n1;
    VectorType c0, c1;
    n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
    n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
    if(n0>n1) res = c0/std::sqrt(n0);
    else      res = c1/std::sqrt(n1);

    return true;
  }

625
  EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
626 627 628 629 630 631 632 633
  static inline void run(SolverType& solver, const MatrixType& mat, int options)
  {
    eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
    eigen_assert((options&~(EigVecMask|GenEigMask))==0
            && (options&EigVecMask)!=EigVecMask
            && "invalid option parameter");
    bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
    
634
    EigenvectorsType& eivecs = solver.m_eivec;
Don Gagne's avatar
Don Gagne committed
635 636
    VectorType& eivals = solver.m_eivalues;
  
637 638 639 640 641 642 643
    // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
    Scalar shift = mat.trace() / Scalar(3);
    // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
    MatrixType scaledMat = mat.template selfadjointView<Lower>();
    scaledMat.diagonal().array() -= shift;
    Scalar scale = scaledMat.cwiseAbs().maxCoeff();
    if(scale > 0) scaledMat /= scale;   // TODO for scale==0 we could save the remaining operations
Don Gagne's avatar
Don Gagne committed
644 645 646 647

    // compute the eigenvalues
    computeRoots(scaledMat,eivals);

648
    // compute the eigenvectors
Don Gagne's avatar
Don Gagne committed
649 650 651 652
    if(computeEigenvectors)
    {
      if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
      {
653
        // All three eigenvalues are numerically the same
Don Gagne's avatar
Don Gagne committed
654 655 656 657 658 659 660
        eivecs.setIdentity();
      }
      else
      {
        MatrixType tmp;
        tmp = scaledMat;

661
        // Compute the eigenvector of the most distinct eigenvalue
Don Gagne's avatar
Don Gagne committed
662 663
        Scalar d0 = eivals(2) - eivals(1);
        Scalar d1 = eivals(1) - eivals(0);
664 665
        Index k(0), l(2);
        if(d0 > d1)
Don Gagne's avatar
Don Gagne committed
666
        {
667
          numext::swap(k,l);
668
          d0 = d1;
Don Gagne's avatar
Don Gagne committed
669 670
        }

671 672 673 674 675 676
        // Compute the eigenvector of index k
        {
          tmp.diagonal().array () -= eivals(k);
          // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
          extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
        }
Don Gagne's avatar
Don Gagne committed
677

678 679 680 681 682 683 684 685
        // Compute eigenvector of index l
        if(d0<=2*Eigen::NumTraits<Scalar>::epsilon()*d1)
        {
          // If d0 is too small, then the two other eigenvalues are numerically the same,
          // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
          eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
          eivecs.col(l).normalize();
        }
Don Gagne's avatar
Don Gagne committed
686 687
        else
        {
688 689 690 691 692
          tmp = scaledMat;
          tmp.diagonal().array () -= eivals(l);

          VectorType dummy;
          extract_kernel(tmp, eivecs.col(l), dummy);
Don Gagne's avatar
Don Gagne committed
693 694
        }

695 696
        // Compute last eigenvector from the other two
        eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
Don Gagne's avatar
Don Gagne committed
697 698
      }
    }
699

Don Gagne's avatar
Don Gagne committed
700 701
    // Rescale back to the original size.
    eivals *= scale;
702
    eivals.array() += shift;
Don Gagne's avatar
Don Gagne committed
703 704 705 706 707 708 709 710
    
    solver.m_info = Success;
    solver.m_isInitialized = true;
    solver.m_eigenvectorsOk = computeEigenvectors;
  }
};

// 2x2 direct eigenvalues decomposition, code from Hauke Heibel
711 712
template<typename SolverType> 
struct direct_selfadjoint_eigenvalues<SolverType,2,false>
Don Gagne's avatar
Don Gagne committed
713 714 715 716
{
  typedef typename SolverType::MatrixType MatrixType;
  typedef typename SolverType::RealVectorType VectorType;
  typedef typename SolverType::Scalar Scalar;
717
  typedef typename SolverType::EigenvectorsType EigenvectorsType;
Don Gagne's avatar
Don Gagne committed
718
  
719
  EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
720 721 722
  static inline void computeRoots(const MatrixType& m, VectorType& roots)
  {
    using std::sqrt;
723
    const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
Don Gagne's avatar
Don Gagne committed
724 725 726 727 728
    const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
    roots(0) = t1 - t0;
    roots(1) = t1 + t0;
  }
  
729
  EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
730 731
  static inline void run(SolverType& solver, const MatrixType& mat, int options)
  {
732 733 734
    EIGEN_USING_STD_MATH(sqrt);
    EIGEN_USING_STD_MATH(abs);
    
Don Gagne's avatar
Don Gagne committed
735 736 737 738 739 740
    eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
    eigen_assert((options&~(EigVecMask|GenEigMask))==0
            && (options&EigVecMask)!=EigVecMask
            && "invalid option parameter");
    bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
    
741
    EigenvectorsType& eivecs = solver.m_eivec;
Don Gagne's avatar
Don Gagne committed
742 743
    VectorType& eivals = solver.m_eivalues;
  
744 745 746 747 748 749 750 751 752
    // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
    Scalar shift = mat.trace() / Scalar(2);
    MatrixType scaledMat = mat;
    scaledMat.coeffRef(0,1) = mat.coeff(1,0);
    scaledMat.diagonal().array() -= shift;
    Scalar scale = scaledMat.cwiseAbs().maxCoeff();
    if(scale > Scalar(0))
      scaledMat /= scale;

Don Gagne's avatar
Don Gagne committed
753 754
    // Compute the eigenvalues
    computeRoots(scaledMat,eivals);
755

Don Gagne's avatar
Don Gagne committed
756 757 758
    // compute the eigen vectors
    if(computeEigenvectors)
    {
759
      if((eivals(1)-eivals(0))<=abs(eivals(1))*Eigen::NumTraits<Scalar>::epsilon())
Don Gagne's avatar
Don Gagne committed
760
      {
761
        eivecs.setIdentity();
Don Gagne's avatar
Don Gagne committed
762 763 764
      }
      else
      {
765 766 767 768 769 770 771 772 773 774 775 776 777 778
        scaledMat.diagonal().array () -= eivals(1);
        Scalar a2 = numext::abs2(scaledMat(0,0));
        Scalar c2 = numext::abs2(scaledMat(1,1));
        Scalar b2 = numext::abs2(scaledMat(1,0));
        if(a2>c2)
        {
          eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
          eivecs.col(1) /= sqrt(a2+b2);
        }
        else
        {
          eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
          eivecs.col(1) /= sqrt(c2+b2);
        }
Don Gagne's avatar
Don Gagne committed
779

780 781
        eivecs.col(0) << eivecs.col(1).unitOrthogonal();
      }
Don Gagne's avatar
Don Gagne committed
782
    }
783

Don Gagne's avatar
Don Gagne committed
784 785
    // Rescale back to the original size.
    eivals *= scale;
786 787
    eivals.array() += shift;

Don Gagne's avatar
Don Gagne committed
788 789 790 791 792 793 794 795 796
    solver.m_info = Success;
    solver.m_isInitialized = true;
    solver.m_eigenvectorsOk = computeEigenvectors;
  }
};

}

template<typename MatrixType>
797
EIGEN_DEVICE_FUNC
Don Gagne's avatar
Don Gagne committed
798 799 800 801 802 803 804
SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
::computeDirect(const MatrixType& matrix, int options)
{
  internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
  return *this;
}

LM's avatar
LM committed
805
namespace internal {
806 807
template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
EIGEN_DEVICE_FUNC
LM's avatar
LM committed
808 809
static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
{
Don Gagne's avatar
Don Gagne committed
810
  using std::abs;
LM's avatar
LM committed
811
  RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
Don Gagne's avatar
Don Gagne committed
812 813 814 815 816 817 818
  RealScalar e = subdiag[end-1];
  // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
  // underflow thus leading to inf/NaN values when using the following commented code:
//   RealScalar e2 = numext::abs2(subdiag[end-1]);
//   RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
  // This explain the following, somewhat more complicated, version:
  RealScalar mu = diag[end];
819
  if(td==RealScalar(0))
Don Gagne's avatar
Don Gagne committed
820 821 822 823 824
    mu -= abs(e);
  else
  {
    RealScalar e2 = numext::abs2(subdiag[end-1]);
    RealScalar h = numext::hypot(td,e);
825 826
    if(e2==RealScalar(0)) mu -= (e / (td + (td>RealScalar(0) ? RealScalar(1) : RealScalar(-1)))) * (e / h);
    else                  mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
Don Gagne's avatar
Don Gagne committed
827 828
  }
  
LM's avatar
LM committed
829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858
  RealScalar x = diag[start] - mu;
  RealScalar z = subdiag[start];
  for (Index k = start; k < end; ++k)
  {
    JacobiRotation<RealScalar> rot;
    rot.makeGivens(x, z);

    // do T = G' T G
    RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
    RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];

    diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
    diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
    subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
    

    if (k > start)
      subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;

    x = subdiag[k];

    if (k < end - 1)
    {
      z = -rot.s() * subdiag[k+1];
      subdiag[k + 1] = rot.c() * subdiag[k+1];
    }
    
    // apply the givens rotation to the unit matrix Q = Q * G
    if (matrixQ)
    {
859 860
      // FIXME if StorageOrder == RowMajor this operation is not very efficient
      Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
LM's avatar
LM committed
861 862 863 864
      q.applyOnTheRight(k,k+1,rot);
    }
  }
}
Don Gagne's avatar
Don Gagne committed
865

LM's avatar
LM committed
866 867
} // end namespace internal

Don Gagne's avatar
Don Gagne committed
868 869
} // end namespace Eigen

LM's avatar
LM committed
870
#endif // EIGEN_SELFADJOINTEIGENSOLVER_H