Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD > Class Template Reference

Public Types

typedef MatrixType::Index Index
 
typedef NumTraits< Scalar >::Real RealScalar
 Real scalar type for MatrixType. More...
 
typedef internal::plain_col_type< MatrixType, RealScalar >::type RealVectorType
 Type for vector of eigenvalues as returned by eigenvalues(). More...
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type MatrixType. More...
 

Public Member Functions

 ArpackGeneralizedSelfAdjointEigenSolver ()
 Default constructor. More...
 
 ArpackGeneralizedSelfAdjointEigenSolver (const MatrixType &A, const MatrixType &B, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0)
 Constructor; computes generalized eigenvalues of given matrix with respect to another matrix. More...
 
 ArpackGeneralizedSelfAdjointEigenSolver (const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0)
 Constructor; computes eigenvalues of given matrix. More...
 
ArpackGeneralizedSelfAdjointEigenSolvercompute (const MatrixType &A, const MatrixType &B, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0)
 Computes generalized eigenvalues / eigenvectors of given matrix using the external ARPACK library. More...
 
ArpackGeneralizedSelfAdjointEigenSolvercompute (const MatrixType &A, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0)
 Computes eigenvalues / eigenvectors of given matrix using the external ARPACK library. More...
 
const Matrix< Scalar, Dynamic, 1 > & eigenvalues () const
 Returns the eigenvalues of given matrix. More...
 
const Matrix< Scalar, Dynamic, Dynamic > & eigenvectors () const
 Returns the eigenvectors of given matrix. More...
 
size_t getNbrConvergedEigenValues () const
 
size_t getNbrIterations () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
Matrix< Scalar, Dynamic, DynamicoperatorInverseSqrt () const
 Computes the inverse square root of the matrix. More...
 
Matrix< Scalar, Dynamic, DynamicoperatorSqrt () const
 Computes the positive-definite square root of the matrix. More...
 

Protected Attributes

bool m_eigenvectorsOk
 
Matrix< Scalar, Dynamic, 1 > m_eivalues
 
Matrix< Scalar, Dynamic, Dynamicm_eivec
 
ComputationInfo m_info
 
bool m_isInitialized
 
size_t m_nbrConverged
 
size_t m_nbrIterations
 

Detailed Description

template<typename MatrixType, typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
class Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >

Definition at line 27 of file ArpackSelfAdjointEigenSolver.h.

Member Typedef Documentation

◆ Index

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
typedef MatrixType::Index Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::Index

Definition at line 34 of file ArpackSelfAdjointEigenSolver.h.

◆ RealScalar

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
typedef NumTraits<Scalar>::Real Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::RealScalar

Real scalar type for MatrixType.

This is just Scalar if Scalar is real (e.g., float or Scalar), and the type of the real part of Scalar if Scalar is complex.

Definition at line 42 of file ArpackSelfAdjointEigenSolver.h.

◆ RealVectorType

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
typedef internal::plain_col_type<MatrixType, RealScalar>::type Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::RealVectorType

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

This is a column vector with entries of type RealScalar. The length of the vector is the size of nbrEigenvalues.

Definition at line 49 of file ArpackSelfAdjointEigenSolver.h.

◆ Scalar

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
typedef MatrixType::Scalar Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::Scalar

Scalar type for matrices of type MatrixType.

Definition at line 33 of file ArpackSelfAdjointEigenSolver.h.

Constructor & Destructor Documentation

◆ ArpackGeneralizedSelfAdjointEigenSolver() [1/3]

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::ArpackGeneralizedSelfAdjointEigenSolver ( )
inline

Default constructor.

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

Definition at line 57 of file ArpackSelfAdjointEigenSolver.h.

◆ ArpackGeneralizedSelfAdjointEigenSolver() [2/3]

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::ArpackGeneralizedSelfAdjointEigenSolver ( const MatrixType A,
const MatrixType B,
Index  nbrEigenvalues,
std::string  eigs_sigma = "LM",
int  options = ComputeEigenvectors,
RealScalar  tol = 0.0 
)
inline

