Eigen::GeneralizedEigenSolver< MatrixType_ > Class Template Reference

Computes the generalized eigenvalues and eigenvectors of a pair of general matrices. More...

#include <GeneralizedEigenSolver.h>

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , Options = internal::traits<MatrixType>::Options , MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime ,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef MatrixType_ MatrixType
 Synonym for the template parameter MatrixType_. More...
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type MatrixType. More...
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Eigen::Index Index
 
typedef std::complex< RealScalarComplexScalar
 Complex scalar type for MatrixType. More...
 
typedef Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > VectorType
 Type for vector of real scalar values eigenvalues as returned by betas(). More...
 
typedef Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ComplexVectorType
 Type for vector of complex scalar values eigenvalues as returned by alphas(). More...
 
typedef CwiseBinaryOp< internal::scalar_quotient_op< ComplexScalar, Scalar >, ComplexVectorType, VectorTypeEigenvalueType
 Expression type for the eigenvalues as returned by eigenvalues(). More...
 
typedef Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTimeEigenvectorsType
 Type for matrix of eigenvectors as returned by eigenvectors(). More...
 

Public Member Functions

 GeneralizedEigenSolver ()
 Default constructor. More...
 
 GeneralizedEigenSolver (Index size)
 Default constructor with memory preallocation. More...
 
 GeneralizedEigenSolver (const MatrixType &A, const MatrixType &B, bool computeEigenvectors=true)
 Constructor; computes the generalized eigendecomposition of given matrix pair. More...
 
EigenvectorsType eigenvectors () const
 Returns the computed generalized eigenvectors. More...
 
EigenvalueType eigenvalues () const
 Returns an expression of the computed generalized eigenvalues. More...
 
const ComplexVectorTypealphas () const
 
const VectorTypebetas () const
 
GeneralizedEigenSolvercompute (const MatrixType &A, const MatrixType &B, bool computeEigenvectors=true)
 Computes generalized eigendecomposition of given matrix. More...
 
ComputationInfo info () const
 
GeneralizedEigenSolversetMaxIterations (Index maxIters)
 

Protected Attributes

EigenvectorsType m_eivec
 
ComplexVectorType m_alphas
 
VectorType m_betas
 
bool m_computeEigenvectors
 
bool m_isInitialized
 
RealQZ< MatrixTypem_realQZ
 
ComplexVectorType m_tmp
 

Detailed Description

template<typename MatrixType_>
class Eigen::GeneralizedEigenSolver< MatrixType_ >

Computes the generalized eigenvalues and eigenvectors of a pair of general matrices.

\eigenvalues_module

Template Parameters
MatrixType_the type of the matrices of which we are computing the eigen-decomposition; this is expected to be an instantiation of the Matrix class template. Currently, only real matrices are supported.

The generalized eigenvalues and eigenvectors of a matrix pair \( A \) and \( B \) are scalars \( \lambda \) and vectors \( v \) such that \( Av = \lambda Bv \). If \( D \) is a diagonal matrix with the eigenvalues on the diagonal, and \( V \) is a matrix with the eigenvectors as its columns, then \( A V = B V D \). The matrix \( V \) is almost always invertible, in which case we have \( A = B V D V^{-1} \). This is called the generalized eigen-decomposition.

The generalized eigenvalues and eigenvectors of a matrix pair may be complex, even when the matrices are real. Moreover, the generalized eigenvalue might be infinite if the matrix B is singular. To workaround this difficulty, the eigenvalues are provided as a pair of complex \( \alpha \) and real \( \beta \) such that: \( \lambda_i = \alpha_i / \beta_i \). If \( \beta_i \) is (nearly) zero, then one can consider the well defined left eigenvalue \( \mu = \beta_i / \alpha_i\) such that: \( \mu_i A v_i = B v_i \), or even \( \mu_i u_i^T A = u_i^T B \) where \( u_i \) is called the left eigenvector.

Call the function compute() to compute the generalized eigenvalues and eigenvectors of a given matrix pair. Alternatively, you can use the GeneralizedEigenSolver(const MatrixType&, const MatrixType&, bool) constructor which computes the eigenvalues and eigenvectors at construction time. Once the eigenvalue and eigenvectors are computed, they can be retrieved with the eigenvalues() and eigenvectors() functions.

Here is an usage example of this class: Example:

