QUDA  0.9.0
multigrid_invert_test.cpp
Go to the documentation of this file.
1 #include <stdlib.h>
2 #include <stdio.h>
3 #include <time.h>
4 #include <math.h>
5 #include <string.h>
6 
7 #include <util_quda.h>
8 #include <test_util.h>
9 #include <dslash_util.h>
10 #include <blas_reference.h>
13 #include "misc.h"
14 
15 #if defined(QMP_COMMS)
16 #include <qmp.h>
17 #elif defined(MPI_COMMS)
18 #include <mpi.h>
19 #endif
20 
21 #include <qio_field.h>
22 
23 #define MAX(a,b) ((a)>(b)?(a):(b))
24 
25 // In a typical application, quda.h is the only QUDA header required.
26 #include <quda.h>
27 
28 // Wilson, clover-improved Wilson, twisted mass, and domain wall are supported.
30 extern int device;
31 extern int xdim;
32 extern int ydim;
33 extern int zdim;
34 extern int tdim;
35 extern int Lsdim;
36 extern int gridsize_from_cmdline[];
38 extern QudaPrecision prec;
43 extern double mass;
44 extern double mu;
45 extern double anisotropy;
46 extern double tol; // tolerance for inverter
47 extern double tol_hq; // heavy-quark tolerance for inverter
48 extern char latfile[];
49 extern int niter;
50 extern int gcrNkrylov; // number of inner iterations for GCR, or l for BiCGstab-l
51 extern int pipeline; // length of pipeline for fused operations in GCR or BiCGstab-l
52 extern int nvec[];
53 extern int mg_levels;
54 
55 extern bool generate_nullspace;
56 extern bool generate_all_levels;
57 extern int nu_pre;
58 extern int nu_post;
60 extern double mu_factor[QUDA_MAX_MG_LEVEL];
61 
63 
65 extern double setup_tol;
66 extern double omega;
68 
71 
72 extern char vec_infile[];
73 extern char vec_outfile[];
74 
75 //Twisted mass flavor type
77 
78 extern void usage(char** );
79 
80 extern double clover_coeff;
81 extern bool compute_clover;
82 
83 namespace quda {
84  extern void setTransferGPU(bool);
85 }
86 
87 void
89 {
90  printfQuda("running the following test:\n");
91 
92  printfQuda("prec sloppy_prec link_recon sloppy_link_recon S_dimension T_dimension Ls_dimension\n");
93  printfQuda("%s %s %s %s %d/%d/%d %d %d\n",
97 
98  printfQuda("MG parameters\n");
99  printfQuda(" - number of levels %d\n", mg_levels);
100  for (int i=0; i<mg_levels-1; i++) printfQuda(" - level %d number of null-space vectors %d\n", i+1, nvec[i]);
101  printfQuda(" - number of pre-smoother applications %d\n", nu_pre);
102  printfQuda(" - number of post-smoother applications %d\n", nu_post);
103 
104  printfQuda("Grid partition info: X Y Z T\n");
105  printfQuda(" %d %d %d %d\n",
106  dimPartitioned(0),
107  dimPartitioned(1),
108  dimPartitioned(2),
109  dimPartitioned(3));
110 
111  return ;
112 
113 }
114 
119 
121  gauge_param.X[0] = xdim;
122  gauge_param.X[1] = ydim;
123  gauge_param.X[2] = zdim;
124  gauge_param.X[3] = tdim;
125 
130 
132 
135 
138 
141 
143 
144  gauge_param.ga_pad = 0;
145  // For multi-GPU, ga_pad must be large enough to store a time-slice
146 #ifdef MULTI_GPU
147  int x_face_size = gauge_param.X[1]*gauge_param.X[2]*gauge_param.X[3]/2;
148  int y_face_size = gauge_param.X[0]*gauge_param.X[2]*gauge_param.X[3]/2;
149  int z_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[3]/2;
150  int t_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[2]/2;
151  int pad_size =MAX(x_face_size, y_face_size);
152  pad_size = MAX(pad_size, z_face_size);
153  pad_size = MAX(pad_size, t_face_size);
154  gauge_param.ga_pad = pad_size;
155 #endif
156 }
157 
160 
161  inv_param.Ls = 1;
162 
163  inv_param.sp_pad = 0;
164  inv_param.cl_pad = 0;
165 
173 
181  }
182 
185 
187 
188  //Free field!
189  inv_param.mass = mass;
190  inv_param.kappa = 1.0 / (2.0 * (1 + 3/anisotropy + mass));
191 
193  inv_param.mu = mu;
196 
198  printfQuda("Twisted-mass doublet non supported (yet)\n");
199  exit(0);
200  }
201  }
202 
205 
208 
210 
211  mg_param.invert_param = &inv_param;
212  mg_param.n_level = mg_levels;
213  for (int i=0; i<mg_param.n_level; i++) {
214  for (int j=0; j<QUDA_MAX_DIM; j++) {
215  // if not defined use 4
216  mg_param.geo_block_size[i][j] = geo_block_size[i][j] ? geo_block_size[i][j] : 4;
217  }
218  mg_param.verbosity[i] = mg_verbosity[i];
219  mg_param.setup_inv_type[i] = setup_inv[i];
220  mg_param.setup_tol[i] = setup_tol;
221  mg_param.spin_block_size[i] = 1;
222  mg_param.n_vec[i] = nvec[i] == 0 ? 24 : nvec[i]; // default to 24 vectors if not set
223  mg_param.nu_pre[i] = nu_pre;
224  mg_param.nu_post[i] = nu_post;
225  mg_param.mu_factor[i] = mu_factor[i];
226 
228 
229  mg_param.smoother[i] = smoother_type;
230 
231  // set the smoother / bottom solver tolerance (for MR smoothing this will be ignored)
232  mg_param.smoother_tol[i] = tol_hq; // repurpose heavy-quark tolerance for now
233 
234  mg_param.global_reduction[i] = QUDA_BOOLEAN_YES;
235 
236  // set to QUDA_DIRECT_SOLVE for no even/odd preconditioning on the smoother
237  // set to QUDA_DIRECT_PC_SOLVE for to enable even/odd preconditioning on the smoother
238  mg_param.smoother_solve_type[i] = QUDA_DIRECT_PC_SOLVE; // EVEN-ODD
239 
240  // set to QUDA_MAT_SOLUTION to inject a full field into coarse grid
241  // set to QUDA_MATPC_SOLUTION to inject single parity field into coarse grid
242 
243  // if we are using an outer even-odd preconditioned solve, then we
244  // use single parity injection into the coarse grid
246 
247  mg_param.omega[i] = omega; // over/under relaxation factor
248 
249  mg_param.location[i] = QUDA_CUDA_FIELD_LOCATION;
250  }
251 
252  // only coarsen the spin on the first restriction
253  mg_param.spin_block_size[0] = 2;
254 
255  // coarse grid solver is GCR
256  mg_param.smoother[mg_levels-1] = QUDA_GCR_INVERTER;
257 
260 
262 
263  mg_param.run_verify = QUDA_BOOLEAN_YES;
264 
265  // set file i/o parameters
266  strcpy(mg_param.vec_infile, vec_infile);
267  strcpy(mg_param.vec_outfile, vec_outfile);
268 
269  // these need to tbe set for now but are actually ignored by the MG setup
270  // needed to make it pass the initialization test
272  inv_param.tol = 1e-10;
273  inv_param.maxiter = 1000;
274  inv_param.reliable_delta = 1e-10;
275  inv_param.gcrNkrylov = 10;
276 
279 }
280 
282  inv_param.Ls = 1;
283 
284  inv_param.sp_pad = 0;
285  inv_param.cl_pad = 0;
286 
290 
295 
302  }
303 
306 
308 
309  //Free field!
310  inv_param.mass = mass;
311  inv_param.kappa = 1.0 / (2.0 * (1 + 3/anisotropy + mass));
312 
314  inv_param.mu = mu;
317 
319  printfQuda("Twisted-mass doublet non supported (yet)\n");
320  exit(0);
321  }
322  }
323 
325 
328 
329  // do we want full solution or single-parity solution
331 
332  // do we want to use an even-odd preconditioned solve or not
335 
337 
340 
341 
345  inv_param.tol = tol;
346 
347  // require both L2 relative and heavy quark residual to determine convergence
349  inv_param.tol_hq = tol_hq; // specify a tolerance for the residual for heavy quark residual
350 
351  // these can be set individually
352  for (int i=0; i<inv_param.num_offset; i++) {
355  }
358 
359  // domain decomposition preconditioner parameters
364  inv_param.omega = 1.0;
365 }
366 
367 int main(int argc, char **argv)
368 {
369  // We give here the default value to some of the array
370  for(int i =0; i<QUDA_MAX_MG_LEVEL; i++) {
373  mu_factor[i] = 1.;
374  }
375 
376  for (int i = 1; i < argc; i++){
377  if(process_command_line_option(argc, argv, &i) == 0){
378  continue;
379  }
380  printf("ERROR: Invalid option:%s\n", argv[i]);
381  usage(argv);
382  }
383 
388 
389  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
390  initComms(argc, argv, gridsize_from_cmdline);
391 
392  // call srand() with a rank-dependent seed
393  initRand();
394 
396 
397  // *** QUDA parameters begin here.
398 
403  printfQuda("dslash_type %d not supported\n", dslash_type);
404  exit(0);
405  }
406 
409 
410  QudaInvertParam mg_inv_param = newQudaInvertParam();
412  mg_param.invert_param = &mg_inv_param;
413 
414  setMultigridParam(mg_param);
415 
416 
419 
420  // *** Everything between here and the call to initQuda() is
421  // *** application-specific.
422 
424 
425  setSpinorSiteSize(24);
426 
427  size_t gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
428  size_t sSize = (inv_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
429 
430  void *gauge[4], *clover=0, *clover_inv=0;
431 
432  for (int dir = 0; dir < 4; dir++) {
433  gauge[dir] = malloc(V*gaugeSiteSize*gSize);
434  }
435 
436  if (strcmp(latfile,"")) { // load in the command line supplied gauge field
439  } else { // else generate a random SU(3) field
440  //generate a random SU(3) field
441  //construct_gauge_field(gauge, 1, gauge_param.cpu_prec, &gauge_param);
442  //generate a unit SU(3) field
444 
445  }
446 
448  double norm = 0.1; // clover components are random numbers in the range (-norm, norm)
449  double diag = 1.0; // constant added to the diagonal
450 
451  size_t cSize = inv_param.clover_cpu_prec;
452  clover = malloc(V*cloverSiteSize*cSize);
453  clover_inv = malloc(V*cloverSiteSize*cSize);
455 
460  }
461 
462  void *spinorIn = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
463  void *spinorCheck = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
464 
465  void *spinorOut = NULL;
467 
468  // start the timer
469  double time0 = -((double)clock());
470 
471  // initialize the QUDA library
472  initQuda(device);
473 
474  // load the gauge field
475  loadGaugeQuda((void*)gauge, &gauge_param);
476 
477  // this line ensure that if we need to construct the clover inverse (in either the smoother or the solver) we do so
480 
481  inv_param.solve_type = solve_type; // restore actual solve_type we want to do
482 
483  // setup the multigrid solver
484  void *mg_preconditioner = newMultigridQuda(&mg_param);
485  inv_param.preconditioner = mg_preconditioner;
486 
487  const int nSrc = 1;
488  for (int i=0; i<nSrc; i++) {
489  // create a point source at 0 (in each subvolume... FIXME)
490  memset(spinorIn, 0, inv_param.Ls*V*spinorSiteSize*sSize);
491  memset(spinorCheck, 0, inv_param.Ls*V*spinorSiteSize*sSize);
493 
495  //((float*)spinorIn)[i] = 1.0;
496  for (int i=0; i<inv_param.Ls*V*spinorSiteSize; i++) ((float*)spinorIn)[i] = rand() / (float)RAND_MAX;
497  } else {
498  //((double*)spinorIn)[i] = 1.0;
499  for (int i=0; i<inv_param.Ls*V*spinorSiteSize; i++) ((double*)spinorIn)[i] = rand() / (double)RAND_MAX;
500  }
501 
502  invertQuda(spinorOut, spinorIn, &inv_param);
503  }
504 
505  // free the multigrid solver
506  destroyMultigridQuda(mg_preconditioner);
507 
508 
509 
510  // stop the timer
511  time0 += clock();
512  time0 /= CLOCKS_PER_SEC;
513 
514  printfQuda("Device memory used:\n Spinor: %f GiB\n Gauge: %f GiB\n",
517  //printfQuda("\nDone: %i iter / %g secs = %g Gflops, total time = %g secs\n",
518  //inv_param.iter, inv_param.secs, inv_param.gflops/inv_param.secs, time0);
519  printfQuda("\nDone: %i iter / %g secs = %g Gflops, total time = %g secs\n",
521 
523 
525  wil_mat(spinorCheck, gauge, spinorOut, inv_param.kappa, 0, inv_param.cpu_prec, gauge_param);
526  } else {
530  } else {
531  printfQuda("Unsupported dslash_type\n");
532  exit(-1);
533  }
534  }
535  }
537  ax(0.5/inv_param.kappa, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
538  }
539 
541 
543  wil_matpc(spinorCheck, gauge, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
545  } else {
550  } else {
551  printfQuda("Unsupported dslash_type\n");
552  exit(-1);
553  }
554  }
555  }
556 
559  }
560 
561  }
562 
563  int vol = inv_param.solution_type == QUDA_MAT_SOLUTION ? V : Vh;
564  mxpy(spinorIn, spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
565  double nrm2 = norm_2(spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
566  double src2 = norm_2(spinorIn, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
567  double l2r = sqrt(nrm2 / src2);
568 
569  printfQuda("Residuals: (L2 relative) tol %g, QUDA = %g, host = %g; (heavy-quark) tol %g, QUDA = %g\n",
571 
572 
573  freeGaugeQuda();
575 
576  // finalize the QUDA library
577  endQuda();
578 
579  // finalize the communications layer
580  finalizeComms();
581 
583  if (clover) free(clover);
584  if (clover_inv) free(clover_inv);
585  }
586 
587  for (int dir = 0; dir<4; dir++) free(gauge[dir]);
588 
589  return 0;
590 }
int maxiter_precondition
Definition: quda.h:267
QudaInverterType setup_inv[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1661
double secs
Definition: quda.h:228
int dimPartitioned(int dim)
Definition: test_util.cpp:1686
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
QudaReconstructType reconstruct_sloppy
Definition: quda.h:46
double anisotropy
Definition: quda.h:31
void freeCloverQuda(void)
int mg_levels
Definition: test_util.cpp:1655
void display_test_info()
void endQuda(void)
void free(void *)
void construct_gauge_field(void **gauge, int type, QudaPrecision precision, QudaGaugeParam *param)
Definition: test_util.cpp:1054
QudaPrecision prec_precondition
Definition: test_util.cpp:1617
QudaSolveType solve_type
Definition: quda.h:182
QudaVerbosity verbosity_precondition
Definition: quda.h:261
int ydim
Definition: test_util.cpp:1621
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:53
QudaMultigridParam newQudaMultigridParam(void)
QudaPrecision & cuda_prec
int nvec[]
Definition: test_util.cpp:1635
char vec_outfile[]
Definition: test_util.cpp:1637
double mu
Definition: quda.h:105
QudaGaugeFixed gauge_fix
Definition: quda.h:51
QudaSchwarzType schwarz_type
Definition: quda.h:276
__host__ __device__ ValueType norm(const complex< ValueType > &z)
Returns the magnitude of z squared.
Definition: complex_quda.h:896
QudaReconstructType link_recon_precondition
Definition: test_util.cpp:1614
double anisotropy
Definition: test_util.cpp:1644
void tm_mat(void *out, void **gauge, void *in, double kappa, double mu, QudaTwistFlavorType flavor, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
enum QudaResidualType_s QudaResidualType
int main(int argc, char **argv)
QudaInverterType inv_type_precondition
Definition: quda.h:248
QudaLinkType type
Definition: quda.h:35
bool generate_all_levels
Definition: test_util.cpp:1666
double kappa
Definition: quda.h:97
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
#define cloverSiteSize
Definition: test_util.h:8
int return_clover_inverse
Definition: quda.h:217
QudaSolveType smoother_solve_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:426
int niter
Definition: test_util.cpp:1630
double mu_factor[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1659
enum QudaSolveType_s QudaSolveType
double setup_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:416
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
int nu_post
Definition: test_util.cpp:1658
QudaPrecision cpu_prec
Definition: quda.h:190
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:1795
bool generate_nullspace
Definition: test_util.cpp:1665
void setMultigridParam(QudaMultigridParam &mg_param)
void destroyMultigridQuda(void *mg_instance)
Free resources allocated by the multigrid solver.
char * strcpy(char *__dst, const char *__src)
QudaDagType dagger
Definition: quda.h:184
void finalizeComms()
Definition: test_util.cpp:107
void ax(const double &a, ColorSpinorField &x)
Definition: blas_quda.cu:209
QudaGaugeParam gauge_param
int tdim
Definition: test_util.cpp:1623
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
double true_res
Definition: quda.h:115
double tol
Definition: test_util.cpp:1647
QudaTwistFlavorType twist_flavor
Definition: test_util.cpp:1649
void tm_matpc(void *outEven, void **gauge, void *inEven, double kappa, double mu, QudaTwistFlavorType flavor, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:704
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
int return_clover
Definition: quda.h:216
#define spinorSiteSize
double mu_factor[QUDA_MAX_MG_LEVEL]
Definition: quda.h:471
int geo_block_size[QUDA_MAX_MG_LEVEL][QUDA_MAX_DIM]
Definition: test_util.cpp:1668
void exit(int) __attribute__((noreturn))
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:202
void setDims(int *)
Definition: test_util.cpp:130
QudaFieldLocation input_location
Definition: quda.h:90
QudaPrecision & cuda_prec_sloppy
void freeGaugeQuda(void)
int nu_pre
Definition: test_util.cpp:1657
double reliable_delta
Definition: quda.h:118
static size_t gSize
Definition: llfat_test.cpp:36
int n_vec[QUDA_MAX_MG_LEVEL]
Definition: quda.h:407
QudaSolutionType solution_type
Definition: quda.h:181
else return(__swbuf(_c, _p))
int strcmp(const char *__s1, const char *__s2)
QudaPrecision clover_cuda_prec
Definition: quda.h:201
int precondition_cycle
Definition: quda.h:273
void * preconditioner
Definition: quda.h:251
char vec_infile[]
Definition: test_util.cpp:1636
QudaPrecision prec
Definition: test_util.cpp:1615
QudaBoolean global_reduction[QUDA_MAX_MG_LEVEL]
Definition: quda.h:444
QudaPrecision & cuda_prec_precondition
void initQuda(int device)
double spinorGiB
Definition: quda.h:225
QudaFieldLocation output_location
Definition: quda.h:91
void * malloc(size_t __size) __attribute__((__warn_unused_result__)) __attribute__((alloc_size(1)))
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:203
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
int printf(const char *,...) __attribute__((__format__(__printf__
QudaMultigridCycleType cycle_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:429
QudaSolveType solve_type
Definition: test_util.cpp:1653
void setTransferGPU(bool)
QudaReconstructType link_recon
Definition: test_util.cpp:1612
QudaPrecision cuda_prec_sloppy
Definition: quda.h:192
QudaVerbosity verbosity
Definition: quda.h:219
void setSpinorSiteSize(int n)
Definition: test_util.cpp:192
double tol_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:156
double mass
Definition: test_util.cpp:1642
QudaInvertParam newQudaInvertParam(void)
double gflops
Definition: quda.h:227
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:770
QudaPrecision cuda_prec_precondition
Definition: quda.h:48
QudaCloverFieldOrder clover_order
Definition: quda.h:205
double omega[QUDA_MAX_MG_LEVEL]
Definition: quda.h:441
int V
Definition: test_util.cpp:28
double tol_hq
Definition: quda.h:112
enum QudaMatPCType_s QudaMatPCType
cpuColorSpinorField * spinorOut
Definition: covdev_test.cpp:41
#define gaugeSiteSize
Definition: test_util.h:6
void * newMultigridQuda(QudaMultigridParam *param)
double sqrt(double)
double true_res_hq
Definition: quda.h:116
double smoother_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:438
int gridsize_from_cmdline[]
Definition: test_util.cpp:50
QudaGammaBasis gamma_basis
Definition: quda.h:197
QudaPrecision cuda_prec_sloppy
Definition: quda.h:45
QudaVerbosity verbosity[QUDA_MAX_MG_LEVEL]
Definition: quda.h:410
double tol_precondition
Definition: quda.h:264
QudaBoolean run_verify
Definition: quda.h:456
QudaReconstructType reconstruct
Definition: quda.h:43
QudaDslashType dslash_type
Definition: test_util.cpp:1626
QudaPrecision cuda_prec
Definition: quda.h:42
QudaPrecision & cpu_prec
int X[4]
Definition: quda.h:29
double mass
Definition: quda.h:96
QudaInverterType smoother_type
Definition: test_util.cpp:1664
char vec_outfile[256]
Definition: quda.h:462
int gcrNkrylov
Definition: quda.h:237
double norm_2(void *v, int len, QudaPrecision precision)
int rand(void) __attribute__((__availability__(swift
int compute_clover_inverse
Definition: quda.h:215
QudaMatPCType matpc_type
Definition: test_util.cpp:1652
int pipeline
Definition: test_util.cpp:1632
void usage(char **)
Definition: test_util.cpp:1693
QudaFieldLocation location[QUDA_MAX_MG_LEVEL]
Definition: quda.h:447
int Lsdim
Definition: test_util.cpp:1624
void construct_clover_field(void *clover, double norm, double diag, QudaPrecision precision)
Definition: test_util.cpp:1166
void * memset(void *__b, int __c, size_t __len)
int spin_block_size[QUDA_MAX_MG_LEVEL]
Definition: quda.h:404
void wil_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
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 geo_block_size[QUDA_MAX_MG_LEVEL][QUDA_MAX_DIM]
Definition: quda.h:401
QudaBoolean generate_all_levels
Definition: quda.h:453
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
void initRand()
Definition: test_util.cpp:117
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.
int xdim
Definition: test_util.cpp:1620
#define printfQuda(...)
Definition: util_quda.h:84
QudaTboundary t_boundary
Definition: quda.h:38
QudaTwistFlavorType twist_flavor
Definition: quda.h:108
int Vh
Definition: test_util.cpp:29
char latfile[]
Definition: test_util.cpp:1627
enum QudaDslashType_s QudaDslashType
void mxpy(ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:192
QudaResidualType residual_type
Definition: quda.h:286
void wil_matpc(void *outEven, void **gauge, void *inEven, double kappa, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
int num_offset
Definition: quda.h:146
enum QudaVerbosity_s QudaVerbosity
double cloverGiB
Definition: quda.h:226
clock_t clock(void) __asm("_" "clock")
int gcrNkrylov
Definition: test_util.cpp:1631
int compute_clover
Definition: quda.h:214
#define QUDA_MAX_DIM
Maximum number of dimensions supported by QUDA. In practice, no routines make use of more than 5...
int nu_post[QUDA_MAX_MG_LEVEL]
Definition: quda.h:435
double omega
Definition: quda.h:270
double mu
Definition: test_util.cpp:1643
QudaComputeNullVector compute_null_vector
Definition: quda.h:450
double omega
Definition: test_util.cpp:1663
void read_gauge_field(const char *filename, void *gauge[], QudaPrecision prec, const int *X, int argc, char *argv[])
Definition: qio_field.h:12
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1613
void setGaugeParam(QudaGaugeParam &gauge_param)
QudaPrecision clover_cpu_prec
Definition: quda.h:200
char vec_infile[256]
Definition: quda.h:459
bool compute_clover
Definition: test_util.cpp:1646
double clover_coeff
Definition: test_util.cpp:1645
QudaPrecision prec_sloppy
Definition: test_util.cpp:1616
void initComms(int argc, char **argv, const int *commDims)
Definition: test_util.cpp:72
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
int zdim
Definition: test_util.cpp:1622
enum QudaInverterType_s QudaInverterType
double setup_tol
Definition: test_util.cpp:1662
#define MAX(a, b)
QudaPrecision cpu_prec
Definition: quda.h:40
double tol_hq
Definition: test_util.cpp:1648
QudaGaugeParam newQudaGaugeParam(void)
void setInvertParam(QudaInvertParam &inv_param)
QudaPreserveSource preserve_source
Definition: quda.h:188
QudaInvertParam * invert_param
Definition: quda.h:395
QudaVerbosity mg_verbosity[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1660
double clover_coeff
Definition: quda.h:208
enum QudaTwistFlavorType_s QudaTwistFlavorType