QUDA  v1.1.0
A library for QCD on GPUs
quda_milc_interface.h
Go to the documentation of this file.
1 #ifndef _QUDA_MILC_INTERFACE_H
2 #define _QUDA_MILC_INTERFACE_H
3 
4 #include <enum_quda.h>
5 #include <quda.h>
6 
15 #if __COMPUTE_CAPABILITY__ >= 600
16 #define USE_QUDA_MANAGED 1
17 #endif
18 
19 #ifdef __cplusplus
20 extern "C" {
21 #endif
22 
26  typedef struct {
27  void *site;
28  void *link;
29  size_t link_offset;
30  void *mom;
31  size_t mom_offset;
32  size_t size;
34 
38  typedef struct {
39  int max_iter;
42  double boundary_phase[4];
46  double tadpole;
49  double naik_epsilon;
51 
56  typedef struct {
58  int nev;
61  double tol_restart;
62 
65  double inc_tol;
66  double eigenval_tol;
67 
70 
73 
74  char *vec_infile;
75  char *vec_outfile;
76 
77  } QudaEigArgs_t;
78 
79 
83  typedef struct {
84  const int* latsize;
85  const int* machsize;
86  int device;
87  } QudaLayout_t;
88 
89 
93  typedef struct {
96  } QudaInitArgs_t; // passed to the initialization struct
97 
98 
102  typedef struct {
107  double force_filter;
109 
110 
114  typedef struct {
118 
124  void qudaSetMPICommHandle(void *mycomm);
125 
131  void qudaInit(QudaInitArgs_t input);
132 
138  void qudaSetLayout(QudaLayout_t layout);
139 
143  void qudaFinalize();
144 
150  void* qudaAllocatePinned(size_t bytes);
151 
156  void qudaFreePinned(void *ptr);
157 
163  void *qudaAllocateManaged(size_t bytes);
164 
169  void qudaFreeManaged(void *ptr);
170 
177  void qudaHisqParamsInit(QudaHisqParams_t hisq_params);
178 
191  void qudaLoadKSLink(int precision,
192  QudaFatLinkArgs_t fatlink_args,
193  const double act_path_coeff[6],
194  void* inlink,
195  void* fatlink,
196  void* longlink);
197 
210  void qudaLoadUnitarizedLink(int precision,
211  QudaFatLinkArgs_t fatlink_args,
212  const double path_coeff[6],
213  void* inlink,
214  void* fatlink,
215  void* ulink);
216 
217 
230  void qudaDslash(int external_precision,
231  int quda_precision,
232  QudaInvertArgs_t inv_args,
233  const void* const milc_fatlink,
234  const void* const milc_longlink,
235  void* source,
236  void* solution,
237  int* num_iters);
238 
261  void qudaDDInvert(int external_precision,
262  int quda_precision,
263  double mass,
264  QudaInvertArgs_t inv_args,
265  double target_residual,
266  double target_fermilab_residual,
267  const int * const domain_overlap,
268  const void* const fatlink,
269  const void* const longlink,
270  void* source,
271  void* solution,
272  double* const final_residual,
273  double* const final_fermilab_residual,
274  int* num_iters);
275 
276 
277 
298  void qudaInvert(int external_precision,
299  int quda_precision,
300  double mass,
301  QudaInvertArgs_t inv_args,
302  double target_residual,
303  double target_fermilab_residual,
304  const void* const milc_fatlink,
305  const void* const milc_longlink,
306  void* source,
307  void* solution,
308  double* const final_resid,
309  double* const final_rel_resid,
310  int* num_iters);
311 
327  void *qudaMultigridCreate(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args,
328  const void *const milc_fatlink, const void *const milc_longlink,
329  const char *const mg_param_file);
330 
355  void qudaInvertMG(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args,
356  double target_residual, double target_fermilab_residual, const void *const milc_fatlink,
357  const void *const milc_longlink, void *mg_pack_ptr, int mg_rebuild_type, void *source,
358  void *solution, double *const final_residual, double *const final_fermilab_residual, int *num_iters);
359 
366  void qudaMultigridDestroy(void *mg_pack_ptr);
367 
389  void qudaInvertMsrc(int external_precision,
390  int quda_precision,
391  double mass,
392  QudaInvertArgs_t inv_args,
393  double target_residual,
394  double target_fermilab_residual,
395  const void* const fatlink,
396  const void* const longlink,
397  void** sourceArray,
398  void** solutionArray,
399  double* const final_residual,
400  double* const final_fermilab_residual,
401  int* num_iters,
402  int num_src);
403 
429  int external_precision,
430  int precision,
431  int num_offsets,
432  double* const offset,
433  QudaInvertArgs_t inv_args,
434  const double* target_residual,
435  const double* target_fermilab_residual,
436  const void* const milc_fatlink,
437  const void* const milc_longlink,
438  void* source,
439  void** solutionArray,
440  double* const final_residual,
441  double* const final_fermilab_residual,
442  int* num_iters);
443 
472  void qudaEigCGInvert(
473  int external_precision,
474  int quda_precision,
475  double mass,
476  QudaInvertArgs_t inv_args,
477  double target_residual,
478  double target_fermilab_residual,
479  const void* const fatlink,
480  const void* const longlink,
481  void* source,
482  void* solution,
483  QudaEigArgs_t eig_args,
484  const int rhs_idx,//current rhs
485  const int last_rhs_flag,//is this the last rhs to solve?
486  double* const final_residual,
487  double* const final_fermilab_residual,
488  int *num_iters);
489 
512  void qudaCloverInvert(int external_precision,
513  int quda_precision,
514  double kappa,
515  double clover_coeff,
516  QudaInvertArgs_t inv_args,
517  double target_residual,
518  double target_fermilab_residual,
519  const void* milc_link,
520  void* milc_clover,
521  void* milc_clover_inv,
522  void* source,
523  void* solution,
524  double* const final_residual,
525  double* const final_fermilab_residual,
526  int* num_iters);
527 
557  int external_precision,
558  int quda_precision,
559  double kappa,
560  double clover_coeff,
561  QudaInvertArgs_t inv_args,
562  double target_residual,
563  double target_fermilab_residual,
564  const void* milc_link,
565  void* milc_clover,
566  void* milc_clover_inv,
567  void* source,
568  void* solution,
569  QudaEigArgs_t eig_args,
570  const int rhs_idx,//current rhs
571  const int last_rhs_flag,//is this the last rhs to solve?
572  double* const final_residual,
573  double* const final_fermilab_residual,
574  int *num_iters);
575 
576 
585  void qudaLoadGaugeField(int external_precision,
586  int quda_precision,
587  QudaInvertArgs_t inv_args,
588  const void* milc_link) ;
589 
593  void qudaFreeGaugeField();
594 
613  void qudaLoadCloverField(int external_precision,
614  int quda_precision,
615  QudaInvertArgs_t inv_args,
616  void* milc_clover,
617  void* milc_clover_inv,
620  double clover_coeff,
621  int compute_trlog,
622  double *trlog) ;
623 
627  void qudaFreeCloverField();
628 
655  void qudaCloverMultishiftInvert(int external_precision,
656  int quda_precision,
657  int num_offsets,
658  double* const offset,
659  double kappa,
660  double clover_coeff,
661  QudaInvertArgs_t inv_args,
662  const double* target_residual,
663  const void* milc_link,
664  void* milc_clover,
665  void* milc_clover_inv,
666  void* source,
667  void** solutionArray,
668  double* const final_residual,
669  int* num_iters
670  );
671 
690  void qudaHisqForce(int precision,
691  int num_terms,
692  int num_naik_terms,
693  double dt,
694  double** coeff,
695  void** quark_field,
696  const double level2_coeff[6],
697  const double fat7_coeff[6],
698  const void* const w_link,
699  const void* const v_link,
700  const void* const u_link,
701  void* const milc_momentum);
702 
714  void qudaGaugeForce(int precision,
715  int num_loop_types,
716  double milc_loop_coeff[3],
717  double eb3,
719 
732  void qudaGaugeForcePhased(int precision, int num_loop_types, double milc_loop_coeff[3], double eb3,
733  QudaMILCSiteArg_t *arg, int phase_in);
734 
743  void qudaUpdateU(int precision,
744  double eps,
746 
756  void qudaUpdateUPhased(int precision, double eps, QudaMILCSiteArg_t *arg, int phase_in);
757 
768  void qudaUpdateUPhasedPipeline(int precision, double eps, QudaMILCSiteArg_t *arg, int phase_in, int want_gaugepipe);
769 
779  void qudaMomLoad(int precision, QudaMILCSiteArg_t *arg);
780 
790  void qudaMomSave(int precision, QudaMILCSiteArg_t *arg);
791 
801  double qudaMomAction(int precision, QudaMILCSiteArg_t *arg);
802 
814  void qudaRephase(int prec, void *gauge, int flag, double i_mu);
815 
824  void qudaUnitarizeSU3(int prec, double tol, QudaMILCSiteArg_t *arg);
825 
835  void qudaUnitarizeSU3Phased(int prec, double tol, QudaMILCSiteArg_t *arg, int phase_in);
836 
855  void qudaCloverForce(void *mom, double dt, void **x, void **p, double *coeff, double kappa,
856  double ck, int nvec, double multiplicity, void *gauge, int precision,
857  QudaInvertArgs_t inv_args);
858 
870  void qudaCloverTrace(void* out,
871  void* dummy,
872  int mu,
873  int nu);
874 
875 
891  void qudaCloverDerivative(void* out,
892  void* gauge,
893  void* oprod,
894  int mu,
895  int nu,
896  double coeff,
897  int precision,
898  int parity,
899  int conjugate);
900 
901 
911  void* qudaCreateExtendedGaugeField(void* gauge,
912  int geometry,
913  int precision);
914 
925  int geometry,
926  int precision);
927 
936  void* qudaCreateGaugeField(void* gauge,
937  int geometry,
938  int precision);
939 
946  void qudaSaveGaugeField(void* gauge,
947  void* inGauge);
948 
954  void qudaDestroyGaugeField(void* gauge);
955 
956 
969  void qudaGaugeFixingOVR( const int precision,
970  const unsigned int gauge_dir,
971  const int Nsteps,
972  const int verbose_interval,
973  const double relax_boost,
974  const double tolerance,
975  const unsigned int reunit_interval,
976  const unsigned int stopWtheta,
977  void* milc_sitelink
978  );
979 
980 
993  void qudaGaugeFixingFFT( int precision,
994  unsigned int gauge_dir,
995  int Nsteps,
996  int verbose_interval,
997  double alpha,
998  unsigned int autotune,
999  double tolerance,
1000  unsigned int stopWtheta,
1001  void* milc_sitelink
1002  );
1003 
1004  /* The below declarations are for removed functions from prior versions of QUDA. */
1005 
1010  void qudaAsqtadForce(int precision,
1011  const double act_path_coeff[6],
1012  const void* const one_link_src[4],
1013  const void* const naik_src[4],
1014  const void* const link,
1015  void* const milc_momentum);
1016 
1021  void qudaComputeOprod(int precision,
1022  int num_terms,
1023  int num_naik_terms,
1024  double** coeff,
1025  double scale,
1026  void** quark_field,
1027  void* oprod[3]);
1028 
1029 #ifdef __cplusplus
1030 }
1031 #endif
1032 
1033 
1034 #endif // _QUDA_MILC_INTERFACE_H
double kappa
QudaSolveType solve_type
double mass
double tol
quda::mgarray< int > nvec
double mu
QudaSolutionType solution_type
double clover_coeff
QudaPrecision prec
QudaParity parity
Definition: covdev_test.cpp:40
enum QudaSolveType_s QudaSolveType
enum QudaPrecision_s QudaPrecision
enum QudaSolutionType_s QudaSolutionType
enum QudaInverterType_s QudaInverterType
enum QudaFieldLocation_s QudaFieldLocation
enum QudaExtLibType_s QudaExtLibType
enum QudaMemoryType_s QudaMemoryType
enum QudaVerbosity_s QudaVerbosity
enum QudaParity_s QudaParity
unsigned long long bytes
__host__ __device__ ValueType arg(const complex< ValueType > &z)
Returns the phase angle of z.
Main header file for the QUDA library.
void qudaLoadGaugeField(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, const void *milc_link)
void qudaDDInvert(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const int *const domain_overlap, const void *const fatlink, const void *const longlink, void *source, void *solution, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaCloverInvert(int external_precision, int quda_precision, double kappa, double clover_coeff, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *milc_link, void *milc_clover, void *milc_clover_inv, void *source, void *solution, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaGaugeFixingOVR(const int precision, const unsigned int gauge_dir, const int Nsteps, const int verbose_interval, const double relax_boost, const double tolerance, const unsigned int reunit_interval, const unsigned int stopWtheta, void *milc_sitelink)
Gauge fixing with overrelaxation with support for single and multi GPU.
void qudaCloverTrace(void *out, void *dummy, int mu, int nu)
void qudaComputeOprod(int precision, int num_terms, int num_naik_terms, double **coeff, double scale, void **quark_field, void *oprod[3])
void qudaUnitarizeSU3Phased(int prec, double tol, QudaMILCSiteArg_t *arg, int phase_in)
void qudaMultigridDestroy(void *mg_pack_ptr)
void qudaInit(QudaInitArgs_t input)
void qudaInvertMsrc(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *const fatlink, const void *const longlink, void **sourceArray, void **solutionArray, double *const final_residual, double *const final_fermilab_residual, int *num_iters, int num_src)
void qudaLoadCloverField(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, void *milc_clover, void *milc_clover_inv, QudaSolutionType solution_type, QudaSolveType solve_type, double clover_coeff, int compute_trlog, double *trlog)
void qudaCloverMultishiftInvert(int external_precision, int quda_precision, int num_offsets, double *const offset, double kappa, double clover_coeff, QudaInvertArgs_t inv_args, const double *target_residual, const void *milc_link, void *milc_clover, void *milc_clover_inv, void *source, void **solutionArray, double *const final_residual, int *num_iters)
void qudaUnitarizeSU3(int prec, double tol, QudaMILCSiteArg_t *arg)
void qudaGaugeForce(int precision, int num_loop_types, double milc_loop_coeff[3], double eb3, QudaMILCSiteArg_t *arg)
void qudaRephase(int prec, void *gauge, int flag, double i_mu)
void * qudaAllocateManaged(size_t bytes)
void qudaFreeManaged(void *ptr)
void qudaHisqForce(int precision, int num_terms, int num_naik_terms, double dt, double **coeff, void **quark_field, 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 *const milc_momentum)
void qudaFinalize()
void qudaInvertMG(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *const milc_fatlink, const void *const milc_longlink, void *mg_pack_ptr, int mg_rebuild_type, void *source, void *solution, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaMultishiftInvert(int external_precision, int precision, int num_offsets, double *const offset, QudaInvertArgs_t inv_args, const double *target_residual, const double *target_fermilab_residual, const void *const milc_fatlink, const void *const milc_longlink, void *source, void **solutionArray, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaSetLayout(QudaLayout_t layout)
void qudaSetMPICommHandle(void *mycomm)
void qudaFreeGaugeField()
double qudaMomAction(int precision, QudaMILCSiteArg_t *arg)
void * qudaCreateExtendedGaugeField(void *gauge, int geometry, int precision)
void qudaGaugeFixingFFT(int precision, unsigned int gauge_dir, int Nsteps, int verbose_interval, double alpha, unsigned int autotune, double tolerance, unsigned int stopWtheta, void *milc_sitelink)
Gauge fixing with Steepest descent method with FFTs with support for single GPU only.
void * qudaMultigridCreate(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, const void *const milc_fatlink, const void *const milc_longlink, const char *const mg_param_file)
void qudaEigCGCloverInvert(int external_precision, int quda_precision, double kappa, double clover_coeff, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *milc_link, void *milc_clover, void *milc_clover_inv, void *source, void *solution, QudaEigArgs_t eig_args, const int rhs_idx, const int last_rhs_flag, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaFreePinned(void *ptr)
void qudaInvert(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *const milc_fatlink, const void *const milc_longlink, void *source, void *solution, double *const final_resid, double *const final_rel_resid, int *num_iters)
void * qudaResidentExtendedGaugeField(void *gauge, int geometry, int precision)
void * qudaAllocatePinned(size_t bytes)
void qudaHisqParamsInit(QudaHisqParams_t hisq_params)
void qudaMomLoad(int precision, QudaMILCSiteArg_t *arg)
void qudaUpdateU(int precision, double eps, QudaMILCSiteArg_t *arg)
void qudaCloverDerivative(void *out, void *gauge, void *oprod, int mu, int nu, double coeff, int precision, int parity, int conjugate)
void qudaAsqtadForce(int precision, const double act_path_coeff[6], const void *const one_link_src[4], const void *const naik_src[4], const void *const link, void *const milc_momentum)
void qudaEigCGInvert(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_residual, double target_fermilab_residual, const void *const fatlink, const void *const longlink, void *source, void *solution, QudaEigArgs_t eig_args, const int rhs_idx, const int last_rhs_flag, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void qudaLoadKSLink(int precision, QudaFatLinkArgs_t fatlink_args, const double act_path_coeff[6], void *inlink, void *fatlink, void *longlink)
void qudaUpdateUPhasedPipeline(int precision, double eps, QudaMILCSiteArg_t *arg, int phase_in, int want_gaugepipe)
void qudaGaugeForcePhased(int precision, int num_loop_types, double milc_loop_coeff[3], double eb3, QudaMILCSiteArg_t *arg, int phase_in)
void qudaDslash(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, const void *const milc_fatlink, const void *const milc_longlink, void *source, void *solution, int *num_iters)
void qudaCloverForce(void *mom, double dt, void **x, void **p, double *coeff, double kappa, double ck, int nvec, double multiplicity, void *gauge, int precision, QudaInvertArgs_t inv_args)
void qudaDestroyGaugeField(void *gauge)
void qudaSaveGaugeField(void *gauge, void *inGauge)
void qudaFreeCloverField()
void qudaLoadUnitarizedLink(int precision, QudaFatLinkArgs_t fatlink_args, const double path_coeff[6], void *inlink, void *fatlink, void *ulink)
void qudaUpdateUPhased(int precision, double eps, QudaMILCSiteArg_t *arg, int phase_in)
void qudaMomSave(int precision, QudaMILCSiteArg_t *arg)
void * qudaCreateGaugeField(void *gauge, int geometry, int precision)
QudaExtLibType deflation_ext_lib
QudaMemoryType mem_type_ritz
QudaFieldLocation location_ritz
QudaPrecision prec_ritz
QudaExtLibType solver_ext_lib
QudaLayout_t layout
QudaVerbosity verbosity
QudaInverterType solver_type
const int * machsize
const int * latsize