QUDA  0.9.0
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;
65  double i_mu;
67  int overlap;
78  size_t gauge_offset;
79  size_t mom_offset;
80  size_t site_size;
83 
84 
88  typedef struct QudaInvertParam_s {
89 
96  double mass;
97  double kappa;
99  double m5;
100  int Ls;
105  double mu;
106  double epsilon;
110  double tol;
111  double tol_restart;
112  double tol_hq;
115  double true_res;
116  double true_res_hq;
117  int maxiter;
118  double reliable_delta;
130 
135 
140 
143 
144  int pipeline;
148  int num_src;
150  int overlap;
154 
157 
160 
163 
166 
169 
172 
175 
179  double action[2];
180 
208  double clover_coeff;
209  double clover_rho;
212  double trlogA[2];
221  int sp_pad;
222  int cl_pad;
224  int iter;
225  double spinorGiB;
226  double cloverGiB;
227  double gflops;
228  double secs;
234  int Nsteps;
235 
238 
239  /*
240  * The following parameters are related to the solver
241  * preconditioner, if enabled.
242  */
243 
249 
252 
255 
262 
265 
268 
270  double omega;
271 
274 
277 
287 
294  int nev;
300  int rhs_idx;
304  double eigenval_tol;
310  double inc_tol;
311 
314 
317 
320 
323 
326 
329 
332 
333  } QudaInvertParam;
334 
335 
336  // Parameter set for solving the eigenvalue problems.
337  // Eigen problems are tightly related with Ritz algorithm.
338  // And the Lanczos algorithm use the Ritz operator.
339  // For Ritz matrix operation,
340  // we need to know about the solution type of dirac operator.
341  // For acceleration, we are also using chevisov polynomial method.
342  // And nk, np values are needed Implicit Restart Lanczos method
343  // which is optimized form of Lanczos algorithm
344  typedef struct QudaEigParam_s {
345 
347 //specific for Lanczos method:
351 
352  double *MatPoly_param;
353  int NPoly;
354  double Stp_residual;
355  int nk;
356  int np;
357  int f_size;
358  double eigen_shift;
359 //more general stuff:
362 
365 
368 
371 
374 
376  char vec_infile[256];
377 
379  char vec_outfile[256];
380 
382  double gflops;
383 
385  double secs;
386 
389 
390  } QudaEigParam;
391 
392 
393  typedef struct QudaMultigridParam_s {
394 
396 
398  int n_level;
399 
402 
405 
408 
411 
414 
417 
420 
424 
427 
430 
433 
436 
439 
442 
445 
448 
451 
454 
457 
459  char vec_infile[256];
460 
462  char vec_outfile[256];
463 
465  double gflops;
466 
468  double secs;
469 
472 
474 
475 
476 
477  /*
478  * Interface functions, found in interface_quda.cpp
479  */
480 
506  void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[],
507  FILE *outfile);
508 
519  typedef int (*QudaCommsMap)(const int *coords, void *fdata);
520 
547  void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata);
548 
559  void initQudaDevice(int device);
560 
567  void initQudaMemory();
568 
578  void initQuda(int device);
579 
583  void endQuda(void);
584 
590  void updateR();
591 
600 
609 
618 
627 
633 
639 
645 
651 
657  void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param);
658 
662  void freeGaugeQuda(void);
663 
669  void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param);
670 
678  void loadCloverQuda(void *h_clover, void *h_clovinv,
680 
684  void freeCloverQuda(void);
685 
695  void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V,
696  void *hp_alpha, void *hp_beta, QudaEigParam *eig_param);
697 
707  void invertQuda(void *h_x, void *h_b, QudaInvertParam *param);
708 
718  void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param);
719 
720 
728  void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param);
729 
738 
743  void destroyMultigridQuda(void *mg_instance);
744 
749  void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param);
750 
759  void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
761 
771  void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param,
773 
783  void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param,
785 
795  void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
796  QudaParity *parity, int inverse);
797 
805  void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
806 
814  void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
815 
816 
817  /*
818  * The following routines are temporary additions used by the HISQ
819  * link-fattening code.
820  */
821 
822  void set_dim(int *);
823  void pack_ghost(void **cpuLink, void **cpuGhost, int nFace,
824  QudaPrecision precision);
825 
826  void computeKSLinkQuda(void* fatlink, void* longlink, void* ulink, void* inlink,
827  double *path_coeff, QudaGaugeParam *param);
828 
829 
830 
844  int computeGaugeForceQuda(void* mom, void* sitelink, int*** input_path_buf, int* path_length,
845  double* loop_coeff, int num_paths, int max_length, double dt,
847 
859  void updateGaugeFieldQuda(void* gauge, void* momentum, double dt,
860  int conj_mom, int exact, QudaGaugeParam* param);
861 
871  void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param);
872 
881  void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param);
882 
891  double momActionQuda(void* momentum, QudaGaugeParam* param);
892 
901  void* createGaugeFieldQuda(void* gauge, int geometry, QudaGaugeParam* param);
902 
910  void saveGaugeFieldQuda(void* outGauge, void* inGauge, QudaGaugeParam* param);
911 
917  void destroyGaugeFieldQuda(void* gauge);
918 
925 
944  void computeCloverForceQuda(void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck,
945  int nvector, double multiplicity, void *gauge,
947 
959  void computeStaggeredForceQuda(void* mom, double dt, double delta, void **x, void *gauge,
960  QudaGaugeParam *gauge_param, QudaInvertParam *invert_param);
961 
976  void computeHISQForceQuda(void* momentum,
977  long long* flops,
978  const double level2_coeff[6],
979  const double fat7_coeff[6],
980  const void* const w_link,
981  const void* const v_link,
982  const void* const u_link,
983  void** quark,
984  int num,
985  int num_naik,
986  double** coeff,
988 
993  void gaussGaugeQuda(long seed);
994 
999  void plaqQuda(double plaq[3]);
1000 
1011  void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param,
1012  unsigned int nSteps, double alpha);
1013 
1019  void performAPEnStep(unsigned int nSteps, double alpha);
1020 
1026  void performSTOUTnStep(unsigned int nSteps, double rho);
1027 
1034  void performOvrImpSTOUTnStep(unsigned int nSteps, double rho, double epsilon);
1035 
1039  double qChargeCuda();
1040 
1054  int computeGaugeFixingOVRQuda(void* gauge,
1055  const unsigned int gauge_dir,
1056  const unsigned int Nsteps,
1057  const unsigned int verbose_interval,
1058  const double relax_boost,
1059  const double tolerance,
1060  const unsigned int reunit_interval,
1061  const unsigned int stopWtheta,
1063  double* timeinfo);
1077  int computeGaugeFixingFFTQuda(void* gauge,
1078  const unsigned int gauge_dir,
1079  const unsigned int Nsteps,
1080  const unsigned int verbose_interval,
1081  const double alpha,
1082  const unsigned int autotune,
1083  const double tolerance,
1084  const unsigned int stopWtheta,
1086  double* timeinfo);
1087 
1092  void flushChronoQuda(int index);
1093 
1094 
1099  void openMagma();
1100 
1101  void closeMagma();
1102 
1109 
1113  void destroyDeflationQuda(void *df_instance);
1114 
1115 #ifdef __cplusplus
1116 }
1117 #endif
1118 
1119 /* #include <quda_new_interface.h> */
1120 
1121 #endif /* _QUDA_H */
void computeCloverForceQuda(void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck, int nvector, double multiplicity, void *gauge, QudaGaugeParam *gauge_param, QudaInvertParam *inv_param)
double action[2]
Definition: quda.h:179
int maxiter_precondition
Definition: quda.h:267
static QudaGaugeParam qudaGaugeParam
double secs
Definition: quda.h:228
void performOvrImpSTOUTnStep(unsigned int nSteps, double rho, double epsilon)
double iter_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:165
enum QudaPreserveSource_s QudaPreserveSource
QudaDiracFieldOrder dirac_order
Definition: quda.h:195
QudaMassNormalization mass_normalization
Definition: quda.h:185
double tol_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:159
enum QudaMassNormalization_s QudaMassNormalization
QudaReconstructType reconstruct_sloppy
Definition: quda.h:46
double anisotropy
Definition: quda.h:31
void freeCloverQuda(void)
double b_5[QUDA_MAX_DWF_LS]
Definition: quda.h:102
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:199
QudaSolveType solve_type
Definition: quda.h:182
QudaVerbosity verbosity_precondition
Definition: quda.h:261
QudaVerbosity verbosity
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:53
void destroyDeflationQuda(void *df_instance)
QudaMultigridParam newQudaMultigridParam(void)
int make_resident_mom
Definition: quda.h:74
int NPoly
Definition: quda.h:353
QudaExtLibType extlib_type
Definition: quda.h:331
size_t gauge_offset
Definition: quda.h:78
double mu
Definition: quda.h:105
QudaGaugeFixed gauge_fix
Definition: quda.h:51
QudaTune tune
Definition: quda.h:230
QudaSchwarzType schwarz_type
Definition: quda.h:276
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:248
void printQudaGaugeParam(QudaGaugeParam *param)
Definition: check_params.h:40
const void * func
QudaLinkType type
Definition: quda.h:35
int max_res_increase_total
Definition: quda.h:139
double kappa
Definition: quda.h:97
void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
QudaPrecision cuda_prec_ritz
Definition: quda.h:290
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
double tol
Definition: quda.h:110
QudaDslashType dslash_type
Definition: quda.h:93
QudaReconstructType reconstruct_precondition
Definition: quda.h:49
QudaInverterType inv_type
Definition: quda.h:94
QudaPrecision cuda_prec
Definition: quda.h:191
double c_5[QUDA_MAX_DWF_LS]
Definition: quda.h:103
int return_clover_inverse
Definition: quda.h:217
QudaSolveType smoother_solve_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:426
void performAPEnStep(unsigned int nSteps, double alpha)
enum QudaSolveType_s QudaSolveType
double setup_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:416
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
double Stp_residual
Definition: quda.h:354
QudaPrecision cpu_prec
Definition: quda.h:190
double momActionQuda(void *momentum, QudaGaugeParam *param)
QudaStaggeredPhase staggered_phase_type
Definition: quda.h:62
void destroyMultigridQuda(void *mg_instance)
Free resources allocated by the multigrid solver.
void plaqQuda(double plaq[3])
double trlogA[2]
Definition: quda.h:212
QudaDagType dagger
Definition: quda.h:184
QudaGaugeParam gauge_param
enum QudaComputeNullVector_s QudaComputeNullVector
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
double true_res
Definition: quda.h:115
void computeKSLinkQuda(void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param)
size_t mom_offset
Definition: quda.h:79
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
int return_clover
Definition: quda.h:216
enum QudaEigType_s QudaEigType
enum QudaTboundary_s QudaTboundary
int make_resident_solution
Definition: quda.h:313
void openMagma()
int overwrite_mom
Definition: quda.h:69
double mu_factor[QUDA_MAX_MG_LEVEL]
Definition: quda.h:471
int compute_clover_trlog
Definition: quda.h:211
int chrono_index
Definition: quda.h:328
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:202
struct QudaInvertParam_s QudaInvertParam
int compute_action
Definition: quda.h:174
void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
QudaFieldLocation input_location
Definition: quda.h:90
void freeGaugeQuda(void)
void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata)
int staggered_phase_applied
Definition: quda.h:63
double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:168
double reliable_delta
Definition: quda.h:118
size_t site_size
Definition: quda.h:80
int solution_accumulator_pipeline
Definition: quda.h:129
QudaUseInitGuess use_init_guess
Definition: quda.h:206
int n_vec[QUDA_MAX_MG_LEVEL]
Definition: quda.h:407
int computeGaugeFixingOVRQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double relax_boost, const double tolerance, const unsigned int reunit_interval, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
Gauge fixing with overrelaxation with support for single and multi GPU.
char * index(const char *, int)
QudaGaugeParam param
Definition: pack_test.cpp:17
void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param, unsigned int nSteps, double alpha)
int llfat_ga_pad
Definition: quda.h:58
QudaSolutionType solution_type
Definition: quda.h:181
void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param)
QudaMemoryType mem_type_ritz
Definition: quda.h:367
QudaSolverNormalization solver_normalization
Definition: quda.h:186
QudaEigType eig_type
Definition: quda.h:350
void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaPrecision clover_cuda_prec
Definition: quda.h:201
int precondition_cycle
Definition: quda.h:273
void * preconditioner
Definition: quda.h:251
void * longlink[4]
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
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)
QudaBoolean global_reduction[QUDA_MAX_MG_LEVEL]
Definition: quda.h:444
double * MatPoly_param
Definition: quda.h:352
QudaInvertParam * invert_param
Definition: quda.h:346
double scale
Definition: quda.h:33
void initQuda(int device)
double spinorGiB
Definition: quda.h:225
void initQudaMemory()
double secs
Definition: quda.h:385
double tol
Definition: test_util.cpp:1647
static unsigned int delta
QudaFieldLocation output_location
Definition: quda.h:91
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:203
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
int site_ga_pad
Definition: quda.h:55
QudaMultigridCycleType cycle_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:429
void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param)
QudaBoolean run_verify
Definition: quda.h:373
double i_mu
Definition: quda.h:65
double m5
Definition: quda.h:99
enum QudaStaggeredPhase_s QudaStaggeredPhase
void createCloverQuda(QudaInvertParam *param)
void * newDeflationQuda(QudaEigParam *param)
QudaPrecision cuda_prec_sloppy
Definition: quda.h:192
QudaVerbosity verbosity
Definition: quda.h:219
double tol_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:156
double true_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:162
QudaInvertParam newQudaInvertParam(void)
double eigen_shift
Definition: quda.h:358
double gflops
Definition: quda.h:227
int eigcg_max_restarts
Definition: quda.h:306
int use_resident_chrono
Definition: quda.h:322
double gflops
Definition: quda.h:382
QudaPrecision cuda_prec_precondition
Definition: quda.h:48
int staple_pad
Definition: quda.h:57
QudaCloverFieldOrder clover_order
Definition: quda.h:205
static __inline__ size_t p
double omega[QUDA_MAX_MG_LEVEL]
Definition: quda.h:441
void saveGaugeFieldQuda(void *outGauge, void *inGauge, QudaGaugeParam *param)
double tol_hq
Definition: quda.h:112
enum QudaMatPCType_s QudaMatPCType
struct QudaEigParam_s QudaEigParam
void * newMultigridQuda(QudaMultigridParam *param)
double true_res_hq
Definition: quda.h:116
double smoother_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:438
enum QudaGaugeFixed_s QudaGaugeFixed
enum QudaSolutionType_s QudaSolutionType
QudaGammaBasis gamma_basis
Definition: quda.h:197
enum QudaSchwarzType_s QudaSchwarzType
QudaPrecision cuda_prec_sloppy
Definition: quda.h:45
int max_search_dim
Definition: quda.h:298
QudaVerbosity verbosity[QUDA_MAX_MG_LEVEL]
Definition: quda.h:410
double tol_precondition
Definition: quda.h:264
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
double qChargeCuda()
double offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:153
int use_sloppy_partial_accumulator
Definition: quda.h:119
enum QudaDagType_s QudaDagType
QudaBoolean run_verify
Definition: quda.h:456
QudaSolutionType RitzMat_Convcheck
Definition: quda.h:349
QudaSolutionType RitzMat_lanczos
Definition: quda.h:348
int heavy_quark_check
Definition: quda.h:142
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
void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param)
Updates the multigrid preconditioner for the new gauge / clover field.
double mass
Definition: quda.h:96
QudaBoolean import_vectors
Definition: quda.h:361
char vec_outfile[256]
Definition: quda.h:462
QudaFieldLocation location
Definition: quda.h:370
int gcrNkrylov
Definition: quda.h:237
void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta, QudaEigParam *eig_param)
int make_resident_chrono
Definition: quda.h:319
int compute_clover_inverse
Definition: quda.h:215
void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param)
enum QudaBoolean_s QudaBoolean
QudaFieldLocation location[QUDA_MAX_MG_LEVEL]
Definition: quda.h:447
void performSTOUTnStep(unsigned int nSteps, double rho)
void printQudaInvertParam(QudaInvertParam *param)
Definition: check_params.h:191
int max_res_increase
Definition: quda.h:134
int computeGaugeFixingFFTQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double alpha, const unsigned int autotune, const double tolerance, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
Gauge fixing with Steepest descent method with FFTs with support for single GPU only.
char vec_outfile[256]
Definition: quda.h:379
void destroyGaugeFieldQuda(void *gauge)
enum QudaFieldLocation_s QudaFieldLocation
void closeMagma()
QudaFieldLocation location
Definition: quda.h:27
int spin_block_size[QUDA_MAX_MG_LEVEL]
Definition: quda.h:404
struct QudaGaugeParam_s QudaGaugeParam
double clover_rho
Definition: quda.h:209
double tadpole_coeff
Definition: quda.h:32
double gaugeGiB
Definition: quda.h:60
QudaPrecision cuda_prec_precondition
Definition: quda.h:193
QudaSolutionType coarse_grid_solution_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:423
int deflation_grid
Definition: quda.h:302
void printQudaEigParam(QudaEigParam *param)
Definition: check_params.h:145
double tol_restart
Definition: quda.h:111
void updateGaugeFieldQuda(void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
void * fatlink[4]
int geo_block_size[QUDA_MAX_MG_LEVEL][QUDA_MAX_DIM]
Definition: quda.h:401
enum QudaGammaBasis_s QudaGammaBasis
enum QudaMultigridCycleType_s QudaMultigridCycleType
QudaBoolean generate_all_levels
Definition: quda.h:453
enum QudaReconstructType_s QudaReconstructType
void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param)
int nu_pre[QUDA_MAX_MG_LEVEL]
Definition: quda.h:432
#define QUDA_MAX_MG_LEVEL
Maximum number of multi-grid levels. This number may be increased if needed.
struct QudaMultigridParam_s QudaMultigridParam
#define QUDA_MAX_DWF_LS
Maximum length of the Ls dimension for domain-wall fermions.
void printQudaMultigridParam(QudaMultigridParam *param)
Definition: check_params.h:504
int mom_ga_pad
Definition: quda.h:59
QudaExtLibType extlib_type
Definition: quda.h:388
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:71
QudaTboundary t_boundary
Definition: quda.h:38
enum QudaTune_s QudaTune
QudaTwistFlavorType twist_flavor
Definition: quda.h:108
int return_result_gauge
Definition: quda.h:75
double residue[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:171
unsigned long long flops
Definition: blas_quda.cu:42
int max_restart_num
Definition: quda.h:308
void * createGaugeFieldQuda(void *gauge, int geometry, QudaGaugeParam *param)
int use_resident_mom
Definition: quda.h:72
enum QudaDslashType_s QudaDslashType
int(* QudaCommsMap)(const int *coords, void *fdata)
Definition: quda.h:519
int f_size
Definition: quda.h:357
void gaussGaugeQuda(long seed)
int compute_true_res
Definition: quda.h:114
QudaResidualType residual_type
Definition: quda.h:286
double inc_tol
Definition: quda.h:310
int num_offset
Definition: quda.h:146
enum QudaVerbosity_s QudaVerbosity
double cloverGiB
Definition: quda.h:226
void computeStaggeredForceQuda(void *mom, double dt, double delta, void **x, void *gauge, QudaGaugeParam *gauge_param, QudaInvertParam *invert_param)
cpuGaugeField * cpuLink
Definition: covdev_test.cpp:39
int return_result_mom
Definition: quda.h:76
int test_type
Definition: test_util.cpp:1634
int compute_clover
Definition: quda.h:214
double epsilon
Definition: quda.h:106
#define QUDA_MAX_DIM
Maximum number of dimensions supported by QUDA. In practice, no routines make use of more than 5...
void flushChronoQuda(int index)
Flush the chronological history for the given index.
int nu_post[QUDA_MAX_MG_LEVEL]
Definition: quda.h:435
void set_dim(int *)
int max_chrono_dim
Definition: quda.h:325
int use_resident_solution
Definition: quda.h:316
double omega
Definition: quda.h:270
void computeHISQForceQuda(void *momentum, long long *flops, const double level2_coeff[6], const double fat7_coeff[6], const void *const w_link, const void *const v_link, const void *const u_link, void **quark, int num, int num_naik, double **coeff, QudaGaugeParam *param)
QudaComputeNullVector compute_null_vector
Definition: quda.h:450
int make_resident_gauge
Definition: quda.h:73
void * deflation_op
Definition: quda.h:254
double eigenval_tol
Definition: quda.h:304
QudaDslashType dslash_type_precondition
Definition: quda.h:259
QudaPrecision clover_cpu_prec
Definition: quda.h:200
QudaPrecision cuda_prec_ritz
Definition: quda.h:364
QudaParity parity
Definition: covdev_test.cpp:53
enum QudaUseInitGuess_s QudaUseInitGuess
char vec_infile[256]
Definition: quda.h:459
enum QudaSolverNormalization_s QudaSolverNormalization
void initQudaDevice(int device)
void pack_ghost(void **cpuLink, void **cpuGhost, int nFace, QudaPrecision precision)
QudaInverterType setup_inv_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:413
QudaInverterType smoother[QUDA_MAX_MG_LEVEL]
Definition: quda.h:419
QudaMatPCType matpc_type
Definition: quda.h:183
void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaEigParam newQudaEigParam(void)
char vec_infile[256]
Definition: quda.h:376
enum QudaInverterType_s QudaInverterType
enum QudaMemoryType_s QudaMemoryType
QudaPrecision cpu_prec
Definition: quda.h:40
enum QudaExtLibType_s QudaExtLibType
int overlap
Definition: quda.h:67
void updateR()
update the radius for halos.
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:188
QudaInvertParam * invert_param
Definition: quda.h:395
double clover_coeff
Definition: quda.h:208
enum QudaTwistFlavorType_s QudaTwistFlavorType
enum cudaDeviceAttr attr int device
enum QudaDiracFieldOrder_s QudaDiracFieldOrder