QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
quda.h
Go to the documentation of this file.
1 #ifndef _QUDA_H
2 #define _QUDA_H
3 
13 #include <enum_quda.h>
14 #include <stdio.h> /* for FILE */
15 #include <quda_define.h>
16 #include <quda_constants.h>
17 
18 #ifndef __CUDACC_RTC__
19 #define double_complex double _Complex
20 #else // keep NVRTC happy since it can't handle C types
21 #define double_complex double2
22 #endif
23 
24 #ifdef __cplusplus
25 extern "C" {
26 #endif
27 
32  typedef struct QudaGaugeParam_s {
33 
36  int X[4];
38  double anisotropy;
39  double tadpole_coeff;
40  double scale;
63  int ga_pad;
67  int staple_pad;
69  int mom_ga_pad;
74  double i_mu;
76  int overlap;
87  size_t gauge_offset;
88  size_t mom_offset;
89  size_t site_size;
92 
93 
97  typedef struct QudaInvertParam_s {
98 
105  double mass;
106  double kappa;
108  double m5;
109  int Ls;
114  double mu;
115  double epsilon;
119  int laplace3D;
121  double tol;
122  double tol_restart;
123  double tol_hq;
126  double true_res;
127  double true_res_hq;
128  int maxiter;
129  double reliable_delta;
143 
148 
153 
158 
163 
166 
167  int pipeline;
171  int num_src;
173  int overlap;
176  double offset[QUDA_MAX_MULTI_SHIFT];
177 
179  double tol_offset[QUDA_MAX_MULTI_SHIFT];
180 
182  double tol_hq_offset[QUDA_MAX_MULTI_SHIFT];
183 
185  double true_res_offset[QUDA_MAX_MULTI_SHIFT];
186 
188  double iter_res_offset[QUDA_MAX_MULTI_SHIFT];
189 
191  double true_res_hq_offset[QUDA_MAX_MULTI_SHIFT];
192 
194  double residue[QUDA_MAX_MULTI_SHIFT];
195 
198 
202  double action[2];
203 
233  double clover_coeff;
234  double clover_rho;
237  double trlogA[2];
246  int sp_pad;
247  int cl_pad;
249  int iter;
250  double gflops;
251  double secs;
256  int Nsteps;
257 
260 
261  /*
262  * The following parameters are related to the solver
263  * preconditioner, if enabled.
264  */
265 
271 
274 
277 
279  void *eig_param;
280 
287 
290 
293 
295  double omega;
296 
299 
302 
305 
308 
311 
321 
328  int nev;
334  int rhs_idx;
338  double eigenval_tol;
344  double inc_tol;
345 
348 
351 
354 
357 
360 
363 
366 
369 
372 
373  } QudaInvertParam;
374 
375  // Parameter set for solving eigenvalue problems.
376  typedef struct QudaEigParam_s {
377 
378  // EIGENSOLVER PARAMS
379  //-------------------------------------------------
382 
385 
388 
390  int poly_deg;
391 
393  double a_min;
394  double a_max;
395 
403 
406 
409 
412 
414  int nEv;
416  int nKr;
420  int nConv;
422  double tol;
427 
431  char arpack_logfile[512];
432 
434  char QUDA_logfile[512];
435 
436  //-------------------------------------------------
437 
438  // EIG-CG PARAMS
439  //-------------------------------------------------
440  int nk;
441  int np;
442 
445 
448 
451 
454 
457 
459  char vec_infile[256];
460 
462  char vec_outfile[256];
463 
465  double gflops;
466 
468  double secs;
469 
472  //-------------------------------------------------
473 
474  } QudaEigParam;
475 
476  typedef struct QudaMultigridParam_s {
477 
479 
481 
483  int n_level;
484 
487 
489  int spin_block_size[QUDA_MAX_MG_LEVEL];
490 
492  int n_vec[QUDA_MAX_MG_LEVEL];
493 
496 
499 
502 
505 
508 
511 
514 
517 
520 
523 
526 
529 
532 
535 
538 
541 
544 
547 
550 
553 
556 
559 
562 
565 
568 
571 
574 
576  QudaPrecision smoother_halo_precision[QUDA_MAX_MG_LEVEL];
577 
579  QudaSchwarzType smoother_schwarz_type[QUDA_MAX_MG_LEVEL];
580 
582  int smoother_schwarz_cycle[QUDA_MAX_MG_LEVEL];
583 
586  QudaSolutionType coarse_grid_solution_type[QUDA_MAX_MG_LEVEL];
587 
590 
593 
595  QudaBoolean global_reduction[QUDA_MAX_MG_LEVEL];
596 
599 
602 
605 
610 
613 
616 
619 
622 
625 
628 
630  char vec_infile[QUDA_MAX_MG_LEVEL][256];
631 
634 
636  char vec_outfile[QUDA_MAX_MG_LEVEL][256];
637 
640 
642  double gflops;
643 
645  double secs;
646 
649 
651 
652 
653 
654  /*
655  * Interface functions, found in interface_quda.cpp
656  */
657 
683  void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[],
684  FILE *outfile);
685 
696  typedef int (*QudaCommsMap)(const int *coords, void *fdata);
697 
702  void qudaSetCommHandle(void *mycomm);
703 
731  void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata);
732 
743  void initQudaDevice(int device);
744 
751  void initQudaMemory();
752 
762  void initQuda(int device);
763 
767  void endQuda(void);
768 
774  void updateR();
775 
784 
793 
802 
811 
817 
823 
829 
835 
841  void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param);
842 
846  void freeGaugeQuda(void);
847 
853  void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param);
854 
862  void loadCloverQuda(void *h_clover, void *h_clovinv,
864 
868  void freeCloverQuda(void);
869 
879  void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta,
880  QudaEigParam *eig_param);
881 
890  void eigensolveQuda(void **h_evecs, double_complex *h_evals, QudaEigParam *param);
891 
901  void invertQuda(void *h_x, void *h_b, QudaInvertParam *param);
902 
912  void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param);
913 
921  void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param);
922 
931 
938  void destroyMultigridQuda(void *mg_instance);
939 
946  void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param);
947 
955  void dumpMultigridQuda(void *mg_instance, QudaMultigridParam *param);
956 
965  void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
967 
977  void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param,
979 
989  void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param,
991 
1001  void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param,
1002  QudaParity *parity, int inverse);
1003 
1011  void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
1012 
1020  void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param);
1021 
1022 
1023  /*
1024  * The following routines are temporary additions used by the HISQ
1025  * link-fattening code.
1026  */
1027 
1028  void set_dim(int *);
1029  void pack_ghost(void **cpuLink, void **cpuGhost, int nFace,
1030  QudaPrecision precision);
1031 
1032  void computeKSLinkQuda(void* fatlink, void* longlink, void* ulink, void* inlink,
1033  double *path_coeff, QudaGaugeParam *param);
1034 
1035 
1036 
1050  int computeGaugeForceQuda(void* mom, void* sitelink, int*** input_path_buf, int* path_length,
1051  double* loop_coeff, int num_paths, int max_length, double dt,
1053 
1065  void updateGaugeFieldQuda(void* gauge, void* momentum, double dt,
1066  int conj_mom, int exact, QudaGaugeParam* param);
1067 
1077  void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param);
1078 
1087  void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param);
1088 
1097  double momActionQuda(void* momentum, QudaGaugeParam* param);
1098 
1107  void* createGaugeFieldQuda(void* gauge, int geometry, QudaGaugeParam* param);
1108 
1116  void saveGaugeFieldQuda(void* outGauge, void* inGauge, QudaGaugeParam* param);
1117 
1123  void destroyGaugeFieldQuda(void* gauge);
1124 
1131 
1150  void computeCloverForceQuda(void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck,
1151  int nvector, double multiplicity, void *gauge,
1153 
1165  void computeStaggeredForceQuda(void* mom, double dt, double delta, void **x, void *gauge,
1166  QudaGaugeParam *gauge_param, QudaInvertParam *invert_param);
1167 
1183  void computeHISQForceQuda(void* momentum,
1184  double dt,
1185  const double level2_coeff[6],
1186  const double fat7_coeff[6],
1187  const void* const w_link,
1188  const void* const v_link,
1189  const void* const u_link,
1190  void** quark,
1191  int num,
1192  int num_naik,
1193  double** coeff,
1195 
1207  void gaussGaugeQuda(unsigned long long seed, double sigma);
1208 
1213  void plaqQuda(double plaq[3]);
1214 
1215  /*
1216  * Performs a deep copy from the internal extendedGaugeResident field.
1217  * @param Pointer to externalGaugeResident cudaGaugeField
1218  * @param Location of gauge field
1219  */
1220  void copyExtendedResidentGaugeQuda(void* resident_gauge, QudaFieldLocation loc);
1221 
1232  void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param, unsigned int nSteps, double alpha);
1233 
1239  void performAPEnStep(unsigned int nSteps, double alpha);
1240 
1246  void performSTOUTnStep(unsigned int nSteps, double rho);
1247 
1254  void performOvrImpSTOUTnStep(unsigned int nSteps, double rho, double epsilon);
1255 
1259  double qChargeQuda();
1260 
1270  void contractQuda(const void *x, const void *y, void *result, const QudaContractType cType, QudaInvertParam *param,
1271  const int *X);
1272 
1278  double qChargeDensityQuda(void *qDensity);
1279 
1293  int computeGaugeFixingOVRQuda(void* gauge,
1294  const unsigned int gauge_dir,
1295  const unsigned int Nsteps,
1296  const unsigned int verbose_interval,
1297  const double relax_boost,
1298  const double tolerance,
1299  const unsigned int reunit_interval,
1300  const unsigned int stopWtheta,
1302  double* timeinfo);
1316  int computeGaugeFixingFFTQuda(void* gauge,
1317  const unsigned int gauge_dir,
1318  const unsigned int Nsteps,
1319  const unsigned int verbose_interval,
1320  const double alpha,
1321  const unsigned int autotune,
1322  const double tolerance,
1323  const unsigned int stopWtheta,
1325  double* timeinfo);
1326 
1331  void flushChronoQuda(int index);
1332 
1333 
1338  void openMagma();
1339 
1340  void closeMagma();
1341 
1348 
1352  void destroyDeflationQuda(void *df_instance);
1353 
1354  void setMPICommHandleQuda(void *mycomm);
1355 
1356 #ifdef __cplusplus
1357 }
1358 #endif
1359 
1360 // remove NVRTC WAR
1361 #undef double_complex
1362 
1363 /* #include <quda_new_interface.h> */
1364 
1365 #endif /* _QUDA_H */
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)
int maxiter_precondition
Definition: quda.h:292
static QudaGaugeParam qudaGaugeParam
double secs
Definition: quda.h:251
void performOvrImpSTOUTnStep(unsigned int nSteps, double rho, double epsilon)
enum QudaPreserveSource_s QudaPreserveSource
QudaDiracFieldOrder dirac_order
Definition: quda.h:219
QudaMassNormalization mass_normalization
Definition: quda.h:208
enum QudaMassNormalization_s QudaMassNormalization
QudaReconstructType reconstruct_sloppy
Definition: quda.h:53
double anisotropy
Definition: quda.h:38
void freeCloverQuda(void)
QudaBoolean use_dagger
Definition: quda.h:401
void invertMultiShiftQuda(void **_hp_x, void *_hp_b, QudaInvertParam *param)
void setVerbosityQuda(QudaVerbosity verbosity, const char prefix[], FILE *outfile)
void endQuda(void)
int max_hq_res_increase
Definition: quda.h:157
QudaFieldLocation clover_location
Definition: quda.h:223
QudaCABasis ca_basis
Definition: quda.h:298
QudaSolveType solve_type
Definition: quda.h:205
QudaVerbosity verbosity_precondition
Definition: quda.h:286
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:63
void destroyDeflationQuda(void *df_instance)
QudaMultigridParam newQudaMultigridParam(void)
int make_resident_mom
Definition: quda.h:83
QudaExtLibType extlib_type
Definition: quda.h:371
size_t gauge_offset
Definition: quda.h:87
double mu
Definition: quda.h:114
void setMPICommHandleQuda(void *mycomm)
QudaGaugeFixed gauge_fix
Definition: quda.h:61
QudaTune tune
Definition: quda.h:253
QudaSchwarzType schwarz_type
Definition: quda.h:310
enum QudaResidualType_s QudaResidualType
#define QUDA_MAX_MULTI_SHIFT
Maximum number of shifts supported by the multi-shift solver. This number may be changed if need be...
QudaInverterType inv_type_precondition
Definition: quda.h:270
void printQudaGaugeParam(QudaGaugeParam *param)
Definition: check_params.h:40
QudaLinkType type
Definition: quda.h:42
int max_res_increase_total
Definition: quda.h:152
double kappa
Definition: quda.h:106
void MatDagMatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
QudaPrecision cuda_prec_ritz
Definition: quda.h:324
int setup_maxiter_refresh[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1682
void invertQuda(void *h_x, void *h_b, QudaInvertParam *param)
double tol
Definition: quda.h:121
QudaDslashType dslash_type
Definition: quda.h:102
QudaReconstructType reconstruct_precondition
Definition: quda.h:59
QudaInverterType inv_type
Definition: quda.h:103
QudaPrecision cuda_prec
Definition: quda.h:214
int return_clover_inverse
Definition: quda.h:242
QudaCABasis coarse_solver_ca_basis[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1697
QudaSetupType setup_type
Definition: quda.h:531
void performAPEnStep(unsigned int nSteps, double alpha)
enum QudaSolveType_s QudaSolveType
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
double epsilon
Definition: test_util.cpp:1649
QudaPrecision cpu_prec
Definition: quda.h:213
double momActionQuda(void *momentum, QudaGaugeParam *param)
QudaStaggeredPhase staggered_phase_type
Definition: quda.h:71
QudaBoolean pre_orthonormalize
Definition: quda.h:534
void destroyMultigridQuda(void *mg_instance)
Free resources allocated by the multigrid solver.
void plaqQuda(double plaq[3])
QudaDagType dagger
Definition: quda.h:207
QudaGaugeParam gauge_param
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:216
enum QudaComputeNullVector_s QudaComputeNullVector
int chrono_replace_last
Definition: quda.h:356
QudaPrecision clover_cuda_prec_refinement_sloppy
Definition: quda.h:227
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
double true_res
Definition: quda.h:126
void computeKSLinkQuda(void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param)
size_t mom_offset
Definition: quda.h:88
void * eig_param
Definition: quda.h:279
void loadCloverQuda(void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
int return_clover
Definition: quda.h:241
enum QudaEigType_s QudaEigType
enum QudaTboundary_s QudaTboundary
int make_resident_solution
Definition: quda.h:347
void openMagma()
double ca_lambda_max
Definition: quda.h:304
QudaBoolean setup_minimize_memory
Definition: quda.h:609
int overwrite_mom
Definition: quda.h:78
QudaBoolean coarse_guess
Definition: quda.h:639
int compute_clover_trlog
Definition: quda.h:236
int chrono_index
Definition: quda.h:365
QudaPrecision clover_cuda_prec_sloppy
Definition: quda.h:226
struct QudaInvertParam_s QudaInvertParam
int compute_action
Definition: quda.h:197
void dslashQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
QudaPrecision chrono_precision
Definition: quda.h:368
QudaFieldLocation input_location
Definition: quda.h:99
void freeGaugeQuda(void)
void initCommsGridQuda(int nDim, const int *dims, QudaCommsMap func, void *fdata)
QudaBoolean use_poly_acc
Definition: quda.h:387
int staggered_phase_applied
Definition: quda.h:72
double reliable_delta
Definition: quda.h:129
size_t site_size
Definition: quda.h:89
int solution_accumulator_pipeline
Definition: quda.h:142
QudaUseInitGuess use_init_guess
Definition: quda.h:231
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.
double ca_lambda_min
Definition: quda.h:301
QudaGaugeParam param
Definition: pack_test.cpp:17
void performWuppertalnStep(void *h_out, void *h_in, QudaInvertParam *param, unsigned int nSteps, double alpha)
int llfat_ga_pad
Definition: quda.h:68
int coarse_solver_maxiter[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1696
QudaSolutionType solution_type
Definition: quda.h:204
void qudaSetCommHandle(void *mycomm)
int nConv
Definition: quda.h:420
int use_alternative_reliable
Definition: quda.h:131
void projectSU3Quda(void *gauge_h, double tol, QudaGaugeParam *param)
QudaMemoryType mem_type_ritz
Definition: quda.h:450
QudaSolverNormalization solver_normalization
Definition: quda.h:209
int geo_block_size[QUDA_MAX_MG_LEVEL][QUDA_MAX_DIM]
Definition: test_util.cpp:1706
QudaEigType eig_type
Definition: quda.h:384
void dslashQuda_4dpc(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaPrecision clover_cuda_prec
Definition: quda.h:225
int precondition_cycle
Definition: quda.h:307
void * preconditioner
Definition: quda.h:273
int num_setup_iter[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1679
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
int chrono_use_resident
Definition: quda.h:359
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)
QudaInvertParam * invert_param
Definition: quda.h:381
double scale
Definition: quda.h:40
QudaBoolean require_convergence
Definition: quda.h:408
void initQuda(int device)
void initQudaMemory()
double secs
Definition: quda.h:468
void gaussGaugeQuda(unsigned long long seed, double sigma)
Generate Gaussian distributed fields and store in the resident gauge field. We create a Gaussian-dist...
double tol
Definition: test_util.cpp:1656
QudaFieldLocation output_location
Definition: quda.h:100
QudaPrecision clover_cuda_prec_precondition
Definition: quda.h:228
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
QudaBoolean post_orthonormalize
Definition: quda.h:537
int site_ga_pad
Definition: quda.h:65
void saveGaugeQuda(void *h_gauge, QudaGaugeParam *param)
QudaBoolean run_verify
Definition: quda.h:456
double i_mu
Definition: quda.h:74
int coarse_solver_ca_basis_size[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1698
double m5
Definition: quda.h:108
enum QudaStaggeredPhase_s QudaStaggeredPhase
void createCloverQuda(QudaInvertParam *param)
void * newDeflationQuda(QudaEigParam *param)
QudaSolveType smoother_solve_type[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1678
QudaPrecision cuda_prec_sloppy
Definition: quda.h:215
QudaVerbosity verbosity
Definition: quda.h:244
void dumpMultigridQuda(void *mg_instance, QudaMultigridParam *param)
Dump the null-space vectors to disk.
QudaInvertParam newQudaInvertParam(void)
double gflops
Definition: quda.h:250
int eigcg_max_restarts
Definition: quda.h:340
double gflops
Definition: quda.h:465
QudaPrecision cuda_prec_precondition
Definition: quda.h:58
double coarse_solver_ca_lambda_max[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1700
int staple_pad
Definition: quda.h:67
QudaCloverFieldOrder clover_order
Definition: quda.h:230
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 saveGaugeFieldQuda(void *outGauge, void *inGauge, QudaGaugeParam *param)
double tol_hq
Definition: quda.h:123
enum QudaMatPCType_s QudaMatPCType
double setup_ca_lambda_min[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1685
struct QudaEigParam_s QudaEigParam
void * newMultigridQuda(QudaMultigridParam *param)
double omega
Definition: test_util.cpp:1690
double true_res_hq
Definition: quda.h:127
int poly_deg
Definition: quda.h:390
enum QudaGaugeFixed_s QudaGaugeFixed
enum QudaSolutionType_s QudaSolutionType
void eigensolveQuda(void **h_evecs, double_complex *h_evals, QudaEigParam *param)
QudaGammaBasis gamma_basis
Definition: quda.h:221
QudaBoolean use_norm_op
Definition: quda.h:402
double a_min
Definition: quda.h:393
double setup_ca_lambda_max[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1686
enum QudaSchwarzType_s QudaSchwarzType
QudaBoolean compute_svd
Definition: quda.h:405
QudaPrecision cuda_prec_sloppy
Definition: quda.h:52
int max_search_dim
Definition: quda.h:332
int chrono_make_resident
Definition: quda.h:353
QudaBoolean arpack_check
Definition: quda.h:429
double tol_precondition
Definition: quda.h:289
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
int use_sloppy_partial_accumulator
Definition: quda.h:132
enum QudaDagType_s QudaDagType
QudaBoolean run_verify
Definition: quda.h:618
int heavy_quark_check
Definition: quda.h:165
enum QudaParity_s QudaParity
QudaFieldLocation setup_location[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1669
QudaReconstructType reconstruct
Definition: quda.h:50
enum QudaLinkType_s QudaLinkType
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
void updateMultigridQuda(void *mg_instance, QudaMultigridParam *param)
Updates the multigrid preconditioner for the new gauge / clover field.
double mass
Definition: quda.h:105
QudaBoolean import_vectors
Definition: quda.h:444
void contractQuda(const void *x, const void *y, void *result, const QudaContractType cType, QudaInvertParam *param, const int *X)
QudaInverterType coarse_solver[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1691
QudaFieldLocation location
Definition: quda.h:453
int gcrNkrylov
Definition: quda.h:259
void lanczosQuda(int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta, QudaEigParam *eig_param)
double smoother_tol[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1695
int max_hq_res_restart_total
Definition: quda.h:162
QudaPrecision cuda_prec_refinement_sloppy
Definition: quda.h:55
int compute_clover_inverse
Definition: quda.h:240
static int dims[4]
Definition: face_gauge.cpp:41
void staggeredPhaseQuda(void *gauge_h, QudaGaugeParam *param)
QudaBoolean run_oblique_proj_check
Definition: quda.h:624
enum QudaBoolean_s QudaBoolean
double coarse_solver_tol[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1692
void performSTOUTnStep(unsigned int nSteps, double rho)
static int index(int ndim, const int *dims, const int *x)
Definition: comm_common.cpp:32
void printQudaInvertParam(QudaInvertParam *param)
Definition: check_params.h:277
QudaBoolean run_low_mode_check
Definition: quda.h:621
int max_res_increase
Definition: quda.h:147
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 destroyGaugeFieldQuda(void *gauge)
void copyExtendedResidentGaugeQuda(void *resident_gauge, QudaFieldLocation loc)
enum QudaFieldLocation_s QudaFieldLocation
void closeMagma()
QudaFieldLocation location
Definition: quda.h:34
int nu_pre[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1671
struct QudaGaugeParam_s QudaGaugeParam
double clover_rho
Definition: quda.h:234
__device__ __host__ Matrix< T, 3 > inverse(const Matrix< T, 3 > &u)
Definition: quda_matrix.h:611
double tadpole_coeff
Definition: quda.h:39
double tol
Definition: quda.h:422
QudaPrecision cuda_prec_precondition
Definition: quda.h:217
int deflation_grid
Definition: quda.h:336
enum QudaEigSpectrumType_s QudaEigSpectrumType
void printQudaEigParam(QudaEigParam *param)
Definition: check_params.h:143
double tol_restart
Definition: quda.h:122
int nLockedMax
Definition: quda.h:418
void updateGaugeFieldQuda(void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
enum QudaGammaBasis_s QudaGammaBasis
enum QudaMultigridCycleType_s QudaMultigridCycleType
QudaBoolean generate_all_levels
Definition: quda.h:615
enum QudaReconstructType_s QudaReconstructType
void invertMultiSrcQuda(void **_hp_x, void **_hp_b, QudaInvertParam *param)
enum QudaCABasis_s QudaCABasis
#define QUDA_MAX_MG_LEVEL
Maximum number of multi-grid levels. This number may be increased if needed.
enum QudaSetupType_s QudaSetupType
struct QudaMultigridParam_s QudaMultigridParam
#define QUDA_MAX_DWF_LS
Maximum length of the Ls dimension for domain-wall fermions.
int setup_ca_basis_size[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1684
double setup_tol[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1680
#define double_complex
Definition: quda.h:19
void printQudaMultigridParam(QudaMultigridParam *param)
Definition: check_params.h:598
int mom_ga_pad
Definition: quda.h:69
QudaExtLibType extlib_type
Definition: quda.h:471
void cloverQuda(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity *parity, int inverse)
void MatQuda(void *h_out, void *h_in, QudaInvertParam *inv_param)
int use_resident_gauge
Definition: quda.h:80
QudaTboundary t_boundary
Definition: quda.h:45
enum QudaTune_s QudaTune
QudaTwistFlavorType twist_flavor
Definition: quda.h:117
int return_result_gauge
Definition: quda.h:84
int chrono_max_dim
Definition: quda.h:362
int max_restart_num
Definition: quda.h:342
double qChargeQuda()
double mu_factor[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1674
void * createGaugeFieldQuda(void *gauge, int geometry, QudaGaugeParam *param)
int use_resident_mom
Definition: quda.h:81
enum QudaDslashType_s QudaDslashType
int(* QudaCommsMap)(const int *coords, void *fdata)
Definition: quda.h:696
int check_interval
Definition: quda.h:424
int max_restarts
Definition: quda.h:426
enum QudaContractType_s QudaContractType
int device
Definition: test_util.cpp:1602
QudaEigSpectrumType spectrum
Definition: quda.h:411
int setup_maxiter[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1681
void * longlink
int compute_true_res
Definition: quda.h:125
QudaResidualType residual_type
Definition: quda.h:320
double inc_tol
Definition: quda.h:344
int num_offset
Definition: quda.h:169
enum QudaVerbosity_s QudaVerbosity
void * fatlink
void computeStaggeredForceQuda(void *mom, double dt, double delta, void **x, void *gauge, QudaGaugeParam *gauge_param, QudaInvertParam *invert_param)
cpuGaugeField * cpuLink
Definition: covdev_test.cpp:39
int return_result_mom
Definition: quda.h:85
int test_type
Definition: test_util.cpp:1636
int compute_clover
Definition: quda.h:239
double epsilon
Definition: quda.h:115
#define QUDA_MAX_DIM
Maximum number of dimensions supported by QUDA. In practice, no routines make use of more than 5...
void flushChronoQuda(int index)
Flush the chronological history for the given index.
void set_dim(int *)
int use_resident_solution
Definition: quda.h:350
double omega
Definition: quda.h:295
QudaComputeNullVector compute_null_vector
Definition: quda.h:612
double coarse_solver_ca_lambda_min[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1699
int make_resident_gauge
Definition: quda.h:82
QudaCABasis setup_ca_basis[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1683
int n_block_ortho[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1673
double qChargeDensityQuda(void *qDensity)
Calculates the topological charge from gaugeSmeared, if it exist, or from gaugePrecise if no smeared ...
void * deflation_op
Definition: quda.h:276
double eigenval_tol
Definition: quda.h:338
double a_max
Definition: quda.h:394
int nu_post[QUDA_MAX_MG_LEVEL]
Definition: test_util.cpp:1672
QudaDslashType dslash_type_precondition
Definition: quda.h:284
QudaPrecision clover_cpu_prec
Definition: quda.h:224
QudaPrecision cuda_prec_ritz
Definition: quda.h:447
QudaParity parity
Definition: covdev_test.cpp:54
enum QudaUseInitGuess_s QudaUseInitGuess
enum QudaSolverNormalization_s QudaSolverNormalization
void initQudaDevice(int device)
void pack_ghost(void **cpuLink, void **cpuGhost, int nFace, QudaPrecision precision)
QudaMatPCType matpc_type
Definition: quda.h:206
void dslashQuda_mdwf(void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
QudaEigParam newQudaEigParam(void)
enum QudaInverterType_s QudaInverterType
enum QudaMemoryType_s QudaMemoryType
QudaReconstructType reconstruct_refinement_sloppy
Definition: quda.h:56
QudaPrecision cpu_prec
Definition: quda.h:47
enum QudaExtLibType_s QudaExtLibType
int overlap
Definition: quda.h:76
void updateR()
update the radius for halos.
QudaGaugeParam newQudaGaugeParam(void)
QudaPreserveSource preserve_source
Definition: quda.h:211
QudaInvertParam * invert_param
Definition: quda.h:478
double reliable_delta_refinement
Definition: quda.h:130
double clover_coeff
Definition: quda.h:233
enum QudaTwistFlavorType_s QudaTwistFlavorType
QudaVerbosity verbosity
Definition: test_util.cpp:1614
enum QudaDiracFieldOrder_s QudaDiracFieldOrder