GeneralizedEigenSolver<MatrixXf> ges;
MatrixXf A = MatrixXf::Random(4, 4);
MatrixXf B = MatrixXf::Random(4, 4);
ges.compute(A, B);
cout << "The (complex) numerators of the generalzied eigenvalues are: " << ges.alphas().transpose() << endl;
cout << "The (real) denominatore of the generalzied eigenvalues are: " << ges.betas().transpose() << endl;
cout << "The (complex) generalzied eigenvalues are (alphas./beta): " << ges.eigenvalues().transpose() << endl;
GeneralizedEigenSolver< MatrixXf > ges
Definition: GeneralizedEigenSolver.cpp:1
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Definition: matrices.h:74

Output:

See also
MatrixBase::eigenvalues(), class ComplexEigenSolver, class SelfAdjointEigenSolver

Member Typedef Documentation

◆ ComplexScalar

template<typename MatrixType_ >
typedef std::complex<RealScalar> Eigen::GeneralizedEigenSolver< MatrixType_ >::ComplexScalar

Complex scalar type for MatrixType.

This is std::complex<Scalar> if Scalar is real (e.g., float or double) and just Scalar if Scalar is complex.

◆ ComplexVectorType

template<typename MatrixType_ >
typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> Eigen::GeneralizedEigenSolver< MatrixType_ >::ComplexVectorType

Type for vector of complex scalar values eigenvalues as returned by alphas().

This is a column vector with entries of type ComplexScalar. The length of the vector is the size of MatrixType.

◆ EigenvalueType

Expression type for the eigenvalues as returned by eigenvalues().

◆ EigenvectorsType

Type for matrix of eigenvectors as returned by eigenvectors().

This is a square matrix with entries of type ComplexScalar. The size is the same as the size of MatrixType.

◆ Index

template<typename MatrixType_ >
typedef Eigen::Index Eigen::GeneralizedEigenSolver< MatrixType_ >::Index
Deprecated:
since Eigen 3.3

◆ MatrixType

template<typename MatrixType_ >
typedef MatrixType_ Eigen::GeneralizedEigenSolver< MatrixType_ >::MatrixType

Synonym for the template parameter MatrixType_.

◆ RealScalar

template<typename MatrixType_ >
typedef NumTraits<Scalar>::Real Eigen::GeneralizedEigenSolver< MatrixType_ >::RealScalar

◆ Scalar

template<typename MatrixType_ >
typedef MatrixType::Scalar Eigen::GeneralizedEigenSolver< MatrixType_ >::Scalar

Scalar type for matrices of type MatrixType.

◆ VectorType

template<typename MatrixType_ >
typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> Eigen::GeneralizedEigenSolver< MatrixType_ >::VectorType

Type for vector of real scalar values eigenvalues as returned by betas().

This is a column vector with entries of type Scalar. The length of the vector is the size of MatrixType.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
67  {
68  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
69  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
70  Options = internal::traits<MatrixType>::Options,
71  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
73  };
@ ColsAtCompileTime
Definition: GeneralizedEigenSolver.h:69
@ MaxColsAtCompileTime
Definition: GeneralizedEigenSolver.h:72
@ MaxRowsAtCompileTime
Definition: GeneralizedEigenSolver.h:71
@ Options
Definition: GeneralizedEigenSolver.h:70
@ RowsAtCompileTime
Definition: GeneralizedEigenSolver.h:68

Constructor & Destructor Documentation

◆ GeneralizedEigenSolver() [1/3]

template<typename MatrixType_ >
Eigen::GeneralizedEigenSolver< MatrixType_ >::GeneralizedEigenSolver ( )
inline

Default constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via EigenSolver::compute(const MatrixType&, bool).

See also
compute() for an example.
124  : m_eivec(), m_alphas(), m_betas(), m_computeEigenvectors(false), m_isInitialized(false), m_realQZ() {}
bool m_isInitialized
Definition: GeneralizedEigenSolver.h:270
EigenvectorsType m_eivec
Definition: GeneralizedEigenSolver.h:266
ComplexVectorType m_alphas
Definition: GeneralizedEigenSolver.h:267
bool m_computeEigenvectors
Definition: GeneralizedEigenSolver.h:269
VectorType m_betas
Definition: GeneralizedEigenSolver.h:268
RealQZ< MatrixType > m_realQZ
Definition: GeneralizedEigenSolver.h:271

◆ GeneralizedEigenSolver() [2/3]

template<typename MatrixType_ >
Eigen::GeneralizedEigenSolver< MatrixType_ >::GeneralizedEigenSolver ( Index  size)
inlineexplicit

