QUDA  v1.1.0
A library for QCD on GPUs
Public Member Functions | Public Attributes | List of all members
quda::MGParam Struct Reference

#include <multigrid.h>

+ Inheritance diagram for quda::MGParam:

Public Member Functions

 MGParam (QudaMultigridParam &param, std::vector< ColorSpinorField * > &B, DiracMatrix *matResidual, DiracMatrix *matSmooth, DiracMatrix *matSmoothSloppy, int level=0)
 
 MGParam (const MGParam &param, std::vector< ColorSpinorField * > &B, DiracMatrix *matResidual, DiracMatrix *matSmooth, DiracMatrix *matSmoothSloppy, int level=0)
 
- Public Member Functions inherited from quda::SolverParam
 SolverParam ()
 
 SolverParam (const QudaInvertParam &param)
 
 SolverParam (const SolverParam &param)
 
 ~SolverParam ()
 
void updateInvertParam (QudaInvertParam &param, int offset=-1)
 
void updateRhsIndex (QudaInvertParam &param)
 

Public Attributes

QudaMultigridParammg_global
 
int level
 
int Nlevel
 
int geoBlockSize [QUDA_MAX_DIM]
 
int spinBlockSize
 
int Nvec
 
int NblockOrtho
 
MGcoarse
 
MGfine
 
std::vector< ColorSpinorField * > & B
 
int nu_pre
 
int nu_post
 
double smoother_tol
 
QudaMultigridCycleType cycle_type
 
QudaBoolean global_reduction
 
DiracMatrixmatResidual
 
DiracMatrixmatSmooth
 
DiracMatrixmatSmoothSloppy
 
QudaInverterType smoother
 
QudaSolutionType coarse_grid_solution_type
 
QudaSolveType smoother_solve_type
 
QudaFieldLocation location
 
QudaFieldLocation setup_location
 
char filename [100]
 
QudaTransferType transfer_type
 
bool use_mma
 
- Public Attributes inherited from quda::SolverParam
QudaInverterType inv_type
 
QudaInverterType inv_type_precondition
 
void * preconditioner
 
void * deflation_op
 
QudaResidualType residual_type
 
bool deflate
 
QudaEigParam eig_param
 
QudaUseInitGuess use_init_guess
 
QudaComputeNullVector compute_null_vector
 
double delta
 
bool use_alternative_reliable
 
bool use_sloppy_partial_accumulator
 
int solution_accumulator_pipeline
 
int max_res_increase
 
int max_res_increase_total
 
int max_hq_res_increase
 
int max_hq_res_restart_total
 
int heavy_quark_check
 
int pipeline
 
double tol
 
double tol_restart
 
double tol_hq
 
bool compute_true_res
 
bool sloppy_converge
 
double true_res
 
double true_res_hq
 
int maxiter
 
int iter
 
QudaPrecision precision
 
QudaPrecision precision_sloppy
 
QudaPrecision precision_refinement_sloppy
 
QudaPrecision precision_precondition
 
QudaPrecision precision_eigensolver
 
QudaPreserveSource preserve_source
 
bool return_residual
 
int overlap_precondition
 
int num_src
 
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 iter_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
 
QudaCABasis ca_basis
 
double ca_lambda_min
 
double ca_lambda_max
 
QudaSchwarzType schwarz_type
 
double secs
 
double gflops
 
QudaPrecision precision_ritz
 
int n_ev
 
int m
 
int deflation_grid
 
int rhs_idx
 
int eigcg_max_restarts
 
int max_restart_num
 
double inc_tol
 
double eigenval_tol
 
QudaVerbosity verbosity_precondition
 
bool is_preconditioner
 verbosity to use for preconditioner More...
 
bool global_reduction
 whether the solver acting as a preconditioner for another solver More...
 
bool mg_instance
 whether to use a global or local (node) reduction for this solver More...
 
QudaExtLibType extlib_type
 

Detailed Description

This struct contains all the metadata required to define the multigrid solver. For each level of multigrid we will have an instance of MGParam describing all the meta data appropriate for given level.

Definition at line 25 of file multigrid.h.

Constructor & Destructor Documentation

◆ MGParam() [1/2]

quda::MGParam::MGParam ( QudaMultigridParam param,
std::vector< ColorSpinorField * > &  B,
DiracMatrix matResidual,
DiracMatrix matSmooth,
DiracMatrix matSmoothSloppy,
int  level = 0 
)
inline

This is top level instantiation done when we start creating the multigrid operator.

