00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #ifndef EIGEN_GENERALIZEDEIGENSOLVER_H
00012 #define EIGEN_GENERALIZEDEIGENSOLVER_H
00013
00014 #include "./RealQZ.h"
00015
00016 namespace Eigen {
00017
00057 template<typename _MatrixType> class GeneralizedEigenSolver
00058 {
00059 public:
00060
00062 typedef _MatrixType MatrixType;
00063
00064 enum {
00065 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00066 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00067 Options = MatrixType::Options,
00068 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00069 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00070 };
00071
00073 typedef typename MatrixType::Scalar Scalar;
00074 typedef typename NumTraits<Scalar>::Real RealScalar;
00075 typedef typename MatrixType::Index Index;
00076
00083 typedef std::complex<RealScalar> ComplexScalar;
00084
00090 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> VectorType;
00091
00097 typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ComplexVectorType;
00098
00101 typedef CwiseBinaryOp<internal::scalar_quotient_op<ComplexScalar,Scalar>,ComplexVectorType,VectorType> EigenvalueType;
00102
00108 typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
00109
00117 GeneralizedEigenSolver() : m_eivec(), m_alphas(), m_betas(), m_isInitialized(false), m_realQZ(), m_matS(), m_tmp() {}
00118
00125 GeneralizedEigenSolver(Index size)
00126 : m_eivec(size, size),
00127 m_alphas(size),
00128 m_betas(size),
00129 m_isInitialized(false),
00130 m_eigenvectorsOk(false),
00131 m_realQZ(size),
00132 m_matS(size, size),
00133 m_tmp(size)
00134 {}
00135
00148 GeneralizedEigenSolver(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true)
00149 : m_eivec(A.rows(), A.cols()),
00150 m_alphas(A.cols()),
00151 m_betas(A.cols()),
00152 m_isInitialized(false),
00153 m_eigenvectorsOk(false),
00154 m_realQZ(A.cols()),
00155 m_matS(A.rows(), A.cols()),
00156 m_tmp(A.cols())
00157 {
00158 compute(A, B, computeEigenvectors);
00159 }
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00198 EigenvalueType eigenvalues() const
00199 {
00200 eigen_assert(m_isInitialized && "GeneralizedEigenSolver is not initialized.");
00201 return EigenvalueType(m_alphas,m_betas);
00202 }
00203
00209 ComplexVectorType alphas() const
00210 {
00211 eigen_assert(m_isInitialized && "GeneralizedEigenSolver is not initialized.");
00212 return m_alphas;
00213 }
00214
00220 VectorType betas() const
00221 {
00222 eigen_assert(m_isInitialized && "GeneralizedEigenSolver is not initialized.");
00223 return m_betas;
00224 }
00225
00249 GeneralizedEigenSolver& compute(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true);
00250
00251 ComputationInfo info() const
00252 {
00253 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00254 return m_realQZ.info();
00255 }
00256
00259 GeneralizedEigenSolver& setMaxIterations(Index maxIters)
00260 {
00261 m_realQZ.setMaxIterations(maxIters);
00262 return *this;
00263 }
00264
00265 protected:
00266 MatrixType m_eivec;
00267 ComplexVectorType m_alphas;
00268 VectorType m_betas;
00269 bool m_isInitialized;
00270 bool m_eigenvectorsOk;
00271 RealQZ<MatrixType> m_realQZ;
00272 MatrixType m_matS;
00273
00274 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
00275 ColumnVectorType m_tmp;
00276 };
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289 template<typename MatrixType>
00290 GeneralizedEigenSolver<MatrixType>&
00291 GeneralizedEigenSolver<MatrixType>::compute(const MatrixType& A, const MatrixType& B, bool computeEigenvectors)
00292 {
00293 using std::sqrt;
00294 using std::abs;
00295 eigen_assert(A.cols() == A.rows() && B.cols() == A.rows() && B.cols() == B.rows());
00296
00297
00298
00299 m_realQZ.compute(A, B, computeEigenvectors);
00300
00301 if (m_realQZ.info() == Success)
00302 {
00303 m_matS = m_realQZ.matrixS();
00304 if (computeEigenvectors)
00305 m_eivec = m_realQZ.matrixZ().transpose();
00306
00307
00308 m_alphas.resize(A.cols());
00309 m_betas.resize(A.cols());
00310 Index i = 0;
00311 while (i < A.cols())
00312 {
00313 if (i == A.cols() - 1 || m_matS.coeff(i+1, i) == Scalar(0))
00314 {
00315 m_alphas.coeffRef(i) = m_matS.coeff(i, i);
00316 m_betas.coeffRef(i) = m_realQZ.matrixT().coeff(i,i);
00317 ++i;
00318 }
00319 else
00320 {
00321 Scalar p = Scalar(0.5) * (m_matS.coeff(i, i) - m_matS.coeff(i+1, i+1));
00322 Scalar z = sqrt(abs(p * p + m_matS.coeff(i+1, i) * m_matS.coeff(i, i+1)));
00323 m_alphas.coeffRef(i) = ComplexScalar(m_matS.coeff(i+1, i+1) + p, z);
00324 m_alphas.coeffRef(i+1) = ComplexScalar(m_matS.coeff(i+1, i+1) + p, -z);
00325
00326 m_betas.coeffRef(i) = m_realQZ.matrixT().coeff(i,i);
00327 m_betas.coeffRef(i+1) = m_realQZ.matrixT().coeff(i,i);
00328 i += 2;
00329 }
00330 }
00331 }
00332
00333 m_isInitialized = true;
00334 m_eigenvectorsOk = false;
00335
00336 return *this;
00337 }
00338
00339 }
00340
00341 #endif // EIGEN_GENERALIZEDEIGENSOLVER_H