QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
quda.h
Go to the documentation of this file.
1 #ifndef _QUDA_H
2 #define _QUDA_H
3 
13 #include <enum_quda.h>
14 #include <stdio.h> /* for FILE */
15 #include <quda_constants.h>
16 
17 #ifdef __cplusplus
18 extern "C" {
19 #endif
20 
25  typedef struct QudaGaugeParam_s {
26 
29  int X[4];
31  double anisotropy;
32  double tadpole_coeff;
33  double scale;
53  int ga_pad;
57  int staple_pad;
59  int mom_ga_pad;
60  double gaugeGiB;
67  int overlap;
75 
76 
80  typedef struct QudaInvertParam_s {
81 
88  double mass;
89  double kappa;
91  double m5;
92  int Ls;
94  double b_5[QUDA_MAX_DWF_LS];
95  double c_5[QUDA_MAX_DWF_LS];
97  double mu;
98  double epsilon;
102  double tol;
103  double tol_restart;
104  double tol_hq;
105  double true_res;
106  double true_res_hq;
107  int maxiter;
108  double reliable_delta;
115 
120 
121  int pipeline;
125  int overlap;
129 
132 
135 
138 
141 
169  double clover_coeff;
172  double trlogA[2];
176  int sp_pad;
177  int cl_pad;
179  int iter;
180  double spinorGiB;
181  double cloverGiB;
182  double gflops;
183  double secs;
189  int Nsteps;
190 
193 
194  /*
195  * The following parameters are related to the domain-decomposed
196  * preconditioner, if enabled.
197  */
198 
204 
211 
214 
217 
219  double omega;
220 
223 
226 
236 
240  int nev;
241 
242  int max_search_dim;//for magma library this parameter must be multiple 16?
243 
244  int rhs_idx;
245 
246  int deflation_grid;//total deflation space is nev*deflation_grid
247 
248  } QudaInvertParam;
249 
250 
251  // Parameter set for solving the eigenvalue problems.
252  // Eigen problems are tightly related with Ritz algorithm.
253  // And the Lanczos algorithm use the Ritz operator.
254  // For Ritz matrix operation,
255  // we need to know about the solution type of dirac operator.
256  // For acceleration, we are also using chevisov polynomial method.
257  // And nk, np values are needed Implicit Restart Lanczos method
258  // which is optimized form of Lanczos algorithm
259  typedef struct QudaEigParam_s {
260 
265 
266  double *MatPoly_param;
267  int NPoly;
268  double Stp_residual;
269  int nk;
270  int np;
271  int f_size;
272  double eigen_shift;
273 
274  } QudaEigParam;
275 
276 
277 
278  /*
279  * Interface functions, found in interface_quda.cpp
280  */
281 
307  void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[],
308  FILE *outfile);
309 
320  typedef int (*QudaCommsMap)(const int *coords, void *fdata);
321 
348  void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata);
349 
360  void initQudaDevice(int device);
361 
368  void initQudaMemory();
369 
379  void initQuda(int device);
380 
384  void endQuda(void);
385 
394 
403 
412 
418 
424 
430 
436  void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param);
437 
441  void freeGaugeQuda(void);
442 
448  void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param);
449 
457  void loadCloverQuda(void *h_clover, void *h_clovinv,
459 
463  void freeCloverQuda(void);
464 
474  void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V,
475  void *hp_alpha, void *hp_beta, QudaEigParam *eig_param);
476 
486  void invertQuda(void *h_x, void *h_b, QudaInvertParam *param);
487 
500  void invertMultiShiftMDQuda(void **_hp_xe, void **_hp_xo, void **_hp_ye,
501  void **_hp_yo, void *_hp_b, QudaInvertParam *param);
502 
510  void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param);
511 
521  void incrementalEigQuda(void *_h_x, void *_h_b, QudaInvertParam *param, void *_h_u, double *inv_eigenvals, int last_rhs);
522 
531  void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
533 
543  void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param,
545 
555  void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param,
557 
567  void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
568  QudaParity *parity, int inverse);
569 
577  void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
578 
586  void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
587 
588 
589  /*
590  * The following routines are temporary additions used by the HISQ
591  * link-fattening code.
592  */
593 
594  void set_dim(int *);
595  void pack_ghost(void **cpuLink, void **cpuGhost, int nFace,
596  QudaPrecision precision);
598 
599  void computeKSLinkQuda(void* fatlink, void* longlink, void* ulink, void* inlink,
600  double *path_coeff, QudaGaugeParam *param, QudaComputeFatMethod method);
601 
602 
603 
618  int computeGaugeForceQuda(void* mom, void* sitelink, int*** input_path_buf, int* path_length,
619  double* loop_coeff, int num_paths, int max_length, double dt,
620  QudaGaugeParam* qudaGaugeParam, double* timeinfo);
621 
633  void updateGaugeFieldQuda(void* gauge, void* momentum, double dt,
634  int conj_mom, int exact, QudaGaugeParam* param);
635 
636 
641  void* createExtendedGaugeField(void* gauge, int geometry, QudaGaugeParam* param);
642 
643  void* createGaugeField(void* gauge, int geometry, QudaGaugeParam* param);
644 
645  void saveGaugeField(void* outGauge, void* inGauge, QudaGaugeParam* param);
646 
647  void extendGaugeField(void* outGauge, void* inGauge);
648 
649 
653  void destroyQudaGaugeField(void* gauge);
654 
656 
657 
658  void computeCloverTraceQuda(void* out, void* clover, int mu, int nu, int dim[4]);
659 
660  void computeCloverDerivativeQuda(void* out, void* gauge, void* oprod, int mu, int nu,
661  double coeff,
662  QudaParity parity, QudaGaugeParam* param, int conjugate);
663 
673  void computeStaggeredOprodQuda(void** oprod, void** quark, int num, double** coeff, QudaGaugeParam* param);
674 
675  void computeStaggeredForceQuda(void* mom, void* quark, double* coeff);
676 
686  void computeAsqtadForceQuda(void* const momentum,
687  long long* flops,
688  const double act_path_coeff[6],
689  const void* const one_link_src[4],
690  const void* const naik_src[4],
691  const void* const link,
692  const QudaGaugeParam* param);
693 
694 
709  void computeHISQForceQuda(void* momentum,
710  long long* flops,
711  const double level2_coeff[6],
712  const double fat7_coeff[6],
713  const void* const staple_src[4],
714  const void* const one_link_src[4],
715  const void* const naik_src[4],
716  const void* const w_link,
717  const void* const v_link,
718  const void* const u_link,
719  const QudaGaugeParam* param);
720 
721 
722 
723  void computeHISQForceCompleteQuda(void* momentum,
724  const double level2_coeff[6],
725  const double fat7_coeff[6],
726  void** quark_array,
727  int num_terms,
728  double** quark_coeff,
729  const void* const w_link,
730  const void* const v_link,
731  const void* const u_link,
732  const QudaGaugeParam* param);
733 
734  void performAPEnStep(unsigned int nSteps, double alpha);
735  double plaqCuda();
740  void openMagma();
741 
742  void closeMagma();
743 
744 
745 #ifdef __cplusplus
746 }
747 #endif
748 
749 /* #include <quda_new_interface.h> */
750 
751 #endif /* _QUDA_H */
752 
int maxiter_precondition
Definition: quda.h:216
double secs
Definition: quda.h:183
enum QudaPreserveSource_s QudaPreserveSource
QudaDiracFieldOrder dirac_order
Definition: quda.h:156
QudaMassNormalization mass_normalization
Definition: quda.h:146
double tol_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:134
enum QudaMassNormalization_s QudaMassNormalization
QudaReconstructType reconstruct_sloppy
Definition: quda.h:46
double anisotropy
Definition: quda.h:31
void computeHISQForceCompleteQuda(void *momentum, const double level2_coeff[6], const double fat7_coeff[6], void **quark_array, int num_terms, double **quark_coeff, const void *const w_link, const void *const v_link, const void *const u_link, const QudaGaugeParam *param)
void freeCloverQuda(void)
void destroyQudaGaugeField(void *gauge)
void computeCloverTraceQuda(void *out, void *clover, int mu, int nu, int dim[4])
double plaqCuda()
double b_5[QUDA_MAX_DWF_LS]
Definition: quda.h:94
void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param)
void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[], FILE *outfile)
void endQuda(void)
QudaFieldLocation clover_location
Definition: quda.h:160
QudaSolveType solve_type
Definition: quda.h:143
QudaVerbosity verbosity_precondition
Definition: quda.h:210
void computeHISQForceQuda(void *momentum, long long *flops, const double level2_coeff[6], const double fat7_coeff[6], const void *const staple_src[4], const void *const one_link_src[4], const void *const naik_src[4], const void *const w_link, const void *const v_link, const void *const u_link, const QudaGaugeParam *param)
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:53
int make_resident_mom
Definition: quda.h:72
int NPoly
Definition: quda.h:267
double mu
Definition: quda.h:97
QudaGaugeFixed gauge_fix
Definition: quda.h:51
QudaTune tune
Definition: quda.h:185
QudaSchwarzType schwarz_type
Definition: quda.h:225
enum QudaResidualType_s QudaResidualType
#define QUDA_MAX_MULTI_SHIFT
Maximum number of shifts supported by the multi-shift solver. This number may be changed if need be...
QudaInverterType inv_type_precondition
Definition: quda.h:203
void printQudaGaugeParam(QudaGaugeParam *param)
Definition: check_params.h:40
QudaLinkType type
Definition: quda.h:35
int max_res_increase_total
Definition: quda.h:119
double kappa
Definition: quda.h:89
void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
QudaPrecision cuda_prec_ritz
Definition: quda.h:238
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
double tol
Definition: quda.h:102
QudaDslashType dslash_type
Definition: quda.h:85
QudaReconstructType reconstruct_precondition
Definition: quda.h:49
QudaInverterType inv_type
Definition: quda.h:86
QudaPrecision cuda_prec
Definition: quda.h:152
double c_5[QUDA_MAX_DWF_LS]
Definition: quda.h:95
void computeKSLinkQuda(void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param, QudaComputeFatMethod method)
void performAPEnStep(unsigned int nSteps, double alpha)
enum QudaSolveType_s QudaSolveType
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
int computeGaugeForceQuda(void *mom, void *sitelink, int ***input_path_buf, int *path_length, double *loop_coeff, int num_paths, int max_length, double dt, QudaGaugeParam *qudaGaugeParam, double *timeinfo)
__global__ void const RealA *const const RealA *const const RealA *const const RealB *const const RealB *const int int mu
double Stp_residual
Definition: quda.h:268
QudaPrecision cpu_prec
Definition: quda.h:151
QudaStaggeredPhase staggered_phase_type
Definition: quda.h:64
double trlogA[2]
Definition: quda.h:172
void * createGaugeField(void *gauge, int geometry, QudaGaugeParam *param)
QudaDagType dagger
Definition: quda.h:145
int test_type
Definition: test_util.cpp:1564
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
double true_res
Definition: quda.h:105
void computeAsqtadForceQuda(void *const momentum, long long *flops, const double act_path_coeff[6], const void *const one_link_src[4], const void *const naik_src[4], const void *const link, const QudaGaugeParam *param)
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
void setFatLinkPadding(QudaComputeFatMethod method, QudaGaugeParam *param)
enum QudaEigType_s QudaEigType
enum QudaTboundary_s QudaTboundary
void openMagma()
int compute_clover_trlog
Definition: quda.h:171
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:163
struct QudaInvertParam_s QudaInvertParam
void invertMultiShiftMDQuda(void **_hp_xe, void **_hp_xo, void **_hp_ye, void **_hp_yo, void *_hp_b, QudaInvertParam *param)
void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
QudaFieldLocation input_location
Definition: quda.h:82
void freeGaugeQuda(void)
void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata)
int staggered_phase_applied
Definition: quda.h:65
double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:140
double reliable_delta
Definition: quda.h:108
void * longlink[4]
QudaUseInitGuess use_init_guess
Definition: quda.h:167
QudaGaugeParam param
Definition: pack_test.cpp:17
int llfat_ga_pad
Definition: quda.h:58
QudaSolutionType solution_type
Definition: quda.h:142
QudaSolverNormalization solver_normalization
Definition: quda.h:147
QudaEigType eig_type
Definition: quda.h:264
void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaPrecision clover_cuda_prec
Definition: quda.h:162
int precondition_cycle
Definition: quda.h:222
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
double * MatPoly_param
Definition: quda.h:266
QudaInvertParam * invert_param
Definition: quda.h:261
double scale
Definition: quda.h:33
void initQuda(int device)
double spinorGiB
Definition: quda.h:180
void initQudaMemory()
QudaFieldLocation output_location
Definition: quda.h:83
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:164
int site_ga_pad
Definition: quda.h:55
void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param)
double m5
Definition: quda.h:91
enum QudaStaggeredPhase_s QudaStaggeredPhase
void createCloverQuda(QudaInvertParam *param)
QudaPrecision cuda_prec_sloppy
Definition: quda.h:153
QudaVerbosity verbosity
Definition: quda.h:174
double tol_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:131
double true_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:137
QudaInvertParam newQudaInvertParam(void)
double eigen_shift
Definition: quda.h:272
double gflops
Definition: quda.h:182
void computeCloverDerivativeQuda(void *out, void *gauge, void *oprod, int mu, int nu, double coeff, QudaParity parity, QudaGaugeParam *param, int conjugate)
QudaPrecision cuda_prec_precondition
Definition: quda.h:48
int staple_pad
Definition: quda.h:57
QudaCloverFieldOrder clover_order
Definition: quda.h:166
double tol_hq
Definition: quda.h:104
enum QudaMatPCType_s QudaMatPCType
int preserve_gauge
Definition: quda.h:62
struct QudaEigParam_s QudaEigParam
void extendGaugeField(void *outGauge, void *inGauge)
double true_res_hq
Definition: quda.h:106
enum QudaGaugeFixed_s QudaGaugeFixed
enum QudaSolutionType_s QudaSolutionType
QudaGammaBasis gamma_basis
Definition: quda.h:158
__constant__ double coeff
enum QudaSchwarzType_s QudaSchwarzType
QudaPrecision cuda_prec_sloppy
Definition: quda.h:45
int max_search_dim
Definition: quda.h:242
double tol_precondition
Definition: quda.h:213
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
double offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:128
void incrementalEigQuda(void *_h_x, void *_h_b, QudaInvertParam *param, void *_h_u, double *inv_eigenvals, int last_rhs)
int use_sloppy_partial_accumulator
Definition: quda.h:109
enum QudaDagType_s QudaDagType
QudaSolutionType RitzMat_Convcheck
Definition: quda.h:263
QudaSolutionType RitzMat_lanczos
Definition: quda.h:262
enum QudaParity_s QudaParity
QudaReconstructType reconstruct
Definition: quda.h:43
enum QudaLinkType_s QudaLinkType
QudaPrecision cuda_prec
Definition: quda.h:42
int X[4]
Definition: quda.h:29
double mass
Definition: quda.h:88
int gcrNkrylov
Definition: quda.h:192
void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta, QudaEigParam *eig_param)
void printQudaInvertParam(QudaInvertParam *param)
Definition: check_params.h:162
int max_res_increase
Definition: quda.h:114
enum QudaFieldLocation_s QudaFieldLocation
void closeMagma()
QudaFieldLocation location
Definition: quda.h:27
QudaInvertParam inv_param
Definition: dslash_test.cpp:38
struct QudaGaugeParam_s QudaGaugeParam
void computeStaggeredOprodQuda(void **oprod, void **quark, int num, double **coeff, QudaGaugeParam *param)
double tadpole_coeff
Definition: quda.h:32
cpuColorSpinorField * out
double gaugeGiB
Definition: quda.h:60
QudaPrecision cuda_prec_precondition
Definition: quda.h:154
int deflation_grid
Definition: quda.h:246
void * createExtendedGaugeField(void *gauge, int geometry, QudaGaugeParam *param)
void printQudaEigParam(QudaEigParam *param)
Definition: check_params.h:126
double tol_restart
Definition: quda.h:103
void updateGaugeFieldQuda(void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
enum QudaGammaBasis_s QudaGammaBasis
enum QudaReconstructType_s QudaReconstructType
#define QUDA_MAX_DWF_LS
Maximum length of the Ls dimension for domain-wall fermions.
int mom_ga_pad
Definition: quda.h:59
void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity *parity, int inverse)
void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
int use_resident_gauge
Definition: quda.h:69
QudaTboundary t_boundary
Definition: quda.h:38
enum QudaTune_s QudaTune
void saveGaugeField(void *outGauge, void *inGauge, QudaGaugeParam *param)
QudaTwistFlavorType twist_flavor
Definition: quda.h:100
int use_resident_mom
Definition: quda.h:70
enum QudaDslashType_s QudaDslashType
int f_size
Definition: quda.h:271
int device
Definition: test_util.cpp:1546
QudaResidualType residual_type
Definition: quda.h:235
int num_offset
Definition: quda.h:123
int(* QudaCommsMap)(const int *coords, void *fdata)
Definition: quda.h:320
enum QudaVerbosity_s QudaVerbosity
double cloverGiB
Definition: quda.h:181
double epsilon
Definition: quda.h:98
void set_dim(int *)
double omega
Definition: quda.h:219
int make_resident_gauge
Definition: quda.h:71
void computeStaggeredForceQuda(void *mom, void *quark, double *coeff)
QudaDslashType dslash_type_precondition
Definition: quda.h:208
QudaPrecision clover_cpu_prec
Definition: quda.h:161
enum QudaComputeFatMethod_s QudaComputeFatMethod
enum QudaUseInitGuess_s QudaUseInitGuess
enum QudaSolverNormalization_s QudaSolverNormalization
void initQudaDevice(int device)
void pack_ghost(void **cpuLink, void **cpuGhost, int nFace, QudaPrecision precision)
const QudaParity parity
Definition: dslash_test.cpp:29
void * gauge[4]
Definition: su3_test.cpp:15
QudaMatPCType matpc_type
Definition: quda.h:144
void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaEigParam newQudaEigParam(void)
enum QudaInverterType_s QudaInverterType
QudaPrecision cpu_prec
Definition: quda.h:40
int overlap
Definition: quda.h:67
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:149
double clover_coeff
Definition: quda.h:169
void * fatlink[4]
enum QudaTwistFlavorType_s QudaTwistFlavorType
enum QudaDiracFieldOrder_s QudaDiracFieldOrder