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

#include <invert_quda.h>

+ Inheritance diagram for quda::Solver:

Public Member Functions

 Solver (const DiracMatrix &mat, const DiracMatrix &matSloppy, const DiracMatrix &matPrecon, const DiracMatrix &matEig, SolverParam &param, TimeProfile &profile)
 
virtual ~Solver ()
 
virtual void operator() (ColorSpinorField &out, ColorSpinorField &in)=0
 
virtual void blocksolve (ColorSpinorField &out, ColorSpinorField &in)
 
const DiracMatrixM ()
 
const DiracMatrixMsloppy ()
 
const DiracMatrixMprecon ()
 
const DiracMatrixMeig ()
 
virtual bool hermitian ()=0
 
bool convergence (double r2, double hq2, double r2_tol, double hq_tol)
 
bool convergenceHQ (double r2, double hq2, double r2_tol, double hq_tol)
 Test for HQ solver convergence – ignore L2 residual. More...
 
bool convergenceL2 (double r2, double hq2, double r2_tol, double hq_tol)
 Test for L2 solver convergence – ignore HQ residual. More...
 
void PrintStats (const char *name, int k, double r2, double b2, double hq2)
 Prints out the running statistics of the solver (requires a verbosity of QUDA_VERBOSE) More...
 
void PrintSummary (const char *name, int k, double r2, double b2, double r2_tol, double hq_tol)
 Prints out the summary of the solver convergence (requires a verbosity of QUDA_SUMMARIZE). Assumes SolverParam.true_res and SolverParam.true_res_hq has been set. More...
 
double precisionEpsilon (QudaPrecision prec=QUDA_INVALID_PRECISION) const
 Returns the epsilon tolerance for a given precision, by default returns the solver precision. More...
 
void constructDeflationSpace (const ColorSpinorField &meta, const DiracMatrix &mat)
 Constructs the deflation space and eigensolver. More...
 
void destroyDeflationSpace ()
 Destroy the allocated deflation space. More...
 
void extendSVDDeflationSpace ()
 Extends the deflation space to twice its size for SVD deflation. More...
 
void injectDeflationSpace (std::vector< ColorSpinorField * > &defl_space)
 Injects a deflation space into the solver from the vector argument. Note the input space is reduced to zero size as a result of calling this function, with responsibility for the space transferred to the solver. More...
 
void extractDeflationSpace (std::vector< ColorSpinorField * > &defl_space)
 Extracts the deflation space from the solver to the vector argument. Note the solver deflation space is reduced to zero size as a result of calling this function, with responsibility for the space transferred to the argument. More...
 
int deflationSpaceSize () const
 Returns the size of deflation space. More...
 
void setDeflateCompute (bool flag)
 Sets the deflation compute boolean. More...
 
void setRecomputeEvals (bool flag)
 Sets the recompute evals boolean. More...
 
virtual double flops () const
 Return flops. More...
 

Static Public Member Functions

static Solvercreate (SolverParam &param, const DiracMatrix &mat, const DiracMatrix &matSloppy, const DiracMatrix &matPrecon, const DiracMatrix &matEig, TimeProfile &profile)
 Solver factory. More...
 
static double stopping (double tol, double b2, QudaResidualType residual_type)
 Set the solver L2 stopping condition. More...
 

Protected Attributes

const DiracMatrixmat
 
const DiracMatrixmatSloppy
 
const DiracMatrixmatPrecon
 
const DiracMatrixmatEig
 
SolverParamparam
 
TimeProfileprofile
 
int node_parity
 
EigenSolvereig_solve
 
bool deflate_init
 
bool deflate_compute
 
bool recompute_evals
 
std::vector< ColorSpinorField * > evecs
 
std::vector< Complexevals
 

Detailed Description

Definition at line 462 of file invert_quda.h.

Constructor & Destructor Documentation

◆ Solver()

quda::Solver::Solver ( const DiracMatrix mat,
const DiracMatrix matSloppy,
const DiracMatrix matPrecon,
const DiracMatrix matEig,
SolverParam param,
TimeProfile profile 
)

Holds the eigenvalues.

Definition at line 14 of file solver.cpp.

◆ ~Solver()

quda::Solver::~Solver ( )
virtual

Definition at line 33 of file solver.cpp.

Member Function Documentation

◆ blocksolve()