Constructor; computes generalized eigenvalues of given matrix with respect to another matrix.

Parameters
[in]ASelf-adjoint matrix whose eigenvalues / eigenvectors will computed. By default, the upper triangular part is used, but can be changed through the template parameter.
[in]BSelf-adjoint matrix for the generalized eigenvalue problem.
[in]nbrEigenvaluesThe number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned.
[in]eigs_sigmaString containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
[in]tolWhat tolerance to find the eigenvalues to. Default is 0, which means machine precision.

This constructor calls compute(const MatrixType&, const MatrixType&, Index, string, int, RealScalar) to compute the eigenvalues of the matrix A with respect to B. The eigenvectors are computed if options equals ComputeEigenvectors.

Definition at line 88 of file ArpackSelfAdjointEigenSolver.h.

91  : m_eivec(),
92  m_eivalues(),
93  m_isInitialized(false),
94  m_eigenvectorsOk(false),
95  m_nbrConverged(0),
97  {
98  compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
99  }
ArpackGeneralizedSelfAdjointEigenSolver & compute(const MatrixType &A, const MatrixType &B, Index nbrEigenvalues, std::string eigs_sigma="LM", int options=ComputeEigenvectors, RealScalar tol=0.0)
Computes generalized eigenvalues / eigenvectors of given matrix using the external ARPACK library.

◆ ArpackGeneralizedSelfAdjointEigenSolver() [3/3]

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::ArpackGeneralizedSelfAdjointEigenSolver ( const MatrixType A,
Index  nbrEigenvalues,
std::string  eigs_sigma = "LM",
int  options = ComputeEigenvectors,
RealScalar  tol = 0.0 
)
inline

Constructor; computes eigenvalues of given matrix.

Parameters
[in]ASelf-adjoint matrix whose eigenvalues / eigenvectors will computed. By default, the upper triangular part is used, but can be changed through the template parameter.
[in]nbrEigenvaluesThe number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned.
[in]eigs_sigmaString containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
[in]tolWhat tolerance to find the eigenvalues to. Default is 0, which means machine precision.

This constructor calls compute(const MatrixType&, Index, string, int, RealScalar) to compute the eigenvalues of the matrix A. The eigenvectors are computed if options equals ComputeEigenvectors.

Definition at line 123 of file ArpackSelfAdjointEigenSolver.h.

126  : m_eivec(),
127  m_eivalues(),
128  m_isInitialized(false),
129  m_eigenvectorsOk(false),
130  m_nbrConverged(0),
131  m_nbrIterations(0)
132  {
133  compute(A, nbrEigenvalues, eigs_sigma, options, tol);
134  }

Member Function Documentation

◆ compute() [1/2]

template<typename MatrixType , typename MatrixSolver , bool BisSPD>
ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD > & Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute ( const MatrixType A,
const MatrixType B,
Index  nbrEigenvalues,
std::string  eigs_sigma = "LM",
int  options = ComputeEigenvectors,
RealScalar  tol = 0.0 
)

Computes generalized eigenvalues / eigenvectors of given matrix using the external ARPACK library.

Parameters
[in]ASelfadjoint matrix whose eigendecomposition is to be computed.
[in]BSelfadjoint matrix for generalized eigenvalues.
[in]nbrEigenvaluesThe number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned.
[in]eigs_sigmaString containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
[in]tolWhat tolerance to find the eigenvalues to. Default is 0, which means machine precision.
Returns
Reference to *this

This function computes the generalized eigenvalues of A with respect to B using ARPACK. The eigenvalues() function can be used to retrieve them. If options equals ComputeEigenvectors, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().

Definition at line 335 of file ArpackSelfAdjointEigenSolver.h.

