GeneralizedSelfAdjointEigenSolver.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
12 #define EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
13 
14 #include "./Tridiagonalization.h"
15 
16 // IWYU pragma: private
17 #include "./InternalHeaderCheck.h"
18 
19 namespace Eigen {
20 
50 template <typename MatrixType_>
53 
54  public:
55  typedef MatrixType_ MatrixType;
56 
65 
79 
107  int options = ComputeEigenvectors | Ax_lBx)
108  : Base(matA.cols()) {
109  compute(matA, matB, options);
110  }
111 
153  int options = ComputeEigenvectors | Ax_lBx);
154 
155  protected:
156 };
157 
158 template <typename MatrixType>
160  const MatrixType& matA, const MatrixType& matB, int options) {
161  eigen_assert(matA.cols() == matA.rows() && matB.rows() == matA.rows() && matB.cols() == matB.rows());
162  eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
163  ((options & GenEigMask) == 0 || (options & GenEigMask) == Ax_lBx || (options & GenEigMask) == ABx_lx ||
164  (options & GenEigMask) == BAx_lx) &&
165  "invalid option parameter");
166 
167  bool computeEigVecs = ((options & EigVecMask) == 0) || ((options & EigVecMask) == ComputeEigenvectors);
168 
169  // Compute the cholesky decomposition of matB = L L' = U'U
170  LLT<MatrixType> cholB(matB);
171 
172  int type = (options & GenEigMask);
173  if (type == 0) type = Ax_lBx;
174 
175  if (type == Ax_lBx) {
176  // compute C = inv(L) A inv(L')
177  MatrixType matC = matA.template selfadjointView<Lower>();
178  cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
179  cholB.matrixU().template solveInPlace<OnTheRight>(matC);
180 
181  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
182 
183  // transform back the eigen vectors: evecs = inv(U) * evecs
184  if (computeEigVecs) cholB.matrixU().solveInPlace(Base::m_eivec);
185  } else if (type == ABx_lx) {
186  // compute C = L' A L
187  MatrixType matC = matA.template selfadjointView<Lower>();
188  matC = matC * cholB.matrixL();
189  matC = cholB.matrixU() * matC;
190 
191  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
192 
193  // transform back the eigen vectors: evecs = inv(U) * evecs
194  if (computeEigVecs) cholB.matrixU().solveInPlace(Base::m_eivec);
195  } else if (type == BAx_lx) {
196  // compute C = L' A L
197  MatrixType matC = matA.template selfadjointView<Lower>();
198  matC = matC * cholB.matrixL();
199  matC = cholB.matrixU() * matC;
200 
201  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
202 
203  // transform back the eigen vectors: evecs = L * evecs
204  if (computeEigVecs) Base::m_eivec = cholB.matrixL() * Base::m_eivec;
205  }
206 
207  return *this;
208 }
209 
210 } // end namespace Eigen
211 
212 #endif // EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
#define eigen_assert(x)
Definition: Macros.h:910
MatrixXf matB(2, 2)
int cols
Definition: Tutorial_commainit_02.cpp:1
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition: GeneralizedSelfAdjointEigenSolver.h:51
GeneralizedSelfAdjointEigenSolver()
Default constructor for fixed-size matrices.
Definition: GeneralizedSelfAdjointEigenSolver.h:64
GeneralizedSelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: GeneralizedSelfAdjointEigenSolver.h:78
SelfAdjointEigenSolver< MatrixType_ > Base
Definition: GeneralizedSelfAdjointEigenSolver.h:52
GeneralizedSelfAdjointEigenSolver & compute(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Computes generalized eigendecomposition of given matrix pencil.
Definition: GeneralizedSelfAdjointEigenSolver.h:159
MatrixType_ MatrixType
Definition: GeneralizedSelfAdjointEigenSolver.h:55
GeneralizedSelfAdjointEigenSolver(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Constructor; computes generalized eigendecomposition of given matrix pencil.
Definition: GeneralizedSelfAdjointEigenSolver.h:106
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
Traits::MatrixU matrixU() const
Definition: LLT.h:117
Traits::MatrixL matrixL() const
Definition: LLT.h:123
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:192
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:191
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: SelfAdjointEigenSolver.h:82
Eigen::Index Index
Definition: SelfAdjointEigenSolver.h:94
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Definition: dense_solvers.cpp:23
@ GenEigMask
Definition: Constants.h:414
@ EigVecMask
Definition: Constants.h:403
@ Ax_lBx
Definition: Constants.h:406
@ ComputeEigenvectors
Definition: Constants.h:401
@ BAx_lx
Definition: Constants.h:412
@ ABx_lx
Definition: Constants.h:409
@ EigenvaluesOnly
Definition: Constants.h:398
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > matA(size, size)
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
type
Definition: compute_granudrum_aor.py:141