Definition at line 110 of file multigrid.h.

◆ MGParam() [2/2]

quda::MGParam::MGParam ( const MGParam param,
std::vector< ColorSpinorField * > &  B,
DiracMatrix matResidual,
DiracMatrix matSmooth,
DiracMatrix matSmoothSloppy,
int  level = 0 
)
inline

Definition at line 143 of file multigrid.h.

Member Data Documentation

◆ B

std::vector<ColorSpinorField*>& quda::MGParam::B

The null space vectors

Definition at line 56 of file multigrid.h.

◆ coarse

MG* quda::MGParam::coarse

This is the next lower level

Definition at line 50 of file multigrid.h.

◆ coarse_grid_solution_type

QudaSolutionType quda::MGParam::coarse_grid_solution_type

The type of residual to send to the next coarse grid, and thus the type of solution to receive back from this coarse grid

Definition at line 87 of file multigrid.h.

◆ cycle_type

QudaMultigridCycleType quda::MGParam::cycle_type

Multigrid cycle type

Definition at line 68 of file multigrid.h.

◆ filename

char quda::MGParam::filename[100]

Filename for where to load/store the null space

Definition at line 99 of file multigrid.h.

◆ fine

MG* quda::MGParam::fine

This is the immediate finer level

Definition at line 53 of file multigrid.h.

◆ geoBlockSize

int quda::MGParam::geoBlockSize[QUDA_MAX_DIM]

Geometric block size

Definition at line 38 of file multigrid.h.

◆ global_reduction

QudaBoolean quda::MGParam::global_reduction

Whether to use global or local (node) reductions

Definition at line 71 of file multigrid.h.

◆ level

int quda::MGParam::level

What is the level of this instance

Definition at line 32 of file multigrid.h.

◆ location

QudaFieldLocation quda::MGParam::location

Where to compute this level of multigrid

Definition at line 93 of file multigrid.h.

◆ matResidual

DiracMatrix* quda::MGParam::matResidual

The Dirac operator to use for residual computation

Definition at line 74 of file multigrid.h.

◆ matSmooth

DiracMatrix* quda::MGParam::matSmooth

The Dirac operator to use for smoothing

Definition at line 77 of file multigrid.h.

◆ matSmoothSloppy

DiracMatrix* quda::MGParam::matSmoothSloppy

The sloppy Dirac operator to use for smoothing

Definition at line 80 of file multigrid.h.

◆ mg_global

QudaMultigridParam& quda::MGParam::mg_global

This points to the parameter struct that is passed into QUDA. We use this to set per-level parameters

Definition at line 29 of file multigrid.h.

◆ NblockOrtho

int quda::MGParam::NblockOrtho

Number of times to apply Gram-Schmidt within a block

Definition at line 47 of file multigrid.h.

◆ Nlevel

int quda::MGParam::Nlevel

Number of levels in the solver

Definition at line 35 of file multigrid.h.

◆ nu_post

int quda::MGParam::nu_post

Number of pre-smoothing applications to perform

Definition at line 62 of file multigrid.h.

◆ nu_pre

int quda::MGParam::nu_pre

Number of pre-smoothing applications to perform

Definition at line 59 of file multigrid.h.

◆ Nvec

int quda::MGParam::Nvec

Number of vectors used to define coarse space

Definition at line 44 of file multigrid.h.

◆ setup_location

QudaFieldLocation quda::MGParam::setup_location

Where to compute this level of the multigrid setup

Definition at line 96 of file multigrid.h.

◆ smoother

QudaInverterType quda::MGParam::smoother

What type of smoother to use

Definition at line 83 of file multigrid.h.

◆ smoother_solve_type

QudaSolveType quda::MGParam::smoother_solve_type

The type of smoother solve to do on each grid (e/o preconditioning or not)

Definition at line 90 of file multigrid.h.

◆ smoother_tol

double quda::MGParam::smoother_tol

Tolerance to use for the solver / smoother (if applicable)

Definition at line 65 of file multigrid.h.

◆ spinBlockSize

int quda::MGParam::spinBlockSize

Spin block size

Definition at line 41 of file multigrid.h.

◆ transfer_type

QudaTransferType quda::MGParam::transfer_type

Whether or not this is a staggered solve or not

Definition at line 102 of file multigrid.h.

◆ use_mma

bool quda::MGParam::use_mma

Whether to use tensor cores (if available)

Definition at line 105 of file multigrid.h.


The documentation for this struct was generated from the following file: