QUDA
v0.7.0
A library for QCD on GPUs
|
#include <invert_quda.h>
Public Member Functions | |
SolverParam (QudaInvertParam ¶m) | |
~SolverParam () | |
void | updateInvertParam (QudaInvertParam ¶m, int offset=-1) |
void | updateRhsIndex (QudaInvertParam ¶m) |
Public Attributes | |
QudaInverterType | inv_type |
QudaInverterType | inv_type_precondition |
QudaResidualType | residual_type |
QudaUseInitGuess | use_init_guess |
double | delta |
bool | use_sloppy_partial_accumulator |
int | max_res_increase |
int | max_res_increase_total |
int | pipeline |
double | tol |
double | tol_restart |
double | tol_hq |
double | true_res |
double | true_res_hq |
int | maxiter |
int | iter |
QudaPrecision | precision |
QudaPrecision | precision_sloppy |
QudaPrecision | precision_precondition |
QudaPreserveSource | preserve_source |
int | overlap_precondition |
int | num_offset |
double | offset [QUDA_MAX_MULTI_SHIFT] |
double | tol_offset [QUDA_MAX_MULTI_SHIFT] |
double | tol_hq_offset [QUDA_MAX_MULTI_SHIFT] |
double | true_res_offset [QUDA_MAX_MULTI_SHIFT] |
double | true_res_hq_offset [QUDA_MAX_MULTI_SHIFT] |
int | Nsteps |
int | Nkrylov |
int | precondition_cycle |
double | tol_precondition |
int | maxiter_precondition |
double | omega |
QudaSchwarzType | schwarz_type |
double | secs |
double | gflops |
QudaPrecision | precision_ritz |
int | nev |
int | m |
int | deflation_grid |
int | rhs_idx |
SolverParam is the meta data used to define linear solvers.
Definition at line 14 of file invert_quda.h.
|
inline |
Constructor that matches the initial values to that of the QudaInvertParam instance
param | The QudaInvertParam instance from which the values are copied |
Definition at line 159 of file invert_quda.h.
|
inline |
Definition at line 188 of file invert_quda.h.
|
inline |
Update the QudaInvertParam with the data from this
param | the QudaInvertParam to be updated |
Definition at line 194 of file invert_quda.h.
|
inline |
Definition at line 213 of file invert_quda.h.
int quda::SolverParam::deflation_grid |
Definition at line 151 of file invert_quda.h.
double quda::SolverParam::delta |
Whether to keep the partial solution accumulator in sloppy precision
Definition at line 41 of file invert_quda.h.
double quda::SolverParam::gflops |
Definition at line 143 of file invert_quda.h.
QudaInverterType quda::SolverParam::inv_type |
Which linear solver to use
Definition at line 18 of file invert_quda.h.
QudaInverterType quda::SolverParam::inv_type_precondition |
The inner Krylov solver used in the preconditioner. Set to QUDA_INVALID_INVERTER to disable the preconditioner entirely.
Definition at line 24 of file invert_quda.h.
int quda::SolverParam::iter |
The precision used by the QUDA solver
Definition at line 78 of file invert_quda.h.
int quda::SolverParam::m |
Definition at line 150 of file invert_quda.h.
int quda::SolverParam::max_res_increase |
This parameter determines how many total reliable update residual increases we tolerate before terminating the solver, i.e., how long do we want to keep trying to converge
Definition at line 49 of file invert_quda.h.
int quda::SolverParam::max_res_increase_total |
Enable pipeline solver
Definition at line 54 of file invert_quda.h.
int quda::SolverParam::maxiter |
The number of iterations performed by the solver
Definition at line 75 of file invert_quda.h.
int quda::SolverParam::maxiter_precondition |
Maximum number of iterations allowed in the inner solver
Definition at line 129 of file invert_quda.h.
int quda::SolverParam::nev |
Definition at line 149 of file invert_quda.h.
int quda::SolverParam::Nkrylov |
Maximum size of Krylov space used by solver
Definition at line 120 of file invert_quda.h.
int quda::SolverParam::Nsteps |
Number of steps in s-step algorithms
Definition at line 117 of file invert_quda.h.
int quda::SolverParam::num_offset |
< Number of offsets in the multi-shift solver
Definition at line 98 of file invert_quda.h.
double quda::SolverParam::offset[QUDA_MAX_MULTI_SHIFT] |
Offsets for multi-shift solver
Definition at line 101 of file invert_quda.h.
double quda::SolverParam::omega |
Relaxation parameter used in GCR-DD (default = 1.0)
Definition at line 132 of file invert_quda.h.
int quda::SolverParam::overlap_precondition |
Definition at line 93 of file invert_quda.h.
int quda::SolverParam::pipeline |
Solver tolerance in the L2 residual norm
Definition at line 57 of file invert_quda.h.
QudaPrecision quda::SolverParam::precision |
The precision used by the QUDA sloppy operator
Definition at line 81 of file invert_quda.h.
QudaPrecision quda::SolverParam::precision_precondition |
Preserve the source or not in the linear solver (deprecated?)
Definition at line 87 of file invert_quda.h.
QudaPrecision quda::SolverParam::precision_ritz |
< The precision of the Ritz vectors
Definition at line 147 of file invert_quda.h.
QudaPrecision quda::SolverParam::precision_sloppy |
The precision used by the QUDA preconditioner
Definition at line 84 of file invert_quda.h.
int quda::SolverParam::precondition_cycle |
Number of preconditioner cycles to perform per iteration
Definition at line 123 of file invert_quda.h.
QudaPreserveSource quda::SolverParam::preserve_source |
Domain overlap to use in the preconditioning
Definition at line 90 of file invert_quda.h.
QudaResidualType quda::SolverParam::residual_type |
Whether to use the L2 relative residual, L2 absolute residual or Fermilab heavy-quark residual, or combinations therein to determine convergence. To require that multiple stopping conditions are satisfied, use a bitwise OR as follows:
p.residual_type = (QudaResidualType) (QUDA_L2_RELATIVE_RESIDUAL | QUDA_HEAVY_QUARK_RESIDUAL);Whether to use an initial guess in the solver or not
Definition at line 35 of file invert_quda.h.
int quda::SolverParam::rhs_idx |
Definition at line 152 of file invert_quda.h.
QudaSchwarzType quda::SolverParam::schwarz_type |
Whether to use additive or multiplicative Schwarz preconditioning The time taken by the solver
Definition at line 137 of file invert_quda.h.
double quda::SolverParam::secs |
The Gflops rate of the solver
Definition at line 140 of file invert_quda.h.
double quda::SolverParam::tol |
Solver tolerance in the L2 residual norm
Definition at line 60 of file invert_quda.h.
double quda::SolverParam::tol_hq |
Actual L2 residual norm achieved in solver
Definition at line 66 of file invert_quda.h.
double quda::SolverParam::tol_hq_offset[QUDA_MAX_MULTI_SHIFT] |
Solver tolerance for each shift when refinement is applied using the heavy-quark residual
Definition at line 107 of file invert_quda.h.
double quda::SolverParam::tol_offset[QUDA_MAX_MULTI_SHIFT] |
Solver tolerance for each offset
Definition at line 104 of file invert_quda.h.
double quda::SolverParam::tol_precondition |
Tolerance in the inner solver
Definition at line 126 of file invert_quda.h.
double quda::SolverParam::tol_restart |
Solver tolerance in the heavy quark residual norm
Definition at line 63 of file invert_quda.h.
double quda::SolverParam::true_res |
Actual heavy quark residual norm achieved in solver
Definition at line 69 of file invert_quda.h.
double quda::SolverParam::true_res_hq |
Maximum number of iterations in the linear solver
Definition at line 72 of file invert_quda.h.
double quda::SolverParam::true_res_hq_offset[QUDA_MAX_MULTI_SHIFT] |
Actual heavy quark residual norm achieved in solver for each offset
Definition at line 113 of file invert_quda.h.
double quda::SolverParam::true_res_offset[QUDA_MAX_MULTI_SHIFT] |
Actual L2 residual norm achieved in solver for each offset
Definition at line 110 of file invert_quda.h.
QudaUseInitGuess quda::SolverParam::use_init_guess |
Reliable update tolerance
Definition at line 38 of file invert_quda.h.
bool quda::SolverParam::use_sloppy_partial_accumulator |
This parameter determines how many consective reliable update residual increases we tolerate before terminating the solver, i.e., how long do we want to keep trying to converge
Definition at line 44 of file invert_quda.h.