QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
staggered_eigensolve_test.cpp
Go to the documentation of this file.
1 #include <iostream>
2 #include <stdio.h>
3 #include <stdlib.h>
4 #include <string.h>
5 
6 #include <quda_internal.h>
7 #include <dirac_quda.h>
8 #include <dslash_quda.h>
9 #include <invert_quda.h>
10 #include <util_quda.h>
11 #include <blas_quda.h>
12 
13 #include <misc.h>
14 #include <test_util.h>
15 #include <dslash_util.h>
16 #include <staggered_gauge_utils.h>
17 #include <unitarization_links.h>
18 #include <llfat_reference.h>
19 #include <gauge_field.h>
20 #include <unitarization_links.h>
21 #include <blas_reference.h>
22 
23 #if defined(QMP_COMMS)
24 #include <qmp.h>
25 #elif defined(MPI_COMMS)
26 #include <mpi.h>
27 #endif
28 
29 #include <qio_field.h>
30 
31 #define MAX(a, b) ((a) > (b) ? (a) : (b))
32 
33 // In a typical application, quda.h is the only QUDA header required.
34 #include <quda.h>
35 
36 #define mySpinorSiteSize 6
37 
38 extern void usage(char **argv);
39 
41 
42 extern int device;
43 
45 size_t gSize = sizeof(double);
46 
47 extern double reliable_delta;
48 extern bool alternative_reliable;
49 extern int test_type;
50 extern int xdim;
51 extern int ydim;
52 extern int zdim;
53 extern int tdim;
54 extern int gridsize_from_cmdline[];
56 extern QudaPrecision prec;
60 extern double mass; // the mass of the Dirac operator
61 extern double kappa;
62 extern int laplace3D;
63 extern double tol; // tolerance for inverter
64 extern double tol_hq; // heavy-quark tolerance for inverter
65 extern char latfile[];
66 extern int Nsrc; // number of spinors to apply to simultaneously
67 extern int niter;
68 extern int gcrNkrylov;
69 extern int pipeline; // length of pipeline for fused operations in GCR or BiCGstab-l
70 extern int solution_accumulator_pipeline; // length of pipeline for fused solution update from the direction vectors
71 extern QudaCABasis ca_basis; // basis for CA-CG solves
72 extern double ca_lambda_min; // minimum eigenvalue for scaling Chebyshev CA-CG solves
73 extern double ca_lambda_max; // maximum eigenvalue for scaling Chebyshev CA-CG solves
74 
75 // Dirac operator type
77 extern QudaMatPCType matpc_type; // preconditioning type
78 extern QudaSolutionType solution_type; // solution type
80 
81 extern bool compute_fatlong; // build the true fat/long links or use random numbers
82 // relativistic correction for naik term
83 extern double eps_naik;
84 // Number of naiks. If eps_naik is 0.0, we only need
85 // to construct one naik.
86 static int n_naiks = 1;
87 
88 // For loading the gauge fields
90 char **argv_copy;
91 
92 // eigensolver params
93 extern int eig_nEv;
94 extern int eig_nKr;
95 extern int eig_nConv;
96 extern bool eig_require_convergence;
97 extern int eig_check_interval;
98 extern int eig_max_restarts;
99 extern double eig_tol;
100 extern int eig_maxiter;
101 extern bool eig_use_poly_acc;
102 extern int eig_poly_deg;
103 extern double eig_amin;
104 extern double eig_amax;
105 extern bool eig_use_normop;
106 extern bool eig_use_dagger;
107 extern bool eig_compute_svd;
109 extern QudaEigType eig_type;
110 extern bool eig_arpack_check;
111 extern char eig_arpack_logfile[];
112 extern char eig_QUDA_logfile[];
113 extern char eig_vec_infile[];
114 extern char eig_vec_outfile[];
115 
116 extern bool verify_results;
117 
118 int X[4];
119 
121 {
122  printfQuda("running the following test:\n");
123  printfQuda("prec sloppy_prec link_recon sloppy_link_recon test_type S_dimension T_dimension\n");
124  printfQuda("%s %s %s %s %s %d/%d/%d %d \n", get_prec_str(prec),
127 
128  printfQuda("Grid partition info: X Y Z T\n");
129  printfQuda(" %d %d %d %d\n", dimPartitioned(0), dimPartitioned(1), dimPartitioned(2),
130  dimPartitioned(3));
131 
132  return;
133 }
134 
135 void usage_extra(char **argv)
136 {
137  printfQuda("Extra options:\n");
138  printfQuda(" --test <0/3/4> # Test method\n");
139  printfQuda(" 0: Full parity inverter\n");
140  printfQuda(" 3: Even even spinor CG inverter\n");
141  printfQuda(" 4: Odd odd spinor CG inverter\n");
142  return;
143 }
144 
146 {
147  gauge_param.X[0] = xdim;
148  gauge_param.X[1] = ydim;
149  gauge_param.X[2] = zdim;
150  gauge_param.X[3] = tdim;
151 
152  gauge_param.cpu_prec = cpu_prec;
153  gauge_param.cuda_prec = prec;
154  gauge_param.reconstruct = link_recon;
155  gauge_param.cuda_prec_sloppy = prec_sloppy;
158 
161 
162  gauge_param.anisotropy = 1.0;
163 
164  // For HISQ, this must always be set to 1.0, since the tadpole
165  // correction is baked into the coefficients for the first fattening.
166  // The tadpole doesn't mean anything for the second fattening
167  // since the input fields are unitarized.
168  gauge_param.tadpole_coeff = 1.0;
169 
171  gauge_param.scale = -1.0 / 24.0;
172  if (eps_naik != 0) { gauge_param.scale *= (1.0 + eps_naik); }
173  } else {
174  gauge_param.scale = 1.0;
175  }
176  gauge_param.gauge_order = QUDA_MILC_GAUGE_ORDER;
177  gauge_param.t_boundary = QUDA_ANTI_PERIODIC_T;
179  gauge_param.gauge_fix = QUDA_GAUGE_FIXED_NO;
180  gauge_param.type = QUDA_WILSON_LINKS;
181 
182  gauge_param.ga_pad = 0;
183 
184 #ifdef MULTI_GPU
185  int x_face_size = gauge_param.X[1] * gauge_param.X[2] * gauge_param.X[3] / 2;
186  int y_face_size = gauge_param.X[0] * gauge_param.X[2] * gauge_param.X[3] / 2;
187  int z_face_size = gauge_param.X[0] * gauge_param.X[1] * gauge_param.X[3] / 2;
188  int t_face_size = gauge_param.X[0] * gauge_param.X[1] * gauge_param.X[2] / 2;
189  int pad_size = MAX(x_face_size, y_face_size);
190  pad_size = MAX(pad_size, z_face_size);
191  pad_size = MAX(pad_size, t_face_size);
192  gauge_param.ga_pad = pad_size;
193 #endif
194 }
195 
197 {
198  // Solver params
199  inv_param.verbosity = QUDA_VERBOSE;
200  inv_param.mass = mass;
201  inv_param.kappa = kappa = 1.0 / (8.0 + mass); // for Laplace operator
202  inv_param.laplace3D = laplace3D; // for Laplace operator
203 
204  // outer solver parameters
205  inv_param.inv_type = QUDA_BICGSTAB_INVERTER; // Dummy setting
206  inv_param.tol = tol;
207  inv_param.tol_restart = 1e-3;
208  inv_param.maxiter = niter;
209  inv_param.reliable_delta = reliable_delta;
211  inv_param.use_sloppy_partial_accumulator = false;
213  inv_param.pipeline = pipeline;
214 
215  inv_param.Ls = Nsrc;
216 
217  if (tol_hq == 0 && tol == 0) {
218  errorQuda("qudaInvert: requesting zero residual\n");
219  exit(1);
220  }
221  // require both L2 relative and heavy quark residual to determine convergence
222  inv_param.residual_type = static_cast<QudaResidualType_s>(0);
223  inv_param.residual_type = (tol != 0) ?
224  static_cast<QudaResidualType_s>(inv_param.residual_type | QUDA_L2_RELATIVE_RESIDUAL) :
225  inv_param.residual_type;
226  inv_param.residual_type = (tol_hq != 0) ?
227  static_cast<QudaResidualType_s>(inv_param.residual_type | QUDA_HEAVY_QUARK_RESIDUAL) :
228  inv_param.residual_type;
229  inv_param.heavy_quark_check = (inv_param.residual_type & QUDA_HEAVY_QUARK_RESIDUAL ? 5 : 0);
230 
231  inv_param.tol_hq = tol_hq; // specify a tolerance for the residual for heavy quark residual
232 
233  inv_param.Nsteps = 2;
234 
235  // domain decomposition preconditioner parameters
237  inv_param.tol_precondition = 1e-1;
238  inv_param.maxiter_precondition = 10;
240  inv_param.cuda_prec_precondition = inv_param.cuda_prec_sloppy;
241 
242  // Specify Krylov sub-size for GCR, BICGSTAB(L), basis size for CA-CG, CA-GCR
243  inv_param.gcrNkrylov = gcrNkrylov;
244 
245  // Specify basis for CA-CG, lambda min/max for Chebyshev basis
246  // lambda_max < lambda_max . use power iters to generate
247  inv_param.ca_basis = ca_basis;
248  inv_param.ca_lambda_min = ca_lambda_min;
249  inv_param.ca_lambda_max = ca_lambda_max;
250 
251  inv_param.solution_type = solution_type;
252  inv_param.solve_type = solve_type;
253  inv_param.matpc_type = matpc_type;
254  inv_param.dagger = QUDA_DAG_NO;
256 
257  inv_param.cpu_prec = cpu_prec;
258  inv_param.cuda_prec = prec;
259  inv_param.cuda_prec_sloppy = prec_sloppy;
262  inv_param.gamma_basis = QUDA_DEGRAND_ROSSI_GAMMA_BASIS; // this is meaningless, but must be thus set
263  inv_param.dirac_order = QUDA_DIRAC_ORDER;
264 
265  inv_param.dslash_type = dslash_type;
266 
269 
270  int tmpint = MAX(X[1] * X[2] * X[3], X[0] * X[2] * X[3]);
271  tmpint = MAX(tmpint, X[0] * X[1] * X[3]);
272  tmpint = MAX(tmpint, X[0] * X[1] * X[2]);
273 
274  inv_param.sp_pad = tmpint;
275 }
276 
277 // Parameters defining the eigensolver
278 void setEigParam(QudaEigParam &eig_param)
279 {
280  eig_param.eig_type = eig_type;
281  eig_param.spectrum = eig_spectrum;
284  errorQuda("Only real spectrum type (LR or SR) can be passed to Lanczos type solver");
285  }
286 
287  // The solver will exit when nConv extremal eigenpairs have converged
288  if (eig_nConv < 0) {
289  eig_param.nConv = eig_nEv;
290  eig_nConv = eig_nEv;
291  } else {
292  eig_param.nConv = eig_nConv;
293  }
294 
295  eig_param.nEv = eig_nEv;
296  eig_param.nKr = eig_nKr;
297  eig_param.tol = eig_tol;
300  eig_param.max_restarts = eig_max_restarts;
301  eig_param.cuda_prec_ritz = prec;
302 
306  if (eig_compute_svd) {
307  eig_param.use_dagger = QUDA_BOOLEAN_FALSE;
308  eig_param.use_norm_op = QUDA_BOOLEAN_TRUE;
309  }
310 
312  eig_param.poly_deg = eig_poly_deg;
313  eig_param.a_min = eig_amin;
314  eig_param.a_max = eig_amax;
315 
317  strcpy(eig_param.arpack_logfile, eig_arpack_logfile);
318  strcpy(eig_param.QUDA_logfile, eig_QUDA_logfile);
319 
320  strcpy(eig_param.vec_infile, eig_vec_infile);
321  strcpy(eig_param.vec_outfile, eig_vec_outfile);
322 }
323 
325 {
326 
328  setGaugeParam(gauge_param);
329 
330  QudaEigParam eig_param = newQudaEigParam();
331  // Though no inversions are performed, the inv_param
332  // structure contains all the information we need to
333  // construct the dirac operator. We encapsualte the
334  // inv_param structure inside the eig_param structure
335  // to avoid any confusion
336  QudaInvertParam eig_inv_param = newQudaInvertParam();
337  setInvertParam(eig_inv_param);
338  eig_param.invert_param = &eig_inv_param;
339  setEigParam(eig_param);
340 
341  // this must be before the FaceBuffer is created (this is because it allocates pinned memory - FIXME)
342  initQuda(device);
343 
344  setDims(gauge_param.X);
345  dw_setDims(gauge_param.X, Nsrc); // so we can use 5-d indexing from dwf
347 
348  void *qdp_inlink[4] = {nullptr, nullptr, nullptr, nullptr};
349  void *qdp_fatlink[4] = {nullptr, nullptr, nullptr, nullptr};
350  void *qdp_longlink[4] = {nullptr, nullptr, nullptr, nullptr};
351  void *milc_fatlink = nullptr;
352  void *milc_longlink = nullptr;
353 
354  gSize = (gauge_param.cpu_prec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
355 
356  for (int dir = 0; dir < 4; dir++) {
357  qdp_inlink[dir] = malloc(V * gaugeSiteSize * gSize);
358  qdp_fatlink[dir] = malloc(V * gaugeSiteSize * gSize);
359  qdp_longlink[dir] = malloc(V * gaugeSiteSize * gSize);
360  }
361  milc_fatlink = malloc(4 * V * gaugeSiteSize * gSize);
362  milc_longlink = malloc(4 * V * gaugeSiteSize * gSize);
363 
364  // for load, etc
365  gauge_param.reconstruct = QUDA_RECONSTRUCT_NO;
366 
367  // load a field WITHOUT PHASES
368  if (strcmp(latfile, "")) {
369  read_gauge_field(latfile, qdp_inlink, gauge_param.cpu_prec, gauge_param.X, argc_copy, argv_copy);
371  applyGaugeFieldScaling_long(qdp_inlink, Vh, &gauge_param, QUDA_STAGGERED_DSLASH, gauge_param.cpu_prec);
372  }
373  } else {
375  construct_gauge_field(qdp_inlink, 1, gauge_param.cpu_prec, &gauge_param);
376  } else {
377  construct_fat_long_gauge_field(qdp_inlink, qdp_longlink, 1, gauge_param.cpu_prec, &gauge_param,
379  }
380  }
381 
382  // Compute plaquette. Routine is aware that the gauge fields already have the phases on them.
383  double plaq[3];
384  computeStaggeredPlaquetteQDPOrder(qdp_inlink, plaq, gauge_param, dslash_type);
385 
386  printfQuda("Computed plaquette is %e (spatial = %e, temporal = %e)\n", plaq[0], plaq[1], plaq[2]);
387 
388  // QUDA_STAGGERED_DSLASH follows the same codepath whether or not you
389  // "compute" the fat/long links or not.
391  for (int dir = 0; dir < 4; dir++) {
392  memcpy(qdp_fatlink[dir], qdp_inlink[dir], V * gaugeSiteSize * gSize);
393  memset(qdp_longlink[dir], 0, V * gaugeSiteSize * gSize);
394  }
395  } else { // QUDA_ASQTAD_DSLASH
396 
397  if (compute_fatlong) {
398  computeFatLongGPU(qdp_fatlink, qdp_longlink, qdp_inlink, gauge_param, gSize, n_naiks, eps_naik);
399  } else {
400  for (int dir = 0; dir < 4; dir++) { memcpy(qdp_fatlink[dir], qdp_inlink[dir], V * gaugeSiteSize * gSize); }
401  }
402 
403  // Compute fat link plaquette.
404  computeStaggeredPlaquetteQDPOrder(qdp_fatlink, plaq, gauge_param, dslash_type);
405 
406  printfQuda("Computed fat link plaquette is %e (spatial = %e, temporal = %e)\n", plaq[0], plaq[1], plaq[2]);
407  }
408 
409  reorderQDPtoMILC(milc_fatlink, qdp_fatlink, V, gaugeSiteSize, gauge_param.cpu_prec, gauge_param.cpu_prec);
410  reorderQDPtoMILC(milc_longlink, qdp_longlink, V, gaugeSiteSize, gauge_param.cpu_prec, gauge_param.cpu_prec);
411 
412 #ifdef MULTI_GPU
413  int tmp_value = MAX(ydim * zdim * tdim / 2, xdim * zdim * tdim / 2);
414  tmp_value = MAX(tmp_value, xdim * ydim * tdim / 2);
415  tmp_value = MAX(tmp_value, xdim * ydim * zdim / 2);
416 
417  int fat_pad = tmp_value;
418  int link_pad = 3 * tmp_value;
419 
420  // FIXME: currently assume staggered is SU(3)
424  gauge_param.reconstruct = QUDA_RECONSTRUCT_NO;
425  GaugeFieldParam cpuFatParam(milc_fatlink, gauge_param);
427  cpuGaugeField *cpuFat = new cpuGaugeField(cpuFatParam);
428  ghost_fatlink = (void **)cpuFat->Ghost();
429 
430  gauge_param.type = QUDA_ASQTAD_LONG_LINKS;
431  GaugeFieldParam cpuLongParam(milc_longlink, gauge_param);
432  cpuLongParam.ghostExchange = QUDA_GHOST_EXCHANGE_PAD;
433  cpuGaugeField *cpuLong = new cpuGaugeField(cpuLongParam);
434  ghost_longlink = (void **)cpuLong->Ghost();
435 
436 #else
437  int fat_pad = 0;
438  int link_pad = 0;
439 #endif
440 
444  gauge_param.ga_pad = fat_pad;
446  gauge_param.reconstruct = link_recon;
448  } else {
449  gauge_param.reconstruct = gauge_param.reconstruct_sloppy = QUDA_RECONSTRUCT_NO;
450  }
451  gauge_param.cuda_prec_precondition = gauge_param.cuda_prec_sloppy;
452  gauge_param.reconstruct_precondition = gauge_param.reconstruct_sloppy;
453  loadGaugeQuda(milc_fatlink, &gauge_param);
454 
456  gauge_param.type = QUDA_ASQTAD_LONG_LINKS;
457  gauge_param.ga_pad = link_pad;
459  gauge_param.reconstruct = link_recon;
461  gauge_param.cuda_prec_precondition = gauge_param.cuda_prec_sloppy;
462  gauge_param.reconstruct_precondition = gauge_param.reconstruct_sloppy;
463  loadGaugeQuda(milc_longlink, &gauge_param);
464  }
465 
466  switch (test_type) {
467  case 0: // full parity solution
468  case 3: // even
469  case 4: {
470  // QUDA eigensolver test
471  //----------------------------------------------------------------------------
472 
473  // Host side arrays to store the eigenpairs computed by QUDA
474  void **host_evecs = (void **)malloc(eig_nConv * sizeof(void *));
475  for (int i = 0; i < eig_nConv; i++) {
476  host_evecs[i] = (void *)malloc(V * mySpinorSiteSize * eig_inv_param.cpu_prec);
477  }
478  double _Complex *host_evals = (double _Complex *)malloc(eig_param.nEv * sizeof(double _Complex));
479 
480  // This function returns the host_evecs and host_evals pointers, populated with
481  // the requested data, at the requested prec. All the information needed to
482  // perfom the solve is in the eig_param container.
483  // If eig_param.arpack_check == true and precision is double, the routine will
484  // use ARPACK rather than the GPU.
485  double time = -((double)clock());
486  if (eig_param.arpack_check && !(prec == QUDA_DOUBLE_PRECISION)) {
487  errorQuda("ARPACK check only available in double precision");
488  }
489 
490  eigensolveQuda(host_evecs, host_evals, &eig_param);
491  time += (double)clock();
492  printfQuda("Time for %s solution = %f\n", eig_param.arpack_check ? "ARPACK" : "QUDA", time / CLOCKS_PER_SEC);
493 
494  // Deallocate host memory
495  for (int i = 0; i < eig_nConv; i++) free(host_evecs[i]);
496  free(host_evecs);
497  free(host_evals);
498 
499  } break;
500  default: errorQuda("Unsupported test type");
501  } // switch
502 
503  // Clean up gauge fields.
504  for (int dir = 0; dir < 4; dir++) {
505  if (qdp_inlink[dir] != nullptr) {
506  free(qdp_inlink[dir]);
507  qdp_inlink[dir] = nullptr;
508  }
509  if (qdp_fatlink[dir] != nullptr) {
510  free(qdp_fatlink[dir]);
511  qdp_fatlink[dir] = nullptr;
512  }
513  if (qdp_longlink[dir] != nullptr) {
514  free(qdp_longlink[dir]);
515  qdp_longlink[dir] = nullptr;
516  }
517  }
518  if (milc_fatlink != nullptr) {
519  free(milc_fatlink);
520  milc_fatlink = nullptr;
521  }
522  if (milc_longlink != nullptr) {
523  free(milc_longlink);
524  milc_longlink = nullptr;
525  }
526 
527 #ifdef MULTI_GPU
528  if (cpuFat != nullptr) {
529  delete cpuFat;
530  cpuFat = nullptr;
531  }
532  if (cpuLong != nullptr) {
533  delete cpuLong;
534  cpuLong = nullptr;
535  }
536 #endif
537 }
538 
539 int main(int argc, char **argv)
540 {
541 
542  // Set a default
544 
545  for (int i = 1; i < argc; i++) {
546 
547  if (process_command_line_option(argc, argv, &i) == 0) { continue; }
548 
549  printf("ERROR: Invalid option:%s\n", argv[i]);
550  usage(argv);
551  }
552 
553  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
554  initComms(argc, argv, gridsize_from_cmdline);
555 
556  if (test_type != 0 && test_type != 3 && test_type != 4) { errorQuda("Test type %d is outside the valid range.\n", test_type); }
557 
558  // Ensure a reasonable default
559  // ensure that the default is improved staggered
561  warningQuda("The dslash_type %d isn't staggered, asqtad, or laplace. Defaulting to asqtad.\n", dslash_type);
563  }
564 
566  // LAPLACE operator path
567  if (test_type != 0) { errorQuda("Test type %d is not supported for the Laplace operator.\n", test_type); }
570  matpc_type = QUDA_MATPC_EVEN_EVEN; // doesn't matter
571  } else {
572  // STAGGERED operator path
574  if (test_type == 0) {
576  } else {
578  }
579  }
580 
581  if (test_type == 3) {
583  } else if (test_type == 4) {
585  } else if (test_type == 0) {
586  matpc_type = QUDA_MATPC_EVEN_EVEN; // it doesn't matter
587  }
588 
589  if (test_type == 0) {
591  } else {
593  }
594  }
595 
597 
600 
601  // Set n_naiks to 2 if eps_naik != 0.0
603  if (eps_naik != 0.0) {
604  if (compute_fatlong) {
605  n_naiks = 2;
606  printfQuda("Note: epsilon-naik != 0, testing epsilon correction links.\n");
607  } else {
608  eps_naik = 0.0;
609  printfQuda("Not computing fat-long, ignoring epsilon correction.\n");
610  }
611  } else {
612  printfQuda("Note: epsilon-naik = 0, testing original HISQ links.\n");
613  }
614  }
615 
617 
618  printfQuda("dslash_type = %d\n", dslash_type);
619 
620  argc_copy = argc;
621  argv_copy = argv;
622 
623  eigensolve_test();
624 
625  endQuda();
626 
627  // finalize the communications layer
628  finalizeComms();
629 }
int maxiter_precondition
Definition: quda.h:292
int device
Definition: test_util.cpp:1602
int gcrNkrylov
Definition: test_util.cpp:1630
void * qdp_longlink[4]
int dimPartitioned(int dim)
Definition: test_util.cpp:1776
QudaDiracFieldOrder dirac_order
Definition: quda.h:219
void applyGaugeFieldScaling_long(Float **gauge, int Vh, QudaGaugeParam *param, QudaDslashType dslash_type)
Definition: test_util.cpp:747
QudaMassNormalization mass_normalization
Definition: quda.h:208
QudaReconstructType reconstruct_sloppy
Definition: quda.h:53
double anisotropy
Definition: quda.h:38
QudaMatPCType matpc_type
Definition: test_util.cpp:1662
QudaGhostExchange ghostExchange
Definition: lattice_field.h:76
QudaBoolean use_dagger
Definition: quda.h:401
double ca_lambda_max
Definition: test_util.cpp:1633
void setGaugeParam(QudaGaugeParam &gauge_param)
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
int niter
Definition: test_util.cpp:1629
enum QudaPrecision_s QudaPrecision
int xdim
Definition: test_util.cpp:1615
int ga_pad
Definition: quda.h:63
void dw_setDims(int *X, const int L5)
Definition: test_util.cpp:187
void setInvertParam(QudaInvertParam &inv_param)
QudaGaugeFixed gauge_fix
Definition: quda.h:61
QudaInverterType inv_type_precondition
Definition: quda.h:270
int eig_nEv
Definition: test_util.cpp:1723
double ca_lambda_min
Definition: test_util.cpp:1632
QudaLinkType type
Definition: quda.h:42
double kappa
Definition: quda.h:106
QudaPrecision prec_sloppy
Definition: test_util.cpp:1609
char eig_QUDA_logfile[]
Definition: test_util.cpp:1741
#define errorQuda(...)
Definition: util_quda.h:121
double tol
Definition: quda.h:121
QudaDslashType dslash_type
Definition: quda.h:102
QudaReconstructType reconstruct_precondition
Definition: quda.h:59
QudaInverterType inv_type
Definition: quda.h:103
void ** ghost_fatlink
QudaPrecision cuda_prec
Definition: quda.h:214
enum QudaSolveType_s QudaSolveType
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
QudaPrecision cpu_prec
Definition: quda.h:213
QudaReconstructType link_recon
Definition: test_util.cpp:1605
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
static int n_naiks
bool eig_use_dagger
Definition: test_util.cpp:1735
QudaStaggeredPhase staggered_phase_type
Definition: quda.h:71
bool alternative_reliable
Definition: test_util.cpp:1659
QudaDagType dagger
Definition: quda.h:207
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1606
void finalizeComms()
Definition: test_util.cpp:128
QudaGaugeParam gauge_param
QudaSolveType solve_type
Definition: test_util.cpp:1663
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:216
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
const char * get_staggered_test_type(int t)
Definition: misc.cpp:827
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:701
enum QudaEigType_s QudaEigType
bool verify_results
Definition: test_util.cpp:1643
QudaPrecision cpu_prec
double ca_lambda_max
Definition: quda.h:304
void usage_extra(char **argv)
void * qdp_fatlink[4]
void setDims(int *)
Definition: test_util.cpp:151
QudaFieldLocation input_location
Definition: quda.h:99
QudaBoolean use_poly_acc
Definition: quda.h:387
double reliable_delta
Definition: quda.h:129
int solution_accumulator_pipeline
Definition: quda.h:142
double ca_lambda_min
Definition: quda.h:301
bool compute_fatlong
Definition: test_util.cpp:1655
QudaSolutionType solution_type
Definition: quda.h:204
int nConv
Definition: quda.h:420
int use_alternative_reliable
Definition: quda.h:131
QudaEigType eig_type
Definition: quda.h:384
QudaEigType eig_type
Definition: test_util.cpp:1738
void computeStaggeredPlaquetteQDPOrder(void **qdp_link, double plaq[3], const QudaGaugeParam &gauge_param_in, const QudaDslashType dslash_type)
QudaInvertParam * invert_param
Definition: quda.h:381
double scale
Definition: quda.h:40
QudaBoolean require_convergence
Definition: quda.h:408
void initQuda(int device)
int test_type
Definition: test_util.cpp:1636
QudaSolutionType solution_type
Definition: test_util.cpp:1664
QudaFieldLocation output_location
Definition: quda.h:100
double mass
Definition: test_util.cpp:1646
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
int ydim
Definition: test_util.cpp:1616
void ** ghost_longlink
void usage(char **argv)
Definition: test_util.cpp:1783
QudaPrecision cuda_prec_sloppy
Definition: quda.h:215
void reorderQDPtoMILC(Out *milc_out, In **qdp_in, int V, int siteSize)
double kappa
Definition: test_util.cpp:1647
QudaVerbosity verbosity
Definition: quda.h:244
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
void setSpinorSiteSize(int n)
Definition: test_util.cpp:211
QudaInvertParam newQudaInvertParam(void)
int main(int argc, char **argv)
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:768
QudaPrecision cuda_prec_precondition
Definition: quda.h:58
double tol_hq
Definition: quda.h:123
enum QudaMatPCType_s QudaMatPCType
int laplace3D
Definition: test_util.cpp:1622
#define warningQuda(...)
Definition: util_quda.h:133
int poly_deg
Definition: quda.h:390
enum QudaSolutionType_s QudaSolutionType
int Nsrc
Definition: test_util.cpp:1627
void eigensolveQuda(void **h_evecs, double_complex *h_evals, QudaEigParam *param)
QudaGammaBasis gamma_basis
Definition: quda.h:221
QudaBoolean use_norm_op
Definition: quda.h:402
bool eig_use_poly_acc
Definition: test_util.cpp:1730
double a_min
Definition: quda.h:393
int eig_check_interval
Definition: test_util.cpp:1727
QudaPrecision prec
Definition: test_util.cpp:1608
QudaBoolean compute_svd
Definition: quda.h:405
QudaPrecision cuda_prec_sloppy
Definition: quda.h:52
double reliable_delta
Definition: test_util.cpp:1658
const void ** Ghost() const
Definition: gauge_field.h:323
QudaBoolean arpack_check
Definition: quda.h:429
double tol_precondition
Definition: quda.h:289
void eigensolve_test()
int use_sloppy_partial_accumulator
Definition: quda.h:132
int heavy_quark_check
Definition: quda.h:165
QudaReconstructType reconstruct
Definition: quda.h:50
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
char latfile[]
Definition: test_util.cpp:1623
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
#define mySpinorSiteSize
char eig_vec_infile[]
Definition: test_util.cpp:1742
QudaPrecision prec_refinement_sloppy
Definition: test_util.cpp:1610
void computeFatLongGPU(void **qdp_fatlink, void **qdp_longlink, void **qdp_inlink, QudaGaugeParam &gauge_param, size_t gSize, int n_naiks, double eps_naik)
void * memset(void *s, int c, size_t n)
cpuGaugeField * cpuLong
int eig_maxiter
void construct_fat_long_gauge_field(void **fatlink, void **longlink, int type, QudaPrecision precision, QudaGaugeParam *param, QudaDslashType dslash_type)
Definition: test_util.cpp:1062
void setEigParam(QudaEigParam &eig_param)
int eig_nConv
Definition: test_util.cpp:1725
char eig_arpack_logfile[]
Definition: test_util.cpp:1740
char vec_outfile[256]
Definition: quda.h:462
QudaResidualType_s
Definition: enum_quda.h:186
bool eig_require_convergence
Definition: test_util.cpp:1726
double eig_tol
Definition: test_util.cpp:1729
double tadpole_coeff
Definition: quda.h:39
double tol_hq
Definition: test_util.cpp:1657
int pipeline
Definition: test_util.cpp:1634
double tol
Definition: quda.h:422
QudaPrecision cuda_prec_precondition
Definition: quda.h:217
enum QudaEigSpectrumType_s QudaEigSpectrumType
double tol_restart
Definition: quda.h:122
int tdim
Definition: test_util.cpp:1618
int eig_poly_deg
Definition: test_util.cpp:1731
QudaEigSpectrumType eig_spectrum
Definition: test_util.cpp:1737
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
double eps_naik
Definition: test_util.cpp:1652
enum QudaCABasis_s QudaCABasis
int eig_max_restarts
Definition: test_util.cpp:1728
void * qdp_inlink[4]
char arpack_logfile[512]
Definition: quda.h:431
bool eig_compute_svd
Definition: test_util.cpp:1736
double eig_amin
Definition: test_util.cpp:1732
QudaDslashType dslash_type
Definition: test_util.cpp:1621
QudaCABasis ca_basis
Definition: test_util.cpp:1631
cpuGaugeField * cpuFat
char eig_vec_outfile[]
Definition: test_util.cpp:1743
#define printfQuda(...)
Definition: util_quda.h:115
QudaTboundary t_boundary
Definition: quda.h:45
char QUDA_logfile[512]
Definition: quda.h:434
double eig_amax
Definition: test_util.cpp:1733
int solution_accumulator_pipeline
Definition: test_util.cpp:1635
bool eig_arpack_check
Definition: test_util.cpp:1739
enum QudaDslashType_s QudaDslashType
int check_interval
Definition: quda.h:424
int max_restarts
Definition: quda.h:426
QudaEigSpectrumType spectrum
Definition: quda.h:411
QudaResidualType residual_type
Definition: quda.h:320
void display_test_info()
int eig_nKr
Definition: test_util.cpp:1724
char ** argv_copy
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
double a_max
Definition: quda.h:394
QudaPrecision cuda_prec_ritz
Definition: quda.h:447
double tol
Definition: test_util.cpp:1656
QudaMatPCType matpc_type
Definition: quda.h:206
QudaEigParam newQudaEigParam(void)
char vec_infile[256]
Definition: quda.h:459
QudaPrecision cpu_prec
Definition: quda.h:47
bool eig_use_normop
Definition: test_util.cpp:1734
#define gaugeSiteSize
Definition: face_gauge.cpp:34
#define MAX(a, b)
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:211
int zdim
Definition: test_util.cpp:1617
int Vh
Definition: test_util.cpp:28