QUDA  v1.1.0
A library for QCD on GPUs
Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
quda::EigenSolver Class Referenceabstract

#include <eigensolve_quda.h>

+ Inheritance diagram for quda::EigenSolver:

Public Member Functions

 EigenSolver (const DiracMatrix &mat, QudaEigParam *eig_param, TimeProfile &profile)
 Constructor for base Eigensolver class. More...
 
virtual ~EigenSolver ()
 
virtual bool hermitian ()=0
 
virtual void operator() (std::vector< ColorSpinorField * > &kSpace, std::vector< Complex > &evals)=0
 Computes the eigen decomposition for the operator passed to create. More...
 
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...
 

Static Public Member Functions

static EigenSolvercreate (QudaEigParam *eig_param, const DiracMatrix &mat, TimeProfile &profile)
 Creates the eigensolver using the parameters given and the matrix. More...
 

Protected Attributes

const DiracMatrixmat
 
QudaEigParameig_param
 
TimeProfileprofile
 
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
 
ColorSpinorFieldtmp1
 
ColorSpinorFieldtmp2
 
QudaPrecision save_prec
 

Detailed Description

Definition at line 14 of file eigensolve_quda.h.

Constructor & Destructor Documentation

◆ EigenSolver()

quda::EigenSolver::EigenSolver ( const DiracMatrix mat,
QudaEigParam eig_param,
TimeProfile profile 
)

Constructor for base Eigensolver class.

Parameters
eig_paramMGParam struct that defines all meta data
profileTimeprofile instance used to profile

Definition at line 22 of file eigensolve_quda.cpp.

◆ ~EigenSolver()

quda::EigenSolver::~EigenSolver ( )
virtual

Destructor for EigenSolver class.

Definition at line 1203 of file eigensolve_quda.cpp.

Member Function Documentation

◆ blockOrthogonalize()

void quda::EigenSolver::blockOrthogonalize ( std::vector< ColorSpinorField * >  v,
std::vector< ColorSpinorField * > &  r,
int  j 
)

Orthogonalise input vectors r against vector space v using block-BLAS.

Parameters
[in]vVector space
[in]rVectors to be orthogonalised
[in]jUse vectors v[0:j]
[in]sarray of

Definition at line 469 of file eigensolve_quda.cpp.

◆ blockReset()

void quda::EigenSolver::blockReset ( std::vector< ColorSpinorField * > &  kSpace,
int  js,
int  je,
int  offset 
)

Copy temp part of kSpace, zero out for next use.

Parameters
[in/out]kSpace The current Krylov space
[in]jsStart of j index
[in]jeEnd of j index
[in]offsetPosition of extra vectors in kSpace

Definition at line 571 of file eigensolve_quda.cpp.

◆ blockRotate()

void quda::EigenSolver::blockRotate ( std::vector< ColorSpinorField * > &  kSpace,
double *  array,
int  rank,
const range &  i,
const range &  j,
blockType  b_type 
)

Rotate part of kSpace.

Parameters
[in/out]kSpace The current Krylov space
[in]arrayThe real rotation matrix
[in]rankrow rank of array
[in]isStart of i index
[in]ieEnd of i index
[in]jsStart of j index
[in]jeEnd of j index
[in]blockTypeType of caxpy(_U/L) to perform
[in]jeEnd of j index
[in]offsetPosition of extra vectors in kSpace

Definition at line 527 of file eigensolve_quda.cpp.

◆ blockRotateComplex()

void quda::EigenSolver::blockRotateComplex ( std::vector< ColorSpinorField * > &  kSpace,
Complex array,
int  rank,
const range &  i,
const range &  j,
blockType  b_type,
int  offset 
)

Rotate part of kSpace.

Parameters
[in/out]kSpace The current Krylov space
[in]arrayThe complex rotation matrix
[in]rankrow rank of array
[in]isStart of i index
[in]ieEnd of i index
[in]jsStart of j index
[in]jeEnd of j index
[in]blockTypeType of caxpy(_U/L) to perform
[in]offsetPosition of extra vectors in kSpace

