QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 #include <limits>
7 
8 #include <util_quda.h>
9 #include <random_quda.h>
10 #include <test_util.h>
11 #include <dslash_util.h>
12 #include <blas_reference.h>
15 #include "misc.h"
16 
17 #include <qio_field.h>
18 
19 #define MAX(a,b) ((a)>(b)?(a):(b))
20 
21 // In a typical application, quda.h is the only QUDA header required.
22 #include <quda.h>
23 
24 // Wilson, clover-improved Wilson, twisted mass, and domain wall are supported.
26 
27 // Twisted mass flavor type
29 
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[];
37 extern QudaPrecision prec;
45 extern double reliable_delta; // reliable update parameter
46 extern bool alternative_reliable;
48 extern int multishift; // whether to test multi-shift or standard solver
49 extern double mass; // mass of Dirac operator
50 extern double kappa; // kappa of Dirac operator
51 extern double mu;
52 extern double epsilon;
53 extern double anisotropy; // temporal anisotropy
54 extern double tol; // tolerance for inverter
55 extern double tol_hq; // heavy-quark tolerance for inverter
56 extern QudaMassNormalization normalization; // mass normalization of Dirac operators
57 extern QudaMatPCType matpc_type; // preconditioning type
58 extern QudaSolutionType solution_type; // the solution we desire
59 extern QudaSolveType solve_type; // the solve type we want to find the solution
60 
61 extern double clover_coeff;
62 extern bool compute_clover;
63 
65 extern QudaVerbosity mg_verbosity[QUDA_MAX_MG_LEVEL]; // use this for preconditioner verbosity
66 
67 extern int Nsrc; // number of spinors to apply to simultaneously
68 extern int niter; // max solver iterations
69 extern int gcrNkrylov; // number of inner iterations for GCR, or l for BiCGstab-l
70 extern QudaCABasis ca_basis; // basis for CA-CG solves
71 extern double ca_lambda_min; // minimum eigenvalue for scaling Chebyshev CA-CG solves
72 extern double ca_lambda_max; // maximum eigenvalue for scaling Chebyshev CA-CG solves
73 extern int pipeline; // length of pipeline for fused operations in GCR or BiCGstab-l
74 extern int solution_accumulator_pipeline; // length of pipeline for fused solution update from the direction vectors
75 extern char latfile[];
76 extern bool unit_gauge;
77 
78 extern void usage(char** );
79 
80 
81 
82 void
84 {
85  printfQuda("running the following test:\n");
86 
87  printfQuda("prec prec_sloppy multishift matpc_type recon recon_sloppy S_dimension T_dimension Ls_dimension dslash_type normalization\n");
88  printfQuda("%6s %6s %d %12s %2s %2s %3d/%3d/%3d %3d %2d %14s %8s\n",
92  xdim, ydim, zdim, tdim, Lsdim,
95 
96  printfQuda("Grid partition info: X Y Z T\n");
97  printfQuda(" %d %d %d %d\n",
98  dimPartitioned(0),
99  dimPartitioned(1),
100  dimPartitioned(2),
101  dimPartitioned(3));
102 
103  return ;
104 
105 }
106 
107 int main(int argc, char **argv)
108 {
109 
110  mg_verbosity[0] = QUDA_SILENT; // set default preconditioner verbosity
111 
112  if (multishift) solution_type = QUDA_MATPCDAG_MATPC_SOLUTION; // set a correct default for the multi-shift solver
113 
114  for (int i = 1; i < argc; i++){
115  if(process_command_line_option(argc, argv, &i) == 0){
116  continue;
117  }
118  printfQuda("ERROR: Invalid option:%s\n", argv[i]);
119  usage(argv);
120  }
121 
127 
128  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
129  initComms(argc, argv, gridsize_from_cmdline);
130 
132 
133  // *** QUDA parameters begin here.
134 
142  printfQuda("dslash_type %d not supported\n", dslash_type);
143  exit(0);
144  }
145 
151 
154 
155  double kappa5;
156 
157  gauge_param.X[0] = xdim;
158  gauge_param.X[1] = ydim;
159  gauge_param.X[2] = zdim;
160  gauge_param.X[3] = tdim;
161  inv_param.Ls = 1;
162 
163  gauge_param.anisotropy = anisotropy;
164  gauge_param.type = QUDA_WILSON_LINKS;
165  gauge_param.gauge_order = QUDA_QDP_GAUGE_ORDER;
166  gauge_param.t_boundary = QUDA_PERIODIC_T;
167 
168  gauge_param.cpu_prec = cpu_prec;
169  gauge_param.cuda_prec = cuda_prec;
170  gauge_param.reconstruct = link_recon;
171  gauge_param.cuda_prec_sloppy = cuda_prec_sloppy;
177 
178  gauge_param.gauge_fix = QUDA_GAUGE_FIXED_NO;
179 
180  inv_param.dslash_type = dslash_type;
181 
182  if (kappa == -1.0) {
183  inv_param.mass = mass;
184  inv_param.kappa = 1.0 / (2.0 * (1 + 3/gauge_param.anisotropy + mass));
185  } else {
186  inv_param.kappa = kappa;
187  inv_param.mass = 0.5/kappa - (1 + 3/gauge_param.anisotropy);
188  }
189  inv_param.mu = mu;
190 
192  inv_param.epsilon = epsilon;
193  inv_param.twist_flavor = twist_flavor;
194  inv_param.Ls = (inv_param.twist_flavor == QUDA_TWIST_NONDEG_DOUBLET) ? 2 : 1;
195  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH ||
198  inv_param.m5 = -1.8;
199  kappa5 = 0.5/(5 + inv_param.m5);
200  inv_param.Ls = Lsdim;
201  for(int k = 0; k < Lsdim; k++) // for mobius only
202  {
203  // b5[k], c[k] values are chosen for arbitrary values,
204  // but the difference of them are same as 1.0
205  inv_param.b_5[k] = 1.452;
206  inv_param.c_5[k] = 0.452;
207  }
208  }
209 
210  // offsets used only by multi-shift solver
211  inv_param.num_offset = 12;
212  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};
213  for (int i=0; i<inv_param.num_offset; i++) inv_param.offset[i] = offset[i];
214 
215  inv_param.inv_type = inv_type;
216  inv_param.solution_type = solution_type;
217  inv_param.solve_type = solve_type;
218  inv_param.matpc_type = matpc_type;
219 
220  inv_param.dagger = QUDA_DAG_NO;
221  inv_param.mass_normalization = normalization;
223 
224 
225  inv_param.pipeline = pipeline;
226 
227  inv_param.Nsteps = 2;
228  inv_param.gcrNkrylov = gcrNkrylov;
229  inv_param.ca_basis = ca_basis;
230  inv_param.ca_lambda_min = ca_lambda_min;
231  inv_param.ca_lambda_max = ca_lambda_max;
232  inv_param.tol = tol;
233  inv_param.tol_restart = 1e-3; //now theoretical background for this parameter...
234  if(tol_hq == 0 && tol == 0){
235  errorQuda("qudaInvert: requesting zero residual\n");
236  exit(1);
237  }
238  // require both L2 relative and heavy quark residual to determine convergence
239  inv_param.residual_type = static_cast<QudaResidualType_s>(0);
240  inv_param.residual_type = (tol != 0) ? static_cast<QudaResidualType_s> ( inv_param.residual_type | QUDA_L2_RELATIVE_RESIDUAL) : inv_param.residual_type;
241  inv_param.residual_type = (tol_hq != 0) ? static_cast<QudaResidualType_s> (inv_param.residual_type | QUDA_HEAVY_QUARK_RESIDUAL) : inv_param.residual_type;
242 
243  inv_param.tol_hq = tol_hq; // specify a tolerance for the residual for heavy quark residual
244 
245  // these can be set individually
246  for (int i=0; i<inv_param.num_offset; i++) {
247  inv_param.tol_offset[i] = inv_param.tol;
248  inv_param.tol_hq_offset[i] = inv_param.tol_hq;
249  }
250  inv_param.maxiter = niter;
251  inv_param.reliable_delta = reliable_delta;
253  inv_param.use_sloppy_partial_accumulator = 0;
255  inv_param.max_res_increase = 1;
256 
257  // domain decomposition preconditioner parameters
259 
261  inv_param.precondition_cycle = 1;
262  inv_param.tol_precondition = 1e-1;
263  inv_param.maxiter_precondition = 10;
264  inv_param.verbosity_precondition = mg_verbosity[0];
266  inv_param.omega = 1.0;
267 
268  inv_param.cpu_prec = cpu_prec;
269  inv_param.cuda_prec = cuda_prec;
274  inv_param.dirac_order = QUDA_DIRAC_ORDER;
275 
278 
279  gauge_param.ga_pad = 0; // 24*24*24/2;
280  inv_param.sp_pad = 0; // 24*24*24/2;
281  inv_param.cl_pad = 0; // 24*24*24/2;
282 
283  // For multi-GPU, ga_pad must be large enough to store a time-slice
284 #ifdef MULTI_GPU
285  int x_face_size = gauge_param.X[1]*gauge_param.X[2]*gauge_param.X[3]/2;
286  int y_face_size = gauge_param.X[0]*gauge_param.X[2]*gauge_param.X[3]/2;
287  int z_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[3]/2;
288  int t_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[2]/2;
289  int pad_size =MAX(x_face_size, y_face_size);
290  pad_size = MAX(pad_size, z_face_size);
291  pad_size = MAX(pad_size, t_face_size);
292  gauge_param.ga_pad = pad_size;
293 #endif
294 
296  inv_param.clover_cpu_prec = cpu_prec;
297  inv_param.clover_cuda_prec = cuda_prec;
302  inv_param.clover_coeff = clover_coeff;
303  }
304 
305  inv_param.verbosity = verbosity;
306 
307  // *** Everything between here and the call to initQuda() is
308  // *** application-specific.
309 
310  // set parameters for the reference Dslash, and prepare fields to be loaded
314  dw_setDims(gauge_param.X, inv_param.Ls);
315  } else {
316  setDims(gauge_param.X);
317  }
318 
319  setSpinorSiteSize(24);
320 
321  size_t gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
322  size_t sSize = (inv_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
323 
324  void *gauge[4], *clover=0, *clover_inv=0;
325 
326  for (int dir = 0; dir < 4; dir++) {
327  gauge[dir] = malloc(V*gaugeSiteSize*gSize);
328  }
329 
330  if (strcmp(latfile,"")) { // load in the command line supplied gauge field
331  read_gauge_field(latfile, gauge, gauge_param.cpu_prec, gauge_param.X, argc, argv);
332  construct_gauge_field(gauge, 2, gauge_param.cpu_prec, &gauge_param);
333  } else { // else generate an SU(3) field
334  if (unit_gauge) {
335  // unit SU(3) field
336  construct_gauge_field(gauge, 0, gauge_param.cpu_prec, &gauge_param);
337  } else {
338  // random SU(3) field
339  construct_gauge_field(gauge, 1, gauge_param.cpu_prec, &gauge_param);
340  }
341  }
342 
344  double norm = 0.01; // clover components are random numbers in the range (-norm, norm)
345  double diag = 1.0; // constant added to the diagonal
346 
347  size_t cSize = inv_param.clover_cpu_prec;
348  clover = malloc(V*cloverSiteSize*cSize);
349  clover_inv = malloc(V*cloverSiteSize*cSize);
350  if (!compute_clover) construct_clover_field(clover, norm, diag, inv_param.clover_cpu_prec);
351 
352  inv_param.compute_clover = compute_clover;
353  if (compute_clover) inv_param.return_clover = 1;
354  inv_param.compute_clover_inverse = 1;
355  inv_param.return_clover_inverse = 1;
356  }
357 
358  void *spinorIn = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
359  void *spinorCheck = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
360 
361  void *spinorOut = NULL, **spinorOutMulti = NULL;
362  if (multishift) {
363  spinorOutMulti = (void**)malloc(inv_param.num_offset*sizeof(void *));
364  for (int i=0; i<inv_param.num_offset; i++) {
365  spinorOutMulti[i] = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
366  }
367  } else {
368  spinorOut = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
369  }
370 
371  // initialize the QUDA library
372  initQuda(device);
373 
374  // load the gauge field
375  loadGaugeQuda((void*)gauge, &gauge_param);
376 
377  double plaq[3];
378  plaqQuda(plaq);
379  printfQuda("Computed plaquette is %e (spatial = %e, temporal = %e)\n", plaq[0], plaq[1], plaq[2]);
380 
381  // load the clover term, if desired
383  loadCloverQuda(clover, clover_inv, &inv_param);
384 
385  double *time = new double[Nsrc];
386  double *gflops = new double[Nsrc];
387  auto *rng = new quda::RNG(quda::LatticeFieldParam(gauge_param), 1234);
388  rng->Init();
389 
390  for (int i = 0; i < Nsrc; i++) {
391 
392  construct_spinor_source(spinorIn, 4, 3, inv_param.cpu_prec, gauge_param.X, *rng);
393 
394  if (multishift) {
395  invertMultiShiftQuda(spinorOutMulti, spinorIn, &inv_param);
396  } else {
397  invertQuda(spinorOut, spinorIn, &inv_param);
398  }
399 
400  time[i] = inv_param.secs;
401  gflops[i] = inv_param.gflops / inv_param.secs;
402  printfQuda("Done: %i iter / %g secs = %g Gflops\n\n", inv_param.iter, inv_param.secs,
403  inv_param.gflops / inv_param.secs);
404  }
405 
406  rng->Release();
407  delete rng;
408 
409  auto mean_time = 0.0;
410  auto mean_time2 = 0.0;
411  auto mean_gflops = 0.0;
412  auto mean_gflops2 = 0.0;
413  for (int i = 0; i < Nsrc; i++) {
414  mean_time += time[i];
415  mean_time2 += time[i] * time[i];
416  mean_gflops += gflops[i];
417  mean_gflops2 += gflops[i] * gflops[i];
418  }
419 
420  mean_time /= Nsrc;
421  mean_time2 /= Nsrc;
422  auto stddev_time = Nsrc > 1 ? sqrt((Nsrc / ((double)Nsrc - 1.0)) * (mean_time2 - mean_time * mean_time)) : std::numeric_limits<double>::infinity();
423  mean_gflops /= Nsrc;
424  mean_gflops2 /= Nsrc;
425  auto stddev_gflops = Nsrc > 1 ? sqrt((Nsrc / ((double)Nsrc - 1.0)) * (mean_gflops2 - mean_gflops * mean_gflops)) : std::numeric_limits<double>::infinity();
426  printfQuda("%d solves, with mean solve time %g (stddev = %g), mean GFLOPS %g (stddev = %g)\n", Nsrc, mean_time,
427  stddev_time, mean_gflops, stddev_gflops);
428 
429  delete[] time;
430  delete[] gflops;
431 
432  if (multishift) {
433  if (inv_param.mass_normalization == QUDA_MASS_NORMALIZATION) {
434  errorQuda("Mass normalization not supported for multi-shift solver in invert_test");
435  }
436 
437  void *spinorTmp = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
438 
439  printfQuda("Host residuum checks: \n");
440  for(int i=0; i < inv_param.num_offset; i++) {
441  ax(0, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
442 
444  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET) {
445  int tm_offset = Vh*spinorSiteSize;
446  void *out0 = spinorCheck;
447  void *out1 = (char*)out0 + tm_offset*cpu_prec;
448 
449  void *tmp0 = spinorTmp;
450  void *tmp1 = (char*)tmp0 + tm_offset*cpu_prec;
451 
452  void *in0 = spinorOutMulti[i];
453  void *in1 = (char*)in0 + tm_offset*cpu_prec;
454 
455  tm_ndeg_matpc(tmp0, tmp1, gauge, in0, in1, inv_param.kappa, inv_param.mu, inv_param.epsilon, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
456  tm_ndeg_matpc(out0, out1, gauge, tmp0, tmp1, inv_param.kappa, inv_param.mu, inv_param.epsilon, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
457  } else {
458  tm_matpc(spinorTmp, gauge, spinorOutMulti[i], inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
459  inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
460  tm_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
461  inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
462  }
463  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
464  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET)
465  errorQuda("Twisted mass solution type not supported");
466  tmc_matpc(spinorTmp, gauge, spinorOutMulti[i], clover, clover_inv, inv_param.kappa, inv_param.mu,
467  inv_param.twist_flavor, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
468  tmc_matpc(spinorCheck, gauge, spinorTmp, clover, clover_inv, inv_param.kappa, inv_param.mu,
469  inv_param.twist_flavor, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
470  } else if (dslash_type == QUDA_WILSON_DSLASH) {
471  wil_matpc(spinorTmp, gauge, spinorOutMulti[i], inv_param.kappa, inv_param.matpc_type, 0,
472  inv_param.cpu_prec, gauge_param);
473  wil_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
474  inv_param.cpu_prec, gauge_param);
475  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
476  clover_matpc(spinorTmp, gauge, clover, clover_inv, spinorOutMulti[i], inv_param.kappa, inv_param.matpc_type, 0,
477  inv_param.cpu_prec, gauge_param);
478  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
479  inv_param.cpu_prec, gauge_param);
480  } else {
481  printfQuda("Domain wall not supported for multi-shift\n");
482  exit(-1);
483  }
484 
485  axpy(inv_param.offset[i], spinorOutMulti[i], spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
486  mxpy(spinorIn, spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
487  double nrm2 = norm_2(spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
488  double src2 = norm_2(spinorIn, Vh*spinorSiteSize, inv_param.cpu_prec);
489  double l2r = sqrt(nrm2 / src2);
490 
491  printfQuda("Shift %d residuals: (L2 relative) tol %g, QUDA = %g, host = %g; (heavy-quark) tol %g, QUDA = %g\n",
492  i, inv_param.tol_offset[i], inv_param.true_res_offset[i], l2r,
493  inv_param.tol_hq_offset[i], inv_param.true_res_hq_offset[i]);
494  }
495  free(spinorTmp);
496 
497  } else {
498 
499  if (inv_param.solution_type == QUDA_MAT_SOLUTION) {
500 
502  if(inv_param.twist_flavor == QUDA_TWIST_SINGLET) {
503  tm_mat(spinorCheck, gauge, spinorOut, inv_param.kappa, inv_param.mu, inv_param.twist_flavor, 0, inv_param.cpu_prec, gauge_param);
504  } else {
505  int tm_offset = V*spinorSiteSize;
506  void *evenOut = spinorCheck;
507  void *oddOut = (char*)evenOut + tm_offset*cpu_prec;
508 
509  void *evenIn = spinorOut;
510  void *oddIn = (char*)evenIn + tm_offset*cpu_prec;
511 
512  tm_ndeg_mat(evenOut, oddOut, gauge, evenIn, oddIn, inv_param.kappa, inv_param.mu, inv_param.epsilon, 0, inv_param.cpu_prec, gauge_param);
513  }
514  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
515  tmc_mat(spinorCheck, gauge, clover, spinorOut, inv_param.kappa, inv_param.mu, inv_param.twist_flavor, 0,
516  inv_param.cpu_prec, gauge_param);
517  } else if (dslash_type == QUDA_WILSON_DSLASH) {
518  wil_mat(spinorCheck, gauge, spinorOut, inv_param.kappa, 0, inv_param.cpu_prec, gauge_param);
519  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
520  clover_mat(spinorCheck, gauge, clover, spinorOut, inv_param.kappa, 0, inv_param.cpu_prec, gauge_param);
521  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
522  dw_mat(spinorCheck, gauge, spinorOut, kappa5, inv_param.dagger, inv_param.cpu_prec, gauge_param, inv_param.mass);
523  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
524  dw_4d_mat(spinorCheck, gauge, spinorOut, kappa5, inv_param.dagger, inv_param.cpu_prec, gauge_param, inv_param.mass);
525  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
526  double _Complex *kappa_b = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
527  double _Complex *kappa_c = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
528  for(int xs = 0 ; xs < Lsdim ; xs++)
529  {
530  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
531  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
532  }
533  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);
534  free(kappa_b);
535  free(kappa_c);
536  } else {
537  errorQuda("Unsupported dslash_type");
538  }
539  if (inv_param.mass_normalization == QUDA_MASS_NORMALIZATION) {
543  ax(0.5/kappa5, spinorCheck, V*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
545  ax(0.5/inv_param.kappa, spinorCheck, 2*V*spinorSiteSize, inv_param.cpu_prec);
546  } else {
547  ax(0.5/inv_param.kappa, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
548  }
549  }
550 
551  } else if(inv_param.solution_type == QUDA_MATPC_SOLUTION) {
552 
554  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET) {
555  int tm_offset = Vh*spinorSiteSize;
556  void *out0 = spinorCheck;
557  void *out1 = (char*)out0 + tm_offset*cpu_prec;
558 
559  void *in0 = spinorOut;
560  void *in1 = (char*)in0 + tm_offset*cpu_prec;
561 
562  tm_ndeg_matpc(out0, out1, gauge, in0, in1, inv_param.kappa, inv_param.mu, inv_param.epsilon, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
563  } else {
564  tm_matpc(spinorCheck, gauge, spinorOut, inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
565  inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
566  }
567  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
568  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET)
569  errorQuda("Twisted mass solution type not supported");
570  tmc_matpc(spinorCheck, gauge, spinorOut, clover, clover_inv, inv_param.kappa, inv_param.mu,
571  inv_param.twist_flavor, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
572  } else if (dslash_type == QUDA_WILSON_DSLASH) {
573  wil_matpc(spinorCheck, gauge, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
574  inv_param.cpu_prec, gauge_param);
575  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
576  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
577  inv_param.cpu_prec, gauge_param);
578  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
579  dw_matpc(spinorCheck, gauge, spinorOut, kappa5, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param, inv_param.mass);
580  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
581  dw_4d_matpc(spinorCheck, gauge, spinorOut, kappa5, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param, inv_param.mass);
582  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
583  double _Complex *kappa_b = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
584  double _Complex *kappa_c = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
585  for(int xs = 0 ; xs < Lsdim ; xs++)
586  {
587  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
588  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
589  }
590  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);
591  free(kappa_b);
592  free(kappa_c);
593  } else {
594  errorQuda("Unsupported dslash_type");
595  }
596 
597  if (inv_param.mass_normalization == QUDA_MASS_NORMALIZATION) {
601  ax(0.25/(kappa5*kappa5), spinorCheck, V*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
602  } else {
603  ax(0.25/(inv_param.kappa*inv_param.kappa), spinorCheck, Vh*spinorSiteSize, inv_param.cpu_prec);
604 
605  }
606  }
607 
608  } else if (inv_param.solution_type == QUDA_MATPCDAG_MATPC_SOLUTION) {
609 
610  void *spinorTmp = malloc(V*spinorSiteSize*sSize*inv_param.Ls);
611 
612  ax(0, spinorCheck, V*spinorSiteSize, inv_param.cpu_prec);
613 
615  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET) {
616  int tm_offset = Vh*spinorSiteSize;
617  void *out0 = spinorCheck;
618  void *out1 = (char*)out0 + tm_offset*cpu_prec;
619 
620  void *tmp0 = spinorTmp;
621  void *tmp1 = (char*)tmp0 + tm_offset*cpu_prec;
622 
623  void *in0 = spinorOut;
624  void *in1 = (char*)in0 + tm_offset*cpu_prec;
625 
626  tm_ndeg_matpc(tmp0, tmp1, gauge, in0, in1, inv_param.kappa, inv_param.mu, inv_param.epsilon, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
627  tm_ndeg_matpc(out0, out1, gauge, tmp0, tmp1, inv_param.kappa, inv_param.mu, inv_param.epsilon, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
628  } else {
629  tm_matpc(spinorTmp, gauge, spinorOut, inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
630  inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
631  tm_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.mu, inv_param.twist_flavor,
632  inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
633  }
634  } else if (dslash_type == QUDA_TWISTED_CLOVER_DSLASH) {
635  if (inv_param.twist_flavor != QUDA_TWIST_SINGLET)
636  errorQuda("Twisted mass solution type not supported");
637  tmc_matpc(spinorTmp, gauge, spinorOut, clover, clover_inv, inv_param.kappa, inv_param.mu,
638  inv_param.twist_flavor, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
639  tmc_matpc(spinorCheck, gauge, spinorTmp, clover, clover_inv, inv_param.kappa, inv_param.mu,
640  inv_param.twist_flavor, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
641  } else if (dslash_type == QUDA_WILSON_DSLASH) {
642  wil_matpc(spinorTmp, gauge, spinorOut, inv_param.kappa, inv_param.matpc_type, 0,
643  inv_param.cpu_prec, gauge_param);
644  wil_matpc(spinorCheck, gauge, spinorTmp, inv_param.kappa, inv_param.matpc_type, 1,
645  inv_param.cpu_prec, gauge_param);
646  } else if (dslash_type == QUDA_CLOVER_WILSON_DSLASH) {
647  clover_matpc(spinorTmp, gauge, clover, clover_inv, spinorOut, inv_param.kappa,
648  inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param);
649  clover_matpc(spinorCheck, gauge, clover, clover_inv, spinorTmp, inv_param.kappa,
650  inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param);
651  } else if (dslash_type == QUDA_DOMAIN_WALL_DSLASH) {
652  dw_matpc(spinorTmp, gauge, spinorOut, kappa5, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param, inv_param.mass);
653  dw_matpc(spinorCheck, gauge, spinorTmp, kappa5, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param, inv_param.mass);
654  } else if (dslash_type == QUDA_DOMAIN_WALL_4D_DSLASH) {
655  dw_4d_matpc(spinorTmp, gauge, spinorOut, kappa5, inv_param.matpc_type, 0, inv_param.cpu_prec, gauge_param, inv_param.mass);
656  dw_4d_matpc(spinorCheck, gauge, spinorTmp, kappa5, inv_param.matpc_type, 1, inv_param.cpu_prec, gauge_param, inv_param.mass);
657  } else if (dslash_type == QUDA_MOBIUS_DWF_DSLASH) {
658  double _Complex *kappa_b = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
659  double _Complex *kappa_c = (double _Complex *)malloc(Lsdim * sizeof(double _Complex));
660  for(int xs = 0 ; xs < Lsdim ; xs++)
661  {
662  kappa_b[xs] = 1.0/(2*(inv_param.b_5[xs]*(4.0 + inv_param.m5) + 1.0));
663  kappa_c[xs] = 1.0/(2*(inv_param.c_5[xs]*(4.0 + inv_param.m5) - 1.0));
664  }
665  mdw_matpc(spinorTmp, 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);
666  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);
667  free(kappa_b);
668  free(kappa_c);
669  } else {
670  errorQuda("Unsupported dslash_type");
671  }
672 
673  if (inv_param.mass_normalization == QUDA_MASS_NORMALIZATION) {
674  errorQuda("Mass normalization not implemented");
675  }
676 
677  free(spinorTmp);
678  }
679 
680 
681  int vol = inv_param.solution_type == QUDA_MAT_SOLUTION ? V : Vh;
682  mxpy(spinorIn, spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
683  double nrm2 = norm_2(spinorCheck, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
684  double src2 = norm_2(spinorIn, vol*spinorSiteSize*inv_param.Ls, inv_param.cpu_prec);
685  double l2r = sqrt(nrm2 / src2);
686 
687  printfQuda("Residuals: (L2 relative) tol %g, QUDA = %g, host = %g; (heavy-quark) tol %g, QUDA = %g\n",
688  inv_param.tol, inv_param.true_res, l2r, inv_param.tol_hq, inv_param.true_res_hq);
689 
690  }
691 
692  freeGaugeQuda();
694 
695  // finalize the QUDA library
696  endQuda();
697 
698  finalizeComms();
699 
701  if (clover) free(clover);
702  if (clover_inv) free(clover_inv);
703  }
704 
705  for (int dir = 0; dir<4; dir++) free(gauge[dir]);
706 
707  return 0;
708 }
double epsilon
Definition: test_util.cpp:1649
int maxiter_precondition
Definition: quda.h:292
int device
Definition: test_util.cpp:1602
static size_t gSize
double secs
Definition: quda.h:251
int dimPartitioned(int dim)
Definition: test_util.cpp:1776
void ax(double a, ColorSpinorField &x)
Definition: blas_quda.cu:508
QudaDiracFieldOrder dirac_order
Definition: quda.h:219
QudaMassNormalization mass_normalization
Definition: quda.h:208
double tol_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:182
enum QudaMassNormalization_s QudaMassNormalization
double anisotropy
Definition: test_util.cpp:1650
QudaReconstructType reconstruct_sloppy
Definition: quda.h:53
double anisotropy
Definition: quda.h:38
void freeCloverQuda(void)
QudaTwistFlavorType twist_flavor
Definition: test_util.cpp:1660
void mdw_matpc(void *out, void **gauge, void *in, double _Complex *kappa_b, double _Complex *kappa_c, QudaMatPCType matpc_type, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm, double _Complex *b5, double _Complex *c5)
cudaColorSpinorField * tmp1
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)
double ca_lambda_max
Definition: test_util.cpp:1633
void endQuda(void)
void construct_gauge_field(void **gauge, int type, QudaPrecision precision, QudaGaugeParam *param)
Definition: test_util.cpp:1047
QudaCABasis ca_basis
Definition: quda.h:298
QudaSolveType solve_type
Definition: quda.h:205
QudaVerbosity verbosity_precondition
Definition: quda.h:286
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:63
double_complex c_5[QUDA_MAX_DWF_LS]
Definition: quda.h:112
int multishift
Definition: test_util.cpp:1642
void mdw_mat(void *out, void **gauge, void *in, double _Complex *kappa_b, double _Complex *kappa_c, int dagger, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm, double _Complex *b5, double _Complex *c5)
int main(int argc, char **argv)
void dw_setDims(int *X, const int L5)
Definition: test_util.cpp:187
double mu
Definition: quda.h:114
QudaGaugeFixed gauge_fix
Definition: quda.h:61
QudaSchwarzType schwarz_type
Definition: quda.h:310
__host__ __device__ ValueType norm(const complex< ValueType > &z)
Returns the magnitude of z squared.
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:270
int Lsdim
Definition: test_util.cpp:1619
QudaLinkType type
Definition: quda.h:42
double kappa
Definition: quda.h:106
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
#define errorQuda(...)
Definition: util_quda.h:121
double tol
Definition: quda.h:121
QudaPrecision prec
Definition: test_util.cpp:1608
void usage(char **)
Definition: test_util.cpp:1783
QudaDslashType dslash_type
Definition: quda.h:102
QudaReconstructType reconstruct_precondition
Definition: quda.h:59
QudaInverterType inv_type
Definition: quda.h:103
QudaPrecision cuda_prec
Definition: quda.h:214
#define cloverSiteSize
Definition: test_util.h:9
int return_clover_inverse
Definition: quda.h:242
enum QudaSolveType_s QudaSolveType
__host__ __device__ ValueType sqrt(ValueType x)
Definition: complex_quda.h:120
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
QudaVerbosity mg_verbosity[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1675
bool compute_clover
Definition: test_util.cpp:1654
QudaPrecision cpu_prec
Definition: quda.h:213
QudaPrecision & cuda_prec
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
double tol_hq
Definition: test_util.cpp:1657
QudaReconstructType link_recon_precondition
Definition: test_util.cpp:1607
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 plaqQuda(double plaq[3])
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:207
QudaPrecision & cuda_prec_refinement_sloppy
const char * get_matpc_str(QudaMatPCType type)
Definition: misc.cpp:1121
void finalizeComms()
Definition: test_util.cpp:128
double reliable_delta
Definition: test_util.cpp:1658
QudaGaugeParam gauge_param
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:216
double clover_coeff
Definition: test_util.cpp:1653
QudaPrecision clover_cuda_prec_refinement_sloppy
Definition: quda.h:227
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
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:126
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:1605
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:701
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
int return_clover
Definition: quda.h:241
void construct_spinor_source(void *v, int nSpin, int nColor, QudaPrecision precision, const int *const x, quda::RNG &rng)
Definition: test_util.cpp:1342
void display_test_info()
Definition: invert_test.cpp:83
double ca_lambda_max
Definition: quda.h:304
#define spinorSiteSize
QudaPrecision & cuda_prec_precondition
bool alternative_reliable
Definition: test_util.cpp:1659
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:226
void setDims(int *)
Definition: test_util.cpp:151
QudaFieldLocation input_location
Definition: quda.h:99
void freeGaugeQuda(void)
double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:191
double reliable_delta
Definition: quda.h:129
int solution_accumulator_pipeline
Definition: quda.h:142
double mu
Definition: test_util.cpp:1648
int pipeline
Definition: test_util.cpp:1634
double ca_lambda_min
Definition: quda.h:301
double_complex b_5[QUDA_MAX_DWF_LS]
Definition: quda.h:111
QudaSolutionType solution_type
Definition: quda.h:204
int use_alternative_reliable
Definition: quda.h:131
QudaSolverNormalization solver_normalization
Definition: quda.h:209
QudaPrecision clover_cuda_prec
Definition: quda.h:225
int precondition_cycle
Definition: quda.h:307
QudaPrecision & cuda_prec_sloppy
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)
QudaSolutionType solution_type
Definition: test_util.cpp:1664
QudaFieldLocation output_location
Definition: quda.h:100
QudaPrecision prec_refinement_sloppy
Definition: test_util.cpp:1610
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:228
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
double m5
Definition: quda.h:108
QudaPrecision cuda_prec_sloppy
Definition: quda.h:215
QudaMassNormalization normalization
Definition: test_util.cpp:1661
QudaVerbosity verbosity
Definition: quda.h:244
void setSpinorSiteSize(int n)
Definition: test_util.cpp:211
const char * get_mass_normalization_str(QudaMassNormalization type)
Definition: misc.cpp:1077
int solution_accumulator_pipeline
Definition: test_util.cpp:1635
double tol_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:179
double true_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:185
void axpy(double a, ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.h:35
QudaInvertParam newQudaInvertParam(void)
double gflops
Definition: quda.h:250
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:768
QudaMatPCType matpc_type
Definition: test_util.cpp:1662
QudaPrecision cuda_prec_precondition
Definition: quda.h:58
QudaCloverFieldOrder clover_order
Definition: quda.h:230
Class declaration to initialize and hold CURAND RNG states.
Definition: random_quda.h:23
void dw_4d_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param, double mferm)
double tol_hq
Definition: quda.h:123
enum QudaMatPCType_s QudaMatPCType
cpuColorSpinorField * spinorOut
Definition: covdev_test.cpp:41
double true_res_hq
Definition: quda.h:127
enum QudaSolutionType_s QudaSolutionType
QudaGammaBasis gamma_basis
Definition: quda.h:221
QudaPrecision cuda_prec_sloppy
Definition: quda.h:52
const char * get_dslash_str(QudaDslashType type)
Definition: misc.cpp:910
double tol_precondition
Definition: quda.h:289
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:176
int use_sloppy_partial_accumulator
Definition: quda.h:132
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:50
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
double tol
Definition: test_util.cpp:1656
double mass
Definition: quda.h:105
int gcrNkrylov
Definition: quda.h:259
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:55
int V
Definition: test_util.cpp:27
double norm_2(void *v, int len, QudaPrecision precision)
int compute_clover_inverse
Definition: quda.h:240
bool unit_gauge
Definition: test_util.cpp:1624
void construct_clover_field(void *clover, double norm, double diag, QudaPrecision precision)
Definition: test_util.cpp:1167
int max_res_increase
Definition: quda.h:147
double ca_lambda_min
Definition: test_util.cpp:1632
QudaResidualType_s
Definition: enum_quda.h:186
int xdim
Definition: test_util.cpp:1615
QudaSolveType solve_type
Definition: test_util.cpp:1663
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:1640
QudaPrecision cuda_prec_precondition
Definition: quda.h:217
double tol_restart
Definition: quda.h:122
int zdim
Definition: test_util.cpp:1617
int tdim
Definition: test_util.cpp:1618
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
char latfile[]
Definition: test_util.cpp:1623
enum QudaCABasis_s QudaCABasis
#define QUDA_MAX_MG_LEVEL
Maximum number of multi-grid levels. This number may be increased if needed.
QudaCABasis ca_basis
Definition: test_util.cpp:1631
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)
#define printfQuda(...)
Definition: util_quda.h:115
QudaTboundary t_boundary
Definition: quda.h:45
QudaTwistFlavorType twist_flavor
Definition: quda.h:117
#define MAX(a, b)
Definition: invert_test.cpp:19
QudaVerbosity verbosity
Definition: test_util.cpp:1614
enum QudaDslashType_s QudaDslashType
int niter
Definition: test_util.cpp:1629
void mxpy(ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.h:34
QudaResidualType residual_type
Definition: quda.h:320
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:169
int Nsrc
Definition: test_util.cpp:1627
enum QudaVerbosity_s QudaVerbosity
QudaPrecision & cpu_prec
QudaDslashType dslash_type
Definition: test_util.cpp:1621
int compute_clover
Definition: quda.h:239
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1606
double epsilon
Definition: quda.h:115
int ydim
Definition: test_util.cpp:1616
cpuColorSpinorField * spinorTmp
double omega
Definition: quda.h:295
double mass
Definition: test_util.cpp:1646
QudaPrecision prec_sloppy
Definition: test_util.cpp:1609
QudaInverterType precon_type
Definition: test_util.cpp:1641
void initComms(int argc, char **argv, int *const commDims)
Definition: test_util.cpp:88
void read_gauge_field(const char *filename, void *gauge[], QudaPrecision prec, const int *X, int argc, char *argv[])
Definition: qio_field.h:14
int gcrNkrylov
Definition: test_util.cpp:1630
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
QudaPrecision clover_cpu_prec
Definition: quda.h:224
QudaMatPCType matpc_type
Definition: quda.h:206
QudaPrecision prec_precondition
Definition: test_util.cpp:1611
enum QudaInverterType_s QudaInverterType
double kappa
Definition: test_util.cpp:1647
double kappa5
QudaReconstructType reconstruct_refinement_sloppy
Definition: quda.h:56
QudaPrecision cpu_prec
Definition: quda.h:47
#define gaugeSiteSize
Definition: face_gauge.cpp:34
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:211
double clover_coeff
Definition: quda.h:233
int Vh
Definition: test_util.cpp:28
enum QudaTwistFlavorType_s QudaTwistFlavorType