QUDA  v1.1.0
A library for QCD on GPUs
quda.h
Go to the documentation of this file.
1 #pragma once
2 
12 #include <enum_quda.h>
13 #include <stdio.h> /* for FILE */
14 #include <quda_define.h>
15 #include <quda_constants.h>
16 
17 #ifndef __CUDACC_RTC__
18 #define double_complex double _Complex
19 #else // keep NVRTC happy since it can't handle C types
20 #define double_complex double2
21 #endif
22 
23 #ifdef __cplusplus
24 extern "C" {
25 #endif
26 
31  typedef struct QudaGaugeParam_s {
32  size_t struct_size;
35  int X[4];
37  double anisotropy;
38  double tadpole_coeff;
39  double scale;
65  int ga_pad;
69  int staple_pad;
71  int mom_ga_pad;
76  double i_mu;
78  int overlap;
89  size_t gauge_offset;
90  size_t mom_offset;
91  size_t site_size;
93 
94 
98  typedef struct QudaInvertParam_s {
99 
101  size_t struct_size;
102 
109  double mass;
110  double kappa;
112  double m5;
113  int Ls;
125  double eofa_shift;
126  int eofa_pm;
127  double mq1;
128  double mq2;
129  double mq3;
130 
131  double mu;
132  double epsilon;
136  int laplace3D;
138  double tol;
139  double tol_restart;
140  double tol_hq;
143  double true_res;
144  double true_res_hq;
145  int maxiter;
146  double reliable_delta;
160 
165 
170 
175 
180 
183 
184  int pipeline;
188  int num_src;
196 
197  int overlap;
201 
204 
207 
210 
213 
216 
219 
222 
226  double action[2];
227 
259  double clover_csw;
260  double clover_coeff;
261  double clover_rho;
264  double trlogA[2];
273  int sp_pad;
274  int cl_pad;
276  int iter;
277  double gflops;
278  double secs;
283  int Nsteps;
284 
287 
288  /*
289  * The following parameters are related to the solver
290  * preconditioner, if enabled.
291  */
292 
298 
301 
304 
306  void *eig_param;
307 
310 
315 
318 
321 
323  double omega;
324 
327 
330 
333 
336 
339 
349 
356  int n_ev;
362  int rhs_idx;
366  double eigenval_tol;
372  double inc_tol;
373 
376 
379 
382 
385 
388 
391 
394 
397 
400 
404 
405  // Parameter set for solving eigenvalue problems.
406  typedef struct QudaEigParam_s {
408  size_t struct_size;
409 
410  // EIGENSOLVER PARAMS
411  //-------------------------------------------------
414 
417 
420 
422  int poly_deg;
423 
425  double a_min;
426  double a_max;
427 
432 
436 
443 
451 
454 
457 
461 
464 
467 
469  int n_ev;
471  int n_kr;
475  int n_conv;
479  double tol;
481  double qr_tol;
490 
494  char arpack_logfile[512];
495 
497  char QUDA_logfile[512];
498 
499  //-------------------------------------------------
500 
501  // EIG-CG PARAMS
502  //-------------------------------------------------
503  int nk;
504  int np;
505 
508 
511 
514 
517 
520 
522  char vec_infile[256];
523 
525  char vec_outfile[256];
526 
529 
534 
536  double gflops;
537 
539  double secs;
540 
543  //-------------------------------------------------
545 
546  typedef struct QudaMultigridParam_s {
547 
549  size_t struct_size;
550 
552 
554 
556  int n_level;
557 
560 
563 
566 
569 
572 
575 
578 
581 
584 
587 
590 
593 
596 
599 
602 
605 
608 
611 
614 
617 
620 
623 
626 
629 
632 
635 
638 
641 
644 
647 
650 
653 
656 
660 
663 
666 
669 
672 
675 
678 
683 
686 
689 
692 
695 
698 
701 
704 
707 
710 
713 
716 
718  double gflops;
719 
721  double secs;
722 
725 
728 
731 
735 
736  typedef struct QudaGaugeObservableParam_s {
737  size_t struct_size;
740  double plaquette[3];
742  double qcharge;
743  double energy[3];
747 
748  typedef struct QudaBLASParam_s {
749  size_t struct_size;
753  int m;
754  int n;
755  int k;
756  int lda;
757  int ldb;
758  int ldc;
759  int a_offset;
760  int b_offset;
761  int c_offset;
762  int a_stride;
763  int b_stride;
764  int c_stride;
774 
775  /*
776  * Interface functions, found in interface_quda.cpp
777  */
778 
804  void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[],
805  FILE *outfile);
806 
817  typedef int (*QudaCommsMap)(const int *coords, void *fdata);
818 
823  void qudaSetCommHandle(void *mycomm);
824 
852  void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata);
853 
864  void initQudaDevice(int device);
865 
872  void initQudaMemory();
873 
883  void initQuda(int device);
884 
888  void endQuda(void);
889 
895  void updateR();
896 
905 
914 
923 
932 
941 
950 
956 
962 
968 
974 
980 
986 
992  void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param);
993 
997  void freeGaugeQuda(void);
998 
1004  void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param);
1005 
1013  void loadCloverQuda(void *h_clover, void *h_clovinv,
1015 
1019  void freeCloverQuda(void);
1020 
1030  void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta,
1031  QudaEigParam *eig_param);
1032 
1041  void eigensolveQuda(void **h_evecs, double_complex *h_evals, QudaEigParam *param);
1042 
1052  void invertQuda(void *h_x, void *h_b, QudaInvertParam *param);
1053 
1069  void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge, QudaGaugeParam *gauge_param);
1070 
1081  void invertMultiSrcStaggeredQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *milc_fatlinks,
1082  void *milc_longlinks, QudaGaugeParam *gauge_param);
1083 
1095  void invertMultiSrcCloverQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge,
1096  QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv);
1097 
1105  void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param);
1106 
1115 
1122  void destroyMultigridQuda(void *mg_instance);
1123 
1131  void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param);
1132 
1140  void dumpMultigridQuda(void *mg_instance, QudaMultigridParam *param);
1141 
1150  void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity);
1151 
1165  void dslashMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge,
1179  void dslashMultiSrcStaggeredQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity,
1180  void *milc_fatlinks, void *milc_longlinks, QudaGaugeParam *gauge_param);
1181 
1194  void dslashMultiSrcCloverQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge,
1195  QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv);
1196 
1206  void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int inverse);
1207 
1215  void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
1216 
1224  void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
1225 
1226 
1227  /*
1228  * The following routines are temporary additions used by the HISQ
1229  * link-fattening code.
1230  */
1231 
1232  void set_dim(int *);
1233  void pack_ghost(void **cpuLink, void **cpuGhost, int nFace,
1234  QudaPrecision precision);
1235 
1236  void computeKSLinkQuda(void* fatlink, void* longlink, void* ulink, void* inlink,
1237  double *path_coeff, QudaGaugeParam *param);
1238 
1246  void momResidentQuda(void *mom, QudaGaugeParam *param);
1247 
1261  int computeGaugeForceQuda(void* mom, void* sitelink, int*** input_path_buf, int* path_length,
1262  double* loop_coeff, int num_paths, int max_length, double dt,
1263  QudaGaugeParam* qudaGaugeParam);
1264 
1276  void updateGaugeFieldQuda(void* gauge, void* momentum, double dt,
1277  int conj_mom, int exact, QudaGaugeParam* param);
1278 
1288  void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param);
1289 
1298  void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param);
1299 
1308  double momActionQuda(void* momentum, QudaGaugeParam* param);
1309 
1318  void* createGaugeFieldQuda(void* gauge, int geometry, QudaGaugeParam* param);
1319 
1327  void saveGaugeFieldQuda(void* outGauge, void* inGauge, QudaGaugeParam* param);
1328 
1334  void destroyGaugeFieldQuda(void* gauge);
1335 
1342 
1361  void computeCloverForceQuda(void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck,
1362  int nvector, double multiplicity, void *gauge,
1364 
1376  void computeStaggeredForceQuda(void *mom, double dt, double delta, void *gauge, void **x, QudaGaugeParam *gauge_param,
1377  QudaInvertParam *invert_param);
1378 
1394  void computeHISQForceQuda(void* momentum,
1395  double dt,
1396  const double level2_coeff[6],
1397  const double fat7_coeff[6],
1398  const void* const w_link,
1399  const void* const v_link,
1400  const void* const u_link,
1401  void** quark,
1402  int num,
1403  int num_naik,
1404  double** coeff,
1406 
1418  void gaussGaugeQuda(unsigned long long seed, double sigma);
1419 
1424  void plaqQuda(double plaq[3]);
1425 
1431  void copyExtendedResidentGaugeQuda(void* resident_gauge, QudaFieldLocation loc);
1432 
1443  void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param, unsigned int n_steps, double alpha);
1444 
1451  void performAPEnStep(unsigned int n_steps, double alpha, int meas_interval);
1452 
1459  void performSTOUTnStep(unsigned int n_steps, double rho, int meas_interval);
1460 
1468  void performOvrImpSTOUTnStep(unsigned int n_steps, double rho, double epsilon, int meas_interval);
1469 
1477  void performWFlownStep(unsigned int n_steps, double step_size, int meas_interval, QudaWFlowType wflow_type);
1478 
1488 
1498  void contractQuda(const void *x, const void *y, void *result, const QudaContractType cType, QudaInvertParam *param,
1499  const int *X);
1500 
1514  int computeGaugeFixingOVRQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps,
1515  const unsigned int verbose_interval, const double relax_boost, const double tolerance,
1516  const unsigned int reunit_interval, const unsigned int stopWtheta,
1517  QudaGaugeParam *param, double *timeinfo);
1532  int computeGaugeFixingFFTQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps,
1533  const unsigned int verbose_interval, const double alpha, const unsigned int autotune,
1534  const double tolerance, const unsigned int stopWtheta, QudaGaugeParam *param,
1535  double *timeinfo);
1536 
1545  void blasGEMMQuda(void *arrayA, void *arrayB, void *arrayC, QudaBoolean native, QudaBLASParam *param);
1546 
1551  void flushChronoQuda(int index);
1552 
1557  void openMagma();
1558 
1559  void closeMagma();
1560 
1567 
1571  void destroyDeflationQuda(void *df_instance);
1572 
1573  void setMPICommHandleQuda(void *mycomm);
1574 
1575 #ifdef __cplusplus
1576 }
1577 #endif
1578 
1579 // remove NVRTC WAR
1580 #undef double_complex
1581 
1582 /* #include <quda_new_interface.h> */
1583 
1584 
QudaWFlowType wflow_type
double tol
QudaVerbosity verbosity
double epsilon
QudaParity parity
Definition: covdev_test.cpp:40
QudaGaugeParam gauge_param
Definition: covdev_test.cpp:26
cpuGaugeField * cpuLink
Definition: covdev_test.cpp:29
QudaInvertParam inv_param
Definition: covdev_test.cpp:27
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
enum QudaSolveType_s QudaSolveType
enum QudaWFlowType_s QudaWFlowType
enum QudaBLASOperation_s QudaBLASOperation
enum QudaPrecision_s QudaPrecision
enum QudaUseInitGuess_s QudaUseInitGuess
enum QudaGaugeFixed_s QudaGaugeFixed
enum QudaStaggeredPhase_s QudaStaggeredPhase
enum QudaTwistFlavorType_s QudaTwistFlavorType
enum QudaDagType_s QudaDagType
enum QudaTransferType_s QudaTransferType
enum QudaMultigridCycleType_s QudaMultigridCycleType
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
enum QudaBLASDataOrder_s QudaBLASDataOrder
enum QudaPreserveSource_s QudaPreserveSource
enum QudaTboundary_s QudaTboundary
enum QudaDslashType_s QudaDslashType
enum QudaBoolean_s QudaBoolean
enum QudaSolutionType_s QudaSolutionType
enum QudaEigSpectrumType_s QudaEigSpectrumType
enum QudaInverterType_s QudaInverterType
enum QudaFieldLocation_s QudaFieldLocation
enum QudaMassNormalization_s QudaMassNormalization
enum QudaBLASDataType_s QudaBLASDataType
enum QudaExtLibType_s QudaExtLibType
enum QudaEigType_s QudaEigType
enum QudaMatPCType_s QudaMatPCType
enum QudaSetupType_s QudaSetupType
enum QudaComputeNullVector_s QudaComputeNullVector
enum QudaResidualType_s QudaResidualType
enum QudaTune_s QudaTune
enum QudaMemoryType_s QudaMemoryType
enum QudaReconstructType_s QudaReconstructType
enum QudaCABasis_s QudaCABasis
enum QudaContractType_s QudaContractType
enum QudaSchwarzType_s QudaSchwarzType
enum QudaVerbosity_s QudaVerbosity
enum QudaDiracFieldOrder_s QudaDiracFieldOrder
enum QudaGammaBasis_s QudaGammaBasis
enum QudaSolverNormalization_s QudaSolverNormalization
enum QudaParity_s QudaParity
enum QudaLinkType_s QudaLinkType
__device__ __host__ Matrix< T, 3 > inverse(const Matrix< T, 3 > &u)
Definition: quda_matrix.h:605
QudaGaugeParam param
Definition: pack_test.cpp:18
struct QudaInvertParam_s QudaInvertParam
void initQudaMemory()
double momActionQuda(void *momentum, QudaGaugeParam *param)
void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge, QudaGaugeParam *gauge_param)
Perform the solve like @invertQuda but for multiple rhs by spliting the comm grid into sub-partitions...
void gaussGaugeQuda(unsigned long long seed, double sigma)
Generate Gaussian distributed fields and store in the resident gauge field. We create a Gaussian-dist...
void dumpMultigridQuda(void *mg_instance, QudaMultigridParam *param)
Dump the null-space vectors to disk.
void createCloverQuda(QudaInvertParam *param)
void * createGaugeFieldQuda(void *gauge, int geometry, QudaGaugeParam *param)
void printQudaMultigridParam(QudaMultigridParam *param)
Definition: check_params.h:689
void destroyGaugeFieldQuda(void *gauge)
void computeStaggeredForceQuda(void *mom, double dt, double delta, void *gauge, void **x, QudaGaugeParam *gauge_param, QudaInvertParam *invert_param)
QudaBLASParam newQudaBLASParam(void)
void destroyDeflationQuda(void *df_instance)
void momResidentQuda(void *mom, QudaGaugeParam *param)
int computeGaugeForceQuda(void *mom, void *sitelink, int ***input_path_buf, int *path_length, double *loop_coeff, int num_paths, int max_length, double dt, QudaGaugeParam *qudaGaugeParam)
void printQudaInvertParam(QudaInvertParam *param)
Definition: check_params.h:342
void setMPICommHandleQuda(void *mycomm)
void * newMultigridQuda(QudaMultigridParam *param)
void computeCloverForceQuda(void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck, int nvector, double multiplicity, void *gauge, QudaGaugeParam *gauge_param, QudaInvertParam *inv_param)
void plaqQuda(double plaq[3])
void performAPEnStep(unsigned int n_steps, double alpha, int meas_interval)
QudaGaugeParam newQudaGaugeParam(void)
void blasGEMMQuda(void *arrayA, void *arrayB, void *arrayC, QudaBoolean native, QudaBLASParam *param)
Strided Batched GEMM.
void invertMultiSrcCloverQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge, QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv)
Really the same with @invertMultiSrcQuda but for clover-style fermions, by accepting pointers to dire...
void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param)
struct QudaGaugeObservableParam_s QudaGaugeObservableParam
void eigensolveQuda(void **h_evecs, double_complex *h_evals, QudaEigParam *param)
void printQudaEigParam(QudaEigParam *param)
Definition: check_params.h:158
void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[], FILE *outfile)
void saveGaugeFieldQuda(void *outGauge, void *inGauge, QudaGaugeParam *param)
void initQudaDevice(int device)
struct QudaEigParam_s QudaEigParam
void freeGaugeQuda(void)
void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
void initQuda(int device)
void openMagma()
void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param, unsigned int n_steps, double alpha)
void invertMultiSrcStaggeredQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, void *milc_fatlinks, void *milc_longlinks, QudaGaugeParam *gauge_param)
Really the same with @invertMultiSrcQuda but for staggered-style fermions, by accepting pointers to f...
void dslashMultiSrcCloverQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge, QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv)
Really the same with @dslashMultiSrcQuda but for clover-style fermions, by accepting pointers to dire...
int computeGaugeFixingFFTQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double alpha, const unsigned int autotune, const double tolerance, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
Gauge fixing with Steepest descent method with FFTs with support for single GPU only.
void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param)
void dslashMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge, QudaGaugeParam *gauge_param)
Perform the solve like @dslashQuda but for multiple rhs by spliting the comm grid into sub-partitions...
void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
void qudaSetCommHandle(void *mycomm)
void contractQuda(const void *x, const void *y, void *result, const QudaContractType cType, QudaInvertParam *param, const int *X)
QudaMultigridParam newQudaMultigridParam(void)
QudaGaugeObservableParam newQudaGaugeObservableParam(void)
void set_dim(int *)
void updateGaugeFieldQuda(void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
void pack_ghost(void **cpuLink, void **cpuGhost, int nFace, QudaPrecision precision)
int computeGaugeFixingOVRQuda(void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double relax_boost, const double tolerance, const unsigned int reunit_interval, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
Gauge fixing with overrelaxation with support for single and multi GPU.
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
void freeCloverQuda(void)
struct QudaBLASParam_s QudaBLASParam
int(* QudaCommsMap)(const int *coords, void *fdata)
Definition: quda.h:817
void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param)
void computeKSLinkQuda(void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param)
struct QudaGaugeParam_s QudaGaugeParam
void performWFlownStep(unsigned int n_steps, double step_size, int meas_interval, QudaWFlowType wflow_type)
QudaInvertParam newQudaInvertParam(void)
void dslashMultiSrcStaggeredQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *milc_fatlinks, void *milc_longlinks, QudaGaugeParam *gauge_param)
Really the same with @dslashMultiSrcQuda but for staggered-style fermions, by accepting pointers to f...
void flushChronoQuda(int index)
Flush the chronological history for the given index.
QudaEigParam newQudaEigParam(void)
void endQuda(void)
void copyExtendedResidentGaugeQuda(void *resident_gauge, QudaFieldLocation loc)
struct QudaMultigridParam_s QudaMultigridParam
void gaugeObservablesQuda(QudaGaugeObservableParam *param)
Calculates a variety of gauge-field observables. If a smeared gauge field is presently loaded (in gau...
void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param)
Updates the multigrid preconditioner for the new gauge / clover field.
void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int inverse)
void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata)
void destroyMultigridQuda(void *mg_instance)
Free resources allocated by the multigrid solver.
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
void computeHISQForceQuda(void *momentum, double dt, const double level2_coeff[6], const double fat7_coeff[6], const void *const w_link, const void *const v_link, const void *const u_link, void **quark, int num, int num_naik, double **coeff, QudaGaugeParam *param)
void performSTOUTnStep(unsigned int n_steps, double rho, int meas_interval)
void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta, QudaEigParam *eig_param)
void printQudaBLASParam(QudaBLASParam *param)
Definition: check_params.h:981
void * newDeflationQuda(QudaEigParam *param)
void updateR()
update the radius for halos.
void printQudaGaugeParam(QudaGaugeParam *param)
Definition: check_params.h:40
void performOvrImpSTOUTnStep(unsigned int n_steps, double rho, double epsilon, int meas_interval)
void printQudaGaugeObservableParam(QudaGaugeObservableParam *param)
Definition: check_params.h:943
void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param)
#define double_complex
Definition: quda.h:18
void closeMagma()
void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
#define QUDA_MAX_DWF_LS
Maximum length of the Ls dimension for domain-wall fermions.
#define QUDA_MAX_MG_LEVEL
Maximum number of multi-grid levels. This number may be increased if needed.
#define QUDA_MAX_DIM
Maximum number of dimensions supported by QUDA. In practice, no routines make use of more than 5.
#define QUDA_MAX_MULTI_SHIFT
Maximum number of shifts supported by the multi-shift solver. This number may be changed if need be.
int c_offset
Definition: quda.h:761
double_complex alpha
Definition: quda.h:766
int a_stride
Definition: quda.h:762
int b_stride
Definition: quda.h:763
size_t struct_size
Definition: quda.h:749
QudaBLASDataOrder data_order
Definition: quda.h:772
int c_stride
Definition: quda.h:764
int b_offset
Definition: quda.h:760
QudaBLASOperation trans_a
Definition: quda.h:751
double_complex beta
Definition: quda.h:767
QudaBLASDataType data_type
Definition: quda.h:771
int a_offset
Definition: quda.h:759
int batch_count
Definition: quda.h:769
QudaBLASOperation trans_b
Definition: quda.h:752
void * preserve_deflation_space
Definition: quda.h:435
QudaEigSpectrumType spectrum
Definition: quda.h:466
QudaBoolean use_eigen_qr
Definition: quda.h:453
QudaPrecision save_prec
Definition: quda.h:528
QudaEigType eig_type
Definition: quda.h:416
int check_interval
Definition: quda.h:483
QudaBoolean import_vectors
Definition: quda.h:507
QudaBoolean use_dagger
Definition: quda.h:449
QudaBoolean arpack_check
Definition: quda.h:492
int n_ev_deflate
Definition: quda.h:477
QudaBoolean compute_svd
Definition: quda.h:456
QudaBoolean io_parity_inflate
Definition: quda.h:533
int batched_rotate
Definition: quda.h:487
QudaBoolean use_poly_acc
Definition: quda.h:419
double secs
Definition: quda.h:539
char vec_outfile[256]
Definition: quda.h:525
QudaBoolean use_norm_op
Definition: quda.h:450
QudaPrecision cuda_prec_ritz
Definition: quda.h:510
QudaBoolean compute_gamma5
Definition: quda.h:460
char vec_infile[256]
Definition: quda.h:522
QudaFieldLocation location
Definition: quda.h:516
QudaExtLibType extlib_type
Definition: quda.h:542
int poly_deg
Definition: quda.h:422
double a_min
Definition: quda.h:425
char QUDA_logfile[512]
Definition: quda.h:497
double gflops
Definition: quda.h:536
double tol
Definition: quda.h:479
QudaBoolean preserve_deflation
Definition: quda.h:431
char arpack_logfile[512]
Definition: quda.h:494
double a_max
Definition: quda.h:426
int nLockedMax
Definition: quda.h:473
int block_size
Definition: quda.h:489
QudaBoolean run_verify
Definition: quda.h:519
QudaBoolean require_convergence
Definition: quda.h:463
int n_ev
Definition: quda.h:469
QudaBoolean preserve_evals
Definition: quda.h:442
int max_restarts
Definition: quda.h:485
QudaInvertParam * invert_param
Definition: quda.h:413
QudaMemoryType mem_type_ritz
Definition: quda.h:513
int n_kr
Definition: quda.h:471
int n_conv
Definition: quda.h:475
double qr_tol
Definition: quda.h:481
size_t struct_size
Definition: quda.h:408
QudaBoolean compute_qcharge_density
Definition: quda.h:744
QudaBoolean compute_qcharge
Definition: quda.h:741
QudaBoolean su_project
Definition: quda.h:738
QudaBoolean compute_plaquette
Definition: quda.h:739
double tadpole_coeff
Definition: quda.h:38
QudaReconstructType reconstruct_precondition
Definition: quda.h:58
int llfat_ga_pad
Definition: quda.h:70
double anisotropy
Definition: quda.h:37
size_t mom_offset
Definition: quda.h:90
QudaReconstructType reconstruct
Definition: quda.h:49
QudaPrecision cuda_prec_precondition
Definition: quda.h:57
int ga_pad
Definition: quda.h:65
QudaLinkType type
Definition: quda.h:41
double scale
Definition: quda.h:39
int make_resident_mom
Definition: quda.h:85
int return_result_mom
Definition: quda.h:87
size_t gauge_offset
Definition: quda.h:89
int use_resident_gauge
Definition: quda.h:82
int mom_ga_pad
Definition: quda.h:71
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:54
QudaFieldLocation location
Definition: quda.h:33
QudaPrecision cuda_prec_sloppy
Definition: quda.h:51
size_t struct_size
Definition: quda.h:32
QudaReconstructType reconstruct_sloppy
Definition: quda.h:52
QudaGaugeFixed gauge_fix
Definition: quda.h:63
int staple_pad
Definition: quda.h:69
QudaGaugeFieldOrder gauge_order
Definition: quda.h:42
int make_resident_gauge
Definition: quda.h:84
double i_mu
Definition: quda.h:76
QudaPrecision cuda_prec
Definition: quda.h:48
size_t site_size
Definition: quda.h:91
QudaReconstructType reconstruct_eigensolver
Definition: quda.h:61
QudaStaggeredPhase staggered_phase_type
Definition: quda.h:73
int X[4]
Definition: quda.h:35
QudaPrecision cpu_prec
Definition: quda.h:46
int site_ga_pad
Definition: quda.h:67
int use_resident_mom
Definition: quda.h:83
int overlap
Definition: quda.h:78
QudaReconstructType reconstruct_refinement_sloppy
Definition: quda.h:55
int staggered_phase_applied
Definition: quda.h:74
QudaTboundary t_boundary
Definition: quda.h:44
int return_result_gauge
Definition: quda.h:86
int overwrite_mom
Definition: quda.h:80
QudaPrecision cuda_prec_eigensolver
Definition: quda.h:60
double ca_lambda_max
Definition: quda.h:332
int compute_clover
Definition: quda.h:266
int use_sloppy_partial_accumulator
Definition: quda.h:149
QudaSolveType solve_type
Definition: quda.h:229
double gflops
Definition: quda.h:277
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:240
int max_res_increase
Definition: quda.h:164
int maxiter_precondition
Definition: quda.h:320
QudaSolutionType solution_type
Definition: quda.h:228
QudaBoolean native_blas_lapack
Definition: quda.h:402
double tol_hq
Definition: quda.h:140
double tol_precondition
Definition: quda.h:317
double iter_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:212
QudaCABasis ca_basis
Definition: quda.h:326
double omega
Definition: quda.h:323
QudaCloverFieldOrder clover_order
Definition: quda.h:256
int compute_true_res
Definition: quda.h:142
QudaSchwarzType schwarz_type
Definition: quda.h:338
double eigenval_tol
Definition: quda.h:366
int chrono_make_resident
Definition: quda.h:381
QudaMassNormalization mass_normalization
Definition: quda.h:232
int use_alternative_reliable
Definition: quda.h:148
QudaPreserveSource preserve_source
Definition: quda.h:235
size_t struct_size
Definition: quda.h:101
double tol_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:206
int solution_accumulator_pipeline
Definition: quda.h:159
int chrono_use_resident
Definition: quda.h:387
double tol_restart
Definition: quda.h:139
QudaFieldLocation clover_location
Definition: quda.h:248
QudaResidualType residual_type
Definition: quda.h:348
double mq2
Definition: quda.h:128
int precondition_cycle
Definition: quda.h:335
double true_res
Definition: quda.h:143
int heavy_quark_check
Definition: quda.h:182
double inc_tol
Definition: quda.h:372
QudaPrecision clover_cuda_prec_refinement_sloppy
Definition: quda.h:252
int make_resident_solution
Definition: quda.h:375
int num_offset
Definition: quda.h:186
QudaPrecision cuda_prec_eigensolver
Definition: quda.h:242
QudaExtLibType extlib_type
Definition: quda.h:399
double mass
Definition: quda.h:109
int max_res_increase_total
Definition: quda.h:169
QudaPrecision clover_cuda_prec
Definition: quda.h:250
int max_search_dim
Definition: quda.h:360
double m5
Definition: quda.h:112
int return_clover
Definition: quda.h:268
double reliable_delta_refinement
Definition: quda.h:147
int eigcg_max_restarts
Definition: quda.h:368
int split_grid[QUDA_MAX_DIM]
Definition: quda.h:195
int compute_clover_inverse
Definition: quda.h:267
QudaTwistFlavorType twist_flavor
Definition: quda.h:134
double ca_lambda_min
Definition: quda.h:329
QudaPrecision clover_cpu_prec
Definition: quda.h:249
QudaBoolean deflate
Definition: quda.h:309
QudaPrecision cuda_prec_ritz
Definition: quda.h:352
int max_hq_res_restart_total
Definition: quda.h:179
QudaDslashType dslash_type
Definition: quda.h:106
int return_clover_inverse
Definition: quda.h:269
double mq1
Definition: quda.h:127
int compute_clover_trlog
Definition: quda.h:263
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:253
int max_hq_res_increase
Definition: quda.h:174
QudaPrecision cuda_prec
Definition: quda.h:238
double mu
Definition: quda.h:131
QudaVerbosity verbosity
Definition: quda.h:271
double_complex b_5[QUDA_MAX_DWF_LS]
Definition: quda.h:115
double true_res_hq
Definition: quda.h:144
double clover_rho
Definition: quda.h:261
QudaPrecision chrono_precision
Definition: quda.h:396
int chrono_replace_last
Definition: quda.h:384
QudaDslashType dslash_type_precondition
Definition: quda.h:312
double clover_coeff
Definition: quda.h:260
double trlogA[2]
Definition: quda.h:264
int gcrNkrylov
Definition: quda.h:286
double eofa_shift
Definition: quda.h:125
double reliable_delta
Definition: quda.h:146
int num_src_per_sub_partition
Definition: quda.h:190
double secs
Definition: quda.h:278
QudaVerbosity verbosity_precondition
Definition: quda.h:314
QudaPrecision clover_cuda_prec_eigensolver
Definition: quda.h:254
double offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:200
QudaDagType dagger
Definition: quda.h:231
QudaInverterType inv_type
Definition: quda.h:107
int max_restart_num
Definition: quda.h:370
QudaTune tune
Definition: quda.h:280
QudaSolverNormalization solver_normalization
Definition: quda.h:233
double epsilon
Definition: quda.h:132
double clover_csw
Definition: quda.h:259
int compute_action
Definition: quda.h:221
QudaPrecision cpu_prec
Definition: quda.h:237
double residue[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:218
double true_res_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:209
int chrono_max_dim
Definition: quda.h:390
int deflation_grid
Definition: quda.h:364
QudaMatPCType matpc_type
Definition: quda.h:230
int chrono_index
Definition: quda.h:393
double_complex c_5[QUDA_MAX_DWF_LS]
Definition: quda.h:116
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:251
double tol_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:203
double mq3
Definition: quda.h:129
QudaInverterType inv_type_precondition
Definition: quda.h:297
QudaGammaBasis gamma_basis
Definition: quda.h:246
void * eig_param
Definition: quda.h:306
double tol
Definition: quda.h:138
QudaPrecision cuda_prec_sloppy
Definition: quda.h:239
QudaFieldLocation input_location
Definition: quda.h:103
QudaFieldLocation output_location
Definition: quda.h:104
void * deflation_op
Definition: quda.h:303
QudaPrecision cuda_prec_precondition
Definition: quda.h:241
double action[2]
Definition: quda.h:226
int use_resident_solution
Definition: quda.h:378
QudaDiracFieldOrder dirac_order
Definition: quda.h:244
double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT]
Definition: quda.h:215
QudaUseInitGuess use_init_guess
Definition: quda.h:257
double kappa
Definition: quda.h:110
void * preconditioner
Definition: quda.h:300
QudaBoolean pre_orthonormalize
Definition: quda.h:607
int smoother_schwarz_cycle[QUDA_MAX_MG_LEVEL]
Definition: quda.h:655
QudaBoolean run_verify
Definition: quda.h:691
double coarse_solver_ca_lambda_min[QUDA_MAX_MG_LEVEL]
Definition: quda.h:628
double setup_ca_lambda_max[QUDA_MAX_MG_LEVEL]
Definition: quda.h:601
QudaSolutionType coarse_grid_solution_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:659
double coarse_solver_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:616
QudaMultigridCycleType cycle_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:665
double omega[QUDA_MAX_MG_LEVEL]
Definition: quda.h:646
QudaBoolean thin_update_only
Definition: quda.h:733
int setup_maxiter_refresh[QUDA_MAX_MG_LEVEL]
Definition: quda.h:589
QudaPrecision precision_null[QUDA_MAX_MG_LEVEL]
Definition: quda.h:568
QudaBoolean use_mma
Definition: quda.h:730
char vec_infile[QUDA_MAX_MG_LEVEL][256]
Definition: quda.h:703
int nu_post[QUDA_MAX_MG_LEVEL]
Definition: quda.h:643
QudaBoolean use_eig_solver[QUDA_MAX_MG_LEVEL]
Definition: quda.h:677
int coarse_solver_maxiter[QUDA_MAX_MG_LEVEL]
Definition: quda.h:619
QudaEigParam * eig_param[QUDA_MAX_MG_LEVEL]
Definition: quda.h:553
int n_vec[QUDA_MAX_MG_LEVEL]
Definition: quda.h:565
int geo_block_size[QUDA_MAX_MG_LEVEL][QUDA_MAX_DIM]
Definition: quda.h:559
QudaInverterType coarse_solver[QUDA_MAX_MG_LEVEL]
Definition: quda.h:613
QudaBoolean coarse_guess
Definition: quda.h:712
int coarse_solver_ca_basis_size[QUDA_MAX_MG_LEVEL]
Definition: quda.h:625
QudaTransferType transfer_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:727
QudaFieldLocation setup_location[QUDA_MAX_MG_LEVEL]
Definition: quda.h:674
QudaBoolean post_orthonormalize
Definition: quda.h:610
int num_setup_iter[QUDA_MAX_MG_LEVEL]
Definition: quda.h:580
double mu_factor[QUDA_MAX_MG_LEVEL]
Definition: quda.h:724
QudaCABasis setup_ca_basis[QUDA_MAX_MG_LEVEL]
Definition: quda.h:592
QudaInverterType setup_inv_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:577
QudaInverterType smoother[QUDA_MAX_MG_LEVEL]
Definition: quda.h:634
int setup_maxiter[QUDA_MAX_MG_LEVEL]
Definition: quda.h:586
int setup_ca_basis_size[QUDA_MAX_MG_LEVEL]
Definition: quda.h:595
QudaBoolean preserve_deflation
Definition: quda.h:715
double coarse_solver_ca_lambda_max[QUDA_MAX_MG_LEVEL]
Definition: quda.h:631
size_t struct_size
Definition: quda.h:549
char vec_outfile[QUDA_MAX_MG_LEVEL][256]
Definition: quda.h:709
QudaBoolean run_low_mode_check
Definition: quda.h:694
double setup_ca_lambda_min[QUDA_MAX_MG_LEVEL]
Definition: quda.h:598
double smoother_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:637
QudaCABasis coarse_solver_ca_basis[QUDA_MAX_MG_LEVEL]
Definition: quda.h:622
int nu_pre[QUDA_MAX_MG_LEVEL]
Definition: quda.h:640
QudaSetupType setup_type
Definition: quda.h:604
QudaBoolean vec_store[QUDA_MAX_MG_LEVEL]
Definition: quda.h:706
QudaSolveType smoother_solve_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:662
QudaFieldLocation location[QUDA_MAX_MG_LEVEL]
Definition: quda.h:671
QudaVerbosity verbosity[QUDA_MAX_MG_LEVEL]
Definition: quda.h:574
QudaBoolean setup_minimize_memory
Definition: quda.h:682
int spin_block_size[QUDA_MAX_MG_LEVEL]
Definition: quda.h:562
int n_block_ortho[QUDA_MAX_MG_LEVEL]
Definition: quda.h:571
QudaBoolean generate_all_levels
Definition: quda.h:688
QudaComputeNullVector compute_null_vector
Definition: quda.h:685
QudaSchwarzType smoother_schwarz_type[QUDA_MAX_MG_LEVEL]
Definition: quda.h:652
QudaBoolean run_oblique_proj_check
Definition: quda.h:697
QudaInvertParam * invert_param
Definition: quda.h:551
double setup_tol[QUDA_MAX_MG_LEVEL]
Definition: quda.h:583
QudaPrecision smoother_halo_precision[QUDA_MAX_MG_LEVEL]
Definition: quda.h:649
QudaBoolean global_reduction[QUDA_MAX_MG_LEVEL]
Definition: quda.h:668
QudaBoolean vec_load[QUDA_MAX_MG_LEVEL]
Definition: quda.h:700