337 {
338  eigen_assert(A.cols() == A.rows());
339  eigen_assert(B.cols() == B.rows());
340  eigen_assert(B.rows() == 0 || A.cols() == B.rows());
341  eigen_assert((options &~ (EigVecMask | GenEigMask)) == 0
342  && (options & EigVecMask) != EigVecMask
343  && "invalid option parameter");
344 
345  bool isBempty = (B.rows() == 0) || (B.cols() == 0);
346 
347  // For clarity, all parameters match their ARPACK name
348  //
349  // Always 0 on the first call
350  //
351  int ido = 0;
352 
353  int n = (int)A.cols();
354 
355  // User options: "LA", "SA", "SM", "LM", "BE"
356  //
357  char whch[3] = "LM";
358 
359  // Specifies the shift if iparam[6] = { 3, 4, 5 }, not used if iparam[6] = { 1, 2 }
360  //
361  RealScalar sigma = 0.0;
362 
363  if (eigs_sigma.length() >= 2 && isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1]))
364  {
365  eigs_sigma[0] = toupper(eigs_sigma[0]);
366  eigs_sigma[1] = toupper(eigs_sigma[1]);
367 
368  // In the following special case we're going to invert the problem, since solving
369  // for larger magnitude is much much faster
370  // i.e., if 'SM' is specified, we're going to really use 'LM', the default
371  //
372  if (eigs_sigma.substr(0,2) != "SM")
373  {
374  whch[0] = eigs_sigma[0];
375  whch[1] = eigs_sigma[1];
376  }
377  }
378  else
379  {
380  eigen_assert(false && "Specifying clustered eigenvalues is not yet supported!");
381 
382  // If it's not scalar values, then the user may be explicitly
383  // specifying the sigma value to cluster the evs around
384  //
385  sigma = atof(eigs_sigma.c_str());
386 
387  // If atof fails, it returns 0.0, which is a fine default
388  //
389  }
390 
391  // "I" means normal eigenvalue problem, "G" means generalized
392  //
393  char bmat[2] = "I";
394  if (eigs_sigma.substr(0,2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) || (!isBempty && !BisSPD))
395  bmat[0] = 'G';
396 
397  // Now we determine the mode to use
398  //
399  int mode = (bmat[0] == 'G') + 1;
400  if (eigs_sigma.substr(0,2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])))
401  {
402  // We're going to use shift-and-invert mode, and basically find
403  // the largest eigenvalues of the inverse operator
404  //
405  mode = 3;
406  }
407 
408  // The user-specified number of eigenvalues/vectors to compute
409  //
410  int nev = (int)nbrEigenvalues;
411 
412  // Allocate space for ARPACK to store the residual
413  //
414  Scalar *resid = new Scalar[n];
415 
416  // Number of Lanczos vectors, must satisfy nev < ncv <= n
417  // Note that this indicates that nev != n, and we cannot compute
418  // all eigenvalues of a mtrix
419  //
420  int ncv = std::min(std::max(2*nev, 20), n);
421 
422  // The working n x ncv matrix, also store the final eigenvectors (if computed)
423  //
424  Scalar *v = new Scalar[n*ncv];
425  int ldv = n;
426 
427  // Working space
428  //
429  Scalar *workd = new Scalar[3*n];
430  int lworkl = ncv*ncv+8*ncv; // Must be at least this length
431  Scalar *workl = new Scalar[lworkl];
432 
433  int *iparam= new int[11];
434  iparam[0] = 1; // 1 means we let ARPACK perform the shifts, 0 means we'd have to do it
435  iparam[2] = std::max(300, (int)std::ceil(2*n/std::max(ncv,1)));
436  iparam[6] = mode; // The mode, 1 is standard ev problem, 2 for generalized ev, 3 for shift-and-invert
437 
438  // Used during reverse communicate to notify where arrays start
439  //
440  int *ipntr = new int[11];
441 
442  // Error codes are returned in here, initial value of 0 indicates a random initial
443  // residual vector is used, any other values means resid contains the initial residual
444  // vector, possibly from a previous run
445  //
446  int info = 0;
447 
448  Scalar scale = 1.0;
449  //if (!isBempty)
450  //{
451  //Scalar scale = B.norm() / std::sqrt(n);
452  //scale = std::pow(2, std::floor(std::log(scale+1)));
454  //for (size_t i=0; i<(size_t)B.outerSize(); i++)
455  // for (typename MatrixType::InnerIterator it(B, i); it; ++it)
456  // it.valueRef() /= scale;
457  //}
458 
459  MatrixSolver OP;
460  if (mode == 1 || mode == 2)
461  {
462  if (!isBempty)
463  OP.compute(B);
464  }
465  else if (mode == 3)
466  {
467  if (sigma == 0.0)
468  {
469  OP.compute(A);
470  }
471  else
472  {
473  // Note: We will never enter here because sigma must be 0.0
474  //
475  if (isBempty)
476  {
477  MatrixType AminusSigmaB(A);
478  for (Index i=0; i<A.rows(); ++i)
479  AminusSigmaB.coeffRef(i,i) -= sigma;
480 
481  OP.compute(AminusSigmaB);
482  }
483  else
484  {
485  MatrixType AminusSigmaB = A - sigma * B;
486  OP.compute(AminusSigmaB);
487  }
488  }
489  }
490 
491  if (!(mode == 1 && isBempty) && !(mode == 2 && isBempty) && OP.info() != Success)
492  std::cout << "Error factoring matrix" << std::endl;
493 
494  do
495  {
496  internal::arpack_wrapper<Scalar, RealScalar>::saupd(&ido, bmat, &n, whch, &nev, &tol, resid,
497  &ncv, v, &ldv, iparam, ipntr, workd, workl,
498  &lworkl, &info);
499 
500  if (ido == -1 || ido == 1)
501  {
502  Scalar *in = workd + ipntr[0] - 1;
503  Scalar *out = workd + ipntr[1] - 1;
504 
505  if (ido == 1 && mode != 2)
506  {
507  Scalar *out2 = workd + ipntr[2] - 1;
508  if (isBempty || mode == 1)
510  else
512 
513  in = workd + ipntr[2] - 1;
514  }
515 
516  if (mode == 1)
517  {
518  if (isBempty)
519  {
520  // OP = A
521  //
523  }
524  else
525  {
526  // OP = L^{-1}AL^{-T}
527  //
528  internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::applyOP(OP, A, n, in, out);
529  }
530  }
531  else if (mode == 2)
532  {
533  if (ido == 1)
535 
536  // OP = B^{-1} A
537  //
539  }
540  else if (mode == 3)
541  {
542  // OP = (A-\sigmaB)B (\sigma could be 0, and B could be I)
543  // The B * in is already computed and stored at in if ido == 1
544  //
545  if (ido == 1 || isBempty)
547  else
549  }
550  }
551  else if (ido == 2)
552  {
553  Scalar *in = workd + ipntr[0] - 1;
554  Scalar *out = workd + ipntr[1] - 1;
555 
556  if (isBempty || mode == 1)
558  else
560  }
561  } while (ido != 99);
562 
563  if (info == 1)
565  else if (info == 3)
567  else if (info < 0)
569  else if (info != 0)
570  eigen_assert(false && "Unknown ARPACK return value!");
571  else
572  {
573  // Do we compute eigenvectors or not?
574  //
575  int rvec = (options & ComputeEigenvectors) == ComputeEigenvectors;
576 
577  // "A" means "All", use "S" to choose specific eigenvalues (not yet supported in ARPACK))
578  //
579  char howmny[2] = "A";
580 
581  // if howmny == "S", specifies the eigenvalues to compute (not implemented in ARPACK)
582  //
583  int *select = new int[ncv];
584 
585  // Final eigenvalues
586  //
587  m_eivalues.resize(nev, 1);
588 
589  internal::arpack_wrapper<Scalar, RealScalar>::seupd(&rvec, howmny, select, m_eivalues.data(), v, &ldv,
590  &sigma, bmat, &n, whch, &nev, &tol, resid, &ncv,
591  v, &ldv, iparam, ipntr, workd, workl, &lworkl, &info);
592 
593  if (info == -14)
595  else if (info != 0)
597  else
598  {
599  if (rvec)
600  {
601  m_eivec.resize(A.rows(), nev);
602  for (int i=0; i<nev; i++)
603  for (int j=0; j<n; j++)
604  m_eivec(j,i) = v[i*n+j] / scale;
605 
606  if (mode == 1 && !isBempty && BisSPD)
607  internal::OP<MatrixSolver, MatrixType, Scalar, BisSPD>::project(OP, n, nev, m_eivec.data());
608 
609  m_eigenvectorsOk = true;
610  }
611 
612  m_nbrIterations = iparam[2];
613  m_nbrConverged = iparam[4];
614 
615  m_info = Success;
616  }
617 
618  delete[] select;
619  }
620 
621  delete[] v;
622  delete[] iparam;
623  delete[] ipntr;
624  delete[] workd;
625  delete[] workl;
626  delete[] resid;
627 
628  m_isInitialized = true;
629 
630  return *this;
631 }
Array< int, Dynamic, 1 > v
int n
SparseMatrix< double > A(n, n)
int i
MatrixXf B
#define eigen_assert(x)
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
Real scalar type for MatrixType.
ComputationInfo info() const
Reports whether previous computation was successful.
MatrixType::Scalar Scalar
Scalar type for matrices of type MatrixType.
constexpr void resize(Index rows, Index cols)
NumericalIssue
ComputeEigenvectors
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
CleanedUpDerType< DerType >::type() max(const AutoDiffScalar< DerType > &x, const T &y)
std::ptrdiff_t j