Definition at line 581 of file eigensolve_quda.cpp.

◆ chebyOp()

void quda::EigenSolver::chebyOp ( const DiracMatrix mat,
ColorSpinorField out,
const ColorSpinorField in 
)

Promoted the specified matVec operation: M, Mdag, MMdag, MdagM to a Chebyshev polynomial.

Parameters
[in]matMatrix operator
[in]outOutput spinor
[in]inInput spinor

Definition at line 308 of file eigensolve_quda.cpp.

◆ checkChebyOpMax()

void quda::EigenSolver::checkChebyOpMax ( const DiracMatrix mat,
std::vector< ColorSpinorField * > &  kSpace 
)

Check for a maximum of the Chebyshev operator.

Parameters
[in]matThe problem operator
[in]kSpaceThe Krylov space vectors

Definition at line 170 of file eigensolve_quda.cpp.

◆ cleanUpEigensolver()

void quda::EigenSolver::cleanUpEigensolver ( std::vector< ColorSpinorField * > &  kSpace,
std::vector< Complex > &  evals 
)

Release memory, save eigenvectors, resize the Krylov space to its original dimension.

Parameters
[in]kSpaceThe Krylov space vectors
[in]evalsThe eigenvalue array

Definition at line 239 of file eigensolve_quda.cpp.

◆ computeEvals() [1/2]

void quda::EigenSolver::computeEvals ( const DiracMatrix mat,
std::vector< ColorSpinorField * > &  evecs,
std::vector< Complex > &  evals 
)
inline

Compute eigenvalues and their residiua. This variant compute the number of converged eigenvalues.

Parameters
[in]matMatrix operator
[in]evecsThe eigenvectors
[in]evalsThe eigenvalues

Definition at line 364 of file eigensolve_quda.h.

◆ computeEvals() [2/2]

void quda::EigenSolver::computeEvals ( const DiracMatrix mat,
std::vector< ColorSpinorField * > &  evecs,
std::vector< Complex > &  evals,
int  size 
)

Compute eigenvalues and their residiua.

Parameters
[in]matMatrix operator
[in]evecsThe eigenvectors
[in]evalsThe eigenvalues
[in]sizeThe number of eigenvalues to compute

Definition at line 718 of file eigensolve_quda.cpp.

◆ computeSVD()

void quda::EigenSolver::computeSVD ( const DiracMatrix mat,
std::vector< ColorSpinorField * > &  evecs,
std::vector< Complex > &  evals 
)

Computes Left/Right SVD from pre computed Right/Left.

Parameters
[in]matMatrix operator
[in]evecsComputed eigenvectors of NormOp
[in]evalsComputed eigenvalues of NormOp

Definition at line 625 of file eigensolve_quda.cpp.

◆ create()

EigenSolver * quda::EigenSolver::create ( QudaEigParam eig_param,
const DiracMatrix mat,
TimeProfile profile 
)
static

Creates the eigensolver using the parameters given and the matrix.

Parameters
eig_paramThe eigensolver parameters
matThe operator to solve
profileTime Profile

Definition at line 97 of file eigensolve_quda.cpp.

◆ deflate() [1/2]

void quda::EigenSolver::deflate ( ColorSpinorField sol,
const ColorSpinorField src,
const std::vector< ColorSpinorField * > &  evecs,
const std::vector< Complex > &  evals,
bool  accumulate = false 
)
inline

Deflate a given source vector with a given eigenspace. This is a wrapper variant for a single source vector.

Parameters
[in]solThe resulting deflated vector
[in]srcThe source vector we are deflating
[in]evecsThe eigenvectors to use in deflation
[in]evalsThe eigenvalues to use in deflation
[in]accumulateWhether to preserve the sol vector content prior to accumulating

Definition at line 288 of file eigensolve_quda.h.

◆ deflate() [2/2]

void quda::EigenSolver::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.

Parameters
[in]solThe resulting deflated vector set
[in]srcThe source vector set we are deflating
[in]evecsThe eigenvectors to use in deflation
[in]evalsThe eigenvalues to use in deflation
[in]accumulateWhether to preserve the sol vector content prior to accumulating