Default constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
GeneralizedEigenSolver()
133  : m_eivec(size, size),
134  m_alphas(size),
135  m_betas(size),
136  m_computeEigenvectors(false),
137  m_isInitialized(false),
138  m_realQZ(size),
139  m_tmp(size) {}
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
ComplexVectorType m_tmp
Definition: GeneralizedEigenSolver.h:272

◆ GeneralizedEigenSolver() [3/3]

template<typename MatrixType_ >
Eigen::GeneralizedEigenSolver< MatrixType_ >::GeneralizedEigenSolver ( const MatrixType A,
const MatrixType B,
bool  computeEigenvectors = true 
)
inline

Constructor; computes the generalized eigendecomposition of given matrix pair.

Parameters
[in]ASquare matrix whose eigendecomposition is to be computed.
[in]BSquare matrix whose eigendecomposition is to be computed.
[in]computeEigenvectorsIf true, both the eigenvectors and the eigenvalues are computed; if false, only the eigenvalues are computed.

This constructor calls compute() to compute the generalized eigenvalues and eigenvectors.

See also
compute()
154  : m_eivec(A.rows(), A.cols()),
155  m_alphas(A.cols()),
156  m_betas(A.cols()),
157  m_computeEigenvectors(false),
158  m_isInitialized(false),
159  m_realQZ(A.cols()),
160  m_tmp(A.cols()) {
161  compute(A, B, computeEigenvectors);
162  }
GeneralizedEigenSolver & compute(const MatrixType &A, const MatrixType &B, bool computeEigenvectors=true)
Computes generalized eigendecomposition of given matrix.
Definition: GeneralizedEigenSolver.h:276

References Eigen::GeneralizedEigenSolver< MatrixType_ >::compute().

Member Function Documentation

◆ alphas()

template<typename MatrixType_ >
const ComplexVectorType& Eigen::GeneralizedEigenSolver< MatrixType_ >::alphas ( ) const
inline
Returns
A const reference to the vectors containing the alpha values

This vector permits to reconstruct the j-th eigenvalues as alphas(i)/betas(j).

See also
betas(), eigenvalues()
210  {
211  eigen_assert(info() == Success && "GeneralizedEigenSolver failed to compute alphas.");
212  return m_alphas;
213  }
#define eigen_assert(x)
Definition: Macros.h:910
ComputationInfo info() const
Definition: GeneralizedEigenSolver.h:250
@ Success
Definition: Constants.h:440

References eigen_assert, Eigen::GeneralizedEigenSolver< MatrixType_ >::info(), Eigen::GeneralizedEigenSolver< MatrixType_ >::m_alphas, and Eigen::Success.

Referenced by generalized_eigensolver_assert(), and generalized_eigensolver_real().

◆ betas()

template<typename MatrixType_ >
const VectorType& Eigen::GeneralizedEigenSolver< MatrixType_ >::betas ( ) const
inline
Returns
A const reference to the vectors containing the beta values

This vector permits to reconstruct the j-th eigenvalues as alphas(i)/betas(j).

See also
alphas(), eigenvalues()
220  {
221  eigen_assert(info() == Success && "GeneralizedEigenSolver failed to compute betas.");
222  return m_betas;
223  }

References eigen_assert, Eigen::GeneralizedEigenSolver< MatrixType_ >::info(), Eigen::GeneralizedEigenSolver< MatrixType_ >::m_betas, and Eigen::Success.

Referenced by generalized_eigensolver_assert(), and generalized_eigensolver_real().

◆ compute()

template<typename MatrixType >
GeneralizedEigenSolver< MatrixType > & Eigen::GeneralizedEigenSolver< MatrixType >::compute ( const MatrixType A,
const MatrixType B,
bool  computeEigenvectors = true 
)

Computes generalized eigendecomposition of given matrix.

Parameters
[in]ASquare matrix whose eigendecomposition is to be computed.
[in]BSquare matrix whose eigendecomposition is to be computed.
[in]computeEigenvectorsIf true, both the eigenvectors and the eigenvalues are computed; if false, only the eigenvalues are computed.
Returns
Reference to *this

This function computes the eigenvalues of the real matrix matrix. The eigenvalues() function can be used to retrieve them. If computeEigenvectors is true, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().

The matrix is first reduced to real generalized Schur form using the RealQZ class. The generalized Schur decomposition is then used to compute the eigenvalues and eigenvectors.

The cost of the computation is dominated by the cost of the generalized Schur decomposition.

This method reuses of the allocated data in the GeneralizedEigenSolver object.