◆ compute() [2/2]

template<typename MatrixType , typename MatrixSolver , bool BisSPD>
ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD > & Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute ( const MatrixType A,
Index  nbrEigenvalues,
std::string  eigs_sigma = "LM",
int  options = ComputeEigenvectors,
RealScalar  tol = 0.0 
)

Computes eigenvalues / eigenvectors of given matrix using the external ARPACK library.

Parameters
[in]ASelfadjoint matrix whose eigendecomposition is to be computed.
[in]nbrEigenvaluesThe number of eigenvalues / eigenvectors to compute. Must be less than the size of the input matrix, or an error is returned.
[in]eigs_sigmaString containing either "LM", "SM", "LA", or "SA", with respective meanings to find the largest magnitude , smallest magnitude, largest algebraic, or smallest algebraic eigenvalues. Alternatively, this value can contain floating point value in string form, in which case the eigenvalues closest to this value will be found.
[in]optionsCan be ComputeEigenvectors (default) or EigenvaluesOnly.
[in]tolWhat tolerance to find the eigenvalues to. Default is 0, which means machine precision.
Returns
Reference to *this

This function computes the eigenvalues of A using ARPACK. The eigenvalues() function can be used to retrieve them. If options equals ComputeEigenvectors, then the eigenvectors are also computed and can be retrieved by calling eigenvectors().

Definition at line 322 of file ArpackSelfAdjointEigenSolver.h.

324 {
325  MatrixType B(0,0);
326  compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
327 
328  return *this;
329 }

◆ eigenvalues()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
const Matrix<Scalar, Dynamic, 1>& Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::eigenvalues ( ) const
inline

Returns the eigenvalues of given matrix.

Returns
A const reference to the column vector containing the eigenvalues.
Precondition
The eigenvalues have been computed before.

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

Example:

MatrixXd ones = MatrixXd::Ones(3,3);
SelfAdjointEigenSolver<MatrixXd> es(ones);
cout << "The eigenvalues of the 3x3 matrix of ones are:"
<< endl << es.eigenvalues() << endl;
MatrixXcf ones
EigenSolver< MatrixXf > es
Matrix< double, Dynamic, Dynamic > MatrixXd

Output:

The eigenvalues of the 3x3 matrix of ones are:
-3.09e-16
        0
        3
See also
eigenvectors(), MatrixBase::eigenvalues()

Definition at line 232 of file ArpackSelfAdjointEigenSolver.h.

233  {
234  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
235  return m_eivalues;
236  }

◆ eigenvectors()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
const Matrix<Scalar, Dynamic, Dynamic>& Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::eigenvectors ( ) const
inline

Returns the eigenvectors of given matrix.

Returns
A const reference to the matrix whose columns are the eigenvectors.
Precondition
The eigenvectors have been computed before.

Column \( k \) of the returned matrix is an eigenvector corresponding to eigenvalue number \( k \) as returned by eigenvalues(). The eigenvectors are normalized to have (Euclidean) norm equal to one. If this object was used to solve the eigenproblem for the selfadjoint matrix \( A \), then the matrix returned by this function is the matrix \( V \) in the eigendecomposition \( A V = D V \). For the generalized eigenproblem, the matrix returned is the solution \( A V = D B V \)

Example:

MatrixXd ones = MatrixXd::Ones(3,3);
SelfAdjointEigenSolver<MatrixXd> es(ones);
cout << "The first eigenvector of the 3x3 matrix of ones is:"
<< endl << es.eigenvectors().col(0) << endl;

Output:

The first eigenvector of the 3x3 matrix of ones is:
-0.816
 0.408
 0.408
See also
eigenvalues()

Definition at line 210 of file ArpackSelfAdjointEigenSolver.h.

211  {
212  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
213  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
214  return m_eivec;
215  }

◆ getNbrConvergedEigenValues()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
size_t Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::getNbrConvergedEigenValues ( ) const
inline

Definition at line 298 of file ArpackSelfAdjointEigenSolver.h.

299  { return m_nbrConverged; }

◆ getNbrIterations()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
size_t Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::getNbrIterations ( ) const
inline

Definition at line 301 of file ArpackSelfAdjointEigenSolver.h.

302  { return m_nbrIterations; }

◆ info()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
ComputationInfo Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NoConvergence otherwise.

Definition at line 292 of file ArpackSelfAdjointEigenSolver.h.

293  {
294  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
295  return m_info;
296  }

◆ operatorInverseSqrt()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Matrix<Scalar, Dynamic, Dynamic> Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::operatorInverseSqrt ( ) const
inline

Computes the inverse square root of the matrix.

Returns
the inverse positive-definite square root of the matrix
Precondition
The eigenvalues and eigenvectors of a positive-definite matrix have been computed before.

This function uses the eigendecomposition \( A = V D V^{-1} \) to compute the inverse square root as \( V D^{-1/2} V^{-1} \). This is cheaper than first computing the square root with operatorSqrt() and then its inverse with MatrixBase::inverse().

Example:

MatrixXd X = MatrixXd::Random(4,4);
MatrixXd A = X * X.transpose();
cout << "Here is a random positive-definite matrix, A:" << endl << A << endl << endl;
SelfAdjointEigenSolver<MatrixXd> es(A);
cout << "The inverse square root of A is: " << endl;
cout << es.operatorInverseSqrt() << endl;
cout << "We can also compute it with operatorSqrt() and inverse(). That yields: " << endl;
cout << es.operatorSqrt().inverse() << endl;
MatrixXf X

Output:

Here is a random positive-definite matrix, A:
  1.41 -0.697 -0.111  0.508
-0.697  0.423 0.0991   -0.4
-0.111 0.0991   1.25  0.902
 0.508   -0.4  0.902    1.4