Definition at line 752 of file eigensolve_quda.cpp.

◆ deflateSVD() [1/2]

void quda::EigenSolver::deflateSVD ( ColorSpinorField sol,
const ColorSpinorField src,
const std::vector< ColorSpinorField * > &  evecs,
const std::vector< Complex > &  evals,
bool  accumulate = false 
)
inline

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.

Parameters
[in]solThe resulting deflated vector set
[in]srcThe source vector set we are deflating
[in]evecsThe singular vectors to use in deflation
[in]evalsThe singular values to use in deflation
[in]accumulateWhether to preserve the sol vector content prior to accumulating

Definition at line 325 of file eigensolve_quda.h.

◆ deflateSVD() [2/2]

void quda::EigenSolver::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.

Parameters
[in]solThe resulting deflated vector set
[in]srcThe source vector set we are deflating
[in]evecsThe singular vectors to use in deflation
[in]evalsThe singular values to use in deflation
[in]accumulateWhether to preserve the sol vector content prior to accumulating

Definition at line 672 of file eigensolve_quda.cpp.

◆ estimateChebyOpMax()

double quda::EigenSolver::estimateChebyOpMax ( const DiracMatrix mat,
ColorSpinorField out,
ColorSpinorField in 
)

Estimate the spectral radius of the operator for the max value of the Chebyshev polynomial.

Parameters
[in]matMatrix operator
[in]outOutput spinor
[in]inInput spinor

Definition at line 379 of file eigensolve_quda.cpp.

◆ hermitian()

virtual bool quda::EigenSolver::hermitian ( )
pure virtual
Returns
Whether the solver is only for Hermitian systems

Implemented in quda::IRAM, quda::BLKTRLM, and quda::TRLM.

◆ loadFromFile()

void quda::EigenSolver::loadFromFile ( const DiracMatrix mat,
std::vector< ColorSpinorField * > &  eig_vecs,
std::vector< Complex > &  evals 
)

Load and check eigenpairs from file.

Parameters
[in]matMatrix operator
[in]eig_vecsThe eigenvectors to save
[in]fileThe filename to save

Definition at line 792 of file eigensolve_quda.cpp.

◆ matVec()

void quda::EigenSolver::matVec ( const DiracMatrix mat,
ColorSpinorField out,
const ColorSpinorField in 
)

Applies the specified matVec operation: M, Mdag, MMdag, MdagM.

Parameters
[in]matMatrix operator
[in]outOutput spinor
[in]inInput spinor

Definition at line 295 of file eigensolve_quda.cpp.

◆ operator()()

virtual void quda::EigenSolver::operator() ( std::vector< ColorSpinorField * > &  kSpace,
std::vector< Complex > &  evals 
)
pure virtual

Computes the eigen decomposition for the operator passed to create.

Parameters
kSpaceThe converged eigenvectors
evalsThe converged eigenvalues

Implemented in quda::IRAM, quda::BLKTRLM, and quda::TRLM.

◆ orthoCheck()

bool quda::EigenSolver::orthoCheck ( std::vector< ColorSpinorField * >  v,
int  j 
)

Check orthonormality of input vector space v.

Parameters
[out]boolIf all vectors are orthonormal to 1e-16 returns true, else false.
[in]vVector space
[in]jUse vectors v[0:j-1]

Definition at line 416 of file eigensolve_quda.cpp.

◆ orthonormalizeMGS()

void quda::EigenSolver::orthonormalizeMGS ( std::vector< ColorSpinorField * > &  v,
int  j 
)

Orthonormalise input vector space v using Modified Gram-Schmidt.

Parameters
[in]vVector space
[in]jUse vectors v[0:j-1]

Definition at line 452 of file eigensolve_quda.cpp.

◆ permuteVecs()

void quda::EigenSolver::permuteVecs ( std::vector< ColorSpinorField * > &  kSpace,
int *  mat,
int  size 
)

Permute the vector space using the permutation matrix.