278  {
279  using std::abs;
280  using std::sqrt;
281  eigen_assert(A.cols() == A.rows() && B.cols() == A.rows() && B.cols() == B.rows());
282  Index size = A.cols();
283  // Reduce to generalized real Schur form:
284  // A = Q S Z and B = Q T Z
285  m_realQZ.compute(A, B, computeEigenvectors);
286  if (m_realQZ.info() == Success) {
287  // Resize storage
290  if (computeEigenvectors) {
292  m_tmp.resize(size);
293  }
294 
295  // Aliases:
296  Map<VectorType> v(reinterpret_cast<Scalar*>(m_tmp.data()), size);
297  ComplexVectorType& cv = m_tmp;
298  const MatrixType& mS = m_realQZ.matrixS();
299  const MatrixType& mT = m_realQZ.matrixT();
300 
301  Index i = 0;
302  while (i < size) {
303  if (i == size - 1 || mS.coeff(i + 1, i) == Scalar(0)) {
304  // Real eigenvalue
305  m_alphas.coeffRef(i) = mS.diagonal().coeff(i);
306  m_betas.coeffRef(i) = mT.diagonal().coeff(i);
307  if (computeEigenvectors) {
308  v.setConstant(Scalar(0.0));
309  v.coeffRef(i) = Scalar(1.0);
310  // For singular eigenvalues do nothing more
312  // Non-singular eigenvalue
313  const Scalar alpha = real(m_alphas.coeffRef(i));
314  const Scalar beta = m_betas.coeffRef(i);
315  for (Index j = i - 1; j >= 0; j--) {
316  const Index st = j + 1;
317  const Index sz = i - j;
318  if (j > 0 && mS.coeff(j, j - 1) != Scalar(0)) {
319  // 2x2 block
320  Matrix<Scalar, 2, 1> rhs = (alpha * mT.template block<2, Dynamic>(j - 1, st, 2, sz) -
321  beta * mS.template block<2, Dynamic>(j - 1, st, 2, sz))
322  .lazyProduct(v.segment(st, sz));
323  Matrix<Scalar, 2, 2> lhs =
324  beta * mS.template block<2, 2>(j - 1, j - 1) - alpha * mT.template block<2, 2>(j - 1, j - 1);
325  v.template segment<2>(j - 1) = lhs.partialPivLu().solve(rhs);
326  j--;
327  } else {
328  v.coeffRef(j) = -v.segment(st, sz)
329  .transpose()
330  .cwiseProduct(beta * mS.block(j, st, 1, sz) - alpha * mT.block(j, st, 1, sz))
331  .sum() /
332  (beta * mS.coeffRef(j, j) - alpha * mT.coeffRef(j, j));
333  }
334  }
335  }
336  m_eivec.col(i).real().noalias() = m_realQZ.matrixZ().transpose() * v;
337  m_eivec.col(i).real().normalize();
338  m_eivec.col(i).imag().setConstant(0);
339  }
340  ++i;
341  } else {
342  // We need to extract the generalized eigenvalues of the pair of a general 2x2 block S and a positive diagonal
343  // 2x2 block T Then taking beta=T_00*T_11, we can avoid any division, and alpha is the eigenvalues of A = (U^-1
344  // * S * U) * diag(T_11,T_00):
345 
346  // T = [a 0]
347  // [0 b]
348  RealScalar a = mT.diagonal().coeff(i), b = mT.diagonal().coeff(i + 1);
349  const RealScalar beta = m_betas.coeffRef(i) = m_betas.coeffRef(i + 1) = a * b;
350 
351  // ^^ NOTE: using diagonal()(i) instead of coeff(i,i) workarounds a MSVC bug.
352  Matrix<RealScalar, 2, 2> S2 = mS.template block<2, 2>(i, i) * Matrix<Scalar, 2, 1>(b, a).asDiagonal();
353 
354  Scalar p = Scalar(0.5) * (S2.coeff(0, 0) - S2.coeff(1, 1));
355  Scalar z = sqrt(abs(p * p + S2.coeff(1, 0) * S2.coeff(0, 1)));
356  const ComplexScalar alpha = ComplexScalar(S2.coeff(1, 1) + p, (beta > 0) ? z : -z);
358  m_alphas.coeffRef(i + 1) = alpha;
359 
360  if (computeEigenvectors) {
361  // Compute eigenvector in position (i+1) and then position (i) is just the conjugate
362  cv.setZero();
363  cv.coeffRef(i + 1) = Scalar(1.0);
364  // here, the "static_cast" workaround expression template issues.
365  cv.coeffRef(i) = -(static_cast<Scalar>(beta * mS.coeffRef(i, i + 1)) - alpha * mT.coeffRef(i, i + 1)) /
366  (static_cast<Scalar>(beta * mS.coeffRef(i, i)) - alpha * mT.coeffRef(i, i));
367  for (Index j = i - 1; j >= 0; j--) {
368  const Index st = j + 1;
369  const Index sz = i + 1 - j;
370  if (j > 0 && mS.coeff(j, j - 1) != Scalar(0)) {
371  // 2x2 block
372  Matrix<ComplexScalar, 2, 1> rhs = (alpha * mT.template block<2, Dynamic>(j - 1, st, 2, sz) -
373  beta * mS.template block<2, Dynamic>(j - 1, st, 2, sz))
374  .lazyProduct(cv.segment(st, sz));
375  Matrix<ComplexScalar, 2, 2> lhs =
376  beta * mS.template block<2, 2>(j - 1, j - 1) - alpha * mT.template block<2, 2>(j - 1, j - 1);
377  cv.template segment<2>(j - 1) = lhs.partialPivLu().solve(rhs);
378  j--;
379  } else {
380  cv.coeffRef(j) = cv.segment(st, sz)
381  .transpose()
382  .cwiseProduct(beta * mS.block(j, st, 1, sz) - alpha * mT.block(j, st, 1, sz))
383  .sum() /
384  (alpha * mT.coeffRef(j, j) - static_cast<Scalar>(beta * mS.coeffRef(j, j)));
385  }
386  }
387  m_eivec.col(i + 1).noalias() = (m_realQZ.matrixZ().transpose() * cv);
388  m_eivec.col(i + 1).normalize();
389  m_eivec.col(i) = m_eivec.col(i + 1).conjugate();
390  }
391  i += 2;
392  }
393  }
394  }
395  m_computeEigenvectors = computeEigenvectors;
396  m_isInitialized = true;
397  return *this;
398 }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
m block< 2, Dynamic >(1, 1, 2, 3).setZero()
float * p
Definition: Tutorial_Map_using.cpp:9
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType.
Definition: GeneralizedEigenSolver.h:86
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ComplexVectorType
Type for vector of complex scalar values eigenvalues as returned by alphas().
Definition: GeneralizedEigenSolver.h:100
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType.
Definition: GeneralizedEigenSolver.h:76
Eigen::Index Index
Definition: GeneralizedEigenSolver.h:78
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:217
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
Definition: CwiseNullaryOp.h:365
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
RealQZ & compute(const MatrixType &A, const MatrixType &B, bool computeQZ=true)
Computes QZ decomposition of given matrix.
Definition: RealQZ.h:502
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealQZ.h:170
const MatrixType & matrixZ() const
Returns matrix Z in the QZ decomposition.
Definition: RealQZ.h:133
const MatrixType & matrixT() const
Returns matrix S in the QZ decomposition.
Definition: RealQZ.h:152
const MatrixType & matrixS() const
Returns matrix S in the QZ decomposition.
Definition: RealQZ.h:143
#define min(a, b)
Definition: datatypes.h:22
RealScalar alpha
Definition: level1_cplx_impl.h:151
const Scalar * a
Definition: level2_cplx_impl.h:32
Scalar beta
Definition: level2_cplx_impl.h:36
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:486
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:482
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References a, abs(), alpha, b, beta, block< 2, Dynamic >(), Eigen::PlainObjectBase< Derived >::coeff(), Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >::coeffRef(), Eigen::PlainObjectBase< Derived >::cols(), Eigen::RealQZ< MatrixType_ >::compute(), Eigen::conj(), Eigen::PlainObjectBase< Derived >::data(), eigen_assert, i, Eigen::RealQZ< MatrixType_ >::info(), j, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_alphas, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_betas, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_computeEigenvectors, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_eivec, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_isInitialized, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_realQZ, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_tmp, Eigen::RealQZ< MatrixType_ >::matrixS(), Eigen::RealQZ< MatrixType_ >::matrixT(), Eigen::RealQZ< MatrixType_ >::matrixZ(), min, p, Eigen::real(), Eigen::PlainObjectBase< Derived >::resize(), Eigen::PlainObjectBase< Derived >::rows(), Eigen::PlainObjectBase< Derived >::setConstant(), Eigen::PlainObjectBase< Derived >::setZero(), size, sqrt(), Eigen::Success, and v.

