![]() |
STAG C++
1.3.0
Spectral Toolkit of Algorithms for Graphs
|
Methods for computing eigenvalues and eigenvectors of sparse matrices.
Typedefs | |
typedef std::tuple< Eigen::VectorXd, Eigen::MatrixXd > | stag::EigenSystem |
Functions | |
stag::EigenSystem | stag::compute_eigensystem (const SprsMat *mat, stag_int num, Spectra::SortRule sort) |
stag::EigenSystem | stag::compute_eigensystem (const SprsMat *mat, stag_int num) |
Eigen::VectorXd | stag::compute_eigenvalues (const SprsMat *mat, stag_int num, Spectra::SortRule sort) |
Eigen::VectorXd | stag::compute_eigenvalues (const SprsMat *mat, stag_int num) |
Eigen::MatrixXd | stag::compute_eigenvectors (const SprsMat *mat, stag_int num, Spectra::SortRule sort) |
Eigen::MatrixXd | stag::compute_eigenvectors (const SprsMat *mat, stag_int num) |
Eigen::VectorXd | stag::power_method (const SprsMat *mat, stag_int num_iterations, Eigen::VectorXd initial_vector) |
Eigen::VectorXd | stag::power_method (const SprsMat *mat, stag_int num_iterations) |
Eigen::VectorXd | stag::power_method (const SprsMat *mat, Eigen::VectorXd initial_vector) |
Eigen::VectorXd | stag::power_method (const SprsMat *mat) |
double | stag::rayleigh_quotient (const SprsMat *mat, Eigen::VectorXd &vec) |
typedef std::tuple<Eigen::VectorXd, Eigen::MatrixXd> stag::EigenSystem |
A convenience type for returning eigenvalues and eigenvectors together.
stag::EigenSystem stag::compute_eigensystem | ( | const SprsMat * | mat, |
stag_int | num, | ||
Spectra::SortRule | sort | ||
) |
Compute the eigenvalues and eigenvectors of a given matrix.
By default, this will compute the eigenvalues of smallest magnitude. This default can be overridden by the sort parameter which takes a Spectra::SortRule
object. This is likely to be one of the following.
Spectra::SortRule::SmallestMagn
will return the eigenvalues with smallest magnitudeSpectra::SortRule::LargestMagn
will return the eigenvalues with largest magnitudeThe following example demonstrates how to compute the 3 largest eigenvectors and eigenvalues of a cycle graph.
mat | the matrix on which to operate |
num | the number of eigenvalues and eigenvectors to compute |
sort | (optional) a Spectra::SortRule indicating which eigenvectors to calculate |
std::runtime_error | if the eigenvalue calculation does not converge |
stag::EigenSystem stag::compute_eigensystem | ( | const SprsMat * | mat, |
stag_int | num | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Eigen::VectorXd stag::compute_eigenvalues | ( | const SprsMat * | mat, |
stag_int | num, | ||
Spectra::SortRule | sort | ||
) |
Compute the eigenvalues of a given matrix.
By default, this will compute the eigenvalues of smallest magnitude. This default can be overridden by the sort parameter which takes a Spectra::SortRule
object. This is likely to be one of the following.
Spectra::SortRule::SmallestMagn
will return the eigenvalues with smallest magnitudeSpectra::SortRule::LargestMagn
will return the eigenvalues with largest magnitudeIf you would like to calculate the eigenvectors and eigenvalues together, then you should instead use stag::compute_eigensystem.
mat | the matrix on which to operate |
num | the number of eigenvalues to compute |
sort | (optional) a Spectra::SortRule indicating which eigenvalues to calculate |
std::runtime_error | if the eigenvalue calculation does not converge |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Eigen::MatrixXd stag::compute_eigenvectors | ( | const SprsMat * | mat, |
stag_int | num, | ||
Spectra::SortRule | sort | ||
) |
Compute the eigenvectors of a given matrix.
By default, this will compute the eigenvectors corresponding to the eigenvalues of smallest magnitude. This default can be overridden by the sort parameter which takes a Spectra::SortRule
object. This is likely to be one of the following.
Spectra::SortRule::SmallestMagn
will return the eigenvalues with smallest magnitudeSpectra::SortRule::LargestMagn
will return the eigenvalues with largest magnitudeIf you would like to calculate the eigenvectors and eigenvalues together, then you should instead use stag::compute_eigensystem.
mat | the matrix on which to operate |
num | the number of eigenvectors to compute |
sort | (optional) a Spectra::SortRule indicating which eigenvectors to calculate |
std::runtime_error | if the eigenvector calculation does not converge |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Eigen::VectorXd stag::power_method | ( | const SprsMat * | mat, |
stag_int | num_iterations, | ||
Eigen::VectorXd | initial_vector | ||
) |
Apply the power method to compute the dominant eigenvector of a matrix.
Given a matrix \(M\), an initial vector \(v_0\), and a number of iterations \(t\), the power method calculates the vector
\[ v_t = M^t v_0, \]
which is close to the eigenvector of \(M\) with largest eigenvalue.
The running time of the power method is \(O(t \cdot \mathrm{nnz}(M))\), where \(\mathrm{nnz}(M)\) is the number of non-zero elements in the matrix \(M\).
mat | the matrix \(M\) on which to operate. |
num_iterations | (optional) the number of iterations of the power method to apply. It this argument is omitted, \(O(\log(n))\) iterations are used which results in a vector whose Rayleigh quotient is a \((1 - \epsilon)\) approximation of the dominant eigenvalue. |
initial_vector | (optional) the initial vector to use for the power iteration. If this argument is omitted, a random unit vector will be used. |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Eigen::VectorXd stag::power_method | ( | const SprsMat * | mat, |
Eigen::VectorXd | initial_vector | ||
) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Eigen::VectorXd stag::power_method | ( | const SprsMat * | mat | ) |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
double stag::rayleigh_quotient | ( | const SprsMat * | mat, |
Eigen::VectorXd & | vec | ||
) |
Compute the Rayleigh quotient of the given vector and matrix.
Given a matrix \(M\), the Rayleigh quotient of vector \(v\) is
\[ R(M, v) = \frac{v^\top M v}{v^\top v}. \]
mat | a sparse matrix \(M \in \mathbb{R}^{n \times n}\). |
vec | a vector \(v \in \mathbb{R}^n\). |