Parameters
[in/out]kSpace The current Krylov space
[in]matEigen object storing the pivots
[in]sizeThe size of the (square) permutation matrix

Definition at line 491 of file eigensolve_quda.cpp.

◆ prepareInitialGuess()

void quda::EigenSolver::prepareInitialGuess ( std::vector< ColorSpinorField * > &  kSpace)

Check for an initial guess. If none present, populate with rands, then orthonormalise.

Parameters
[in]kSpaceThe Krylov space vectors

Definition at line 139 of file eigensolve_quda.cpp.

◆ prepareKrylovSpace()

void quda::EigenSolver::prepareKrylovSpace ( std::vector< ColorSpinorField * > &  kSpace,
std::vector< Complex > &  evals 
)

Extend the Krylov space.

Parameters
[in]kSpaceThe Krylov space vectors
[in]evalsThe eigenvalue array

Definition at line 179 of file eigensolve_quda.cpp.

◆ printEigensolverSetup()

void quda::EigenSolver::printEigensolverSetup ( )

Dump the eigensolver parameters to stdout.

Definition at line 192 of file eigensolve_quda.cpp.

◆ queryPrec()

void quda::EigenSolver::queryPrec ( const QudaPrecision  prec)

Query the eigensolver precision to stdout.

Parameters
[in]precPrecision of the solver instance

Definition at line 228 of file eigensolve_quda.cpp.

◆ rotateVecs()

void quda::EigenSolver::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.

Parameters
[in]kSpacethe Krylov space
[in]rot_arrayThe real rotation matrix
[in]offsetThe position of the start of unused vectors in kSpace
[in]dimThe number of rows in the rotation array
[in]keepThe number of columns in the rotation array
[in]lockedThe number of locked vectors in kSpace
[in]profileTime profiler

Definition at line 1062 of file eigensolve_quda.cpp.

◆ rotateVecsComplex()

void quda::EigenSolver::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.

Parameters
[in]kSpacethe Krylov space
[in]rot_arrayThe complex rotation matrix
[in]offsetThe position of the start of unused vector in kSpace
[in]dimThe number of rows in the rotation array
[in]keepThe number of columns in the rotation array
[in]lockedThe number of locked vectors in kSpace
[in]profileTime profiler

Definition at line 918 of file eigensolve_quda.cpp.

◆ setEpsilon()

double quda::EigenSolver::setEpsilon ( const QudaPrecision  prec)

Set the epsilon parameter.

Parameters
[in]precPrecision of the solver instance
[out]epsilonThe deduced epsilon value

Definition at line 215 of file eigensolve_quda.cpp.

◆ sortArrays() [1/4]

void quda::EigenSolver::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.

Parameters
[in]spec_typeThe spectrum type (Largest/Smallest)(Modulus/Imaginary/Real)
[in]nThe number of elements to sort
[in]xThe array to sort
[in]yAn array whose elements will be permuted in tandem with x

Definition at line 821 of file eigensolve_quda.cpp.

◆ sortArrays() [2/4]

void quda::EigenSolver::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.

Parameters
[in]spec_typeThe spectrum type (Largest/Smallest)(Modulus/Imaginary/Real)
[in]nThe number of elements to sort
[in]xThe array to sort
[in]yAn array whose elements will be permuted in tandem with x

Definition at line 882 of file eigensolve_quda.cpp.

◆ sortArrays() [3/4]

void quda::EigenSolver::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.

Parameters
[in]spec_typeThe spectrum type (Largest/Smallest)(Modulus/Imaginary/Real)
[in]nThe number of elements to sort
[in]xThe array to sort
[in]yAn array whose elements will be permuted in tandem with x

Definition at line 892 of file eigensolve_quda.cpp.

◆ sortArrays() [4/4]

void quda::EigenSolver::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.

Parameters
[in]spec_typeThe spectrum type (Largest/Smallest)(Modulus/Imaginary/Real) that determines the sorting condition
[in]nThe number of elements to sort
[in]xThe array to sort
[in]yAn array whose elements will be permuted in tandem with x

Definition at line 902 of file eigensolve_quda.cpp.

