QUDA  0.9.0
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 
31 // Twisted mass flavor type
33 
34 extern int device;
35 extern int xdim;
36 extern int ydim;
37 extern int zdim;
38 extern int tdim;
39 extern int Lsdim;
40 extern int gridsize_from_cmdline[];
41 extern QudaPrecision prec;
49 extern int multishift; // whether to test multi-shift or standard solver
50 extern double mass; // mass of Dirac operator
51 extern double mu;
52 extern double anisotropy; // temporal anisotropy
53 extern double tol; // tolerance for inverter
54 extern double tol_hq; // heavy-quark tolerance for inverter
55 extern QudaMassNormalization normalization; // mass normalization of Dirac operators
56 extern QudaMatPCType matpc_type; // preconditioning type
57 
58 extern double clover_coeff;
59 extern bool compute_clover;
60 
61 extern int niter; // max solver iterations
62 extern int gcrNkrylov; // number of inner iterations for GCR, or l for BiCGstab-l
63 extern int pipeline; // length of pipeline for fused operations in GCR or BiCGstab-l
64 extern int solution_accumulator_pipeline; // length of pipeline for fused solution update from the direction vectors
65 extern char latfile[];
66 
67 extern void usage(char** );
68 
69 
70 
71 void
73 {
74  printfQuda("running the following test:\n");
75 
76  printfQuda("prec prec_sloppy multishift matpc_type recon recon_sloppy S_dimension T_dimension Ls_dimension dslash_type normalization\n");
77  printfQuda("%6s %6s %d %12s %2s %2s %3d/%3d/%3d %3d %2d %14s %8s\n",
81  xdim, ydim, zdim, tdim, Lsdim,
84 
85  printfQuda("Grid partition info: X Y Z T\n");
86  printfQuda(" %d %d %d %d\n",
87  dimPartitioned(0),
88  dimPartitioned(1),
89  dimPartitioned(2),
90  dimPartitioned(3));
91 
92  return ;
93 
94 }
95 
96 int main(int argc, char **argv)
97 {
98 
99  for (int i = 1; i < argc; i++){
100  if(process_command_line_option(argc, argv, &i) == 0){
101  continue;
102  }
103  printfQuda("ERROR: Invalid option:%s\n", argv[i]);
104  usage(argv);
105  }
106 
111 
112  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
113  initComms(argc, argv, gridsize_from_cmdline);
114 
116 
117  // *** QUDA parameters begin here.
118 
126  printfQuda("dslash_type %d not supported\n", dslash_type);
127  exit(0);
128  }
129 
134 
137 
138  double kappa5;
139 
140  gauge_param.X[0] = xdim;
141  gauge_param.X[1] = ydim;
142  gauge_param.X[2] = zdim;
143  gauge_param.X[3] = tdim;
144  inv_param.Ls = 1;
145 
150 
159 
161 
162  inv_param.mass = mass;
163  inv_param.mu = mu;
164  inv_param.kappa = 1.0 / (2.0 * (1 + 3/gauge_param.anisotropy + mass));
165 
167  inv_param.epsilon = 0.1385;
170  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH ||
173  inv_param.m5 = -1.8;
174  kappa5 = 0.5/(5 + inv_param.m5);
175  inv_param.Ls = Lsdim;
176  for(int k = 0; k < Lsdim; k++) // for mobius only
177  {
178  // b5[k], c[k] values are chosen for arbitrary values,
179  // but the difference of them are same as 1.0
180  inv_param.b_5[k] = 1.452;
181  inv_param.c_5[k] = 0.452;
182  }
183  }
184 
185  // offsets used only by multi-shift solver
186  inv_param.num_offset = 12;
187  double offset[12] = {0.01, 0.02, 0.03, 0.04, 0.05, 0.06, 0.07, 0.08, 0.09, 0.1, 0.11, 0.12};
188  for (int i=0; i<inv_param.num_offset; i++) inv_param.offset[i] = offset[i];
189 
191  if (multishift) {
197  } else {
199  }
200 
201  //inv_param.solution_type = QUDA_MAT_SOLUTION;
202 
204 
208 
216  } else {
218  }
219 
220 
222 
223  inv_param.Nsteps = 2;
225  inv_param.tol = tol;
226  inv_param.tol_restart = 1e-3; //now theoretical background for this parameter...
227  if(tol_hq == 0 && tol == 0){
228  errorQuda("qudaInvert: requesting zero residual\n");
229  exit(1);
230  }
231  // require both L2 relative and heavy quark residual to determine convergence
232  inv_param.residual_type = static_cast<QudaResidualType_s>(0);
233  inv_param.residual_type = (tol != 0) ? static_cast<QudaResidualType_s> ( inv_param.residual_type | QUDA_L2_RELATIVE_RESIDUAL) : inv_param.residual_type;
235 
236  inv_param.tol_hq = tol_hq; // specify a tolerance for the residual for heavy quark residual
237 
238  // these can be set individually
239  for (int i=0; i<inv_param.num_offset; i++) {
242  }
248 
249  // domain decomposition preconditioner parameters
251 
258  inv_param.omega = 1.0;
259 
266 
269 
270  gauge_param.ga_pad = 0; // 24*24*24/2;
271  inv_param.sp_pad = 0; // 24*24*24/2;
272  inv_param.cl_pad = 0; // 24*24*24/2;
273 
274  // For multi-GPU, ga_pad must be large enough to store a time-slice
275 #ifdef MULTI_GPU
276  int x_face_size = gauge_param.X[1]*gauge_param.X[2]*gauge_param.X[3]/2;
277  int y_face_size = gauge_param.X[0]*gauge_param.X[2]*gauge_param.X[3]/2;
278  int z_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[3]/2;
279  int t_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[2]/2;
280  int pad_size =MAX(x_face_size, y_face_size);
281  pad_size = MAX(pad_size, z_face_size);
282  pad_size = MAX(pad_size, t_face_size);
283  gauge_param.ga_pad = pad_size;
284 #endif
285 
293  }
294 
296 
297  // *** Everything between here and the call to initQuda() is
298  // *** application-specific.
299 
300  // set parameters for the reference Dslash, and prepare fields to be loaded
305  } else {
307  }
308 
309  setSpinorSiteSize(24);
310 
311  size_t gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
312  size_t sSize = (inv_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
313 
314  void *gauge[4], *clover=0, *clover_inv=0;
315 
316  for (int dir = 0; dir < 4; dir++) {
317  gauge[dir] = malloc(V*gaugeSiteSize*gSize);
318  }
319 
320  if (strcmp(latfile,"")) { // load in the command line supplied gauge field
323  } else { // else generate a random SU(3) field
325  }
326 
328  double norm = 0.01; // clover components are random numbers in the range (-norm, norm)
329  double diag = 1.0; // constant added to the diagonal
330 
331  size_t cSize = inv_param.clover_cpu_prec;
332  clover = malloc(V*cloverSiteSize*cSize);
333  clover_inv = malloc(V*cloverSiteSize*cSize);
335 
340  }
341 
342  void *spinorIn = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
343  void *spinorCheck = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
344 
345  void *spinorOut = NULL, **spinorOutMulti = NULL;
346  if (multishift) {
347  spinorOutMulti = (void**)malloc(inv_param.num_offset*sizeof(void *));
348  for (int i=0; i<inv_param.num_offset; i++) {
349  spinorOutMulti[i] = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
350  }
351  } else {
353  }
354 
355  memset(spinorIn, 0, inv_param.Ls*V*spinorSiteSize*sSize);
356  memset(spinorCheck, 0, inv_param.Ls*V*spinorSiteSize*sSize);
357  if (multishift) {
358  for (int i=0; i<inv_param.num_offset; i++) memset(spinorOutMulti[i], 0, inv_param.Ls*V*spinorSiteSize*sSize);
359  } else {
361  }
362 
363  // create a point source at 0 (in each subvolume... FIXME)
364 
365  // create a point source at 0 (in each subvolume... FIXME)
366 
368  //((float*)spinorIn)[0] = 1.0;
369  for (int i=0; i<inv_param.Ls*V*spinorSiteSize; i++) ((float*)spinorIn)[i] = rand() / (float)RAND_MAX;
370  } else {
371  //((double*)spinorIn)[0] = 1.0;
372  for (int i=0; i<inv_param.Ls*V*spinorSiteSize; i++) ((double*)spinorIn)[i] = rand() / (double)RAND_MAX;
373  }
374 
375  // start the timer
376  double time0 = -((double)clock());
377 
378  // initialize the QUDA library
379  initQuda(device);
380 
381  // load the gauge field
382  loadGaugeQuda((void*)gauge, &gauge_param);
383 
384  // load the clover term, if desired
386  loadCloverQuda(clover, clover_inv, &inv_param);
387 
388  // perform the inversion
389  if (multishift) {
390  invertMultiShiftQuda(spinorOutMulti, spinorIn, &inv_param);
391  } else {
392  invertQuda(spinorOut, spinorIn, &inv_param);
393  }
394 
395  // stop the timer
396  time0 += clock();
397  time0 /= CLOCKS_PER_SEC;
398 
399  printfQuda("Device memory used:\n Spinor: %f GiB\n Gauge: %f GiB\n",
402  printfQuda(" Clover: %f GiB\n", inv_param.cloverGiB);
403  }
404  printfQuda("\nDone: %i iter / %g secs = %g Gflops, total time = %g secs\n",
406 
407  if (multishift) {
409  errorQuda("Mass normalization not supported for multi-shift solver in invert_test");
410  }
411 
412  void *spinorTmp = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
413 
414  printfQuda("Host residuum checks: \n");
415  for(int i=0; i < inv_param.num_offset; i++) {
416  ax(0, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
417 
420  int tm_offset = Vh*spinorSiteSize;
421  void *out0 = spinorCheck;
422  void *out1 = (char*)out0 + tm_offset*cpu_prec;
423 
424  void *tmp0 = spinorTmp;
425  void *tmp1 = (char*)tmp0 + tm_offset*cpu_prec;
426 
427  void *in0 = spinorOutMulti[i];
428  void *in1 = (char*)in0 + tm_offset*cpu_prec;
429 
432  } else {
433  tm_matpc(spinorTmp, gauge, spinorOutMulti[i], inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
437  }
438  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
440  errorQuda("Twisted mass solution type not supported");
441  tmc_matpc(spinorTmp, gauge, spinorOutMulti[i], clover, clover_inv, inv_param.kappa, inv_param.mu,
443  tmc_matpc(spinorCheck, gauge, spinorTmp, clover, clover_inv, inv_param.kappa, inv_param.mu,
445  } else if (dslash_type == QUDA_WILSON_DSLASH) {
446  wil_matpc(spinorTmp, gauge, spinorOutMulti[i], inv_param.kappa, inv_param.matpc_type, 0,
448  wil_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
450  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
451  clover_matpc(spinorTmp, gauge, clover, clover_inv, spinorOutMulti[i], inv_param.kappa, inv_param.matpc_type, 0,
453  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
455  } else {
456  printfQuda("Domain wall not supported for multi-shift\n");
457  exit(-1);
458  }
459 
460  axpy(inv_param.offset[i], spinorOutMulti[i], spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
461  mxpy(spinorIn, spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
462  double nrm2 = norm_2(spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
463  double src2 = norm_2(spinorIn, Vh*spinorSiteSize, inv_param.cpu_prec);
464  double l2r = sqrt(nrm2 / src2);
465 
466  printfQuda("Shift %d residuals: (L2 relative) tol %g, QUDA = %g, host = %g; (heavy-quark) tol %g, QUDA = %g\n",
469  }
470  free(spinorTmp);
471 
472  } else {
473 
475 
479  } else {
480  int tm_offset = V*spinorSiteSize;
481  void *evenOut = spinorCheck;
482  void *oddOut = (char*)evenOut + tm_offset*cpu_prec;
483 
484  void *evenIn = spinorOut;
485  void *oddIn = (char*)evenIn + tm_offset*cpu_prec;
486 
487  tm_ndeg_mat(evenOut, oddOut, gauge, evenIn, oddIn, inv_param.kappa, inv_param.mu, inv_param.epsilon, 0, inv_param.cpu_prec, gauge_param);
488  }
489  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
492  } else if (dslash_type == QUDA_WILSON_DSLASH) {
493  wil_mat(spinorCheck, gauge, spinorOut, inv_param.kappa, 0, inv_param.cpu_prec, gauge_param);
494  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
496  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
498  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
500  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
501  double *kappa_b, *kappa_c;
502  kappa_b = (double*)malloc(Lsdim*sizeof(double));
503  kappa_c = (double*)malloc(Lsdim*sizeof(double));
504  for(int xs = 0 ; xs < Lsdim ; xs++)
505  {
506  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
507  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
508  }
509  mdw_mat(spinorCheck, gauge, spinorOut, kappa_b, kappa_c, inv_param.dagger, inv_param.cpu_prec, gauge_param, inv_param.mass, inv_param.b_5, inv_param.c_5);
510  free(kappa_b);
511  free(kappa_c);
512  } else {
513  errorQuda("Unsupported dslash_type");
514  }
519  ax(0.5/kappa5, spinorCheck, V*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
521  ax(0.5/inv_param.kappa, spinorCheck, 2*V*spinorSiteSize, inv_param.cpu_prec);
522  } else {
523  ax(0.5/inv_param.kappa, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
524  }
525  }
526 
528 
531  int tm_offset = Vh*spinorSiteSize;
532  void *out0 = spinorCheck;
533  void *out1 = (char*)out0 + tm_offset*cpu_prec;
534 
535  void *in0 = spinorOut;
536  void *in1 = (char*)in0 + tm_offset*cpu_prec;
537 
539  } else {
542  }
543  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
545  errorQuda("Twisted mass solution type not supported");
546  tmc_matpc(spinorCheck, gauge, spinorOut, clover, clover_inv, inv_param.kappa, inv_param.mu,
548  } else if (dslash_type == QUDA_WILSON_DSLASH) {
549  wil_matpc(spinorCheck, gauge, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
551  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
552  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
554  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
556  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
558  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
559  double *kappa_b, *kappa_c;
560  kappa_b = (double*)malloc(Lsdim*sizeof(double));
561  kappa_c = (double*)malloc(Lsdim*sizeof(double));
562  for(int xs = 0 ; xs < Lsdim ; xs++)
563  {
564  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
565  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
566  }
567  mdw_matpc(spinorCheck, gauge, spinorOut, kappa_b, kappa_c, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param, inv_param.mass, inv_param.b_5, inv_param.c_5);
568  free(kappa_b);
569  free(kappa_c);
570  } else {
571  errorQuda("Unsupported dslash_type");
572  }
573 
578  ax(0.25/(kappa5*kappa5), spinorCheck, V*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
579  } else {
581 
582  }
583  }
584 
586 
587  void *spinorTmp = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
588 
589  ax(0, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
590 
593  int tm_offset = Vh*spinorSiteSize;
594  void *out0 = spinorCheck;
595  void *out1 = (char*)out0 + tm_offset*cpu_prec;
596 
597  void *tmp0 = spinorTmp;
598  void *tmp1 = (char*)tmp0 + tm_offset*cpu_prec;
599 
600  void *in0 = spinorOut;
601  void *in1 = (char*)in0 + tm_offset*cpu_prec;
602 
605  } else {
610  }
611  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
613  errorQuda("Twisted mass solution type not supported");
616  tmc_matpc(spinorCheck, gauge, spinorTmp, clover, clover_inv, inv_param.kappa, inv_param.mu,
618  } else if (dslash_type == QUDA_WILSON_DSLASH) {
621  wil_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
623  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
624  clover_matpc(spinorTmp, gauge, clover, clover_inv, spinorOut, inv_param.kappa,
626  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorTmp, inv_param.kappa,
628  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
631  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
634  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
635  double *kappa_b, *kappa_c;
636  kappa_b = (double*)malloc(Lsdim*sizeof(double));
637  kappa_c = (double*)malloc(Lsdim*sizeof(double));
638  for(int xs = 0 ; xs < Lsdim ; xs++)
639  {
640  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
641  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
642  }
644  mdw_matpc(spinorCheck, gauge, spinorTmp, kappa_b, kappa_c, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param, inv_param.mass, inv_param.b_5, inv_param.c_5);
645  free(kappa_b);
646  free(kappa_c);
647  } else {
648  errorQuda("Unsupported dslash_type");
649  }
650 
652  errorQuda("Mass normalization not implemented");
653  }
654 
655  free(spinorTmp);
656  }
657 
658 
659  int vol = inv_param.solution_type == QUDA_MAT_SOLUTION ? V : Vh;
660  mxpy(spinorIn, spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
661  double nrm2 = norm_2(spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
662  double src2 = norm_2(spinorIn, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
663  double l2r = sqrt(nrm2 / src2);
664 
665  printfQuda("Residuals: (L2 relative) tol %g, QUDA = %g, host = %g; (heavy-quark) tol %g, QUDA = %g\n",
667 
668  }
669 
670  freeGaugeQuda();
672 
673  // finalize the QUDA library
674  endQuda();
675 
676  finalizeComms();
677 
679  if (clover) free(clover);
680  if (clover_inv) free(clover_inv);
681  }
682 
683  for (int dir = 0; dir<4; dir++) free(gauge[dir]);
684 
685  return 0;
686 }
int maxiter_precondition
Definition: quda.h:267
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
enum QudaMassNormalization_s QudaMassNormalization
double anisotropy
Definition: test_util.cpp:1644
QudaReconstructType reconstruct_sloppy
Definition: quda.h:46
double anisotropy
Definition: quda.h:31
void freeCloverQuda(void)
QudaTwistFlavorType twist_flavor
Definition: test_util.cpp:1649
double b_5[QUDA_MAX_DWF_LS]
Definition: quda.h:102
void dw_4d_matpc(void *out, void **gauge, void *in, double kappa, QudaMatPCType matpc_type, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm)
void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param)
void endQuda(void)
void free(void *)
void construct_gauge_field(void **gauge, int type, QudaPrecision precision, QudaGaugeParam *param)
Definition: test_util.cpp:1054
QudaSolveType solve_type
Definition: quda.h:182
QudaVerbosity verbosity_precondition
Definition: quda.h:261
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:53
int multishift
Definition: test_util.cpp:1640
int main(int argc, char **argv)
Definition: invert_test.cpp:96
void dw_setDims(int *X, const int L5)
Definition: test_util.cpp:167
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
void tm_mat(void *out, void **gauge, void *in, double kappa, double mu, QudaTwistFlavorType flavor, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaInverterType inv_type_precondition
Definition: quda.h:248
int Lsdim
Definition: test_util.cpp:1624
QudaLinkType type
Definition: quda.h:35
double kappa
Definition: quda.h:97
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
#define errorQuda(...)
Definition: util_quda.h:90
double tol
Definition: quda.h:110
QudaPrecision prec
Definition: test_util.cpp:1615
void usage(char **)
Definition: test_util.cpp:1693
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
#define cloverSiteSize
Definition: test_util.h:8
int return_clover_inverse
Definition: quda.h:217
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
bool compute_clover
Definition: test_util.cpp:1646
QudaPrecision cpu_prec
Definition: quda.h:190
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:1795
double tol_hq
Definition: test_util.cpp:1648
QudaReconstructType link_recon_precondition
Definition: test_util.cpp:1614
void tm_ndeg_mat(void *evenOut, void *oddOut, void **gauge, void *evenIn, void *oddIn, double kappa, double mu, double epsilon, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
void clover_matpc(void *out, void **gauge, void *clover, void *clover_inv, void *in, double kappa, QudaMatPCType matpc_type, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaDagType dagger
Definition: quda.h:184
const char * get_matpc_str(QudaMatPCType type)
Definition: misc.cpp:987
void finalizeComms()
Definition: test_util.cpp:107
void ax(const double &a, ColorSpinorField &x)
Definition: blas_quda.cu:209
QudaGaugeParam gauge_param
double clover_coeff
Definition: test_util.cpp:1645
void mdw_mat(void *out, void **gauge, void *in, double *kappa_b, double *kappa_c, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm, double *b5, double *c5)
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
void dw_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm)
double true_res
Definition: quda.h:115
void tmc_mat(void *out, void **gauge, void *clover, void *in, double kappa, double mu, QudaTwistFlavorType flavor, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaReconstructType link_recon
Definition: test_util.cpp:1612
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)
QudaPrecision & cuda_prec_precondition
int return_clover
Definition: quda.h:216
void display_test_info()
Definition: invert_test.cpp:72
#define spinorSiteSize
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
void freeGaugeQuda(void)
double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:168
double reliable_delta
Definition: quda.h:118
QudaPrecision cpu_prec
Definition: covdev_test.cpp:33
int solution_accumulator_pipeline
Definition: quda.h:129
size_t size_t offset
double mu
Definition: test_util.cpp:1643
static size_t gSize
Definition: llfat_test.cpp:36
int pipeline
Definition: test_util.cpp:1632
QudaSolutionType solution_type
Definition: quda.h:181
else return(__swbuf(_c, _p))
QudaSolverNormalization solver_normalization
Definition: quda.h:186
int strcmp(const char *__s1, const char *__s2)
QudaPrecision clover_cuda_prec
Definition: quda.h:201
int precondition_cycle
Definition: quda.h:273
void initQuda(int device)
void dw_matpc(void *out, void **gauge, void *in, double kappa, QudaMatPCType matpc_type, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm)
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
QudaPrecision & cuda_prec_sloppy
double m5
Definition: quda.h:99
QudaPrecision cuda_prec
Definition: covdev_test.cpp:34
QudaPrecision cuda_prec_sloppy
Definition: quda.h:192
QudaMassNormalization normalization
Definition: test_util.cpp:1651
QudaVerbosity verbosity
Definition: quda.h:219
void setSpinorSiteSize(int n)
Definition: test_util.cpp:192
const char * get_mass_normalization_str(QudaMassNormalization type)
Definition: misc.cpp:943
int solution_accumulator_pipeline
Definition: test_util.cpp:1633
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 gflops
Definition: quda.h:227
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:770
QudaMatPCType matpc_type
Definition: test_util.cpp:1652
QudaPrecision cuda_prec_precondition
Definition: quda.h:48
QudaCloverFieldOrder clover_order
Definition: quda.h:205
void dw_4d_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm)
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
double sqrt(double)
void mdw_matpc(void *out, void **gauge, void *in, double *kappa_b, double *kappa_c, QudaMatPCType matpc_type, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm, double *b5, double *c5)
double true_res_hq
Definition: quda.h:116
QudaGammaBasis gamma_basis
Definition: quda.h:197
QudaPrecision cuda_prec_sloppy
Definition: quda.h:45
const char * get_dslash_str(QudaDslashType type)
Definition: misc.cpp:878
double tol_precondition
Definition: quda.h:264
#define tmp1
Definition: tmc_core.h:15
void clover_mat(void *out, void **gauge, void *clover, void *in, double kappa, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param)
double offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:153
int use_sloppy_partial_accumulator
Definition: quda.h:119
void tm_ndeg_matpc(void *outEven1, void *outEven2, void **gauge, void *inEven1, void *inEven2, double kappa, double mu, double epsilon, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaReconstructType reconstruct
Definition: quda.h:43
QudaPrecision cuda_prec
Definition: quda.h:42
int X[4]
Definition: quda.h:29
double tol
Definition: test_util.cpp:1647
double mass
Definition: quda.h:96
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
void axpy(const double &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:150
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 max_res_increase
Definition: quda.h:134
QudaResidualType_s
Definition: enum_quda.h:165
int xdim
Definition: test_util.cpp:1620
void wil_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaInverterType inv_type
Definition: test_util.cpp:1638
double gaugeGiB
Definition: quda.h:60
QudaPrecision cuda_prec_precondition
Definition: quda.h:193
double tol_restart
Definition: quda.h:111
int zdim
Definition: test_util.cpp:1622
int tdim
Definition: test_util.cpp:1623
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
char latfile[]
Definition: test_util.cpp:1627
void tmc_matpc(void *out, void **gauge, void *in, void *clover, void *cInv, double kappa, double mu, QudaTwistFlavorType flavor, QudaMatPCType matpc_type, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param)
if(err !=cudaSuccess)
#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
#define MAX(a, b)
Definition: invert_test.cpp:23
enum QudaDslashType_s QudaDslashType
int niter
Definition: test_util.cpp:1630
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
double cloverGiB
Definition: quda.h:226
clock_t clock(void) __asm("_" "clock")
QudaDslashType dslash_type
Definition: test_util.cpp:1626
int compute_clover
Definition: quda.h:214
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1613
double epsilon
Definition: quda.h:106
int ydim
Definition: test_util.cpp:1621
cpuColorSpinorField * spinorTmp
double omega
Definition: quda.h:270
double mass
Definition: test_util.cpp:1642
QudaPrecision prec_sloppy
Definition: test_util.cpp:1616
QudaInverterType precon_type
Definition: test_util.cpp:1639
void read_gauge_field(const char *filename, void *gauge[], QudaPrecision prec, const int *X, int argc, char *argv[])
Definition: qio_field.h:12
int gcrNkrylov
Definition: test_util.cpp:1631
int gridsize_from_cmdline[]
Definition: test_util.cpp:50
QudaPrecision clover_cpu_prec
Definition: quda.h:200
void initComms(int argc, char **argv, const int *commDims)
Definition: test_util.cpp:72
#define tmp0
Definition: tmc_core.h:14
QudaMatPCType matpc_type
Definition: quda.h:183
QudaPrecision prec_precondition
Definition: test_util.cpp:1617
enum QudaInverterType_s QudaInverterType
double kappa5
QudaPrecision cpu_prec
Definition: quda.h:40
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:188
double clover_coeff
Definition: quda.h:208
enum QudaTwistFlavorType_s QudaTwistFlavorType