Referenced by generalized_eigensolver_assert(), generalized_eigensolver_real(), and Eigen::GeneralizedEigenSolver< MatrixType_ >::GeneralizedEigenSolver().

◆ eigenvalues()

template<typename MatrixType_ >
EigenvalueType Eigen::GeneralizedEigenSolver< MatrixType_ >::eigenvalues ( ) const
inline

Returns an expression of the computed generalized eigenvalues.

Returns
An expression of the column vector containing the eigenvalues.

It is a shortcut for

this->alphas().cwiseQuotient(this->betas());
const VectorType & betas() const
Definition: GeneralizedEigenSolver.h:220
const ComplexVectorType & alphas() const
Definition: GeneralizedEigenSolver.h:210

Not that betas might contain zeros. It is therefore not recommended to use this function, but rather directly deal with the alphas and betas vectors.

Precondition
Either the constructor GeneralizedEigenSolver(const MatrixType&,const MatrixType&,bool) or the member function compute(const MatrixType&,const MatrixType&,bool) has been called before.

The eigenvalues are repeated according to their algebraic multiplicity, so there are as many eigenvalues as rows in the matrix. The eigenvalues are not sorted in any particular order.

See also
alphas(), betas(), eigenvectors()
200  {
201  eigen_assert(info() == Success && "GeneralizedEigenSolver failed to compute eigenvalues.");
203  }
CwiseBinaryOp< internal::scalar_quotient_op< ComplexScalar, Scalar >, ComplexVectorType, VectorType > EigenvalueType
Expression type for the eigenvalues as returned by eigenvalues().
Definition: GeneralizedEigenSolver.h:105