Member Data Documentation

◆ batched_rotate

int quda::EigenSolver::batched_rotate
protected

Definition at line 39 of file eigensolve_quda.h.

◆ block_size

int quda::EigenSolver::block_size
protected

Definition at line 40 of file eigensolve_quda.h.

◆ check_interval

int quda::EigenSolver::check_interval
protected

Definition at line 38 of file eigensolve_quda.h.

◆ converged

bool quda::EigenSolver::converged
protected

Part of the spectrum to be computed

Definition at line 35 of file eigensolve_quda.h.

◆ d_vecs_tmp

std::vector<ColorSpinorField *> quda::EigenSolver::d_vecs_tmp
protected

Definition at line 53 of file eigensolve_quda.h.

◆ eig_param

QudaEigParam* quda::EigenSolver::eig_param
protected

Definition at line 20 of file eigensolve_quda.h.

◆ iter

int quda::EigenSolver::iter
protected

Definition at line 41 of file eigensolve_quda.h.

◆ iter_converged

int quda::EigenSolver::iter_converged
protected

Definition at line 42 of file eigensolve_quda.h.

◆ iter_keep

int quda::EigenSolver::iter_keep
protected

Definition at line 44 of file eigensolve_quda.h.

◆ iter_locked

int quda::EigenSolver::iter_locked
protected

Definition at line 43 of file eigensolve_quda.h.

◆ mat

const DiracMatrix& quda::EigenSolver::mat
protected

Definition at line 19 of file eigensolve_quda.h.

◆ max_restarts

int quda::EigenSolver::max_restarts
protected

Definition at line 37 of file eigensolve_quda.h.

◆ n_conv

int quda::EigenSolver::n_conv
protected

Size of Krylov space after extension

Definition at line 27 of file eigensolve_quda.h.

◆ n_ev

int quda::EigenSolver::n_ev
protected

Definition at line 25 of file eigensolve_quda.h.

◆ n_ev_deflate

int quda::EigenSolver::n_ev_deflate
protected

Number of converged eigenvalues requested

Definition at line 28 of file eigensolve_quda.h.

◆ n_kr

int quda::EigenSolver::n_kr
protected

Size of initial factorisation

Definition at line 26 of file eigensolve_quda.h.

◆ num_converged

int quda::EigenSolver::num_converged
protected

Definition at line 45 of file eigensolve_quda.h.

◆ num_keep

int quda::EigenSolver::num_keep
protected

Definition at line 47 of file eigensolve_quda.h.

◆ num_locked

int quda::EigenSolver::num_locked
protected

Definition at line 46 of file eigensolve_quda.h.

◆ profile

TimeProfile& quda::EigenSolver::profile
protected

Definition at line 21 of file eigensolve_quda.h.

◆ r

std::vector<ColorSpinorField *> quda::EigenSolver::r
protected

Definition at line 52 of file eigensolve_quda.h.

◆ residua

std::vector<double> quda::EigenSolver::residua
protected

Definition at line 49 of file eigensolve_quda.h.

◆ restart_iter

int quda::EigenSolver::restart_iter
protected

Definition at line 36 of file eigensolve_quda.h.

◆ reverse

bool quda::EigenSolver::reverse
protected

Tolerance on eigenvalues

Definition at line 30 of file eigensolve_quda.h.

◆ save_prec

QudaPrecision quda::EigenSolver::save_prec
protected

Definition at line 58 of file eigensolve_quda.h.

◆ spectrum

char quda::EigenSolver::spectrum[3]
protected

True if using polynomial acceleration

Definition at line 31 of file eigensolve_quda.h.

◆ tmp1

ColorSpinorField* quda::EigenSolver::tmp1
protected

Definition at line 55 of file eigensolve_quda.h.

◆ tmp2

ColorSpinorField* quda::EigenSolver::tmp2
protected

Definition at line 56 of file eigensolve_quda.h.

◆ tol

double quda::EigenSolver::tol
protected

Number of converged eigenvalues to use in deflation

Definition at line 29 of file eigensolve_quda.h.


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