QUDA  v1.1.0
A library for QCD on GPUs
command_line_params.cpp
Go to the documentation of this file.
1 #include "command_line_params.h"
2 #include <comm_quda.h>
3 
4 // parameters parsed from the command line
5 
6 #ifdef MULTI_GPU
7 int device_ordinal = -1;
8 #else
9 int device_ordinal = 0;
10 #endif
11 
13 std::array<int, 4> gridsize_from_cmdline = {1, 1, 1, 1};
18 
19 bool native_blas_lapack = true;
20 
21 std::array<int, 4> dim_partitioned = {0, 0, 0, 0};
34 std::array<int, 4> dim = {24, 24, 24, 24};
35 int &xdim = dim[0];
36 int &ydim = dim[1];
37 int &zdim = dim[2];
38 int &tdim = dim[3];
39 int Lsdim = 16;
40 bool dagger = false;
42 int laplace3D = 4;
43 char latfile[256] = "";
44 bool unit_gauge = false;
45 double gaussian_sigma = 0.2;
46 char gauge_outfile[256] = "";
47 int Nsrc = 1;
48 int Msrc = 1;
49 int niter = 100;
51 int gcrNkrylov = 10;
53 double ca_lambda_min = 0.0;
54 double ca_lambda_max = -1.0;
55 int pipeline = 0;
57 int test_type = 0;
62 bool inv_deflate = false;
63 bool inv_multigrid = false;
67 int multishift = 1;
68 bool verify_results = true;
69 bool low_mode_check = false;
70 bool oblique_proj_check = false;
71 double mass = 0.1;
72 double kappa = -1.0;
73 double mu = 0.1;
74 double epsilon = 0.01;
75 double m5 = -1.5;
76 double b5 = 1.5;
77 double c5 = 0.5;
78 double anisotropy = 1.0;
79 double tadpole_factor = 1.0;
80 double eps_naik = 0.0;
81 int n_naiks = 1;
82 double clover_csw = 1.0;
83 double clover_coeff = 0.0;
84 bool compute_clover = false;
85 bool compute_fatlong = false;
86 double tol = 1e-7;
87 double tol_precondition = 1e-1;
88 double tol_hq = 0.;
89 double reliable_delta = 0.1;
90 bool alternative_reliable = false;
97 
98 int mg_levels = 2;
99 
117 bool pre_orthonormalize = false;
119 double omega = 0.85;
130 bool generate_nullspace = true;
135 
136 // Aggregation type for the top level of staggered
137 // FIXME: replace with QUDA_TRANSFER_OPTIMIZED_KD when ready
139 
140 // we only actually support 4 here currently
142 
143 #if (CUDA_VERSION >= 10010 && __COMPUTE_CAPABILITY__ >= 700)
144 bool mg_use_mma = true;
145 #else
146 bool mg_use_mma = false;
147 #endif
148 
149 int n_ev = 8;
150 int max_search_dim = 64;
151 int deflation_grid = 16;
152 double tol_restart = 5e+3 * tol;
153 
156 double inc_tol = 1e-2;
157 double eigenval_tol = 1e-1;
158 
163 
164 // Parameters for the stand alone eigensolver
166 int eig_n_ev = 16;
167 int eig_n_kr = 32;
168 int eig_n_conv = -1; // If unchanged, will be set to n_ev
169 int eig_n_ev_deflate = -1; // If unchanged, will be set to n_conv
170 int eig_batched_rotate = 0; // If unchanged, will be set to maximum
173 int eig_max_restarts = 1000;
174 double eig_tol = 1e-6;
175 double eig_qr_tol = 1e-11;
176 bool eig_use_eigen_qr = true;
177 bool eig_use_poly_acc = true;
178 int eig_poly_deg = 100;
179 double eig_amin = 0.1;
180 double eig_amax = 0.0; // If zero is passed to the solver, an estimate will be computed
181 bool eig_use_normop = true;
182 bool eig_use_dagger = false;
183 bool eig_compute_svd = false;
184 bool eig_compute_gamma5 = false;
187 bool eig_arpack_check = false;
188 char eig_arpack_logfile[256] = "arpack_logfile.log";
189 char eig_vec_infile[256] = "";
190 char eig_vec_outfile[256] = "";
193 
194 // Parameters for the MG eigensolver.
195 // The coarsest grid params are for deflation,
196 // all others are for PR vectors.
218 
219 bool mg_eig_coarse_guess = false;
221 
222 double heatbath_beta_value = 6.2;
227 bool heatbath_coldstart = false;
228 
229 int eofa_pm = 1;
230 double eofa_shift = -1.2345;
231 double eofa_mq1 = 1.0;
232 double eofa_mq2 = 0.085;
233 double eofa_mq3 = 1.0;
234 
235 double stout_smear_rho = 0.1;
236 double stout_smear_epsilon = -0.25;
237 double ape_smear_rho = 0.6;
238 int smear_steps = 50;
239 double wflow_epsilon = 0.01;
240 int wflow_steps = 100;
243 
245 
246 std::array<int, 4> grid_partition = {1, 1, 1, 1};
251 
252 std::array<int, 3> blas_mnk = {64, 64, 64};
253 auto &blas_m = blas_mnk[0];
254 auto &blas_n = blas_mnk[1];
255 auto &blas_k = blas_mnk[2];
256 
257 std::array<int, 3> blas_leading_dims = {128, 128, 128};
261 
262 std::array<int, 3> blas_offsets = {0, 0, 0};
266 
267 std::array<int, 3> blas_strides = {1, 1, 1};
271 
272 std::array<double, 2> blas_alpha_re_im = {1.0, 0.0};
273 std::array<double, 2> blas_beta_re_im = {1.0, 0.0};
274 int blas_batch = 16;
275 
276 namespace
277 {
278  CLI::TransformPairs<QudaCABasis> ca_basis_map {{"power", QUDA_POWER_BASIS}, {"chebyshev", QUDA_CHEBYSHEV_BASIS}};
279 
280  CLI::TransformPairs<QudaBLASDataType> blas_dt_map {
282 
283  CLI::TransformPairs<QudaBLASDataOrder> blas_data_order_map {{"row", QUDA_BLAS_DATAORDER_ROW},
284  {"col", QUDA_BLAS_DATAORDER_COL}};
285 
286  CLI::TransformPairs<QudaBLASOperation> blas_op_map {{"N", QUDA_BLAS_OP_N}, {"T", QUDA_BLAS_OP_T}, {"C", QUDA_BLAS_OP_C}};
287 
288  CLI::TransformPairs<QudaContractType> contract_type_map {{"open", QUDA_CONTRACT_TYPE_OPEN},
289  {"dr", QUDA_CONTRACT_TYPE_DR}};
290 
291  CLI::TransformPairs<QudaDslashType> dslash_type_map {{"wilson", QUDA_WILSON_DSLASH},
292  {"clover", QUDA_CLOVER_WILSON_DSLASH},
293  {"twisted-mass", QUDA_TWISTED_MASS_DSLASH},
294  {"twisted-clover", QUDA_TWISTED_CLOVER_DSLASH},
295  {"clover-hasenbusch-twist", QUDA_CLOVER_HASENBUSCH_TWIST_DSLASH},
296  {"staggered", QUDA_STAGGERED_DSLASH},
297  {"asqtad", QUDA_ASQTAD_DSLASH},
298  {"domain-wall", QUDA_DOMAIN_WALL_DSLASH},
299  {"domain-wall-4d", QUDA_DOMAIN_WALL_4D_DSLASH},
300  {"mobius", QUDA_MOBIUS_DWF_DSLASH},
301  {"mobius-eofa", QUDA_MOBIUS_DWF_EOFA_DSLASH},
302  {"laplace", QUDA_LAPLACE_DSLASH}};
303 
304  CLI::TransformPairs<QudaTwistFlavorType> twist_flavor_type_map {{"singlet", QUDA_TWIST_SINGLET},
305  {"deg-doublet", QUDA_TWIST_DEG_DOUBLET},
306  {"nondeg-doublet", QUDA_TWIST_NONDEG_DOUBLET},
307  {"no", QUDA_TWIST_NO}};
308 
309  CLI::TransformPairs<QudaInverterType> inverter_type_map {{"invalid", QUDA_INVALID_INVERTER},
310  {"cg", QUDA_CG_INVERTER},
311  {"bicgstab", QUDA_BICGSTAB_INVERTER},
312  {"gcr", QUDA_GCR_INVERTER},
313  {"pcg", QUDA_PCG_INVERTER},
314  {"mpcg", QUDA_MPCG_INVERTER},
315  {"mpbicgstab", QUDA_MPBICGSTAB_INVERTER},
316  {"mr", QUDA_MR_INVERTER},
317  {"sd", QUDA_SD_INVERTER},
318  {"eigcg", QUDA_EIGCG_INVERTER},
319  {"inc-eigcg", QUDA_INC_EIGCG_INVERTER},
320  {"gmresdr", QUDA_GMRESDR_INVERTER},
321  {"gmresdr-proj", QUDA_GMRESDR_PROJ_INVERTER},
322  {"gmresdr-sh", QUDA_GMRESDR_SH_INVERTER},
323  {"fgmresdr", QUDA_FGMRESDR_INVERTER},
324  {"mg", QUDA_MG_INVERTER},
325  {"bicgstab-l", QUDA_BICGSTABL_INVERTER},
326  {"cgne", QUDA_CGNE_INVERTER},
327  {"cgnr", QUDA_CGNR_INVERTER},
328  {"cg3", QUDA_CG3_INVERTER},
329  {"cg3ne", QUDA_CG3NE_INVERTER},
330  {"cg3nr", QUDA_CG3NR_INVERTER},
331  {"ca-cg", QUDA_CA_CG_INVERTER},
332  {"ca-cgne", QUDA_CA_CGNE_INVERTER},
333  {"ca-cgnr", QUDA_CA_CGNR_INVERTER},
334  {"ca-gcr", QUDA_CA_GCR_INVERTER}};
335 
336  CLI::TransformPairs<QudaPrecision> precision_map {{"double", QUDA_DOUBLE_PRECISION},
337  {"single", QUDA_SINGLE_PRECISION},
338  {"half", QUDA_HALF_PRECISION},
339  {"quarter", QUDA_QUARTER_PRECISION}};
340 
341  CLI::TransformPairs<QudaSchwarzType> schwarz_type_map {{"invalid", QUDA_INVALID_SCHWARZ},
342  {"additive", QUDA_ADDITIVE_SCHWARZ},
343  {"multiplicative", QUDA_MULTIPLICATIVE_SCHWARZ}};
344 
345  CLI::TransformPairs<QudaSolutionType> solution_type_map {{"mat", QUDA_MAT_SOLUTION},
346  {"mat-dag-mat", QUDA_MATDAG_MAT_SOLUTION},
347  {"mat-pc", QUDA_MATPC_SOLUTION},
348  {"mat-pc-dag", QUDA_MATPC_DAG_SOLUTION},
349  {"mat-pc-dag-mat-pc", QUDA_MATPCDAG_MATPC_SOLUTION}};
350 
351  CLI::TransformPairs<QudaEigType> eig_type_map {{"trlm", QUDA_EIG_TR_LANCZOS},
352  {"blktrlm", QUDA_EIG_BLK_TR_LANCZOS},
353  {"iram", QUDA_EIG_IR_ARNOLDI},
354  {"blkiram", QUDA_EIG_BLK_IR_ARNOLDI}};
355 
356  CLI::TransformPairs<QudaTransferType> transfer_type_map {{"aggregate", QUDA_TRANSFER_AGGREGATE},
357  {"kd-coarse", QUDA_TRANSFER_COARSE_KD},
358  {"kd-optimized", QUDA_TRANSFER_OPTIMIZED_KD}};
359 
360  CLI::TransformPairs<QudaTboundary> fermion_t_boundary_map {{"periodic", QUDA_PERIODIC_T},
361  {"anti-periodic", QUDA_ANTI_PERIODIC_T}};
362 
363  CLI::TransformPairs<QudaSolveType> solve_type_map {
364  {"direct", QUDA_DIRECT_SOLVE}, {"direct-pc", QUDA_DIRECT_PC_SOLVE}, {"normop", QUDA_NORMOP_SOLVE},
365  {"normop-pc", QUDA_NORMOP_PC_SOLVE}, {"normerr", QUDA_NORMERR_SOLVE}, {"normerr-pc", QUDA_NORMERR_PC_SOLVE}};
366 
367  CLI::TransformPairs<QudaEigSpectrumType> seig_pectrum_map {
370 
371  CLI::TransformPairs<QudaFieldLocation> field_location_map {{"cpu", QUDA_CPU_FIELD_LOCATION},
372  {"host", QUDA_CPU_FIELD_LOCATION},
373  {"gpu", QUDA_CUDA_FIELD_LOCATION},
374  {"device", QUDA_CUDA_FIELD_LOCATION}};
375 
376  CLI::TransformPairs<QudaVerbosity> verbosity_map {
377  {"silent", QUDA_SILENT}, {"summarize", QUDA_SUMMARIZE}, {"verbose", QUDA_VERBOSE}, {"debug", QUDA_DEBUG_VERBOSE}};
378 
379  CLI::TransformPairs<QudaMassNormalization> mass_normalization_map {{"kappa", QUDA_KAPPA_NORMALIZATION},
380  {"mass", QUDA_MASS_NORMALIZATION},
381  {"asym-mass", QUDA_ASYMMETRIC_MASS_NORMALIZATION}};
382 
383  CLI::TransformPairs<QudaMatPCType> matpc_type_map {{"even-even", QUDA_MATPC_EVEN_EVEN},
384  {"odd-odd", QUDA_MATPC_ODD_ODD},
385  {"even-even-asym", QUDA_MATPC_EVEN_EVEN_ASYMMETRIC},
386  {"odd-odd-asym", QUDA_MATPC_ODD_ODD_ASYMMETRIC}};
387 
388  CLI::TransformPairs<QudaReconstructType> reconstruct_type_map {{"18", QUDA_RECONSTRUCT_NO},
389  {"13", QUDA_RECONSTRUCT_13},
390  {"12", QUDA_RECONSTRUCT_12},
391  {"9", QUDA_RECONSTRUCT_9},
392  {"8", QUDA_RECONSTRUCT_8}};
393 
394  CLI::TransformPairs<QudaEigSpectrumType> eig_spectrum_map {
397 
398  CLI::TransformPairs<QudaWFlowType> wflow_type_map {{"wilson", QUDA_WFLOW_TYPE_WILSON},
399  {"symanzik", QUDA_WFLOW_TYPE_SYMANZIK}};
400 
401  CLI::TransformPairs<QudaSetupType> setup_type_map {{"test", QUDA_TEST_VECTOR_SETUP}, {"null", QUDA_TEST_VECTOR_SETUP}};
402 
403  CLI::TransformPairs<QudaExtLibType> extlib_map {{"eigen", QUDA_EIGEN_EXTLIB}, {"magma", QUDA_MAGMA_EXTLIB}};
404 
405 } // namespace
406 
407 std::shared_ptr<QUDAApp> make_app(std::string app_description, std::string app_name)
408 {
409  auto quda_app = std::make_shared<QUDAApp>(app_description, app_name);
410  quda_app->option_defaults()->always_capture_default();
411 
412  quda_app->add_option("--alternative-reliable", alternative_reliable, "use alternative reliable updates");
413  quda_app->add_option("--anisotropy", anisotropy, "Temporal anisotropy factor (default 1.0)");
414 
415  quda_app->add_option("--ca-basis-type", ca_basis, "The basis to use for CA-CG (default power)")
416  ->transform(CLI::QUDACheckedTransformer(ca_basis_map));
417  quda_app->add_option(
418  "--cheby-basis-eig-max",
419  ca_lambda_max, "Conservative estimate of largest eigenvalue for Chebyshev basis CA-CG (default is to guess with power iterations)");
420  quda_app->add_option("--cheby-basis-eig-min", ca_lambda_min,
421  "Conservative estimate of smallest eigenvalue for Chebyshev basis CA-CG (default 0)");
422  quda_app->add_option("--clover-csw", clover_csw, "Clover Csw coefficient 1.0")->capture_default_str();
423  quda_app->add_option("--clover-coeff", clover_coeff, "The overall clover coefficient, kappa * Csw. (default 0.0. Will be inferred from clover-csw (default 1.0) and kappa. "
424  "If the user populates this value with anything other than 0.0, the passed value will override the inferred value)")->capture_default_str();
425 
426  quda_app->add_option("--compute-clover", compute_clover,
427  "Compute the clover field or use random numbers (default false)");
428  quda_app->add_option("--compute-fat-long", compute_fatlong,
429  "Compute the fat/long field or use random numbers (default false)");
430 
431  quda_app
432  ->add_option("--contraction-type", contract_type,
433  "Whether to leave spin elemental open, or use a gamma basis and contract on "
434  "spin (default open)")
435  ->transform(CLI::QUDACheckedTransformer(contract_type_map));
436 
437  quda_app
438  ->add_option("--blas-data-type", blas_data_type,
439  "Whether to use single(S), double(D), and/or complex(C/Z) data types (default C)")
440  ->transform(CLI::QUDACheckedTransformer(blas_dt_map));
441 
442  quda_app
443  ->add_option("--blas-data-order", blas_data_order, "Whether data is in row major or column major order (default row)")
444  ->transform(CLI::QUDACheckedTransformer(blas_data_order_map));
445 
446  quda_app
447  ->add_option(
448  "--blas-trans-a", blas_trans_a,
449  "Whether to leave the A GEMM matrix as is (N), to transpose (T) or transpose conjugate (C) (default N) ")
450  ->transform(CLI::QUDACheckedTransformer(blas_op_map));
451 
452  quda_app
453  ->add_option(
454  "--blas-trans-b", blas_trans_b,
455  "Whether to leave the B GEMM matrix as is (N), to transpose (T) or transpose conjugate (C) (default N) ")
456  ->transform(CLI::QUDACheckedTransformer(blas_op_map));
457 
458  quda_app->add_option("--blas-alpha", blas_alpha_re_im, "Set the complex value of alpha for GEMM (default {1.0,0.0}")
459  ->expected(2);
460 
461  quda_app->add_option("--blas-beta", blas_beta_re_im, "Set the complex value of beta for GEMM (default {1.0,0.0}")
462  ->expected(2);
463 
464  quda_app
465  ->add_option("--blas-mnk", blas_mnk, "Set the dimensions of the A, B, and C matrices GEMM (default 128 128 128)")
466  ->expected(3);
467 
468  quda_app
469  ->add_option("--blas-leading-dims", blas_leading_dims,
470  "Set the leading dimensions A, B, and C matrices GEMM (default 128 128 128) ")
471  ->expected(3);
472 
473  quda_app->add_option("--blas-offsets", blas_offsets, "Set the offsets for matrices A, B, and C (default 0 0 0)")
474  ->expected(3);
475 
476  quda_app->add_option("--blas-strides", blas_strides, "Set the strides for matrices A, B, and C (default 1 1 1)")
477  ->expected(3);
478 
479  quda_app->add_option("--blas-batch", blas_batch, "Set the number of batches for GEMM (default 16)");
480 
481  quda_app->add_flag("--dagger", dagger, "Set the dagger to 1 (default 0)");
482  quda_app->add_option("--device", device_ordinal, "Set the CUDA device to use (default 0, single GPU only)")
483  ->check(CLI::Range(0, 16));
484 
485  quda_app->add_option("--dslash-type", dslash_type, "Set the dslash type")
486  ->transform(CLI::QUDACheckedTransformer(dslash_type_map));
487 
488  quda_app->add_option("--epsilon", epsilon, "Twisted-Mass flavor twist of Dirac operator (default 0.01)");
489  quda_app->add_option("--epsilon-naik", eps_naik, "Epsilon factor on Naik term (default 0.0, suggested non-zero -0.1)");
490 
491  quda_app
492  ->add_option("--flavor", twist_flavor,
493  "Set the twisted mass flavor type (singlet (default), deg-doublet, nondeg-doublet)")
494  ->transform(CLI::QUDACheckedTransformer(twist_flavor_type_map));
495  ;
496  quda_app->add_option("--gaussian-sigma", gaussian_sigma,
497  "Width of the Gaussian noise used for random gauge field contruction (default 0.2)");
498 
499  quda_app->add_option("--heatbath-beta", heatbath_beta_value, "Beta value used in heatbath test (default 6.2)");
500  quda_app->add_option("--heatbath-coldstart", heatbath_coldstart,
501  "Whether to use a cold or hot start in heatbath test (default false)");
502  quda_app->add_option("--heatbath-num-hb-per-step", heatbath_num_heatbath_per_step,
503  "Number of heatbath hits per heatbath step (default 5)");
504  quda_app->add_option("--heatbath-num-or-per-step", heatbath_num_overrelax_per_step,
505  "Number of overrelaxation hits per heatbath step (default 5)");
506  quda_app->add_option("--heatbath-num-steps", heatbath_num_steps,
507  "Number of measurement steps in heatbath test (default 10)");
508  quda_app->add_option("--heatbath-warmup-steps", heatbath_warmup_steps,
509  "Number of warmup steps in heatbath test (default 10)");
510 
511  quda_app->add_option("--inv-type", inv_type, "The type of solver to use (default cg)")
512  ->transform(CLI::QUDACheckedTransformer(inverter_type_map));
513  quda_app->add_option("--inv-deflate", inv_deflate, "Deflate the inverter using the eigensolver");
514  quda_app->add_option("--inv-multigrid", inv_multigrid, "Precondition the inverter using multigrid");
515  quda_app->add_option("--kappa", kappa, "Kappa of Dirac operator (default 0.12195122... [equiv to mass])");
516  quda_app->add_option(
517  "--laplace3D", laplace3D,
518  "Restrict laplace operator to omit the t dimension (n=3), or include all dims (n=4) (default 4)");
519  quda_app->add_option("--load-gauge", latfile, "Load gauge field \" file \" for the test (requires QIO)");
520  quda_app->add_option("--Lsdim", Lsdim, "Set Ls dimension size(default 16)");
521  quda_app->add_option("--mass", mass, "Mass of Dirac operator (default 0.1)");
522 
523  quda_app->add_option("--mass-normalization", normalization, "Mass normalization (kappa (default) / mass / asym-mass)")
524  ->transform(CLI::QUDACheckedTransformer(mass_normalization_map));
525 
526  quda_app
527  ->add_option("--matpc", matpc_type, "Matrix preconditioning type (even-even, odd-odd, even-even-asym, odd-odd-asym)")
528  ->transform(CLI::QUDACheckedTransformer(matpc_type_map));
529  quda_app->add_option("--msrc", Msrc,
530  "Used for testing non-square block blas routines where nsrc defines the other dimension");
531  quda_app->add_option("--mu", mu, "Twisted-Mass chiral twist of Dirac operator (default 0.1)");
532  quda_app->add_option("--m5", m5, "Mass of shift of five-dimensional Dirac operators (default -1.5)");
533  quda_app->add_option("--b5", b5, "Mobius b5 parameter (default 1.5)");
534  quda_app->add_option("--c5", c5, "Mobius c5 parameter (default 0.5)");
535  quda_app->add_option(
536  "--multishift", multishift,
537  "Whether to do a multi-shift solver test or not. Default is 1 (single mass)"
538  "If a value N > 1 is passed, heavier masses will be constructed and the multi-shift solver will be called");
539  quda_app->add_option("--ngcrkrylov", gcrNkrylov,
540  "The number of inner iterations to use for GCR, BiCGstab-l, CA-CG (default 10)");
541  quda_app->add_option("--niter", niter, "The number of iterations to perform (default 100)");
542  quda_app->add_option("--native-blas-lapack", native_blas_lapack,
543  "Use the native or generic BLAS LAPACK implementation (default true)");
544  quda_app->add_option("--maxiter-precondition", maxiter_precondition,
545  "The number of iterations to perform for any preconditioner (default 10)");
546  quda_app->add_option("--nsrc", Nsrc,
547  "How many spinors to apply the dslash to simultaneusly (experimental for staggered only)");
548 
549  quda_app->add_option("--pipeline", pipeline,
550  "The pipeline length for fused operations in GCR, BiCGstab-l (default 0, no pipelining)");
551 
552  // precision options
553 
554  CLI::QUDACheckedTransformer prec_transform(precision_map);
555  quda_app->add_option("--prec", prec, "Precision in GPU")->transform(prec_transform);
556  quda_app->add_option("--prec-precondition", prec_precondition, "Preconditioner precision in GPU")->transform(prec_transform);
557 
558  quda_app->add_option("--prec-eigensolver", prec_eigensolver, "Eigensolver precision in GPU")->transform(prec_transform);
559 
560  quda_app->add_option("--prec-refine", prec_refinement_sloppy, "Sloppy precision for refinement in GPU")
561  ->transform(prec_transform);
562 
563  quda_app->add_option("--prec-ritz", prec_ritz, "Eigenvector precision in GPU")->transform(prec_transform);
564 
565  quda_app->add_option("--prec-sloppy", prec_sloppy, "Sloppy precision in GPU")->transform(prec_transform);
566 
567  quda_app->add_option("--prec-null", prec_null, "Precison TODO")->transform(prec_transform);
568 
569  quda_app->add_option("--precon-type", precon_type, "The type of solver to use (default none (=unspecified)).")
570  ->transform(CLI::QUDACheckedTransformer(inverter_type_map));
571  quda_app
572  ->add_option("--precon-schwarz-type", precon_schwarz_type,
573  "The type of Schwarz preconditioning to use (default=invalid)")
574  ->transform(CLI::QUDACheckedTransformer(schwarz_type_map));
575  quda_app->add_option("--precon-schwarz-cycle", precon_schwarz_cycle,
576  "The number of Schwarz cycles to apply per smoother application (default=1)");
577 
578  CLI::TransformPairs<int> rank_order_map {{"col", 0}, {"row", 1}};
579  quda_app
580  ->add_option("--rank-order", rank_order,
581  "Set the [t][z][y][x] rank order as either column major (t fastest, default) or row major (x fastest)")
582  ->transform(CLI::QUDACheckedTransformer(rank_order_map));
583 
584  quda_app->add_option("--recon", link_recon, "Link reconstruction type")
585  ->transform(CLI::QUDACheckedTransformer(reconstruct_type_map));
586  quda_app->add_option("--recon-precondition", link_recon_precondition, "Preconditioner link reconstruction type")
587  ->transform(CLI::QUDACheckedTransformer(reconstruct_type_map));
588  quda_app->add_option("--recon-eigensolver", link_recon_eigensolver, "Eigensolver link reconstruction type")
589  ->transform(CLI::QUDACheckedTransformer(reconstruct_type_map));
590  quda_app->add_option("--recon-sloppy", link_recon_sloppy, "Sloppy link reconstruction type")
591  ->transform(CLI::QUDACheckedTransformer(reconstruct_type_map));
592 
593  quda_app->add_option("--reliable-delta", reliable_delta, "Set reliable update delta factor");
594  quda_app->add_option("--save-gauge", gauge_outfile,
595  "Save gauge field \" file \" for the test (requires QIO, heatbath test only)");
596 
597  quda_app->add_option("--solution-pipeline", solution_accumulator_pipeline,
598  "The pipeline length for fused solution accumulation (default 0, no pipelining)");
599 
600  quda_app
601  ->add_option(
602  "--solution-type", solution_type,
603  "The solution we desire (mat (default), mat-dag-mat, mat-pc, mat-pc-dag-mat-pc (default for multi-shift))")
604  ->transform(CLI::QUDACheckedTransformer(solution_type_map));
605 
606  quda_app
607  ->add_option("--fermion-t-boundary", fermion_t_boundary,
608  "The fermoinic temporal boundary conditions (anti-periodic (default), periodic")
609  ->transform(CLI::QUDACheckedTransformer(fermion_t_boundary_map));
610 
611  quda_app
612  ->add_option("--solve-type", solve_type,
613  "The type of solve to do (direct, direct-pc, normop, normop-pc, normerr, normerr-pc)")
614  ->transform(CLI::QUDACheckedTransformer(solve_type_map));
615  quda_app
616  ->add_option("--solver-ext-lib-type", solver_ext_lib, "Set external library for the solvers (default Eigen library)")
617  ->transform(CLI::QUDACheckedTransformer(extlib_map));
618 
619  quda_app->add_option("--tadpole-coeff", tadpole_factor,
620  "Tadpole coefficient for HISQ fermions (default 1.0, recommended [Plaq]^1/4)");
621 
622  quda_app->add_option("--tol", tol, "Set L2 residual tolerance");
623  quda_app->add_option("--tolhq", tol_hq, "Set heavy-quark residual tolerance");
624  quda_app->add_option("--tol-precondition", tol_precondition, "Set L2 residual tolerance for preconditioner");
625  quda_app->add_option(
626  "--unit-gauge", unit_gauge,
627  "Generate a unit valued gauge field in the tests. If false, a random gauge is generated (default false)");
628 
629  quda_app->add_option("--verbosity", verbosity, "The the verbosity on the top level of QUDA( default summarize)")
630  ->transform(CLI::QUDACheckedTransformer(verbosity_map));
631  quda_app->add_option("--verify", verify_results, "Verify the GPU results using CPU results (default true)");
632 
633  // lattice dimensions
634  auto dimopt = quda_app->add_option("--dim", dim, "Set space-time dimension (X Y Z T)")->check(CLI::Range(1, 512));
635  auto sdimopt = quda_app
636  ->add_option(
637  "--sdim",
638  [](CLI::results_t res) {
639  return CLI::detail::lexical_cast(res[0], xdim) && CLI::detail::lexical_cast(res[0], ydim)
640  && CLI::detail::lexical_cast(res[0], zdim);
641  },
642  "Set space dimension(X/Y/Z) size")
643  ->type_name("INT")
644  ->check(CLI::Range(1, 512));
645 
646  quda_app->add_option("--xdim", xdim, "Set X dimension size(default 24)")
647  ->check(CLI::Range(1, 512))
648  ->excludes(dimopt)
649  ->excludes(sdimopt);
650  quda_app->add_option("--ydim", ydim, "Set X dimension size(default 24)")
651  ->check(CLI::Range(1, 512))
652  ->excludes(dimopt)
653  ->excludes(sdimopt);
654  quda_app->add_option("--zdim", zdim, "Set X dimension size(default 24)")
655  ->check(CLI::Range(1, 512))
656  ->excludes(dimopt)
657  ->excludes(sdimopt);
658  quda_app->add_option("--tdim", tdim, "Set T dimension size(default 24)")->check(CLI::Range(1, 512))->excludes(dimopt);
659 
660  // multi-gpu partitioning
661 
662  quda_app->add_option(
663  "--partition",
664  [](CLI::results_t res) {
665  int p;
666  auto retval = CLI::detail::lexical_cast(res[0], p);
667  for (int j = 0; j < 4; j++) {
668  if (p & (1 << j)) { dim_partitioned[j] = 1; }
669  }
670  return retval;
671  },
672  "Set the communication topology (X=1, Y=2, Z=4, T=8, and combinations of these)");
673 
674  auto gridsizeopt
675  = quda_app
676  ->add_option("--gridsize", gridsize_from_cmdline, "Set the grid size in all four dimension (default 1 1 1 1)")
677  ->expected(4);
678  quda_app->add_option("--xgridsize", grid_x, "Set grid size in X dimension (default 1)")->excludes(gridsizeopt);
679  quda_app->add_option("--ygridsize", grid_y, "Set grid size in Y dimension (default 1)")->excludes(gridsizeopt);
680  quda_app->add_option("--zgridsize", grid_z, "Set grid size in Z dimension (default 1)")->excludes(gridsizeopt);
681  quda_app->add_option("--tgridsize", grid_t, "Set grid size in T dimension (default 1)")->excludes(gridsizeopt);
682 
683  return quda_app;
684 }
685 
686 void add_eigen_option_group(std::shared_ptr<QUDAApp> quda_app)
687 {
688 
689  CLI::QUDACheckedTransformer prec_transform(precision_map);
690  // Option group for Eigensolver related options
691  auto opgroup = quda_app->add_option_group("Eigensolver", "Options controlling eigensolver");
692 
693  opgroup->add_option("--eig-amax", eig_amax, "The maximum in the polynomial acceleration")->check(CLI::PositiveNumber);
694  opgroup->add_option("--eig-amin", eig_amin, "The minimum in the polynomial acceleration")->check(CLI::PositiveNumber);
695 
696  opgroup->add_option("--eig-ARPACK-logfile", eig_arpack_logfile, "The filename storing the log from arpack");
697  opgroup->add_option("--eig-arpack-check", eig_arpack_check,
698  "Cross check the device data against ARPACK (requires ARPACK, default false)");
699  opgroup->add_option("--eig-use-eigen-qr", eig_use_eigen_qr,
700  "Use Eigen to eigensolve the upper Hessenberg in IRAM, else use QUDA's QR code. (default true)");
701  opgroup->add_option("--eig-compute-svd", eig_compute_svd,
702  "Solve the MdagM problem, use to compute SVD of M (default false)");
703 
704  opgroup->add_option("--eig-compute-gamma5", eig_compute_gamma5,
705  "Solve the gamma5 OP problem. Solve for OP then multiply by gamma_5 (default false)");
706 
707  opgroup->add_option("--eig-max-restarts", eig_max_restarts, "Perform n iterations of the restart in the eigensolver");
708  opgroup->add_option("--eig-block-size", eig_block_size, "The block size to use in the block variant eigensolver");
709  opgroup->add_option(
710  "--eig-n-ev-deflate", eig_n_ev_deflate,
711  "The number of converged eigenpairs that will be used in the deflation routines (default eig_n_conv)");
712  opgroup->add_option("--eig-n-conv", eig_n_conv, "The number of converged eigenvalues requested (default eig_n_ev)");
713  opgroup->add_option("--eig-n-ev", eig_n_ev, "The size of eigenvector search space in the eigensolver");
714  opgroup->add_option("--eig-n-kr", eig_n_kr, "The size of the Krylov subspace to use in the eigensolver");
715  opgroup->add_option("--eig-batched-rotate", eig_batched_rotate,
716  "The maximum number of extra eigenvectors the solver may allocate to perform a Ritz rotation.");
717  opgroup->add_option("--eig-poly-deg", eig_poly_deg, "TODO");
718  opgroup->add_option(
719  "--eig-require-convergence",
720  eig_require_convergence, "If true, the solver will error out if convergence is not attained. If false, a warning will be given (default true)");
721  opgroup->add_option("--eig-save-vec", eig_vec_outfile, "Save eigenvectors to <file> (requires QIO)");
722  opgroup->add_option("--eig-load-vec", eig_vec_infile, "Load eigenvectors to <file> (requires QIO)")
723  ->check(CLI::ExistingFile);
724  opgroup
725  ->add_option("--eig-save-prec", eig_save_prec,
726  "If saving eigenvectors, use this precision to save. No-op if eig-save-prec is greater than or equal "
727  "to precision of eigensolver (default = double)")
728  ->transform(prec_transform);
729 
730  opgroup->add_option(
731  "--eig-io-parity-inflate", eig_io_parity_inflate,
732  "Whether to inflate single-parity eigenvectors onto dual parity full fields for file I/O (default = false)");
733 
734  opgroup
735  ->add_option("--eig-spectrum", eig_spectrum,
736  "The spectrum part to be calulated. S=smallest L=largest R=real M=modulus I=imaginary")
737  ->transform(CLI::QUDACheckedTransformer(eig_spectrum_map));
738  opgroup->add_option("--eig-tol", eig_tol, "The tolerance to use in the eigensolver (default 1e-6)");
739  opgroup->add_option("--eig-qr-tol", eig_qr_tol, "The tolerance to use in the qr (default 1e-11)");
740 
741  opgroup->add_option("--eig-type", eig_type, "The type of eigensolver to use (default trlm)")
742  ->transform(CLI::QUDACheckedTransformer(eig_type_map));
743 
744  opgroup->add_option("--eig-use-dagger", eig_use_dagger,
745  "Solve the Mdag problem instead of M (MMdag if eig-use-normop == true) (default false)");
746  opgroup->add_option("--eig-use-normop", eig_use_normop,
747  "Solve the MdagM problem instead of M (MMdag if eig-use-dagger == true) (default false)");
748  opgroup->add_option("--eig-use-poly-acc", eig_use_poly_acc, "Use Chebyshev polynomial acceleration in the eigensolver");
749 }
750 
751 void add_deflation_option_group(std::shared_ptr<QUDAApp> quda_app)
752 {
753  auto opgroup = quda_app->add_option_group("Deflation", "Options controlling deflation");
754 
755  opgroup
756  ->add_option("--df-deflation-grid", deflation_grid,
757  "Set maximum number of cycles needed to compute eigenvectors(default 1)")
758  ->check(CLI::PositiveNumber);
759  opgroup
760  ->add_option(
761  "--df-eigcg-max-restarts",
762  eigcg_max_restarts, "Set how many iterative refinement cycles will be solved with eigCG within a single physical right hand site solve (default 4)")
763  ->check(CLI::PositiveNumber);
764  ;
765  opgroup->add_option("--df-ext-lib-type", deflation_ext_lib,
766  "Set external library for the deflation methods (default Eigen library)");
767  opgroup->add_option("--df-location-ritz", location_ritz,
768  "Set memory location for the ritz vectors (default cuda memory location)");
769  opgroup->add_option("--df-max-restart-num", max_restart_num,
770  "Set maximum number of the initCG restarts in the deflation stage (default 3)");
771  opgroup->add_option("--df-max-search-dim", max_search_dim, "Set the size of eigenvector search space (default 64)");
772  opgroup->add_option("--df-mem-type-ritz", mem_type_ritz,
773  "Set memory type for the ritz vectors (default device memory type)");
774  opgroup->add_option("--df-n-ev", n_ev, "Set number of eigenvectors computed within a single solve cycle (default 8)");
775  opgroup->add_option("--df-tol-eigenval", eigenval_tol, "Set maximum eigenvalue residual norm (default 1e-1)");
776  opgroup->add_option("--df-tol-inc", inc_tol,
777  "Set tolerance for the subsequent restarts in the initCG solver (default 1e-2)");
778  opgroup->add_option("--df-tol-restart", tol_restart,
779  "Set tolerance for the first restart in the initCG solver(default 5e-5)");
780 }
781 
782 void add_multigrid_option_group(std::shared_ptr<QUDAApp> quda_app)
783 {
784  auto opgroup = quda_app->add_option_group("MultiGrid", "Options controlling deflation");
785 
786  // MWTODO: clean this up - code duplication
787 
788  auto solve_type_transform = CLI::QUDACheckedTransformer(solve_type_map);
789 
790  CLI::QUDACheckedTransformer prec_transform(precision_map);
791  // TODO
792  quda_app->add_mgoption(
793  opgroup, "--mg-block-size", geo_block_size, CLI::Validator(),
794  "Set the geometric block size for the each multigrid levels transfer operator (default 4 4 4 4)");
795  quda_app->add_mgoption(opgroup, "--mg-coarse-solve-type", coarse_solve_type, solve_type_transform,
796  "The type of solve to do on each level (direct, direct-pc) (default = solve_type)");
797 
798  auto solver_trans = CLI::QUDACheckedTransformer(inverter_type_map);
799  quda_app->add_mgoption(opgroup, "--mg-coarse-solver", coarse_solver, solver_trans,
800  "The solver to wrap the V cycle on each level (default gcr, only for levels 1+)");
801 
802  quda_app->add_mgoption(opgroup, "--mg-coarse-solver-ca-basis-size", coarse_solver_ca_basis_size, CLI::PositiveNumber,
803  "The basis size to use for CA-CG setup of multigrid (default 4)");
804 
805  quda_app->add_mgoption(opgroup, "--mg-coarse-solver-ca-basis-type", coarse_solver_ca_basis,
806  CLI::QUDACheckedTransformer(ca_basis_map),
807  "The basis to use for CA-CG setup of multigrid(default power)");
808  quda_app->add_mgoption(opgroup, "--mg-coarse-solver-cheby-basis-eig-max", coarse_solver_ca_lambda_max,
809  CLI::PositiveNumber,
810  "Conservative estimate of largest eigenvalue for Chebyshev basis CA-CG in setup of multigrid "
811  "(default is to guess with power iterations)");
812  quda_app->add_mgoption(
813  opgroup, "--mg-coarse-solver-cheby-basis-eig-min", coarse_solver_ca_lambda_min, CLI::PositiveNumber,
814  "Conservative estimate of smallest eigenvalue for Chebyshev basis CA-CG in setup of multigrid (default 0)");
815  quda_app->add_mgoption(opgroup, "--mg-coarse-solver-maxiter", coarse_solver_maxiter, CLI::PositiveNumber,
816  "The coarse solver maxiter for each level (default 100)");
817  quda_app->add_mgoption(opgroup, "--mg-coarse-solver-tol", coarse_solver_tol, CLI::PositiveNumber,
818  "The coarse solver tolerance for each level (default 0.25, only for levels 1+)");
819  quda_app->add_mgoption(opgroup, "--mg-eig", mg_eig, CLI::Validator(),
820  "Use the eigensolver on this level (default false)");
821  quda_app->add_mgoption(opgroup, "--mg-eig-amax", mg_eig_amax, CLI::PositiveNumber,
822  "The maximum in the polynomial acceleration (default 4.0)");
823  quda_app->add_mgoption(opgroup, "--mg-eig-amin", mg_eig_amin, CLI::PositiveNumber,
824  "The minimum in the polynomial acceleration (default 0.1)");
825  quda_app->add_mgoption(
826  opgroup, "--mg-eig-check-interval", mg_eig_check_interval, CLI::Validator(),
827  "Perform a convergence check every nth restart/iteration (only used in Implicit Restart types)");
828  quda_app->add_option("--mg-eig-coarse-guess", mg_eig_coarse_guess,
829  "If deflating on the coarse grid, optionally use an initial guess (default = false)");
830  quda_app->add_option("--mg-eig-preserve-deflation", mg_eig_preserve_deflation,
831  "If the multigrid operator is updated, preserve generated deflation space (default = false)");
832  quda_app->add_mgoption(opgroup, "--mg-eig-max-restarts", mg_eig_max_restarts, CLI::PositiveNumber,
833  "Perform a maximun of n restarts in eigensolver (default 100)");
834  quda_app->add_mgoption(
835  opgroup, "--mg-eig-use-eigen-qr", mg_eig_use_eigen_qr, CLI::Validator(),
836  "Use Eigen to eigensolve the upper Hessenberg in IRAM, else use QUDA's QR code. (default true)");
837  quda_app->add_mgoption(opgroup, "--mg-eig-block-size", mg_eig_block_size, CLI::Validator(),
838  "The block size to use in the block variant eigensolver");
839  quda_app->add_mgoption(opgroup, "--mg-eig-n-ev", mg_eig_n_ev, CLI::Validator(),
840  "The size of eigenvector search space in the eigensolver");
841  quda_app->add_mgoption(opgroup, "--mg-eig-n-kr", mg_eig_n_kr, CLI::Validator(),
842  "The size of the Krylov subspace to use in the eigensolver");
843  quda_app->add_mgoption(opgroup, "--mg-eig-n-ev-deflate", mg_eig_n_ev_deflate, CLI::Validator(),
844  "The number of converged eigenpairs that will be used in the deflation routines");
845  quda_app->add_mgoption(
846  opgroup, "--mg-eig-batched-rotate", mg_eig_batched_rotate, CLI::Validator(),
847  "The maximum number of extra eigenvectors the solver may allocate to perform a Ritz rotation.");
848  quda_app->add_mgoption(opgroup, "--mg-eig-poly-deg", mg_eig_poly_deg, CLI::PositiveNumber,
849  "Set the degree of the Chebyshev polynomial (default 100)");
850  quda_app->add_mgoption(
851  opgroup, "--mg-eig-require-convergence", mg_eig_require_convergence,
852  CLI::Validator(), "If true, the solver will error out if convergence is not attained. If false, a warning will be given (default true)");
853 
854  quda_app->add_mgoption(
855  opgroup, "--mg-eig-spectrum", mg_eig_spectrum, CLI::QUDACheckedTransformer(eig_spectrum_map),
856  "The spectrum part to be calulated. S=smallest L=largest R=real M=modulus I=imaginary (default SR)");
857  quda_app->add_mgoption(opgroup, "--mg-eig-tol", mg_eig_tol, CLI::PositiveNumber,
858  "The tolerance to use in the eigensolver (default 1e-6)");
859  quda_app->add_mgoption(opgroup, "--mg-eig-qr-tol", mg_eig_qr_tol, CLI::PositiveNumber,
860  "The tolerance to use in the QR (default 1e-11)");
861 
862  quda_app->add_mgoption(opgroup, "--mg-eig-type", mg_eig_type, CLI::QUDACheckedTransformer(eig_type_map),
863  "The type of eigensolver to use (default trlm)");
864  quda_app->add_mgoption(opgroup, "--mg-eig-use-dagger", mg_eig_use_dagger, CLI::Validator(),
865  "Solve the MMdag problem instead of M (MMdag if eig-use-normop == true) (default false)");
866  quda_app->add_mgoption(opgroup, "--mg-eig-use-normop", mg_eig_use_normop, CLI::Validator(),
867  "Solve the MdagM problem instead of M (MMdag if eig-use-dagger == true) (default false)");
868  quda_app->add_mgoption(opgroup, "--mg-eig-use-poly-acc", mg_eig_use_poly_acc, CLI::Validator(),
869  "Use Chebyshev polynomial acceleration in the eigensolver (default true)");
870  opgroup->add_option(
871  "--mg-generate-all-levels",
872  generate_all_levels, "true=generate null-space on all levels, false=generate on level 0 and create other levels from that (default true)");
873  opgroup->add_option("--mg-evolve-thin-updates", mg_evolve_thin_updates,
874  "Utilize thin updates for multigrid evolution tests (default false)");
875  opgroup->add_option("--mg-generate-nullspace", generate_nullspace,
876  "Generate the null-space vector dynamically (default true, if set false and mg-load-vec isn't "
877  "set, creates free-field null vectors)");
878  opgroup->add_option("--mg-levels", mg_levels, "The number of multigrid levels to do (default 2)");
879 
880  // TODO
881  quda_app->add_mgoption(opgroup, "--mg-load-vec", mg_vec_infile, CLI::Validator(),
882  "Load the vectors <file> for the multigrid_test (requires QIO)");
883  quda_app->add_mgoption(opgroup, "--mg-save-vec", mg_vec_outfile, CLI::Validator(),
884  "Save the generated null-space vectors <file> from the multigrid_test (requires QIO)");
885 
886  quda_app
887  ->add_mgoption("--mg-eig-save-prec", mg_eig_save_prec, CLI::Validator(),
888  "If saving eigenvectors, use this precision to save. No-op if mg-eig-save-prec is greater than or "
889  "equal to precision of eigensolver (default = double)")
890  ->transform(prec_transform);
891 
892  opgroup->add_option(
893  "--mg-low-mode-check", low_mode_check,
894  "Measure how well the null vector subspace overlaps with the low eigenmode subspace (default false)");
895  quda_app->add_mgoption(opgroup, "--mg-mu-factor", mu_factor, CLI::Validator(),
896  "Set the multiplicative factor for the twisted mass mu parameter on each level (default 1)");
897  quda_app->add_mgoption(opgroup, "--mg-n-block-ortho", n_block_ortho, CLI::PositiveNumber,
898  "The number of times to run Gram-Schmidt during block orthonormalization (default 1)");
899  quda_app->add_mgoption(opgroup, "--mg-nu-post", nu_post, CLI::PositiveNumber,
900  "The number of post-smoother applications to do at a given multigrid level (default 2)");
901  quda_app->add_mgoption(opgroup, "--mg-nu-pre", nu_pre, CLI::PositiveNumber,
902  "The number of pre-smoother applications to do at a given multigrid level (default 2)");
903  quda_app->add_mgoption(opgroup, "--mg-nvec", nvec, CLI::PositiveNumber,
904  "Number of null-space vectors to define the multigrid transfer operator on a given level");
905  opgroup->add_option("--mg-oblique-proj-check", oblique_proj_check,
906  "Measure how well the null vector subspace adjusts the low eigenmode subspace (default false)");
907  opgroup->add_option("--mg-omega", omega,
908  "The over/under relaxation factor for the smoother of multigrid (default 0.85)");
909  opgroup->add_option("--mg-post-orth", post_orthonormalize,
910  "If orthonormalize the vector after inverting in the setup of multigrid (default true)");
911  opgroup->add_option("--mg-pre-orth", pre_orthonormalize,
912  "If orthonormalize the vector before inverting in the setup of multigrid (default false)");
913 
914  quda_app
915  ->add_mgoption(opgroup, "--mg-schwarz-type", mg_schwarz_type, CLI::Validator(),
916  "The type of preconditioning to use (requires MR smoother and GCR setup solver) (default=invalid)")
917  ->transform(CLI::QUDACheckedTransformer(schwarz_type_map));
918  quda_app->add_mgoption(opgroup, "--mg-schwarz-cycle", mg_schwarz_cycle, CLI::PositiveNumber,
919  "The number of Schwarz cycles to apply per smoother application (default=1)");
920  quda_app->add_mgoption(opgroup, "--mg-setup-ca-basis-size", setup_ca_basis_size, CLI::PositiveNumber,
921  "The basis size to use for CA-CG setup of multigrid (default 4)");
922  quda_app->add_mgoption(opgroup, "--mg-setup-ca-basis-type", setup_ca_basis, CLI::QUDACheckedTransformer(ca_basis_map),
923  "The basis to use for CA-CG setup of multigrid(default power)");
924  quda_app->add_mgoption(opgroup, "--mg-setup-cheby-basis-eig-max", setup_ca_lambda_max, CLI::PositiveNumber,
925  "Conservative estimate of largest eigenvalue for Chebyshev basis CA-CG in setup of multigrid "
926  "(default is to guess with power iterations)");
927  quda_app->add_mgoption(
928  opgroup, "--mg-setup-cheby-basis-eig-min", setup_ca_lambda_min, CLI::PositiveNumber,
929  "Conservative estimate of smallest eigenvalue for Chebyshev basis CA-CG in setup of multigrid (default 0)");
930  quda_app->add_mgoption(opgroup, "--mg-setup-inv", setup_inv, solver_trans,
931  "The inverter to use for the setup of multigrid (default bicgstab)");
932  quda_app->add_mgoption(opgroup, "--mg-setup-iters", num_setup_iter, CLI::PositiveNumber,
933  "The number of setup iterations to use for the multigrid (default 1)");
934 
935  quda_app->add_mgoption(
936  opgroup, "--mg-setup-maxiter", setup_maxiter, CLI::Validator(),
937  "The maximum number of solver iterations to use when relaxing on a null space vector (default 500)");
938  quda_app->add_mgoption(
939  opgroup, "--mg-setup-maxiter-refresh", setup_maxiter_refresh, CLI::Validator(),
940  "The maximum number of solver iterations to use when refreshing the pre-existing null space vectors (default 100)");
941  quda_app->add_mgoption(opgroup, "--mg-setup-tol", setup_tol, CLI::Validator(),
942  "The tolerance to use for the setup of multigrid (default 5e-6)");
943 
944  opgroup->add_option("--mg-setup-type", setup_type, "The type of setup to use for the multigrid (default null)")
945  ->transform(CLI::QUDACheckedTransformer(setup_type_map));
946 
947  // FIXME: default should become kd-optimized
948  opgroup
949  ->add_option(
950  "--mg-staggered-coarsen-type",
951  staggered_transfer_type, "The type of coarsening to use for the top level staggered operator (aggregate, kd-coarse (default), kd-optimized)")
952  ->transform(CLI::QUDACheckedTransformer(transfer_type_map));
953 
954  quda_app->add_mgoption(opgroup, "--mg-smoother", smoother_type, solver_trans,
955  "The smoother to use for multigrid (default mr)");
956 
957  opgroup
958  ->add_option("--mg-smoother-halo-prec", smoother_halo_prec,
959  "The smoother halo precision (applies to all levels - defaults to null_precision)")
960  ->transform(prec_transform);
961 
962  quda_app->add_mgoption(opgroup, "--mg-smoother-solve-type", smoother_solve_type, solve_type_transform,
963  "The type of solve to do in smoother (direct, direct-pc (default) )");
964  quda_app->add_mgoption(opgroup, "--mg-smoother-tol", smoother_tol, CLI::Validator(),
965  "The smoother tolerance to use for each multigrid (default 0.25)");
966 
967  quda_app->add_mgoption(opgroup, "--mg-verbosity", mg_verbosity, CLI::QUDACheckedTransformer(verbosity_map),
968  "The verbosity to use on each level of the multigrid (default summarize)");
969 
970  opgroup->add_option(
971  "--mg-use-mma", mg_use_mma,
972  "Use tensor-core to accelerate multigrid (default = true on Volta or later with CUDA >=10.1, otherwise false)");
973 }
974 
975 void add_eofa_option_group(std::shared_ptr<QUDAApp> quda_app)
976 {
977  auto opgroup = quda_app->add_option_group("EOFA", "Options controlling EOFA parameteres");
978 
979  CLI::TransformPairs<int> eofa_pm_map {{"plus", 1}, {"minus", 0}};
980  opgroup->add_option("--eofa-pm", eofa_pm, "Set to evalute \"plus\" or \"minus\" EOFA operator (default plus)")
981  ->transform(CLI::QUDACheckedTransformer(eofa_pm_map));
982  opgroup->add_option("--eofa-shift", eofa_shift, "Set the shift for the EOFA operator (default -0.12345)");
983  opgroup->add_option("--eofa-mq1", eofa_mq1, "Set mq1 for EOFA operator (default 1.0)");
984  opgroup->add_option("--eofa-mq2", eofa_mq1, "Set mq2 for EOFA operator (default 0.085)");
985  opgroup->add_option("--eofa-mq3", eofa_mq1, "Set mq3 for EOFA operator (default 1.0)");
986 }
987 
988 void add_su3_option_group(std::shared_ptr<QUDAApp> quda_app)
989 {
990 
991  // Option group for SU(3) related options
992  auto opgroup = quda_app->add_option_group("SU(3)", "Options controlling SU(3) tests");
993  opgroup->add_option("--su3-ape-rho", ape_smear_rho, "rho coefficient for APE smearing (default 0.6)");
994 
995  opgroup->add_option("--su3-stout-rho", stout_smear_rho,
996  "rho coefficient for Stout and Over-Improved Stout smearing (default 0.08)");
997 
998  opgroup->add_option("--su3-stout-epsilon", stout_smear_epsilon,
999  "epsilon coefficient for Over-Improved Stout smearing (default -0.25)");
1000 
1001  opgroup->add_option("--su3-smear-steps", smear_steps, "The number of smearing steps to perform (default 50)");
1002 
1003  opgroup->add_option("--su3-wflow-epsilon", wflow_epsilon, "The step size in the Runge-Kutta integrator (default 0.01)");
1004 
1005  opgroup->add_option("--su3-wflow-steps", wflow_steps,
1006  "The number of steps in the Runge-Kutta integrator (default 100)");
1007 
1008  opgroup->add_option("--su3-wflow-type", wflow_type, "The type of action to use in the wilson flow (default wilson)")
1009  ->transform(CLI::QUDACheckedTransformer(wflow_type_map));
1010  ;
1011 
1012  opgroup->add_option("--su3-measurement-interval", measurement_interval,
1013  "Measure the field energy and topological charge every Nth step (default 5) ");
1014 }
1015 
1016 void add_comms_option_group(std::shared_ptr<QUDAApp> quda_app)
1017 {
1018  auto opgroup
1019  = quda_app->add_option_group("Communication", "Options controlling communication (split grid) parameteres");
1020  opgroup->add_option("--grid-partition", grid_partition, "Set the grid partition (default 1 1 1 1)")->expected(4);
1021 }
bool inv_multigrid
int heatbath_num_steps
bool mg_use_mma
QudaTransferType staggered_transfer_type
bool generate_all_levels
QudaWFlowType wflow_type
double kappa
int eig_check_interval
double stout_smear_epsilon
int eig_batched_rotate
QudaInverterType precon_type
double reliable_delta
QudaSolveType solve_type
QudaPrecision prec_refinement_sloppy
int laplace3D
quda::mgarray< int > num_setup_iter
int eig_poly_deg
bool eig_io_parity_inflate
double inc_tol
double c5
QudaEigType eig_type
quda::mgarray< int > mg_eig_poly_deg
QudaInverterType inv_type
quda::mgarray< QudaEigType > mg_eig_type
QudaReconstructType link_recon_sloppy
int eig_n_ev
bool mg_eig_preserve_deflation
double clover_csw
quda::mgarray< int > mg_eig_check_interval
quda::mgarray< char[256]> mg_vec_outfile
double eigenval_tol
std::shared_ptr< QUDAApp > make_app(std::string app_description, std::string app_name)
double tadpole_factor
char eig_vec_infile[256]
bool heatbath_coldstart
bool eig_arpack_check
auto & blas_lda
double eig_tol
double eofa_mq1
quda::mgarray< double > setup_ca_lambda_max
double mass
int precon_schwarz_cycle
QudaPrecision prec_null
double tol
QudaExtLibType solver_ext_lib
QudaReconstructType link_recon
int deflation_grid
auto & blas_b_offset
auto & blas_c_offset
int niter
int test_type
quda::mgarray< int > mg_eig_n_ev_deflate
quda::mgarray< bool > mg_eig
std::array< int, 4 > dim
auto & blas_c_stride
double anisotropy
int device_ordinal
bool eig_compute_svd
int blas_batch
auto & blas_ldb
QudaPrecision eig_save_prec
QudaVerbosity verbosity
quda::mgarray< QudaInverterType > coarse_solver
QudaReconstructType link_recon_precondition
quda::mgarray< QudaCABasis > coarse_solver_ca_basis
quda::mgarray< int > n_block_ortho
double tol_hq
int & ydim
void add_multigrid_option_group(std::shared_ptr< QUDAApp > quda_app)
double eofa_shift
bool eig_use_poly_acc
QudaTwistFlavorType twist_flavor
std::array< int, 4 > grid_partition
int eig_n_ev_deflate
double stout_smear_rho
QudaBLASOperation blas_trans_b
QudaPrecision prec_ritz
double gaussian_sigma
auto & blas_n
auto & blas_b_stride
quda::mgarray< int > coarse_solver_maxiter
auto & blas_a_stride
double epsilon
int Msrc
bool compute_clover
bool verify_results
quda::mgarray< char[256]> mg_vec_infile
quda::mgarray< int > nu_post
quda::mgarray< int > nu_pre
char gauge_outfile[256]
quda::mgarray< int > setup_ca_basis_size
bool pre_orthonormalize
int rank_order
quda::mgarray< int > mg_eig_n_kr
quda::mgarray< int > coarse_solver_ca_basis_size
char eig_vec_outfile[256]
quda::mgarray< double > mg_eig_amin
double eofa_mq2
auto & grid_x
bool post_orthonormalize
quda::mgarray< QudaVerbosity > mg_verbosity
quda::mgarray< double > setup_ca_lambda_min
double wflow_epsilon
char latfile[256]
int eig_block_size
int smear_steps
double eig_amax
auto & blas_ldc
quda::mgarray< int > setup_maxiter
quda::mgarray< int > nvec
std::array< double, 2 > blas_alpha_re_im
void add_eofa_option_group(std::shared_ptr< QUDAApp > quda_app)
bool eig_use_dagger
bool inv_deflate
auto & blas_a_offset
QudaBLASDataOrder blas_data_order
auto & grid_y
quda::mgarray< QudaCABasis > setup_ca_basis
double mu
double eig_qr_tol
int max_search_dim
quda::mgarray< QudaSchwarzType > mg_schwarz_type
bool oblique_proj_check
int gcrNkrylov
int wflow_steps
QudaMemoryType mem_type_ritz
quda::mgarray< QudaEigSpectrumType > mg_eig_spectrum
quda::mgarray< int > mg_eig_max_restarts
int & zdim
double ca_lambda_min
QudaSolutionType solution_type
bool mg_eig_coarse_guess
double m5
bool eig_use_eigen_qr
double eps_naik
int heatbath_num_heatbath_per_step
int max_restart_num
auto & grid_t
int measurement_interval
QudaDslashType dslash_type
double ca_lambda_max
QudaTboundary fermion_t_boundary
bool eig_use_normop
quda::mgarray< bool > mg_eig_use_dagger
std::array< int, 3 > blas_leading_dims
quda::mgarray< double > mu_factor
QudaExtLibType deflation_ext_lib
bool generate_nullspace
quda::mgarray< double > coarse_solver_ca_lambda_max
bool alternative_reliable
int solution_accumulator_pipeline
bool low_mode_check
std::array< int, 4 > dim_partitioned
std::array< int, 3 > blas_strides
quda::mgarray< int > mg_eig_n_ev
int eig_n_kr
int mg_levels
quda::mgarray< double > mg_eig_amax
double omega
void add_eigen_option_group(std::shared_ptr< QUDAApp > quda_app)
int eig_n_conv
double ape_smear_rho
quda::mgarray< QudaInverterType > setup_inv
bool mg_evolve_thin_updates
double clover_coeff
bool native_blas_lapack
int pipeline
std::array< double, 2 > blas_beta_re_im
QudaMatPCType matpc_type
quda::mgarray< double > setup_tol
int multishift
QudaFieldLocation location_ritz
quda::mgarray< double > coarse_solver_ca_lambda_min
QudaEigSpectrumType eig_spectrum
quda::mgarray< bool > mg_eig_use_poly_acc
double tol_precondition
quda::mgarray< double > mg_eig_qr_tol
QudaPrecision prec_eigensolver
quda::mgarray< QudaSolveType > coarse_solve_type
double tol_restart
auto & grid_z
bool compute_fatlong
std::array< int, 3 > blas_offsets
int Nsrc
double b5
quda::mgarray< double > mg_eig_tol
QudaSchwarzType precon_schwarz_type
bool unit_gauge
void add_deflation_option_group(std::shared_ptr< QUDAApp > quda_app)
quda::mgarray< int > mg_eig_batched_rotate
auto & blas_m
QudaPrecision prec_precondition
quda::mgarray< bool > mg_eig_use_normop
quda::mgarray< bool > mg_eig_use_eigen_qr
QudaPrecision prec
QudaReconstructType link_recon_eigensolver
quda::mgarray< QudaPrecision > mg_eig_save_prec
char eig_arpack_logfile[256]
auto & blas_k
void add_su3_option_group(std::shared_ptr< QUDAApp > quda_app)
QudaCABasis ca_basis
int n_naiks
quda::mgarray< double > smoother_tol
QudaSetupType setup_type
QudaContractType contract_type
double heatbath_beta_value
int heatbath_num_overrelax_per_step
int Lsdim
quda::mgarray< QudaSolveType > smoother_solve_type
double eofa_mq3
bool eig_require_convergence
int eig_max_restarts
int & tdim
quda::mgarray< bool > mg_eig_require_convergence
quda::mgarray< int > mg_eig_block_size
int & xdim
quda::mgarray< double > coarse_solver_tol
bool eig_compute_gamma5
QudaPrecision smoother_halo_prec
int maxiter_precondition
int eigcg_max_restarts
std::array< int, 3 > blas_mnk
quda::mgarray< QudaInverterType > smoother_type
quda::mgarray< int > setup_maxiter_refresh
void add_comms_option_group(std::shared_ptr< QUDAApp > quda_app)
int eofa_pm
QudaBLASOperation blas_trans_a
double eig_amin
quda::mgarray< int > mg_schwarz_cycle
std::array< int, 4 > gridsize_from_cmdline
int heatbath_warmup_steps
bool dagger
QudaMassNormalization normalization
QudaBLASDataType blas_data_type
quda::mgarray< std::array< int, 4 > > geo_block_size
QudaPrecision prec_sloppy
enum QudaSolveType_s QudaSolveType
enum QudaWFlowType_s QudaWFlowType
enum QudaBLASOperation_s QudaBLASOperation
enum QudaPrecision_s QudaPrecision
@ QUDA_TEST_VECTOR_SETUP
Definition: enum_quda.h:448
@ QUDA_NULL_VECTOR_SETUP
Definition: enum_quda.h:447
@ QUDA_CHEBYSHEV_BASIS
Definition: enum_quda.h:202
@ QUDA_POWER_BASIS
Definition: enum_quda.h:201
@ QUDA_WILSON_DSLASH
Definition: enum_quda.h:90
@ QUDA_TWISTED_CLOVER_DSLASH
Definition: enum_quda.h:100
@ QUDA_STAGGERED_DSLASH
Definition: enum_quda.h:97
@ QUDA_MOBIUS_DWF_DSLASH
Definition: enum_quda.h:95
@ QUDA_CLOVER_WILSON_DSLASH
Definition: enum_quda.h:91
@ QUDA_TWISTED_MASS_DSLASH
Definition: enum_quda.h:99
@ QUDA_DOMAIN_WALL_DSLASH
Definition: enum_quda.h:93
@ QUDA_ASQTAD_DSLASH
Definition: enum_quda.h:98
@ QUDA_MOBIUS_DWF_EOFA_DSLASH
Definition: enum_quda.h:96
@ QUDA_LAPLACE_DSLASH
Definition: enum_quda.h:101
@ QUDA_CLOVER_HASENBUSCH_TWIST_DSLASH
Definition: enum_quda.h:92
@ QUDA_DOMAIN_WALL_4D_DSLASH
Definition: enum_quda.h:94
@ QUDA_CUDA_FIELD_LOCATION
Definition: enum_quda.h:326
@ QUDA_CPU_FIELD_LOCATION
Definition: enum_quda.h:325
enum QudaTwistFlavorType_s QudaTwistFlavorType
@ QUDA_KAPPA_NORMALIZATION
Definition: enum_quda.h:226
@ QUDA_ASYMMETRIC_MASS_NORMALIZATION
Definition: enum_quda.h:228
@ QUDA_MASS_NORMALIZATION
Definition: enum_quda.h:227
@ QUDA_SILENT
Definition: enum_quda.h:265
@ QUDA_DEBUG_VERBOSE
Definition: enum_quda.h:268
@ QUDA_SUMMARIZE
Definition: enum_quda.h:266
@ QUDA_VERBOSE
Definition: enum_quda.h:267
@ QUDA_BLAS_DATATYPE_Z
Definition: enum_quda.h:480
@ QUDA_BLAS_DATATYPE_D
Definition: enum_quda.h:478
@ QUDA_BLAS_DATATYPE_C
Definition: enum_quda.h:479
@ QUDA_BLAS_DATATYPE_S
Definition: enum_quda.h:477
enum QudaTransferType_s QudaTransferType
enum QudaBLASDataOrder_s QudaBLASDataOrder
enum QudaTboundary_s QudaTboundary
@ QUDA_RECONSTRUCT_NO
Definition: enum_quda.h:70
@ QUDA_RECONSTRUCT_INVALID
Definition: enum_quda.h:76
@ QUDA_RECONSTRUCT_12
Definition: enum_quda.h:71
@ QUDA_RECONSTRUCT_13
Definition: enum_quda.h:74
@ QUDA_RECONSTRUCT_8
Definition: enum_quda.h:72
@ QUDA_RECONSTRUCT_9
Definition: enum_quda.h:73
@ QUDA_ANTI_PERIODIC_T
Definition: enum_quda.h:56
@ QUDA_PERIODIC_T
Definition: enum_quda.h:57
@ QUDA_MEMORY_DEVICE
Definition: enum_quda.h:13
enum QudaDslashType_s QudaDslashType
@ QUDA_TRANSFER_COARSE_KD
Definition: enum_quda.h:454
@ QUDA_TRANSFER_AGGREGATE
Definition: enum_quda.h:453
@ QUDA_TRANSFER_OPTIMIZED_KD
Definition: enum_quda.h:455
enum QudaSolutionType_s QudaSolutionType
enum QudaEigSpectrumType_s QudaEigSpectrumType
enum QudaInverterType_s QudaInverterType
@ QUDA_EIG_BLK_IR_ARNOLDI
Definition: enum_quda.h:140
@ QUDA_EIG_IR_ARNOLDI
Definition: enum_quda.h:139
@ QUDA_EIG_BLK_TR_LANCZOS
Definition: enum_quda.h:138
@ QUDA_EIG_TR_LANCZOS
Definition: enum_quda.h:137
enum QudaFieldLocation_s QudaFieldLocation
enum QudaMassNormalization_s QudaMassNormalization
enum QudaBLASDataType_s QudaBLASDataType
enum QudaExtLibType_s QudaExtLibType
enum QudaEigType_s QudaEigType
@ QUDA_MATPC_ODD_ODD_ASYMMETRIC
Definition: enum_quda.h:219
@ QUDA_MATPC_EVEN_EVEN_ASYMMETRIC
Definition: enum_quda.h:218
@ QUDA_MATPC_ODD_ODD
Definition: enum_quda.h:217
@ QUDA_MATPC_EVEN_EVEN
Definition: enum_quda.h:216
enum QudaMatPCType_s QudaMatPCType
enum QudaSetupType_s QudaSetupType
@ QUDA_BLAS_DATAORDER_COL
Definition: enum_quda.h:486
@ QUDA_BLAS_DATAORDER_ROW
Definition: enum_quda.h:485
@ QUDA_CA_CGNE_INVERTER
Definition: enum_quda.h:130
@ QUDA_GCR_INVERTER
Definition: enum_quda.h:109
@ QUDA_CGNE_INVERTER
Definition: enum_quda.h:124
@ QUDA_GMRESDR_PROJ_INVERTER
Definition: enum_quda.h:119
@ QUDA_INC_EIGCG_INVERTER
Definition: enum_quda.h:117
@ QUDA_FGMRESDR_INVERTER
Definition: enum_quda.h:121
@ QUDA_MR_INVERTER
Definition: enum_quda.h:110
@ QUDA_CA_CG_INVERTER
Definition: enum_quda.h:129
@ QUDA_MPCG_INVERTER
Definition: enum_quda.h:115
@ QUDA_CG3NE_INVERTER
Definition: enum_quda.h:127
@ QUDA_BICGSTABL_INVERTER
Definition: enum_quda.h:123
@ QUDA_CGNR_INVERTER
Definition: enum_quda.h:125
@ QUDA_GMRESDR_SH_INVERTER
Definition: enum_quda.h:120
@ QUDA_CG3NR_INVERTER
Definition: enum_quda.h:128
@ QUDA_PCG_INVERTER
Definition: enum_quda.h:114
@ QUDA_CA_GCR_INVERTER
Definition: enum_quda.h:132
@ QUDA_SD_INVERTER
Definition: enum_quda.h:112
@ QUDA_MPBICGSTAB_INVERTER
Definition: enum_quda.h:111
@ QUDA_CG_INVERTER
Definition: enum_quda.h:107
@ QUDA_CA_CGNR_INVERTER
Definition: enum_quda.h:131
@ QUDA_INVALID_INVERTER
Definition: enum_quda.h:133
@ QUDA_CG3_INVERTER
Definition: enum_quda.h:126
@ QUDA_EIGCG_INVERTER
Definition: enum_quda.h:116
@ QUDA_MG_INVERTER
Definition: enum_quda.h:122
@ QUDA_BICGSTAB_INVERTER
Definition: enum_quda.h:108
@ QUDA_GMRESDR_INVERTER
Definition: enum_quda.h:118
enum QudaMemoryType_s QudaMemoryType
@ QUDA_BLAS_OP_C
Definition: enum_quda.h:472
@ QUDA_BLAS_OP_N
Definition: enum_quda.h:470
@ QUDA_BLAS_OP_T
Definition: enum_quda.h:471
enum QudaReconstructType_s QudaReconstructType
@ QUDA_MATPC_DAG_SOLUTION
Definition: enum_quda.h:160
@ QUDA_MATPC_SOLUTION
Definition: enum_quda.h:159
@ QUDA_MATDAG_MAT_SOLUTION
Definition: enum_quda.h:158
@ QUDA_MATPCDAG_MATPC_SOLUTION
Definition: enum_quda.h:161
@ QUDA_MAT_SOLUTION
Definition: enum_quda.h:157
@ QUDA_DOUBLE_PRECISION
Definition: enum_quda.h:65
@ QUDA_SINGLE_PRECISION
Definition: enum_quda.h:64
@ QUDA_INVALID_PRECISION
Definition: enum_quda.h:66
@ QUDA_QUARTER_PRECISION
Definition: enum_quda.h:62
@ QUDA_HALF_PRECISION
Definition: enum_quda.h:63
@ QUDA_ADDITIVE_SCHWARZ
Definition: enum_quda.h:187
@ QUDA_MULTIPLICATIVE_SCHWARZ
Definition: enum_quda.h:188
@ QUDA_INVALID_SCHWARZ
Definition: enum_quda.h:189
enum QudaCABasis_s QudaCABasis
@ QUDA_EIGEN_EXTLIB
Definition: enum_quda.h:557
@ QUDA_MAGMA_EXTLIB
Definition: enum_quda.h:558
enum QudaContractType_s QudaContractType
@ QUDA_SPECTRUM_LM_EIG
Definition: enum_quda.h:147
@ QUDA_SPECTRUM_SM_EIG
Definition: enum_quda.h:148
@ QUDA_SPECTRUM_LR_EIG
Definition: enum_quda.h:149
@ QUDA_SPECTRUM_SR_EIG
Definition: enum_quda.h:150
@ QUDA_SPECTRUM_SI_EIG
Definition: enum_quda.h:152
@ QUDA_SPECTRUM_LI_EIG
Definition: enum_quda.h:151
enum QudaSchwarzType_s QudaSchwarzType
enum QudaVerbosity_s QudaVerbosity
@ QUDA_TWIST_SINGLET
Definition: enum_quda.h:400
@ QUDA_TWIST_NO
Definition: enum_quda.h:403
@ QUDA_TWIST_NONDEG_DOUBLET
Definition: enum_quda.h:401
@ QUDA_TWIST_DEG_DOUBLET
Definition: enum_quda.h:402
@ QUDA_DIRECT_SOLVE
Definition: enum_quda.h:167
@ QUDA_NORMERR_SOLVE
Definition: enum_quda.h:171
@ QUDA_NORMERR_PC_SOLVE
Definition: enum_quda.h:172
@ QUDA_NORMOP_PC_SOLVE
Definition: enum_quda.h:170
@ QUDA_DIRECT_PC_SOLVE
Definition: enum_quda.h:169
@ QUDA_NORMOP_SOLVE
Definition: enum_quda.h:168
@ QUDA_CONTRACT_TYPE_OPEN
Definition: enum_quda.h:523
@ QUDA_CONTRACT_TYPE_DR
Definition: enum_quda.h:524
@ QUDA_WFLOW_TYPE_SYMANZIK
Definition: enum_quda.h:550
@ QUDA_WFLOW_TYPE_WILSON
Definition: enum_quda.h:549
std::array< T, QUDA_MAX_MG_LEVEL > mgarray
::std::string string
Definition: gtest-port.h:891
internal::ParamGenerator< T > Range(T start, T end, IncrementT step)