References eigen_assert, Eigen::GeneralizedEigenSolver< MatrixType_ >::info(), Eigen::GeneralizedEigenSolver< MatrixType_ >::m_alphas, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_betas, and Eigen::Success.

Referenced by generalized_eigensolver_assert(), and generalized_eigensolver_real().

◆ eigenvectors()

template<typename MatrixType_ >
EigenvectorsType Eigen::GeneralizedEigenSolver< MatrixType_ >::eigenvectors ( ) const
inline

Returns the computed generalized eigenvectors.

Returns
Matrix whose columns are the (possibly complex) right eigenvectors. i.e. the eigenvectors that solve (A - l*B)x = 0. The ordering matches the eigenvalues.
Precondition
Either the constructor GeneralizedEigenSolver(const MatrixType&,const MatrixType&, bool) or the member function compute(const MatrixType&, const MatrixType& bool) has been called before, and computeEigenvectors was set to true (the default).
See also
eigenvalues()
176  {
177  eigen_assert(info() == Success && "GeneralizedEigenSolver failed to compute eigenvectors");
178  eigen_assert(m_computeEigenvectors && "Eigenvectors for GeneralizedEigenSolver were not calculated");
179  return m_eivec;
180  }

References eigen_assert, Eigen::GeneralizedEigenSolver< MatrixType_ >::info(), Eigen::GeneralizedEigenSolver< MatrixType_ >::m_computeEigenvectors, Eigen::GeneralizedEigenSolver< MatrixType_ >::m_eivec, and Eigen::Success.

Referenced by generalized_eigensolver_assert(), and generalized_eigensolver_real().

◆ info()

◆ setMaxIterations()

template<typename MatrixType_ >
GeneralizedEigenSolver& Eigen::GeneralizedEigenSolver< MatrixType_ >::setMaxIterations ( Index  maxIters)
inline

Sets the maximal number of iterations allowed.

257  {
258  m_realQZ.setMaxIterations(maxIters);
259  return *this;
260  }
RealQZ & setMaxIterations(Index maxIters)
Definition: RealQZ.h:185

References Eigen::GeneralizedEigenSolver< MatrixType_ >::m_realQZ, and Eigen::RealQZ< MatrixType_ >::setMaxIterations().

Referenced by generalized_eigensolver_assert().

Member Data Documentation

◆ m_alphas

◆ m_betas

◆ m_computeEigenvectors

template<typename MatrixType_ >
bool Eigen::GeneralizedEigenSolver< MatrixType_ >::m_computeEigenvectors
protected

◆ m_eivec

◆ m_isInitialized

template<typename MatrixType_ >
bool Eigen::GeneralizedEigenSolver< MatrixType_ >::m_isInitialized
protected

◆ m_realQZ

◆ m_tmp

template<typename MatrixType_ >
ComplexVectorType Eigen::GeneralizedEigenSolver< MatrixType_ >::m_tmp
protected

The documentation for this class was generated from the following file: