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

#include <ArpackSelfAdjointEigenSolver.h>

Public Types

typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type MatrixType. More...
 
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...
 

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, Dynamic > & eigenvectors () const
 Returns the eigenvectors of given matrix. More...
 
const Matrix< Scalar, Dynamic, 1 > & eigenvalues () const
 Returns the eigenvalues of given matrix. More...
 
Matrix< Scalar, Dynamic, DynamicoperatorSqrt () const
 Computes the positive-definite square root of the matrix. More...
 
Matrix< Scalar, Dynamic, DynamicoperatorInverseSqrt () const
 Computes the inverse square root of the matrix. More...
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
size_t getNbrConvergedEigenValues () const
 
size_t getNbrIterations () const
 

Protected Attributes

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

Member Typedef Documentation

◆ Index

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

◆ 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.

◆ 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.

◆ 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.

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().

58  : m_eivec(),
59  m_eivalues(),
60  m_isInitialized(false),
61  m_eigenvectorsOk(false),
62  m_nbrConverged(0),
63  m_nbrIterations(0) {}
size_t m_nbrIterations
Definition: ArpackSelfAdjointEigenSolver.h:299
bool m_eigenvectorsOk
Definition: ArpackSelfAdjointEigenSolver.h:296
Matrix< Scalar, Dynamic, Dynamic > m_eivec
Definition: ArpackSelfAdjointEigenSolver.h:292
bool m_isInitialized
Definition: ArpackSelfAdjointEigenSolver.h:295
Matrix< Scalar, Dynamic, 1 > m_eivalues
Definition: ArpackSelfAdjointEigenSolver.h:293
size_t m_nbrConverged
Definition: ArpackSelfAdjointEigenSolver.h:298

◆ 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.

90  : m_eivec(),
91  m_eivalues(),
92  m_isInitialized(false),
93  m_eigenvectorsOk(false),
94  m_nbrConverged(0),
95  m_nbrIterations(0) {
96  compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
97  }
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
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.
Definition: ArpackSelfAdjointEigenSolver.h:316
Definition: matrices.h:74

References Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute().

◆ 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.

123  : m_eivec(),
124  m_eivalues(),
125  m_isInitialized(false),
126  m_eigenvectorsOk(false),
127  m_nbrConverged(0),
128  m_nbrIterations(0) {
129  compute(A, nbrEigenvalues, eigs_sigma, options, tol);
130  }

References Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::compute().

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().

320  {
321  eigen_assert(A.cols() == A.rows());
322  eigen_assert(B.cols() == B.rows());
323  eigen_assert(B.rows() == 0 || A.cols() == B.rows());
324  eigen_assert((options & ~(EigVecMask | GenEigMask)) == 0 && (options & EigVecMask) != EigVecMask &&
325  "invalid option parameter");
326 
327  bool isBempty = (B.rows() == 0) || (B.cols() == 0);
328 
329  // For clarity, all parameters match their ARPACK name
330  //
331  // Always 0 on the first call
332  //
333  int ido = 0;
334 
335  int n = (int)A.cols();
336 
337  // User options: "LA", "SA", "SM", "LM", "BE"
338  //
339  char whch[3] = "LM";
340 
341  // Specifies the shift if iparam[6] = { 3, 4, 5 }, not used if iparam[6] = { 1, 2 }
342  //
343  RealScalar sigma = 0.0;
344 
345  if (eigs_sigma.length() >= 2 && isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) {
346  eigs_sigma[0] = toupper(eigs_sigma[0]);
347  eigs_sigma[1] = toupper(eigs_sigma[1]);
348 
349  // In the following special case we're going to invert the problem, since solving
350  // for larger magnitude is much much faster
351  // i.e., if 'SM' is specified, we're going to really use 'LM', the default
352  //
353  if (eigs_sigma.substr(0, 2) != "SM") {
354  whch[0] = eigs_sigma[0];
355  whch[1] = eigs_sigma[1];
356  }
357  } else {
358  eigen_assert(false && "Specifying clustered eigenvalues is not yet supported!");
359 
360  // If it's not scalar values, then the user may be explicitly
361  // specifying the sigma value to cluster the evs around
362  //
363  sigma = atof(eigs_sigma.c_str());
364 
365  // If atof fails, it returns 0.0, which is a fine default
366  //
367  }
368 
369  // "I" means normal eigenvalue problem, "G" means generalized
370  //
371  char bmat[2] = "I";
372  if (eigs_sigma.substr(0, 2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1])) || (!isBempty && !BisSPD))
373  bmat[0] = 'G';
374 
375  // Now we determine the mode to use
376  //
377  int mode = (bmat[0] == 'G') + 1;
378  if (eigs_sigma.substr(0, 2) == "SM" || !(isalpha(eigs_sigma[0]) && isalpha(eigs_sigma[1]))) {
379  // We're going to use shift-and-invert mode, and basically find
380  // the largest eigenvalues of the inverse operator
381  //
382  mode = 3;
383  }
384 
385  // The user-specified number of eigenvalues/vectors to compute
386  //
387  int nev = (int)nbrEigenvalues;
388 
389  // Allocate space for ARPACK to store the residual
390  //
391  Scalar *resid = new Scalar[n];
392 
393  // Number of Lanczos vectors, must satisfy nev < ncv <= n
394  // Note that this indicates that nev != n, and we cannot compute
395  // all eigenvalues of a mtrix
396  //
397  int ncv = std::min(std::max(2 * nev, 20), n);
398 
399  // The working n x ncv matrix, also store the final eigenvectors (if computed)
400  //
401  Scalar *v = new Scalar[n * ncv];
402  int ldv = n;
403 
404  // Working space
405  //
406  Scalar *workd = new Scalar[3 * n];
407  int lworkl = ncv * ncv + 8 * ncv; // Must be at least this length
408  Scalar *workl = new Scalar[lworkl];
409 
410  int *iparam = new int[11];
411  iparam[0] = 1; // 1 means we let ARPACK perform the shifts, 0 means we'd have to do it
412  iparam[2] = std::max(300, (int)std::ceil(2 * n / std::max(ncv, 1)));
413  iparam[6] = mode; // The mode, 1 is standard ev problem, 2 for generalized ev, 3 for shift-and-invert
414 
415  // Used during reverse communicate to notify where arrays start
416  //
417  int *ipntr = new int[11];
418 
419  // Error codes are returned in here, initial value of 0 indicates a random initial
420  // residual vector is used, any other values means resid contains the initial residual
421  // vector, possibly from a previous run
422  //
423  int info = 0;
424 
425  Scalar scale = 1.0;
426  // if (!isBempty)
427  //{
428  // Scalar scale = B.norm() / std::sqrt(n);
429  // scale = std::pow(2, std::floor(std::log(scale+1)));
431  // for (size_t i=0; i<(size_t)B.outerSize(); i++)
432  // for (typename MatrixType::InnerIterator it(B, i); it; ++it)
433  // it.valueRef() /= scale;
434  // }
435 
436  MatrixSolver OP;
437  if (mode == 1 || mode == 2) {
438  if (!isBempty) OP.compute(B);
439  } else if (mode == 3) {
440  if (sigma == 0.0) {
441  OP.compute(A);
442  } else {
443  // Note: We will never enter here because sigma must be 0.0
444  //
445  if (isBempty) {
446  MatrixType AminusSigmaB(A);
447  for (Index i = 0; i < A.rows(); ++i) AminusSigmaB.coeffRef(i, i) -= sigma;
448 
449  OP.compute(AminusSigmaB);
450  } else {
451  MatrixType AminusSigmaB = A - sigma * B;
452  OP.compute(AminusSigmaB);
453  }
454  }
455  }
456 
457  if (!(mode == 1 && isBempty) && !(mode == 2 && isBempty) && OP.info() != Success)
458  std::cout << "Error factoring matrix" << std::endl;
459 
460  do {
461  internal::arpack_wrapper<Scalar, RealScalar>::saupd(&ido, bmat, &n, whch, &nev, &tol, resid, &ncv, v, &ldv, iparam,
462  ipntr, workd, workl, &lworkl, &info);
463 
464  if (ido == -1 || ido == 1) {
465  Scalar *in = workd + ipntr[0] - 1;
466  Scalar *out = workd + ipntr[1] - 1;
467 
468  if (ido == 1 && mode != 2) {
469  Scalar *out2 = workd + ipntr[2] - 1;
470  if (isBempty || mode == 1)
472  else
474 
475  in = workd + ipntr[2] - 1;
476  }
477 
478  if (mode == 1) {
479  if (isBempty) {
480  // OP = A
481  //
483  } else {
484  // OP = L^{-1}AL^{-T}
485  //
487  }
488  } else if (mode == 2) {
490 
491  // OP = B^{-1} A
492  //
494  } else if (mode == 3) {
495  // OP = (A-\sigmaB)B (\sigma could be 0, and B could be I)
496  // The B * in is already computed and stored at in if ido == 1
497  //
498  if (ido == 1 || isBempty)
500  else
502  }
503  } else if (ido == 2) {
504  Scalar *in = workd + ipntr[0] - 1;
505  Scalar *out = workd + ipntr[1] - 1;
506 
507  if (isBempty || mode == 1)
509  else
511  }
512  } while (ido != 99);
513 
514  if (info == 1)
516  else if (info == 3)
518  else if (info < 0)
520  else if (info != 0)
521  eigen_assert(false && "Unknown ARPACK return value!");
522  else {
523  // Do we compute eigenvectors or not?
524  //
525  int rvec = (options & ComputeEigenvectors) == ComputeEigenvectors;
526 
527  // "A" means "All", use "S" to choose specific eigenvalues (not yet supported in ARPACK))
528  //
529  char howmny[2] = "A";
530 
531  // if howmny == "S", specifies the eigenvalues to compute (not implemented in ARPACK)
532  //
533  int *select = new int[ncv];
534 
535  // Final eigenvalues
536  //
537  m_eivalues.resize(nev, 1);
538 
539  internal::arpack_wrapper<Scalar, RealScalar>::seupd(&rvec, howmny, select, m_eivalues.data(), v, &ldv, &sigma, bmat,
540  &n, whch, &nev, &tol, resid, &ncv, v, &ldv, iparam, ipntr,
541  workd, workl, &lworkl, &info);
542 
543  if (info == -14)
545  else if (info != 0)
547  else {
548  if (rvec) {
549  m_eivec.resize(A.rows(), nev);
550  for (int i = 0; i < nev; i++)
551  for (int j = 0; j < n; j++) m_eivec(j, i) = v[i * n + j] / scale;
552 
553  if (mode == 1 && !isBempty && BisSPD)
555 
556  m_eigenvectorsOk = true;
557  }
558 
559  m_nbrIterations = iparam[2];
560  m_nbrConverged = iparam[4];
561 
562  m_info = Success;
563  }
564 
565  delete[] select;
566  }
567 
568  delete[] v;
569  delete[] iparam;
570  delete[] ipntr;
571  delete[] workd;
572  delete[] workl;
573  delete[] resid;
574 
575  m_isInitialized = true;
576 
577  return *this;
578 }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define eigen_assert(x)
Definition: Macros.h:910
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:48
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
MatrixType::Index Index
Definition: ArpackSelfAdjointEigenSolver.h:34
ComputationInfo m_info
Definition: ArpackSelfAdjointEigenSolver.h:294
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: ArpackSelfAdjointEigenSolver.h:282
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:595
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
#define OP(X)
Definition: common.h:54
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
@ NumericalIssue
Definition: Constants.h:442
@ InvalidInput
Definition: Constants.h:447
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
@ GenEigMask
Definition: Constants.h:414
@ EigVecMask
Definition: Constants.h:403
@ ComputeEigenvectors
Definition: Constants.h:401
return int(ret)+1
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 ceil(const bfloat16 &a)
Definition: BFloat16.h:644
int sigma
Definition: calibrate.py:179
static void project(MatrixSolver &OP, int n, int k, Scalar *vecs)
static void applyOP(MatrixSolver &OP, const MatrixType &A, int n, Scalar *in, Scalar *out)
static void saupd(int *ido, char *bmat, int *n, char *which, int *nev, RealScalar *tol, Scalar *resid, int *ncv, Scalar *v, int *ldv, int *iparam, int *ipntr, Scalar *workd, Scalar *workl, int *lworkl, int *info)
Definition: ArpackSelfAdjointEigenSolver.h:604
static void seupd(int *rvec, char *All, int *select, Scalar *d, Scalar *z, int *ldz, RealScalar *sigma, char *bmat, int *n, char *which, int *nev, RealScalar *tol, Scalar *resid, int *ncv, Scalar *v, int *ldv, int *iparam, int *ipntr, Scalar *workd, Scalar *workl, int *lworkl, int *ierr)
Definition: ArpackSelfAdjointEigenSolver.h:610
std::ofstream out("Result.txt")
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Eigen::internal::OP< MatrixSolver, MatrixType, Scalar, BisSPD >::applyOP(), Eigen::bfloat16_impl::ceil(), Eigen::PlainObjectBase< Derived >::cols(), Eigen::ComputeEigenvectors, eigen_assert, Eigen::EigVecMask, Eigen::GenEigMask, i, info, int(), Eigen::InvalidInput, j, Eigen::PlainObjectBase< Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ > >::Map(), max, min, n, Eigen::NoConvergence, Eigen::NumericalIssue, OP, out(), Eigen::internal::OP< MatrixSolver, MatrixType, Scalar, BisSPD >::project(), Eigen::PlainObjectBase< Derived >::rows(), Eigen::internal::arpack_wrapper< Scalar, RealScalar >::saupd(), Eigen::internal::arpack_wrapper< Scalar, RealScalar >::seupd(), calibrate::sigma, Eigen::Success, and v.

Referenced by Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::ArpackGeneralizedSelfAdjointEigenSolver().

◆ 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().

307  {
308  MatrixType B(0, 0);
309  compute(A, B, nbrEigenvalues, eigs_sigma, options, tol);
310 
311  return *this;
312 }

References compute().

◆ 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
Definition: ComplexEigenSolver_eigenvalues.cpp:1
EigenSolver< MatrixXf > es
Definition: EigenSolver_compute.cpp:1

Output:

See also
eigenvectors(), MatrixBase::eigenvalues()
225  {
226  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
227  return m_eivalues;
228  }

References eigen_assert, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivalues, and Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized.

◆ 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:

See also
eigenvalues()
204  {
205  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
206  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
207  return m_eivec;
208  }

References eigen_assert, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eigenvectorsOk, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivec, and Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized.

◆ getNbrConvergedEigenValues()

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

◆ getNbrIterations()

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

◆ 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.
282  {
283  eigen_assert(m_isInitialized && "ArpackGeneralizedSelfAdjointEigenSolver is not initialized.");
284  return m_info;
285  }

References eigen_assert, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_info, and Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized.

◆ 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;
#define X
Definition: icosphere.cpp:20

Output:

See also
operatorSqrt(), MatrixBase::inverse(), MatrixFunctions Module
272  {
273  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
274  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
275  return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
276  }

References eigen_assert, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eigenvectorsOk, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivalues, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivec, and Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized.

◆ 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
Definition: SelfAdjointEigenSolver_operatorSqrt.cpp:6

Output:

See also
operatorInverseSqrt(), MatrixFunctions Module
248  {
249  eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
250  eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
251  return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
252  }

References eigen_assert, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eigenvectorsOk, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivalues, Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_eivec, and Eigen::ArpackGeneralizedSelfAdjointEigenSolver< MatrixType, MatrixSolver, BisSPD >::m_isInitialized.

Member Data Documentation

◆ m_eigenvectorsOk

◆ m_eivalues

◆ m_eivec

◆ m_info

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

◆ m_isInitialized

◆ m_nbrConverged

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

◆ m_nbrIterations

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

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