void quda::Solver::blocksolve ( ColorSpinorField out,
ColorSpinorField in 
)
virtual

Reimplemented in quda::CG.

Definition at line 302 of file solver.cpp.

◆ constructDeflationSpace()

void quda::Solver::constructDeflationSpace ( const ColorSpinorField meta,
const DiracMatrix mat 
)

Constructs the deflation space and eigensolver.

Parameters
[in]metaA sample ColorSpinorField with which to instantiate the eigensolver
[in]matThe operator to eigensolve
[in]Whetherto compute the SVD

Definition at line 168 of file solver.cpp.

◆ convergence()

bool quda::Solver::convergence ( double  r2,
double  hq2,
double  r2_tol,
double  hq_tol 
)

@briefTest for solver convergence

Parameters
[in]r2L2 norm squared of the residual
[in]hq2Heavy quark residual
[in]r2_tolSolver L2 tolerance
[in]hq_tolSolver heavy-quark tolerance
Returns
Whether converged

Definition at line 328 of file solver.cpp.

◆ convergenceHQ()

bool quda::Solver::convergenceHQ ( double  r2,
double  hq2,
double  r2_tol,
double  hq_tol 
)

Test for HQ solver convergence – ignore L2 residual.

Parameters
[in]r2L2 norm squared of the residual
[in]hq2Heavy quark residual
[in]r2_tolSolver L2 tolerance
[in[hq_tol Solver heavy-quark tolerance
Returns
Whether converged

Definition at line 348 of file solver.cpp.

◆ convergenceL2()

bool quda::Solver::convergenceL2 ( double  r2,
double  hq2,
double  r2_tol,
double  hq_tol 
)

Test for L2 solver convergence – ignore HQ residual.

Parameters
[in]r2L2 norm squared of the residual
[in]hq2Heavy quark residual
[in]r2_tolSolver L2 tolerance
[in]hq_tolSolver heavy-quark tolerance

Definition at line 361 of file solver.cpp.

◆ create()

Solver * quda::Solver::create ( SolverParam param,
const DiracMatrix mat,
const DiracMatrix matSloppy,
const DiracMatrix matPrecon,
const DiracMatrix matEig,
TimeProfile profile 
)
static

Solver factory.

Definition at line 42 of file solver.cpp.

◆ deflationSpaceSize()

int quda::Solver::deflationSpaceSize ( ) const
inline

Returns the size of deflation space.

Definition at line 615 of file invert_quda.h.

◆ destroyDeflationSpace()

void quda::Solver::destroyDeflationSpace ( )

Destroy the allocated deflation space.

Definition at line 229 of file solver.cpp.

◆ extendSVDDeflationSpace()

void quda::Solver::extendSVDDeflationSpace ( )

Extends the deflation space to twice its size for SVD deflation.

Definition at line 287 of file solver.cpp.

◆ extractDeflationSpace()

void quda::Solver::extractDeflationSpace ( std::vector< ColorSpinorField * > &  defl_space)

Extracts the deflation space from the solver to the vector argument. Note the solver deflation space is reduced to zero size as a result of calling this function, with responsibility for the space transferred to the argument.

Parameters
[in,out]defl_spacethe extracted deflation space. On input, this vector should have zero size.

Definition at line 276 of file solver.cpp.

◆ flops()

virtual double quda::Solver::flops ( ) const
inlinevirtual

Return flops.

Returns
flops expended by this operator

Reimplemented in quda::MG.

Definition at line 633 of file invert_quda.h.

◆ hermitian()

virtual bool quda::Solver::hermitian ( )
pure virtual

◆ injectDeflationSpace()

void quda::Solver::injectDeflationSpace ( std::vector< ColorSpinorField * > &  defl_space)

Injects a deflation space into the solver from the vector argument. Note the input space is reduced to zero size as a result of calling this function, with responsibility for the space transferred to the solver.

Parameters
[in,out]defl_spacethe deflation space we wish to transfer to the solver.

Definition at line 266 of file solver.cpp.

◆ M()

const DiracMatrix& quda::Solver::M ( )
inline

Definition at line 489 of file invert_quda.h.

◆ Meig()

const DiracMatrix& quda::Solver::Meig ( )
inline

Definition at line 492 of file invert_quda.h.

◆ Mprecon()

const DiracMatrix& quda::Solver::Mprecon ( )
inline

Definition at line 491 of file invert_quda.h.

◆ Msloppy()

const DiracMatrix& quda::Solver::Msloppy ( )
inline

Definition at line 490 of file invert_quda.h.

◆ operator()()

virtual void quda::Solver::operator() ( ColorSpinorField out,
ColorSpinorField in 
)
pure virtual

◆ precisionEpsilon()

double quda::Solver::precisionEpsilon ( QudaPrecision  prec = QUDA_INVALID_PRECISION) const

Returns the epsilon tolerance for a given precision, by default returns the solver precision.

Parameters
[in]precInput precision, default value is solver precision

Definition at line 412 of file solver.cpp.

◆ PrintStats()

void quda::Solver::PrintStats ( const char *  name,
int  k,
double  r2,
double  b2,
double  hq2 
)

Prints out the running statistics of the solver (requires a verbosity of QUDA_VERBOSE)

Parameters
[in]nameName of solver that called this
[in]kiteration count
[in]r2L2 norm squared of the residual
[in]hq2Heavy quark residual

Definition at line 373 of file solver.cpp.

◆ PrintSummary()

void quda::Solver::PrintSummary ( const char *  name,
int  k,
double  r2,
double  b2,
double  r2_tol,
double  hq_tol 
)

Prints out the summary of the solver convergence (requires a verbosity of QUDA_SUMMARIZE). Assumes SolverParam.true_res and SolverParam.true_res_hq has been set.

Parameters
[in]nameName of solver that called this
[in]kiteration count
[in]r2L2 norm squared of the residual
[in]hq2Heavy quark residual
[in]r2_tolSolver L2 tolerance
[in]hq_tolSolver heavy-quark tolerance

Definition at line 386 of file solver.cpp.

◆ setDeflateCompute()

void quda::Solver::setDeflateCompute ( bool  flag)
inline

Sets the deflation compute boolean.

Parameters
[in]flagSet to this boolean value

Definition at line 621 of file invert_quda.h.

◆ setRecomputeEvals()

void quda::Solver::setRecomputeEvals ( bool  flag)
inline

Sets the recompute evals boolean.

Parameters
[in]flagSet to this boolean value

Definition at line 627 of file invert_quda.h.

◆ stopping()

double quda::Solver::stopping ( double  tol,
double  b2,
QudaResidualType  residual_type 
)
static

Set the solver L2 stopping condition.

Parameters
[in]Desiredsolver tolerance
[in]b2L2 norm squared of the source vector
[in]residual_typeThe type of residual we want to solve for
Returns
L2 stopping condition

Definition at line 311 of file solver.cpp.

Member Data Documentation

◆ deflate_compute

bool quda::Solver::deflate_compute
protected

If true, the deflation space has been computed.

Definition at line 475 of file invert_quda.h.

◆ deflate_init

bool quda::Solver::deflate_init
protected

Eigensolver object.

Definition at line 474 of file invert_quda.h.

◆ eig_solve

EigenSolver* quda::Solver::eig_solve
protected

Definition at line 473 of file invert_quda.h.

◆ evals

std::vector<Complex> quda::Solver::evals
protected

Holds the eigenvectors.

Definition at line 478 of file invert_quda.h.

◆ evecs

std::vector<ColorSpinorField *> quda::Solver::evecs
protected

If true, instruct the solver to recompute evals from an existing deflation space.

Definition at line 477 of file invert_quda.h.

◆ mat

const DiracMatrix& quda::Solver::mat
protected

Definition at line 465 of file invert_quda.h.

◆ matEig

const DiracMatrix& quda::Solver::matEig
protected

Definition at line 468 of file invert_quda.h.

◆ matPrecon

const DiracMatrix& quda::Solver::matPrecon
protected

Definition at line 467 of file invert_quda.h.

◆ matSloppy

const DiracMatrix& quda::Solver::matSloppy
protected

Definition at line 466 of file invert_quda.h.

◆ node_parity

int quda::Solver::node_parity
protected

Definition at line 472 of file invert_quda.h.

◆ param

SolverParam& quda::Solver::param
protected

Definition at line 470 of file invert_quda.h.

◆ profile

TimeProfile& quda::Solver::profile
protected

Definition at line 471 of file invert_quda.h.

◆ recompute_evals

bool quda::Solver::recompute_evals
protected

If true, instruct the solver to create a deflation space.

Definition at line 476 of file invert_quda.h.


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