CUSP
Loading...
Searching...
No Matches
EigenSolvers

Iterative methods for computing eigenpairs of hermitian and non-hermitian linear systems. More...

Functions

template<typename Matrix , typename Array2d >
void cusp::eigen::arnoldi (const Matrix &A, Array2d &H, size_t k=10)
 Approximate spectral radius of A using Arnoldi.
 
template<typename LinearOperator , typename Array1d , typename Array2d , typename LanczosOptions >
void cusp::eigen::lanczos (const LinearOperator &A, Array1d &eigVals, Array2d &eigVecs, LanczosOptions &options)
 Lanczos method.
 
template<typename LinearOperator , typename Array1d , typename Array2d , typename Monitor , typename Preconditioner >
void cusp::eigen::lobpcg (LinearOperator &A, Array1d &S, Array2d &X, Monitor &monitor, Preconditioner &M, bool largest=true)
 LOBPCG method.
 
template<typename MatrixType >
double cusp::eigen::disks_spectral_radius (const MatrixType &A)
 Uses Gershgorin disks to approximate spectral radius.
 
template<typename MatrixType >
double cusp::eigen::estimate_rho_Dinv_A (const MatrixType &A)
 Approximate spectral radius of (D^-1)A.
 
template<typename MatrixType >
double cusp::eigen::estimate_spectral_radius (const MatrixType &A, size_t k=20)
 Approximate spectral radius of A using Lanczos.
 
template<typename MatrixType >
double cusp::eigen::ritz_spectral_radius (const MatrixType &A, size_t k=10, bool symmetric=false)
 Approximate spectral radius of A using Lanczos or Arnoldi.
 

Detailed Description

Iterative methods for computing eigenpairs of hermitian and non-hermitian linear systems.

Function Documentation

◆ arnoldi()

template<typename Matrix , typename Array2d >
void cusp::eigen::arnoldi ( const Matrix &  A,
Array2d &  H,
size_t  k = 10 
)

Approximate spectral radius of A using Arnoldi.

Template Parameters
Matrixtype of a sparse or dense matrix
Array2dtype of dense matrix of partials
Parameters
Amatrix of the linear system
Hdense matrix of ritz values
kmaximum number of outer Arnoldi iterations
Overview
Approximates the spectral radius A using a specified number of Arnoldi iterations.
Example
The following code snippet demonstrates how to use arnoldi to compute the spectral radius of a 16x16 Laplacian matrix.
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// compute the largest eigenpair of A using 20 Arnoldi iterations
float rho = cusp::eigen::arnoldi(A, 20);
std::cout << "Spectral radius of A : " << rho << std::endl;
return 0;
}
Arnoldi method.
Compressed sparse row (CSR) representation a sparse matrix.
Definition csr_matrix.h:108
Compressed Sparse Row matrix format.
void arnoldi(const Matrix &A, Array2d &H, size_t k=10)
Approximate spectral radius of A using Arnoldi.
Poisson matrix generators.

◆ disks_spectral_radius()

template<typename MatrixType >
double cusp::eigen::disks_spectral_radius ( const MatrixType &  A)

Uses Gershgorin disks to approximate spectral radius.

Template Parameters
MatrixTypetype of a sparse or dense matrix
Parameters
Amatrix of the linear system
Returns
spectral radius approximation
Overview
Approximates the spectral radius of a matrix using Gershgorin disks.
See also
https://en.wikipedia.org/wiki/Gershgorin_circle_theorem
Example
The following code snippet demonstrates how to use disks_spectral_radius to compute the spectral radius of a 16x16 Laplacian matrix.
#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// compute the largest eigenpair of A
std::cout << "Spectral radius of A : " << rho << std::endl;
return 0;
}
double disks_spectral_radius(const MatrixType &A)
Uses Gershgorin disks to approximate spectral radius.
Monitor iterative solver convergence.
Various methods to compute spectral radius.

◆ estimate_rho_Dinv_A()

template<typename MatrixType >
double cusp::eigen::estimate_rho_Dinv_A ( const MatrixType &  A)

Approximate spectral radius of (D^-1)A.

Template Parameters
MatrixTypetype of a sparse or dense matrix
Parameters
Amatrix of the linear system
Returns
spectral radius approximation
Overview
Approximates the spectral radius (D^-1)A, where D is a diagonal matrix containing the diagonal entries of A. The spectral radius of (D^-1)A is computed using either Lanczos or Arnoldi.
Example
The following code snippet demonstrates how to use estimate_rho_Dinv_A to compute the spectral radius of a 16x16 Laplacian matrix.
#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// compute the largest eigenpair of A
std::cout << "Spectral radius of (D^-1)A : " << rho << std::endl;
return 0;
}
double estimate_rho_Dinv_A(const MatrixType &A)
Approximate spectral radius of (D^-1)A.

◆ estimate_spectral_radius()

template<typename MatrixType >
double cusp::eigen::estimate_spectral_radius ( const MatrixType &  A,
size_t  k = 20 
)

Approximate spectral radius of A using Lanczos.

