QUDA
v1.1.0
A library for QCD on GPUs
|
Thick Restarted Lanczos Method. More...
#include <eigensolve_quda.h>
Public Member Functions | |
TRLM (const DiracMatrix &mat, QudaEigParam *eig_param, TimeProfile &profile) | |
Constructor for Thick Restarted Eigensolver class. More... | |
virtual | ~TRLM () |
Destructor for Thick Restarted Eigensolver class. More... | |
virtual bool | hermitian () |
void | operator() (std::vector< ColorSpinorField * > &kSpace, std::vector< Complex > &evals) |
Compute eigenpairs. More... | |
void | lanczosStep (std::vector< ColorSpinorField * > v, int j) |
Lanczos step: extends the Krylov space. More... | |
void | reorder (std::vector< ColorSpinorField * > &kSpace) |
Reorder the Krylov space by eigenvalue. More... | |
void | eigensolveFromArrowMat () |
Get the eigendecomposition from the arrow matrix. More... | |
void | computeKeptRitz (std::vector< ColorSpinorField * > &kSpace) |
Rotate the Ritz vectors usinng the arrow matrix eigendecomposition. More... | |
Public Member Functions inherited from quda::EigenSolver | |
EigenSolver (const DiracMatrix &mat, QudaEigParam *eig_param, TimeProfile &profile) | |
Constructor for base Eigensolver class. More... | |
virtual | ~EigenSolver () |
void | prepareInitialGuess (std::vector< ColorSpinorField * > &kSpace) |
Check for an initial guess. If none present, populate with rands, then orthonormalise. More... | |
void | checkChebyOpMax (const DiracMatrix &mat, std::vector< ColorSpinorField * > &kSpace) |
Check for a maximum of the Chebyshev operator. More... | |
void | prepareKrylovSpace (std::vector< ColorSpinorField * > &kSpace, std::vector< Complex > &evals) |
Extend the Krylov space. More... | |
double | setEpsilon (const QudaPrecision prec) |
Set the epsilon parameter. More... | |
void | queryPrec (const QudaPrecision prec) |
Query the eigensolver precision to stdout. More... | |
void | printEigensolverSetup () |
Dump the eigensolver parameters to stdout. More... | |
void | cleanUpEigensolver (std::vector< ColorSpinorField * > &kSpace, std::vector< Complex > &evals) |
Release memory, save eigenvectors, resize the Krylov space to its original dimension. More... | |
void | matVec (const DiracMatrix &mat, ColorSpinorField &out, const ColorSpinorField &in) |
Applies the specified matVec operation: M, Mdag, MMdag, MdagM. More... | |
void | chebyOp (const DiracMatrix &mat, ColorSpinorField &out, const ColorSpinorField &in) |
Promoted the specified matVec operation: M, Mdag, MMdag, MdagM to a Chebyshev polynomial. More... | |
double | estimateChebyOpMax (const DiracMatrix &mat, ColorSpinorField &out, ColorSpinorField &in) |
Estimate the spectral radius of the operator for the max value of the Chebyshev polynomial. More... | |
void | blockOrthogonalize (std::vector< ColorSpinorField * > v, std::vector< ColorSpinorField * > &r, int j) |
Orthogonalise input vectors r against vector space v using block-BLAS. More... | |
void | orthonormalizeMGS (std::vector< ColorSpinorField * > &v, int j) |
Orthonormalise input vector space v using Modified Gram-Schmidt. More... | |
bool | orthoCheck (std::vector< ColorSpinorField * > v, int j) |
Check orthonormality of input vector space v. More... | |
void | rotateVecs (std::vector< ColorSpinorField * > &kSpace, const double *rot_array, const int offset, const int dim, const int keep, const int locked, TimeProfile &profile) |
Rotate the Krylov space. More... | |
void | rotateVecsComplex (std::vector< ColorSpinorField * > &kSpace, const Complex *rot_array, const int offset, const int dim, const int keep, const int locked, TimeProfile &profile) |
Rotate the Krylov space. More... | |
void | permuteVecs (std::vector< ColorSpinorField * > &kSpace, int *mat, int size) |
Permute the vector space using the permutation matrix. More... | |
void | blockRotate (std::vector< ColorSpinorField * > &kSpace, double *array, int rank, const range &i, const range &j, blockType b_type) |
Rotate part of kSpace. More... | |
void | blockRotateComplex (std::vector< ColorSpinorField * > &kSpace, Complex *array, int rank, const range &i, const range &j, blockType b_type, int offset) |
Rotate part of kSpace. More... | |
void | blockReset (std::vector< ColorSpinorField * > &kSpace, int js, int je, int offset) |
Copy temp part of kSpace, zero out for next use. More... | |
void | deflate (std::vector< ColorSpinorField * > &sol, const std::vector< ColorSpinorField * > &src, const std::vector< ColorSpinorField * > &evecs, const std::vector< Complex > &evals, bool accumulate=false) const |
Deflate a set of source vectors with a given eigenspace. More... | |
void | deflate (ColorSpinorField &sol, const ColorSpinorField &src, const std::vector< ColorSpinorField * > &evecs, const std::vector< Complex > &evals, bool accumulate=false) |
Deflate a given source vector with a given eigenspace. This is a wrapper variant for a single source vector. More... | |
void | deflateSVD (std::vector< ColorSpinorField * > &sol, const std::vector< ColorSpinorField * > &vec, const std::vector< ColorSpinorField * > &evecs, const std::vector< Complex > &evals, bool accumulate=false) const |
Deflate a set of source vectors with a set of left and right singular vectors. More... | |
void | deflateSVD (ColorSpinorField &sol, const ColorSpinorField &src, const std::vector< ColorSpinorField * > &evecs, const std::vector< Complex > &evals, bool accumulate=false) |
Deflate a a given source vector with a given with a set of left and right singular vectors This is a wrapper variant for a single source vector. More... | |
void | computeSVD (const DiracMatrix &mat, std::vector< ColorSpinorField * > &evecs, std::vector< Complex > &evals) |
Computes Left/Right SVD from pre computed Right/Left. More... | |
void | computeEvals (const DiracMatrix &mat, std::vector< ColorSpinorField * > &evecs, std::vector< Complex > &evals, int size) |
Compute eigenvalues and their residiua. More... | |
void | computeEvals (const DiracMatrix &mat, std::vector< ColorSpinorField * > &evecs, std::vector< Complex > &evals) |
Compute eigenvalues and their residiua. This variant compute the number of converged eigenvalues. More... | |
void | loadFromFile (const DiracMatrix &mat, std::vector< ColorSpinorField * > &eig_vecs, std::vector< Complex > &evals) |
Load and check eigenpairs from file. More... | |
void | sortArrays (QudaEigSpectrumType spec_type, int n, std::vector< Complex > &x, std::vector< Complex > &y) |
Sort array the first n elements of x according to spec_type, y comes along for the ride. More... | |
void | sortArrays (QudaEigSpectrumType spec_type, int n, std::vector< double > &x, std::vector< Complex > &y) |
Sort array the first n elements of x according to spec_type, y comes along for the ride Overloaded version with real x. More... | |
void | sortArrays (QudaEigSpectrumType spec_type, int n, std::vector< Complex > &x, std::vector< double > &y) |
Sort array the first n elements of x according to spec_type, y comes along for the ride Overloaded version with real y. More... | |
void | sortArrays (QudaEigSpectrumType spec_type, int n, std::vector< double > &x, std::vector< double > &y) |
Sort array the first n elements of x according to spec_type, y comes along for the ride Overloaded version with real x and real y. More... | |
Public Attributes | |
std::vector< double > | ritz_mat |
double * | alpha |
double * | beta |
Additional Inherited Members | |
Static Public Member Functions inherited from quda::EigenSolver | |
static EigenSolver * | create (QudaEigParam *eig_param, const DiracMatrix &mat, TimeProfile &profile) |
Creates the eigensolver using the parameters given and the matrix. More... | |
Protected Attributes inherited from quda::EigenSolver | |
const DiracMatrix & | mat |
QudaEigParam * | eig_param |
TimeProfile & | profile |
int | n_ev |
int | n_kr |
int | n_conv |
int | n_ev_deflate |
double | tol |
bool | reverse |
char | spectrum [3] |
bool | converged |
int | restart_iter |
int | max_restarts |
int | check_interval |
int | batched_rotate |
int | block_size |
int | iter |
int | iter_converged |
int | iter_locked |
int | iter_keep |
int | num_converged |
int | num_locked |
int | num_keep |
std::vector< double > | residua |
std::vector< ColorSpinorField * > | r |
std::vector< ColorSpinorField * > | d_vecs_tmp |
ColorSpinorField * | tmp1 |
ColorSpinorField * | tmp2 |
QudaPrecision | save_prec |
Thick Restarted Lanczos Method.
Definition at line 421 of file eigensolve_quda.h.
quda::TRLM::TRLM | ( | const DiracMatrix & | mat, |
QudaEigParam * | eig_param, | ||
TimeProfile & | profile | ||
) |
Constructor for Thick Restarted Eigensolver class.
eig_param | The eigensolver parameters |
mat | The operator to solve |
profile | Time Profile |
Definition at line 19 of file eig_trlm.cpp.
|
virtual |
Destructor for Thick Restarted Eigensolver class.
Definition at line 195 of file eig_trlm.cpp.
void quda::TRLM::computeKeptRitz | ( | std::vector< ColorSpinorField * > & | kSpace | ) |
Rotate the Ritz vectors usinng the arrow matrix eigendecomposition.
[in] | nKspace | current Krylov space |
Definition at line 339 of file eig_trlm.cpp.
void quda::TRLM::eigensolveFromArrowMat | ( | ) |
Get the eigendecomposition from the arrow matrix.
Definition at line 276 of file eig_trlm.cpp.
|
inlinevirtual |
Implements quda::EigenSolver.
Reimplemented in quda::BLKTRLM.
Definition at line 441 of file eigensolve_quda.h.
void quda::TRLM::lanczosStep | ( | std::vector< ColorSpinorField * > | v, |
int | j | ||
) |
Lanczos step: extends the Krylov space.
[in] | v | Vector space |
[in] | j | Index of vector being computed |
Definition at line 205 of file eig_trlm.cpp.
|
virtual |
Compute eigenpairs.
[in] | kSpace | Krylov vector space |
[in] | evals | Computed eigenvalues |
Implements quda::EigenSolver.
Reimplemented in quda::BLKTRLM.
Definition at line 43 of file eig_trlm.cpp.
void quda::TRLM::reorder | ( | std::vector< ColorSpinorField * > & | kSpace | ) |
Reorder the Krylov space by eigenvalue.
[in] | kSpace | the Krylov space |
Definition at line 249 of file eig_trlm.cpp.
double* quda::TRLM::alpha |
Definition at line 447 of file eigensolve_quda.h.
double* quda::TRLM::beta |
Definition at line 448 of file eigensolve_quda.h.
std::vector<double> quda::TRLM::ritz_mat |
TRLM is only for Hermitian systems
Definition at line 444 of file eigensolve_quda.h.