Eigen::ComplexEigenSolver< MatrixType_ > Class Template Reference

Computes eigenvalues and eigenvectors of general complex matrices. More...

Public Types

enum  {
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  Options ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef std::complex< RealScalarComplexScalar
 Complex scalar type for MatrixType. More...
 
typedef Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &(~RowMajor), MaxColsAtCompileTime, 1 > EigenvalueType
 Type for vector of eigenvalues as returned by eigenvalues(). More...
 
typedef Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTimeEigenvectorType
 Type for matrix of eigenvectors as returned by eigenvectors(). More...
 
typedef Eigen::Index Index
 
typedef MatrixType_ MatrixType
 Synonym for the template parameter MatrixType_. More...
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type MatrixType. More...
 

Public Member Functions

 ComplexEigenSolver ()
 Default constructor. More...
 
template<typename InputType >
 ComplexEigenSolver (const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
 Constructor; computes eigendecomposition of given matrix. More...
 
 ComplexEigenSolver (Index size)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
ComplexEigenSolver< MatrixType > & compute (const EigenBase< InputType > &matrix, bool computeEigenvectors)
 
template<typename InputType >
ComplexEigenSolvercompute (const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
 Computes eigendecomposition of given matrix. More...
 
const EigenvalueTypeeigenvalues () const
 Returns the eigenvalues of given matrix. More...
 
const EigenvectorTypeeigenvectors () const
 Returns the eigenvectors of given matrix. More...
 
Index getMaxIterations ()
 Returns the maximum number of iterations. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
ComplexEigenSolversetMaxIterations (Index maxIters)
 Sets the maximum number of iterations allowed. More...
 

Protected Attributes

bool m_eigenvectorsOk
 
EigenvalueType m_eivalues
 
EigenvectorType m_eivec
 
bool m_isInitialized
 
EigenvectorType m_matX
 
ComplexSchur< MatrixTypem_schur
 

Private Member Functions

void doComputeEigenvectors (RealScalar matrixnorm)
 
void sortEigenvalues (bool computeEigenvectors)
 

Detailed Description

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

Computes eigenvalues and eigenvectors of general complex matrices.

This is defined in the Eigenvalues module.

Template Parameters
MatrixType_the type of the matrix of which we are computing the eigendecomposition; this is expected to be an instantiation of the Matrix class template.

The eigenvalues and eigenvectors of a matrix \( A \) are scalars \( \lambda \) and vectors \( v \) such that \( Av = \lambda v \). 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 = V D \). The matrix \( V \) is almost always invertible, in which case we have \( A = V D V^{-1} \). This is called the eigendecomposition.

The main function in this class is compute(), which computes the eigenvalues and eigenvectors of a given function. The documentation for that function contains an example showing the main features of the class.

See also
class EigenSolver, class SelfAdjointEigenSolver

Definition at line 47 of file ComplexEigenSolver.h.

Member Typedef Documentation

◆ ComplexScalar

template<typename MatrixType_ >
typedef std::complex<RealScalar> Eigen::ComplexEigenSolver< 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.

Definition at line 73 of file ComplexEigenSolver.h.

◆ EigenvalueType

template<typename MatrixType_ >
typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options&(~RowMajor), MaxColsAtCompileTime, 1> Eigen::ComplexEigenSolver< MatrixType_ >::EigenvalueType

Type for vector of eigenvalues as returned by eigenvalues().

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

Definition at line 80 of file ComplexEigenSolver.h.

◆ EigenvectorType

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.

Definition at line 87 of file ComplexEigenSolver.h.

◆ Index

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

Definition at line 65 of file ComplexEigenSolver.h.

◆ MatrixType

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

Synonym for the template parameter MatrixType_.

Definition at line 52 of file ComplexEigenSolver.h.

◆ RealScalar

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

Definition at line 64 of file ComplexEigenSolver.h.

◆ Scalar

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

Scalar type for matrices of type MatrixType.

Definition at line 63 of file ComplexEigenSolver.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 54 of file ComplexEigenSolver.h.

54  {
55  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57  Options = MatrixType::Options,
58  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
59  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
60  };

Constructor & Destructor Documentation

◆ ComplexEigenSolver() [1/3]

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

Default constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via compute().

Definition at line 94 of file ComplexEigenSolver.h.

95  : m_eivec(),
96  m_eivalues(),
97  m_schur(),
98  m_isInitialized(false),
99  m_eigenvectorsOk(false),
100  m_matX()
101  {}
ComplexSchur< MatrixType > m_schur

◆ ComplexEigenSolver() [2/3]

template<typename MatrixType_ >
Eigen::ComplexEigenSolver< MatrixType_ >::ComplexEigenSolver ( 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
ComplexEigenSolver()

Definition at line 109 of file ComplexEigenSolver.h.

110  : m_eivec(size, size),
111  m_eivalues(size),
112  m_schur(size),
113  m_isInitialized(false),
114  m_eigenvectorsOk(false),
115  m_matX(size, size)
116  {}

◆ ComplexEigenSolver() [3/3]

template<typename MatrixType_ >
template<typename InputType >
Eigen::ComplexEigenSolver< MatrixType_ >::ComplexEigenSolver ( const EigenBase< InputType > &  matrix,
bool  computeEigenvectors = true 
)
inlineexplicit

Constructor; computes eigendecomposition of given matrix.

Parameters
[in]matrixSquare 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 eigendecomposition.

Definition at line 128 of file ComplexEigenSolver.h.

129  : m_eivec(matrix.rows(),matrix.cols()),
130  m_eivalues(matrix.cols()),
131  m_schur(matrix.rows()),
132  m_isInitialized(false),
133  m_eigenvectorsOk(false),
134  m_matX(matrix.rows(),matrix.cols())
135  {
136  compute(matrix.derived(), computeEigenvectors);
137  }
ComplexEigenSolver & compute(const EigenBase< InputType > &matrix, bool computeEigenvectors=true)
Computes eigendecomposition of given matrix.

Member Function Documentation

◆ compute() [1/2]

template<typename MatrixType_ >
template<typename InputType >
ComplexEigenSolver<MatrixType>& Eigen::ComplexEigenSolver< MatrixType_ >::compute ( const EigenBase< InputType > &  matrix,
bool  computeEigenvectors 
)

Definition at line 260 of file ComplexEigenSolver.h.

261 {
262  // this code is inspired from Jampack
263  eigen_assert(matrix.cols() == matrix.rows());
264 
265  // Do a complex Schur decomposition, A = U T U^*
266  // The eigenvalues are on the diagonal of T.
267  m_schur.compute(matrix.derived(), computeEigenvectors);
268 
269  if(m_schur.info() == Success)
270  {
271  m_eivalues = m_schur.matrixT().diagonal();
272  if(computeEigenvectors)
274  sortEigenvalues(computeEigenvectors);
275  }
276 
277  m_isInitialized = true;
278  m_eigenvectorsOk = computeEigenvectors;
279  return *this;
280 }
#define eigen_assert(x)
Definition: Macros.h:902
void sortEigenvalues(bool computeEigenvectors)
void doComputeEigenvectors(RealScalar matrixnorm)
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: ComplexSchur.h:219
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition: ComplexSchur.h:164
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
@ Success
Definition: Constants.h:446

◆ compute() [2/2]

template<typename MatrixType_ >
template<typename InputType >
ComplexEigenSolver& Eigen::ComplexEigenSolver< MatrixType_ >::compute ( const EigenBase< InputType > &  matrix,
bool  computeEigenvectors = true 
)

Computes eigendecomposition of given matrix.

Parameters
[in]matrixSquare 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 complex 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 Schur form using the ComplexSchur class. The Schur decomposition is then used to compute the eigenvalues and eigenvectors.

The cost of the computation is dominated by the cost of the Schur decomposition, which is \( O(n^3) \) where \( n \) is the size of the matrix.

Example:

cout << "Here is a random 4x4 matrix, A:" << endl << A << endl << endl;
ComplexEigenSolver<MatrixXcf> ces;
ces.compute(A);
cout << "The eigenvalues of A are:" << endl << ces.eigenvalues() << endl;
cout << "The matrix of eigenvectors, V, is:" << endl << ces.eigenvectors() << endl << endl;
complex<float> lambda = ces.eigenvalues()[0];
cout << "Consider the first eigenvalue, lambda = " << lambda << endl;
VectorXcf v = ces.eigenvectors().col(0);
cout << "If v is the corresponding eigenvector, then lambda * v = " << endl << lambda * v << endl;
cout << "... and A * v = " << endl << A * v << endl << endl;
cout << "Finally, V * D * V^(-1) = " << endl
<< ces.eigenvectors() * ces.eigenvalues().asDiagonal() * ces.eigenvectors().inverse() << endl;
Array< int, Dynamic, 1 > v
MatrixXcf A
cout<< "The eigenvalues of A are:"<< endl<< ces.eigenvalues()<< endl;cout<< "The matrix of eigenvectors, V, is:"<< endl<< ces.eigenvectors()<< endl<< endl;complex< float > lambda
cout<< "Here is a random 4x4 matrix, A:"<< endl<< A<< endl<< endl;ComplexEigenSolver< MatrixXcf > ces
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< std::complex< float >, Dynamic, 1 > VectorXcf
DynamicĂ—1 vector of type std::complex<float>.
Definition: Matrix.h:503
Matrix< std::complex< float >, Dynamic, Dynamic > MatrixXcf
DynamicĂ—Dynamic matrix of type std::complex<float>.
Definition: Matrix.h:503

Output:

Here is a random 4x4 matrix, A:
  (-0.211,0.68)  (0.108,-0.444)   (0.435,0.271) (-0.198,-0.687)
  (0.597,0.566) (0.258,-0.0452)  (0.214,-0.717)  (-0.782,-0.74)
 (-0.605,0.823)  (0.0268,-0.27) (-0.514,-0.967)  (-0.563,0.998)
  (0.536,-0.33)   (0.832,0.904)  (0.608,-0.726)  (0.678,0.0259)

The eigenvalues of A are:
 (0.137,0.505)
 (-0.758,1.22)
 (1.52,-0.402)
(-0.691,-1.63)
The matrix of eigenvectors, V, is:
  (-0.246,-0.106)     (0.418,0.263)   (0.0417,-0.296)    (-0.122,0.271)
  (-0.205,-0.629)    (0.466,-0.457)    (0.244,-0.456)      (0.247,0.23)
 (-0.432,-0.0359) (-0.0651,-0.0146)    (-0.191,0.334)   (0.859,-0.0877)
    (-0.301,0.46)    (-0.41,-0.397)     (0.623,0.328)    (-0.116,0.195)

Consider the first eigenvalue, lambda = (0.137,0.505)
If v is the corresponding eigenvector, then lambda * v = 
 (0.0197,-0.139)
    (0.29,-0.19)
(-0.0412,-0.223)
(-0.274,-0.0891)
... and A * v = 
 (0.0197,-0.139)
    (0.29,-0.19)
(-0.0412,-0.223)
(-0.274,-0.0891)

Finally, V * D * V^(-1) = 
  (-0.211,0.68)  (0.108,-0.444)   (0.435,0.271) (-0.198,-0.687)
  (0.597,0.566) (0.258,-0.0452)  (0.214,-0.717)  (-0.782,-0.74)
 (-0.605,0.823)  (0.0268,-0.27) (-0.514,-0.967)  (-0.563,0.998)
  (0.536,-0.33)   (0.832,0.904)  (0.608,-0.726)  (0.678,0.0259)

◆ doComputeEigenvectors()

template<typename MatrixType >
void Eigen::ComplexEigenSolver< MatrixType >::doComputeEigenvectors ( RealScalar  matrixnorm)
private

Definition at line 284 of file ComplexEigenSolver.h.

285 {
286  const Index n = m_eivalues.size();
287 
288  matrixnorm = numext::maxi(matrixnorm,(std::numeric_limits<RealScalar>::min)());
289 
290  // Compute X such that T = X D X^(-1), where D is the diagonal of T.
291  // The matrix X is unit triangular.
292  m_matX = EigenvectorType::Zero(n, n);
293  for(Index k=n-1 ; k>=0 ; k--)
294  {
295  m_matX.coeffRef(k,k) = ComplexScalar(1.0,0.0);
296  // Compute X(i,k) using the (i,k) entry of the equation X T = D X
297  for(Index i=k-1 ; i>=0 ; i--)
298  {
299  m_matX.coeffRef(i,k) = -m_schur.matrixT().coeff(i,k);
300  if(k-i-1>0)
301  m_matX.coeffRef(i,k) -= (m_schur.matrixT().row(i).segment(i+1,k-i-1) * m_matX.col(k).segment(i+1,k-i-1)).value();
303  if(z==ComplexScalar(0))
304  {
305  // If the i-th and k-th eigenvalue are equal, then z equals 0.
306  // Use a small value instead, to prevent division by zero.
307  numext::real_ref(z) = NumTraits<RealScalar>::epsilon() * matrixnorm;
308  }
309  m_matX.coeffRef(i,k) = m_matX.coeff(i,k) / z;
310  }
311  }
312 
313  // Compute V as V = U X; now A = U T U^* = U X D X^(-1) U^* = V D V^(-1)
314  m_eivec.noalias() = m_schur.matrixU() * m_matX;
315  // .. and normalize the eigenvectors
316  for(Index k=0 ; k<n ; k++)
317  {
318  m_eivec.col(k).normalize();
319  }
320 }
int n
std::complex< RealScalar > ComplexScalar
Complex scalar type for MatrixType.
const ComplexMatrixType & matrixU() const
Returns the unitary matrix in the Schur decomposition.
Definition: ComplexSchur.h:140
constexpr const Scalar & coeff(Index rowId, Index colId) const
constexpr Scalar & coeffRef(Index rowId, Index colId)
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
internal::add_const_on_value_type_t< EIGEN_MATHFUNC_RETVAL(real_ref, Scalar) > real_ref(const Scalar &x)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)

◆ eigenvalues()

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

Returns the eigenvalues of given matrix.

Returns
A const reference to the column vector containing the eigenvalues.
Precondition
Either the constructor ComplexEigenSolver(const MatrixType& matrix, bool) or the member function compute(const MatrixType& matrix, bool) has been called before to compute the eigendecomposition of a matrix.

This function returns a column vector containing the eigenvalues. 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.

Example:

ComplexEigenSolver<MatrixXcf> ces(ones, /* computeEigenvectors = */ false);
cout << "The eigenvalues of the 3x3 matrix of ones are:"
<< endl << ces.eigenvalues() << endl;
static const ConstantReturnType Ones()

Output:

The eigenvalues of the 3x3 matrix of ones are:
(0,-0)
 (0,0)
 (3,0)

Definition at line 184 of file ComplexEigenSolver.h.

185  {
186  eigen_assert(m_isInitialized && "ComplexEigenSolver is not initialized.");
187  return m_eivalues;
188  }

◆ eigenvectors()

template<typename MatrixType_ >
const EigenvectorType& Eigen::ComplexEigenSolver< MatrixType_ >::eigenvectors ( ) const
inline

Returns the eigenvectors of given matrix.

Returns
A const reference to the matrix whose columns are the eigenvectors.
Precondition
Either the constructor ComplexEigenSolver(const MatrixType& matrix, bool) or the member function compute(const MatrixType& matrix, bool) has been called before to compute the eigendecomposition of a matrix, and computeEigenvectors was set to true (the default).

This function returns a matrix whose columns are the eigenvectors. Column \( k \) is an eigenvector corresponding to eigenvalue number \( k \) as returned by eigenvalues(). The eigenvectors are normalized to have (Euclidean) norm equal to one. The matrix returned by this function is the matrix \( V \) in the eigendecomposition \( A = V D V^{-1} \), if it exists.

Example:

ComplexEigenSolver<MatrixXcf> ces(ones);
cout << "The first eigenvector of the 3x3 matrix of ones is:"
<< endl << ces.eigenvectors().col(0) << endl;

Output:

The first eigenvector of the 3x3 matrix of ones is:
(-0.816,0)
 (0.408,0)
 (0.408,0)

Definition at line 159 of file ComplexEigenSolver.h.

160  {
161  eigen_assert(m_isInitialized && "ComplexEigenSolver is not initialized.");
162  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
163  return m_eivec;
164  }

◆ getMaxIterations()

template<typename MatrixType_ >
Index Eigen::ComplexEigenSolver< MatrixType_ >::getMaxIterations ( )
inline

Returns the maximum number of iterations.

Definition at line 235 of file ComplexEigenSolver.h.

236  {
237  return m_schur.getMaxIterations();
238  }
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: ComplexSchur.h:237

◆ info()

template<typename MatrixType_ >
ComputationInfo Eigen::ComplexEigenSolver< MatrixType_ >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NoConvergence otherwise.

Definition at line 221 of file ComplexEigenSolver.h.

222  {
223  eigen_assert(m_isInitialized && "ComplexEigenSolver is not initialized.");
224  return m_schur.info();
225  }

◆ setMaxIterations()

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

Sets the maximum number of iterations allowed.

Definition at line 228 of file ComplexEigenSolver.h.

229  {
230  m_schur.setMaxIterations(maxIters);
231  return *this;
232  }
ComplexSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: ComplexSchur.h:230

◆ sortEigenvalues()

template<typename MatrixType >
void Eigen::ComplexEigenSolver< MatrixType >::sortEigenvalues ( bool  computeEigenvectors)
private

Definition at line 324 of file ComplexEigenSolver.h.

325 {
326  const Index n = m_eivalues.size();
327  for (Index i=0; i<n; i++)
328  {
329  Index k;
330  m_eivalues.cwiseAbs().tail(n-i).minCoeff(&k);
331  if (k != 0)
332  {
333  k += i;
335  if(computeEigenvectors)
336  m_eivec.col(i).swap(m_eivec.col(k));
337  }
338  }
339 }
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:788

Member Data Documentation

◆ m_eigenvectorsOk

template<typename MatrixType_ >
bool Eigen::ComplexEigenSolver< MatrixType_ >::m_eigenvectorsOk
protected

Definition at line 248 of file ComplexEigenSolver.h.

◆ m_eivalues

template<typename MatrixType_ >
EigenvalueType Eigen::ComplexEigenSolver< MatrixType_ >::m_eivalues
protected

Definition at line 245 of file ComplexEigenSolver.h.

◆ m_eivec

template<typename MatrixType_ >
EigenvectorType Eigen::ComplexEigenSolver< MatrixType_ >::m_eivec
protected

Definition at line 244 of file ComplexEigenSolver.h.

◆ m_isInitialized

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

Definition at line 247 of file ComplexEigenSolver.h.

◆ m_matX

template<typename MatrixType_ >
EigenvectorType Eigen::ComplexEigenSolver< MatrixType_ >::m_matX
protected

Definition at line 249 of file ComplexEigenSolver.h.

◆ m_schur

template<typename MatrixType_ >
ComplexSchur<MatrixType> Eigen::ComplexEigenSolver< MatrixType_ >::m_schur
protected

Definition at line 246 of file ComplexEigenSolver.h.


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