Template Parameters
MatrixTypetype of a sparse or dense matrix
Parameters
Amatrix of the linear system
knumber of Lanczos iterations
Returns
spectral radius approximation
Overview
Approximates the spectral radius A using a specified number of Lanczos iterations.
Example
The following code snippet demonstrates how to use estimate_spectral_radius to compute the spectral radius of a 16x16 Laplacian matrix.
#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// compute the largest eigenpair of A using 20 Lanczos iterations
std::cout << "Spectral radius of A : " << rho << std::endl;
return 0;
}
double estimate_spectral_radius(const MatrixType &A, size_t k=20)
Approximate spectral radius of A using Lanczos.

◆ lanczos()

template<typename LinearOperator , typename Array1d , typename Array2d , typename LanczosOptions >
void cusp::eigen::lanczos ( const LinearOperator &  A,
Array1d &  eigVals,
Array2d &  eigVecs,
LanczosOptions &  options 
)

Lanczos method.

Template Parameters
LinearOperatoris a matrix or subclass of linear_operator
Array1darray of eigenvalues
Array2dmatrix of eigenvectors
Parameters
Amatrix of the linear system
eigValseigenvalues
eigVecseigenvectors
optionsLanczos solver options
Overview
Computes the extreme eigenpairs of hermitian linear systems A x = s x.
Note
A must be symmetric.
Example
The following code snippet demonstrates how to use lanczos to compute the eigenpairs of a 10x10 Laplacian matrix.
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage and initialize eigenpairs
// initialize Lanzcos option
cusp::eigen::lanczos_options<float> options;
options.tol = 1e-6;
options.maxIter = 100;
options.verbose = true;
options.computeEigVecs = false;
options.reorth = cusp::eigen::Full;
// compute the largest eigenpair of A
cusp::eigen::lanczos(A, eigvals, eigvecs, options);
// print largest eigenvalue
std::cout << "Largest eigenvalue : " << eigvals.back() << std::endl;
return 0;
}
The array1d class is a 1D vector container that may contain elements stored in "host" or "device" mem...
Definition array1d.h:99
The array2d class is a 2D vector container that may contain elements stored in "host" or "device" mem...
Definition array2d.h:94
void lanczos(const LinearOperator &A, Array1d &eigVals, Array2d &eigVecs, LanczosOptions &options)
Lanczos method.
Lanczos method.

◆ lobpcg()

template<typename LinearOperator , typename Array1d , typename Array2d , typename Monitor , typename Preconditioner >
void cusp::eigen::lobpcg ( LinearOperator &  A,
Array1d &  S,
Array2d &  X,
Monitor &  monitor,
Preconditioner &  M,
bool  largest = true 
)

LOBPCG method.

Template Parameters
LinearOperatoris a matrix or subtypename of linear_operator
Vectorvector
Monitoris a monitor
Preconditioneris a matrix or subtypename of linear_operator
Parameters
Amatrix of the linear system
Seigenvalues
Xeigenvectors
monitormonitors iteration and determines stopping conditions
Mpreconditioner for A
largestIf true compute the eigenpair corresponding to the largest eigenvalue otherwise compute the smallest.
Overview
Computes the extreme eigenpairs of hermitian linear systems A x = s x using LOBPCG.
Note
A and M must be symmetric.
See also
https://en.wikipedia.org/wiki/LOBPCG
Example
The following code snippet demonstrates how to use lobpcg to solve a 10x10 Poisson problem.
#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// allocate storage and initialize eigenpairs
cusp::random_array<double> randx(A.num_rows);
// set stopping criteria:
// iteration_limit = 100
// relative_tolerance = 1e-6
// absolute_tolerance = 0
// verbose = true
cusp::monitor<double> monitor(X, 100, 1e-6, 0, true);
// set preconditioner (identity)
// Compute the largest eigenpair of A
cusp::eigen::lobpcg(A, S, X, monitor, M, true);
std::cout << "Largest eigenvalue : " << S[0] << std::endl;
return 0;
}
Simple identity operator.
Implements standard convergence criteria and reporting for iterative solvers.
Definition monitor.h:102
Specialized array1d_view wrapping cusp::random_iterator.
Definition array1d.h:664
void lobpcg(LinearOperator &A, Array1d &S, Array2d &X, Monitor &monitor, Preconditioner &M, bool largest=true)
LOBPCG method.
LOBPCG method.
See also
monitor

◆ ritz_spectral_radius()

template<typename MatrixType >
double cusp::eigen::ritz_spectral_radius ( const MatrixType &  A,
size_t  k = 10,
bool  symmetric = false 
)

Approximate spectral radius of A using Lanczos or Arnoldi.

Template Parameters
MatrixTypetype of a sparse or dense matrix
Parameters
Amatrix of the linear system
knumber of iterations
symmetricif true use Lanczos, otherwise use Arnoldi
Returns
spectral radius approximation
Overview
Approximates the spectral radius A using a specified number of Lanczos or Arnoldi iterations.
Example
The following code snippet demonstrates how to use ritz_spectral_radius to compute the spectral radius of a 16x16 Laplacian matrix.
#include <cusp/monitor.h>
int main(void)
{
// create an empty sparse matrix structure (CSR format)
// initialize matrix
// compute the largest eigenpair of A using 20 Arnoldi iterations
float rho = cusp::eigen::ritz_spectral_radius(A, 20, false);
std::cout << "Spectral radius of A : " << rho << std::endl;
return 0;
}
double ritz_spectral_radius(const MatrixType &A, size_t k=10, bool symmetric=false)
Approximate spectral radius of A using Lanczos or Arnoldi.