GeneralizedSelfAdjointEigenSolver.h 9.49 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_GENERALIZEDSELFADJOINTEIGENSOLVER_H
#define EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H

#include "./Tridiagonalization.h"

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

LM's avatar
LM committed
18 19 20 21 22 23 24 25 26 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
/** \eigenvalues_module \ingroup Eigenvalues_Module
  *
  *
  * \class GeneralizedSelfAdjointEigenSolver
  *
  * \brief Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem
  *
  * \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.
  *
  * This class solves the generalized eigenvalue problem
  * \f$ Av = \lambda Bv \f$. In this case, the matrix \f$ A \f$ should be
  * selfadjoint and the matrix \f$ B \f$ should be positive definite.
  *
  * 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
  * GeneralizedSelfAdjointEigenSolver(const MatrixType&, 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 GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)
  * contains an example of the typical use of this class.
  *
  * \sa class SelfAdjointEigenSolver, class EigenSolver, class ComplexEigenSolver
  */
template<typename _MatrixType>
class GeneralizedSelfAdjointEigenSolver : public SelfAdjointEigenSolver<_MatrixType>
{
    typedef SelfAdjointEigenSolver<_MatrixType> Base;
  public:

    typedef _MatrixType MatrixType;

    /** \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
      * GeneralizedSelfAdjointEigenSolver(Index) for dynamic-size matrices.
      */
    GeneralizedSelfAdjointEigenSolver() : Base() {}

    /** \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
      */
76
    explicit GeneralizedSelfAdjointEigenSolver(Index size)
LM's avatar
LM committed
77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 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 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223
        : Base(size)
    {}

    /** \brief Constructor; computes generalized eigendecomposition of given matrix pencil.
      *
      * \param[in]  matA  Selfadjoint matrix in matrix pencil.
      *                   Only the lower triangular part of the matrix is referenced.
      * \param[in]  matB  Positive-definite matrix in matrix pencil.
      *                   Only the lower triangular part of the matrix is referenced.
      * \param[in]  options A or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}.
      *                     Default is #ComputeEigenvectors|#Ax_lBx.
      *
      * This constructor calls compute(const MatrixType&, const MatrixType&, int)
      * to compute the eigenvalues and (if requested) the eigenvectors of the
      * generalized eigenproblem \f$ Ax = \lambda B x \f$ with \a matA the
      * selfadjoint matrix \f$ A \f$ and \a matB the positive definite matrix
      * \f$ B \f$. Each eigenvector \f$ x \f$ satisfies the property
      * \f$ x^* B x = 1 \f$. The eigenvectors are computed if
      * \a options contains ComputeEigenvectors.
      *
      * In addition, the two following variants can be solved via \p options:
      * - \c ABx_lx: \f$ ABx = \lambda x \f$
      * - \c BAx_lx: \f$ BAx = \lambda x \f$
      *
      * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType2.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType2.out
      *
      * \sa compute(const MatrixType&, const MatrixType&, int)
      */
    GeneralizedSelfAdjointEigenSolver(const MatrixType& matA, const MatrixType& matB,
                                      int options = ComputeEigenvectors|Ax_lBx)
      : Base(matA.cols())
    {
      compute(matA, matB, options);
    }

    /** \brief Computes generalized eigendecomposition of given matrix pencil.
      *
      * \param[in]  matA  Selfadjoint matrix in matrix pencil.
      *                   Only the lower triangular part of the matrix is referenced.
      * \param[in]  matB  Positive-definite matrix in matrix pencil.
      *                   Only the lower triangular part of the matrix is referenced.
      * \param[in]  options A or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}.
      *                     Default is #ComputeEigenvectors|#Ax_lBx.
      *
      * \returns    Reference to \c *this
      *
      * Accoring to \p options, this function computes eigenvalues and (if requested)
      * the eigenvectors of one of the following three generalized eigenproblems:
      * - \c Ax_lBx: \f$ Ax = \lambda B x \f$
      * - \c ABx_lx: \f$ ABx = \lambda x \f$
      * - \c BAx_lx: \f$ BAx = \lambda x \f$
      * with \a matA the selfadjoint matrix \f$ A \f$ and \a matB the positive definite
      * matrix \f$ B \f$.
      * In addition, each eigenvector \f$ x \f$ satisfies the property \f$ x^* B x = 1 \f$.
      *
      * The eigenvalues() function can be used to retrieve
      * the eigenvalues. If \p options contains ComputeEigenvectors, then the
      * eigenvectors are also computed and can be retrieved by calling
      * eigenvectors().
      *
      * The implementation uses LLT to compute the Cholesky decomposition
      * \f$ B = LL^* \f$ and computes the classical eigendecomposition
      * of the selfadjoint matrix \f$ L^{-1} A (L^*)^{-1} \f$ if \p options contains Ax_lBx
      * and of \f$ L^{*} A L \f$ otherwise. This solves the
      * generalized eigenproblem, because any solution of the generalized
      * eigenproblem \f$ Ax = \lambda B x \f$ corresponds to a solution
      * \f$ L^{-1} A (L^*)^{-1} (L^* x) = \lambda (L^* x) \f$ of the
      * eigenproblem for \f$ L^{-1} A (L^*)^{-1} \f$. Similar statements
      * can be made for the two other variants.
      *
      * Example: \include SelfAdjointEigenSolver_compute_MatrixType2.cpp
      * Output: \verbinclude SelfAdjointEigenSolver_compute_MatrixType2.out
      *
      * \sa GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)
      */
    GeneralizedSelfAdjointEigenSolver& compute(const MatrixType& matA, const MatrixType& matB,
                                               int options = ComputeEigenvectors|Ax_lBx);

  protected:

};


template<typename MatrixType>
GeneralizedSelfAdjointEigenSolver<MatrixType>& GeneralizedSelfAdjointEigenSolver<MatrixType>::
compute(const MatrixType& matA, const MatrixType& matB, int options)
{
  eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
  eigen_assert((options&~(EigVecMask|GenEigMask))==0
          && (options&EigVecMask)!=EigVecMask
          && ((options&GenEigMask)==0 || (options&GenEigMask)==Ax_lBx
           || (options&GenEigMask)==ABx_lx || (options&GenEigMask)==BAx_lx)
          && "invalid option parameter");

  bool computeEigVecs = ((options&EigVecMask)==0) || ((options&EigVecMask)==ComputeEigenvectors);

  // Compute the cholesky decomposition of matB = L L' = U'U
  LLT<MatrixType> cholB(matB);

  int type = (options&GenEigMask);
  if(type==0)
    type = Ax_lBx;

  if(type==Ax_lBx)
  {
    // compute C = inv(L) A inv(L')
    MatrixType matC = matA.template selfadjointView<Lower>();
    cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
    cholB.matrixU().template solveInPlace<OnTheRight>(matC);

    Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly );

    // transform back the eigen vectors: evecs = inv(U) * evecs
    if(computeEigVecs)
      cholB.matrixU().solveInPlace(Base::m_eivec);
  }
  else if(type==ABx_lx)
  {
    // compute C = L' A L
    MatrixType matC = matA.template selfadjointView<Lower>();
    matC = matC * cholB.matrixL();
    matC = cholB.matrixU() * matC;

    Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);

    // transform back the eigen vectors: evecs = inv(U) * evecs
    if(computeEigVecs)
      cholB.matrixU().solveInPlace(Base::m_eivec);
  }
  else if(type==BAx_lx)
  {
    // compute C = L' A L
    MatrixType matC = matA.template selfadjointView<Lower>();
    matC = matC * cholB.matrixL();
    matC = cholB.matrixU() * matC;

    Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);

    // transform back the eigen vectors: evecs = L * evecs
    if(computeEigVecs)
      Base::m_eivec = cholB.matrixL() * Base::m_eivec;
  }

  return *this;
}

Don Gagne's avatar
Don Gagne committed
224 225
} // end namespace Eigen

LM's avatar
LM committed
226
#endif // EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H