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 namespace Eigen {
17 
47 template<typename _MatrixType>
49 {
51  public:
52 
53  typedef _MatrixType MatrixType;
54 
63 
77  : Base(size)
78  {}
79 
106  GeneralizedSelfAdjointEigenSolver(const MatrixType& matA, const MatrixType& matB,
107  int options = ComputeEigenvectors|Ax_lBx)
108  : Base(matA.cols())
109  {
110  compute(matA, matB, options);
111  }
112 
153  GeneralizedSelfAdjointEigenSolver& compute(const MatrixType& matA, const MatrixType& matB,
155 
156  protected:
157 
158 };
159 
160 
161 template<typename MatrixType>
163 compute(const MatrixType& matA, const MatrixType& matB, int options)
164 {
165  eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
166  eigen_assert((options&~(EigVecMask|GenEigMask))==0
167  && (options&EigVecMask)!=EigVecMask
168  && ((options&GenEigMask)==0 || (options&GenEigMask)==Ax_lBx
169  || (options&GenEigMask)==ABx_lx || (options&GenEigMask)==BAx_lx)
170  && "invalid option parameter");
171 
172  bool computeEigVecs = ((options&EigVecMask)==0) || ((options&EigVecMask)==ComputeEigenvectors);
173 
174  // Compute the cholesky decomposition of matB = L L' = U'U
175  LLT<MatrixType> cholB(matB);
176 
177  int type = (options&GenEigMask);
178  if(type==0)
179  type = Ax_lBx;
180 
181  if(type==Ax_lBx)
182  {
183  // compute C = inv(L) A inv(L')
184  MatrixType matC = matA.template selfadjointView<Lower>();
185  cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
186  cholB.matrixU().template solveInPlace<OnTheRight>(matC);
187 
188  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly );
189 
190  // transform back the eigen vectors: evecs = inv(U) * evecs
191  if(computeEigVecs)
192  cholB.matrixU().solveInPlace(Base::m_eivec);
193  }
194  else if(type==ABx_lx)
195  {
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 = inv(U) * evecs
204  if(computeEigVecs)
205  cholB.matrixU().solveInPlace(Base::m_eivec);
206  }
207  else if(type==BAx_lx)
208  {
209  // compute C = L' A L
210  MatrixType matC = matA.template selfadjointView<Lower>();
211  matC = matC * cholB.matrixL();
212  matC = cholB.matrixU() * matC;
213 
214  Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
215 
216  // transform back the eigen vectors: evecs = L * evecs
217  if(computeEigVecs)
218  Base::m_eivec = cholB.matrixL() * Base::m_eivec;
219  }
220 
221  return *this;
222 }
223 
224 } // end namespace Eigen
225 
226 #endif // EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
GeneralizedSelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
EIGEN_DEVICE_FUNC SelfAdjointEigenSolver & compute(const EigenBase< InputType > &matrix, int options=ComputeEigenvectors)
Computes eigendecomposition of given matrix.
Traits::MatrixU matrixU() const
Definition: LLT.h:115
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: LDLT.h:16
static constexpr size_t size(Tuple< Args... > &)
Provides access to the number of elements in a tuple as a compile-time constant expression.
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:52
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
options
#define eigen_assert(x)
Definition: Macros.h:577
GeneralizedSelfAdjointEigenSolver & compute(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Computes generalized eigendecomposition of given matrix pencil.
Traits::MatrixL matrixL() const
Definition: LLT.h:122
GeneralizedSelfAdjointEigenSolver(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Constructor; computes generalized eigendecomposition of given matrix pencil.
GeneralizedSelfAdjointEigenSolver()
Default constructor for fixed-size matrices.


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:08:12