Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifndef EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
00027 #define EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
00028
00029 #include "./EigenvaluesCommon.h"
00030 #include "./Tridiagonalization.h"
00031
00061 template<typename _MatrixType>
00062 class GeneralizedSelfAdjointEigenSolver : public SelfAdjointEigenSolver<_MatrixType>
00063 {
00064 typedef SelfAdjointEigenSolver<_MatrixType> Base;
00065 public:
00066
00067 typedef typename Base::Index Index;
00068 typedef _MatrixType MatrixType;
00069
00077 GeneralizedSelfAdjointEigenSolver() : Base() {}
00078
00091 GeneralizedSelfAdjointEigenSolver(Index size)
00092 : Base(size)
00093 {}
00094
00121 GeneralizedSelfAdjointEigenSolver(const MatrixType& matA, const MatrixType& matB,
00122 int options = ComputeEigenvectors|Ax_lBx)
00123 : Base(matA.cols())
00124 {
00125 compute(matA, matB, options);
00126 }
00127
00168 GeneralizedSelfAdjointEigenSolver& compute(const MatrixType& matA, const MatrixType& matB,
00169 int options = ComputeEigenvectors|Ax_lBx);
00170
00171 protected:
00172
00173 };
00174
00175
00176 template<typename MatrixType>
00177 GeneralizedSelfAdjointEigenSolver<MatrixType>& GeneralizedSelfAdjointEigenSolver<MatrixType>::
00178 compute(const MatrixType& matA, const MatrixType& matB, int options)
00179 {
00180 eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
00181 eigen_assert((options&~(EigVecMask|GenEigMask))==0
00182 && (options&EigVecMask)!=EigVecMask
00183 && ((options&GenEigMask)==0 || (options&GenEigMask)==Ax_lBx
00184 || (options&GenEigMask)==ABx_lx || (options&GenEigMask)==BAx_lx)
00185 && "invalid option parameter");
00186
00187 bool computeEigVecs = ((options&EigVecMask)==0) || ((options&EigVecMask)==ComputeEigenvectors);
00188
00189
00190 LLT<MatrixType> cholB(matB);
00191
00192 int type = (options&GenEigMask);
00193 if(type==0)
00194 type = Ax_lBx;
00195
00196 if(type==Ax_lBx)
00197 {
00198
00199 MatrixType matC = matA.template selfadjointView<Lower>();
00200 cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
00201 cholB.matrixU().template solveInPlace<OnTheRight>(matC);
00202
00203 Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly );
00204
00205
00206 if(computeEigVecs)
00207 cholB.matrixU().solveInPlace(Base::m_eivec);
00208 }
00209 else if(type==ABx_lx)
00210 {
00211
00212 MatrixType matC = matA.template selfadjointView<Lower>();
00213 matC = matC * cholB.matrixL();
00214 matC = cholB.matrixU() * matC;
00215
00216 Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
00217
00218
00219 if(computeEigVecs)
00220 cholB.matrixU().solveInPlace(Base::m_eivec);
00221 }
00222 else if(type==BAx_lx)
00223 {
00224
00225 MatrixType matC = matA.template selfadjointView<Lower>();
00226 matC = matC * cholB.matrixL();
00227 matC = cholB.matrixU() * matC;
00228
00229 Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
00230
00231
00232 if(computeEigVecs)
00233 Base::m_eivec = cholB.matrixL() * Base::m_eivec;
00234 }
00235
00236 return *this;
00237 }
00238
00239 #endif // EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H