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 
108  : Base(matA.cols())
109  {
111  }
112 
155 
156  protected:
157 
158 };
159 
160 
161 template<typename MatrixType>
162 GeneralizedSelfAdjointEigenSolver<MatrixType>& GeneralizedSelfAdjointEigenSolver<MatrixType>::
164 {
165  eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
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
matA
MatrixXf matA(2, 2)
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Eigen::BAx_lx
@ BAx_lx
Definition: Constants.h:416
matB
MatrixXf matB(2, 2)
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
eigen_assert
#define eigen_assert(x)
Definition: Macros.h:1037
type
Definition: pytypes.h:1525
Eigen::GeneralizedSelfAdjointEigenSolver::Base
SelfAdjointEigenSolver< _MatrixType > Base
Definition: GeneralizedSelfAdjointEigenSolver.h:50
Eigen::ComputeEigenvectors
@ ComputeEigenvectors
Definition: Constants.h:405
size
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Eigen::GeneralizedSelfAdjointEigenSolver::GeneralizedSelfAdjointEigenSolver
GeneralizedSelfAdjointEigenSolver(Index size)
Constructor, pre-allocates memory for dynamic-size matrices.
Definition: GeneralizedSelfAdjointEigenSolver.h:76
Eigen::GeneralizedSelfAdjointEigenSolver::GeneralizedSelfAdjointEigenSolver
GeneralizedSelfAdjointEigenSolver(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Constructor; computes generalized eigendecomposition of given matrix pencil.
Definition: GeneralizedSelfAdjointEigenSolver.h:106
Eigen::SelfAdjointEigenSolver
Computes eigenvalues and eigenvectors of selfadjoint matrices.
Definition: SelfAdjointEigenSolver.h:76
Eigen::LLT::matrixL
Traits::MatrixL matrixL() const
Definition: LLT.h:135
Eigen::GenEigMask
@ GenEigMask
Definition: Constants.h:418
Eigen::Ax_lBx
@ Ax_lBx
Definition: Constants.h:410
Eigen::GeneralizedSelfAdjointEigenSolver::GeneralizedSelfAdjointEigenSolver
GeneralizedSelfAdjointEigenSolver()
Default constructor for fixed-size matrices.
Definition: GeneralizedSelfAdjointEigenSolver.h:62
compute
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Definition: dense_solvers.cpp:25
Eigen::EigenvaluesOnly
@ EigenvaluesOnly
Definition: Constants.h:402
Eigen::SelfAdjointEigenSolver::Index
Eigen::Index Index
Definition: SelfAdjointEigenSolver.h:90
Eigen::LLT
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:66
Tridiagonalization.h
Eigen::GeneralizedSelfAdjointEigenSolver
Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem.
Definition: GeneralizedSelfAdjointEigenSolver.h:48
Eigen::GeneralizedSelfAdjointEigenSolver::MatrixType
_MatrixType MatrixType
Definition: GeneralizedSelfAdjointEigenSolver.h:53
Eigen::ABx_lx
@ ABx_lx
Definition: Constants.h:413
Eigen::LLT::matrixU
Traits::MatrixU matrixU() const
Definition: LLT.h:128
Eigen::EigVecMask
@ EigVecMask
Definition: Constants.h:407
Eigen::GeneralizedSelfAdjointEigenSolver::compute
GeneralizedSelfAdjointEigenSolver & compute(const MatrixType &matA, const MatrixType &matB, int options=ComputeEigenvectors|Ax_lBx)
Computes generalized eigendecomposition of given matrix pencil.
Definition: GeneralizedSelfAdjointEigenSolver.h:163
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
options
Definition: options.h:16


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:02:20