The inverse square root of A is: 
  1.88   2.78 -0.546  0.605
  2.78   8.61   -2.3   2.74
-0.546   -2.3   1.92  -1.36
 0.605   2.74  -1.36   2.18
We can also compute it with operatorSqrt() and inverse(). That yields: 
  1.88   2.78 -0.546  0.605
  2.78   8.61   -2.3   2.74
-0.546   -2.3   1.92  -1.36
 0.605   2.74  -1.36   2.18
See also
operatorSqrt(), MatrixBase::inverse(), MatrixFunctions Module

Definition at line 281 of file ArpackSelfAdjointEigenSolver.h.

282  {
283  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
284  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
285  return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
286  }

◆ operatorSqrt()

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Matrix<Scalar, Dynamic, Dynamic> Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::operatorSqrt ( ) const
inline

Computes the positive-definite square root of the matrix.

Returns
the positive-definite square root of the matrix
Precondition
The eigenvalues and eigenvectors of a positive-definite matrix have been computed before.

The square root of a positive-definite matrix \( A \) is the positive-definite matrix whose square equals \( A \). This function uses the eigendecomposition \( A = V D V^{-1} \) to compute the square root as \( A^{1/2} = V D^{1/2} V^{-1} \).

Example:

MatrixXd X = MatrixXd::Random(4,4);
MatrixXd A = X * X.transpose();
cout << "Here is a random positive-definite matrix, A:" << endl << A << endl << endl;
SelfAdjointEigenSolver<MatrixXd> es(A);
MatrixXd sqrtA = es.operatorSqrt();
cout << "The square root of A is: " << endl << sqrtA << endl;
cout << "If we square this, we get: " << endl << sqrtA*sqrtA << endl;
MatrixXd sqrtA

Output:

Here is a random positive-definite matrix, A:
  1.41 -0.697 -0.111  0.508
-0.697  0.423 0.0991   -0.4
-0.111 0.0991   1.25  0.902
 0.508   -0.4  0.902    1.4

The square root of A is: 
   1.09  -0.432 -0.0685     0.2
 -0.432   0.379   0.141  -0.269
-0.0685   0.141       1   0.468
    0.2  -0.269   0.468    1.04
If we square this, we get: 
  1.41 -0.697 -0.111  0.508
-0.697  0.423 0.0991   -0.4
-0.111 0.0991   1.25  0.902
 0.508   -0.4  0.902    1.4
See also
operatorInverseSqrt(), MatrixFunctions Module

Definition at line 256 of file ArpackSelfAdjointEigenSolver.h.

257  {
258  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
259  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
260  return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
261  }

Member Data Documentation

◆ m_eigenvectorsOk

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
bool Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eigenvectorsOk
protected

Definition at line 309 of file ArpackSelfAdjointEigenSolver.h.

◆ m_eivalues

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Matrix<Scalar, Dynamic, 1> Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivalues
protected

Definition at line 306 of file ArpackSelfAdjointEigenSolver.h.

◆ m_eivec

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
Matrix<Scalar, Dynamic, Dynamic> Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivec
protected

Definition at line 305 of file ArpackSelfAdjointEigenSolver.h.

◆ m_info

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
ComputationInfo Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_info
protected

Definition at line 307 of file ArpackSelfAdjointEigenSolver.h.

◆ m_isInitialized

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
bool Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized
protected

Definition at line 308 of file ArpackSelfAdjointEigenSolver.h.

◆ m_nbrConverged

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
size_t Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_nbrConverged
protected

Definition at line 311 of file ArpackSelfAdjointEigenSolver.h.

◆ m_nbrIterations

template<typename MatrixType , typename MatrixSolver = SimplicialLLT<MatrixType>, bool BisSPD = false>
size_t Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_nbrIterations
protected

Definition at line 312 of file ArpackSelfAdjointEigenSolver.h.


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