QUDA  0.9.0
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 
16 #ifdef __cplusplus
17 extern "C" {
18 #endif
19 
23  typedef struct {
24  void *site;
25  void *link;
26  size_t link_offset;
27  void *mom;
28  size_t mom_offset;
29  size_t size;
31 
35  typedef struct {
36  int max_iter;
39  double boundary_phase[4];
44 
49  typedef struct {
51  int nev;
54  double tol_restart;
55 
58  double inc_tol;
59  double eigenval_tol;
60 
63 
66 
67  char *vec_infile;
68  char *vec_outfile;
69 
70  } QudaEigArgs_t;
71 
72 
76  typedef struct {
77  const int* latsize;
78  const int* machsize;
79  int device;
80  } QudaLayout_t;
81 
82 
86  typedef struct {
89  } QudaInitArgs_t; // passed to the initialization struct
90 
91 
95  typedef struct {
100  double force_filter;
102 
103 
107  typedef struct {
111 
117  void qudaInit(QudaInitArgs_t input);
118 
125 
129  void qudaFinalize();
130 
136  void* qudaAllocatePinned(size_t bytes);
137 
142  void qudaFreePinned(void *ptr);
143 
150  void qudaHisqParamsInit(QudaHisqParams_t hisq_params);
151 
164  void qudaLoadKSLink(int precision,
165  QudaFatLinkArgs_t fatlink_args,
166  const double act_path_coeff[6],
167  void* inlink,
168  void* fatlink,
169  void* longlink);
170 
183  void qudaLoadUnitarizedLink(int precision,
184  QudaFatLinkArgs_t fatlink_args,
185  const double path_coeff[6],
186  void* inlink,
187  void* fatlink,
188  void* ulink);
189 
190 
204  void qudaDslash(int external_precision,
205  int quda_precision,
206  QudaInvertArgs_t inv_args,
207  const void* const milc_fatlink,
208  const void* const milc_longlink,
209  const double tadpole,
210  void* source,
211  void* solution,
212  int* num_iters);
213 
236  void qudaDDInvert(int external_precision,
237  int quda_precision,
238  double mass,
239  QudaInvertArgs_t inv_args,
240  double target_residual,
241  double target_fermilab_residual,
242  const int * const domain_overlap,
243  const void* const fatlink,
244  const void* const longlink,
245  void* source,
246  void* solution,
247  double* const final_residual,
248  double* const final_fermilab_residual,
249  int* num_iters);
250 
251 
252 
274  void qudaInvert(int external_precision,
275  int quda_precision,
276  double mass,
277  QudaInvertArgs_t inv_args,
278  double target_residual,
279  double target_fermilab_residual,
280  const void* const milc_fatlink,
281  const void* const milc_longlink,
282  const double tadpole,
283  void* source,
284  void* solution,
285  double* const final_resid,
286  double* const final_rel_resid,
287  int* num_iters);
288 
311  void qudaInvertMsrc(int external_precision,
312  int quda_precision,
313  double mass,
314  QudaInvertArgs_t inv_args,
315  double target_residual,
316  double target_fermilab_residual,
317  const void* const fatlink,
318  const void* const longlink,
319  const double tadpole,
320  void** sourceArray,
321  void** solutionArray,
322  double* const final_residual,
323  double* const final_fermilab_residual,
324  int* num_iters,
325  int num_src);
326 
353  int external_precision,
354  int precision,
355  int num_offsets,
356  double* const offset,
357  QudaInvertArgs_t inv_args,
358  const double* target_residual,
359  const double* target_fermilab_residual,
360  const void* const milc_fatlink,
361  const void* const milc_longlink,
362  const double tadpole,
363  void* source,
364  void** solutionArray,
365  double* const final_residual,
366  double* const final_fermilab_residual,
367  int* num_iters);
368 
398  void qudaEigCGInvert(
399  int external_precision,
400  int quda_precision,
401  double mass,
402  QudaInvertArgs_t inv_args,
403  double target_residual,
404  double target_fermilab_residual,
405  const void* const fatlink,
406  const void* const longlink,
407  const double tadpole,
408  void* source,
409  void* solution,
410  QudaEigArgs_t eig_args,
411  const int rhs_idx,//current rhs
412  const int last_rhs_flag,//is this the last rhs to solve?
413  double* const final_residual,
414  double* const final_fermilab_residual,
415  int *num_iters);
416 
439  void qudaCloverInvert(int external_precision,
440  int quda_precision,
441  double kappa,
442  double clover_coeff,
443  QudaInvertArgs_t inv_args,
444  double target_residual,
445  double target_fermilab_residual,
446  const void* milc_link,
447  void* milc_clover,
448  void* milc_clover_inv,
449  void* source,
450  void* solution,
451  double* const final_residual,
452  double* const final_fermilab_residual,
453  int* num_iters);
454 
484  int external_precision,
485  int quda_precision,
486  double kappa,
487  double clover_coeff,
488  QudaInvertArgs_t inv_args,
489  double target_residual,
490  double target_fermilab_residual,
491  const void* milc_link,
492  void* milc_clover,
493  void* milc_clover_inv,
494  void* source,
495  void* solution,
496  QudaEigArgs_t eig_args,
497  const int rhs_idx,//current rhs
498  const int last_rhs_flag,//is this the last rhs to solve?
499  double* const final_residual,
500  double* const final_fermilab_residual,
501  int *num_iters);
502 
503 
512  void qudaLoadGaugeField(int external_precision,
513  int quda_precision,
514  QudaInvertArgs_t inv_args,
515  const void* milc_link) ;
516 
520  void qudaFreeGaugeField();
521 
540  void qudaLoadCloverField(int external_precision,
541  int quda_precision,
542  QudaInvertArgs_t inv_args,
543  void* milc_clover,
544  void* milc_clover_inv,
545  QudaSolutionType solution_type,
547  double clover_coeff,
548  int compute_trlog,
549  double *trlog) ;
550 
554  void qudaFreeCloverField();
555 
582  void qudaCloverMultishiftInvert(int external_precision,
583  int quda_precision,
584  int num_offsets,
585  double* const offset,
586  double kappa,
587  double clover_coeff,
588  QudaInvertArgs_t inv_args,
589  const double* target_residual,
590  const void* milc_link,
591  void* milc_clover,
592  void* milc_clover_inv,
593  void* source,
594  void** solutionArray,
595  double* const final_residual,
596  int* num_iters
597  );
598 
616  void qudaHisqForce(int precision,
617  int num_terms,
618  int num_naik_terms,
619  double** coeff,
620  void** quark_field,
621  const double level2_coeff[6],
622  const double fat7_coeff[6],
623  const void* const w_link,
624  const void* const v_link,
625  const void* const u_link,
626  void* const milc_momentum);
627 
639  void qudaGaugeForce(int precision,
640  int num_loop_types,
641  double milc_loop_coeff[3],
642  double eb3,
644 
653  void qudaUpdateU(int precision,
654  double eps,
656 
667  double qudaMomAction(int precision,
668  void *momentum);
669 
681  void qudaRephase(int prec, void *gauge, int flag, double i_mu);
682 
691  void qudaUnitarizeSU3(int prec, double tol, QudaMILCSiteArg_t *arg);
692 
711  void qudaCloverForce(void *mom, double dt, void **x, void **p, double *coeff, double kappa,
712  double ck, int nvec, double multiplicity, void *gauge, int precision,
713  QudaInvertArgs_t inv_args);
714 
726  void qudaCloverTrace(void* out,
727  void* dummy,
728  int mu,
729  int nu);
730 
731 
747  void qudaCloverDerivative(void* out,
748  void* gauge,
749  void* oprod,
750  int mu,
751  int nu,
752  double coeff,
753  int precision,
754  int parity,
755  int conjugate);
756 
757 
767  void* qudaCreateExtendedGaugeField(void* gauge,
768  int geometry,
769  int precision);
770 
780  void* qudaResidentExtendedGaugeField(void* gauge,
781  int geometry,
782  int precision);
783 
792  void* qudaCreateGaugeField(void* gauge,
793  int geometry,
794  int precision);
795 
802  void qudaSaveGaugeField(void* gauge,
803  void* inGauge);
804 
810  void qudaDestroyGaugeField(void* gauge);
811 
812 
825  void qudaGaugeFixingOVR( const int precision,
826  const unsigned int gauge_dir,
827  const int Nsteps,
828  const int verbose_interval,
829  const double relax_boost,
830  const double tolerance,
831  const unsigned int reunit_interval,
832  const unsigned int stopWtheta,
833  void* milc_sitelink
834  );
835 
836 
849  void qudaGaugeFixingFFT( int precision,
850  unsigned int gauge_dir,
851  int Nsteps,
852  int verbose_interval,
853  double alpha,
854  unsigned int autotune,
855  double tolerance,
856  unsigned int stopWtheta,
857  void* milc_sitelink
858  );
859 
860  /* The below declarations are for removed functions from prior versions of QUDA. */
861 
866  void qudaAsqtadForce(int precision,
867  const double act_path_coeff[6],
868  const void* const one_link_src[4],
869  const void* const naik_src[4],
870  const void* const link,
871  void* const milc_momentum);
872 
877  void qudaComputeOprod(int precision,
878  int num_terms,
879  int num_naik_terms,
880  double** coeff,
881  double scale,
882  void** quark_field,
883  void* oprod[3]);
884 
885 
886 #ifdef __cplusplus
887 }
888 #endif
889 
890 
891 #endif // _QUDA_MILC_INTERFACE_H
void qudaHisqParamsInit(QudaHisqParams_t hisq_params)
double mu
Definition: test_util.cpp:1643
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)
enum QudaPrecision_s QudaPrecision
void qudaUnitarizeSU3(int prec, double tol, 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.
QudaExtLibType deflation_ext_lib
QudaPrecision prec_ritz
enum QudaSolveType_s QudaSolveType
void qudaInit(QudaInitArgs_t input)
QudaMemoryType mem_type_ritz
void qudaLoadGaugeField(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, const void *milc_link)
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, const double tadpole, 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 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 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)
int nvec
Definition: test_util.cpp:1635
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, const double tadpole, void **sourceArray, void **solutionArray, double *const final_residual, double *const final_fermilab_residual, int *num_iters, int num_src)
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 qudaSaveGaugeField(void *gauge, void *inGauge)
double qudaMomAction(int precision, void *momentum)
void qudaSetLayout(QudaLayout_t layout)
size_t size_t offset
void * longlink[4]
double tol
Definition: test_util.cpp:1647
void qudaFreePinned(void *ptr)
void qudaUpdateU(int precision, double eps, QudaMILCSiteArg_t *arg)
VOLATILE spinorFloat kappa
QIO_Layout layout
Definition: qio_field.cpp:6
void qudaFreeCloverField()
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, const double tadpole, void *source, void **solutionArray, double *const final_residual, double *const final_fermilab_residual, int *num_iters)
void * qudaCreateGaugeField(void *gauge, int geometry, int precision)
static __inline__ size_t p
QudaInverterType solver_type
void qudaRephase(int prec, void *gauge, int flag, double i_mu)
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, const double tadpole, void *source, void *solution, double *const final_resid, double *const final_rel_resid, int *num_iters)
enum QudaSolutionType_s QudaSolutionType
void qudaComputeOprod(int precision, int num_terms, int num_naik_terms, double **coeff, double scale, void **quark_field, void *oprod[3])
const int * machsize
void qudaFreeGaugeField()
enum QudaParity_s QudaParity
void qudaDslash(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, const void *const milc_fatlink, const void *const milc_longlink, const double tadpole, void *source, void *solution, int *num_iters)
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 qudaFinalize()
QudaFieldLocation location_ritz
const void * ptr
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)
QudaSolveType solve_type
Definition: test_util.cpp:1653
void qudaCloverDerivative(void *out, void *gauge, void *oprod, int mu, int nu, double coeff, int precision, int parity, int conjugate)
void qudaLoadUnitarizedLink(int precision, QudaFatLinkArgs_t fatlink_args, const double path_coeff[6], void *inlink, void *fatlink, void *ulink)
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
void * fatlink[4]
double clover_coeff
Definition: test_util.cpp:1645
Main header file for the QUDA library.
void * qudaResidentExtendedGaugeField(void *gauge, int geometry, int precision)
void qudaCloverTrace(void *out, void *dummy, int mu, int nu)
const int * latsize
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 * qudaAllocatePinned(size_t bytes)
__host__ __device__ ValueType arg(const complex< ValueType > &z)
Returns the phase angle of z.
Definition: complex_quda.h:880
void qudaHisqForce(int precision, int num_terms, int num_naik_terms, 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)
enum QudaVerbosity_s QudaVerbosity
QudaVerbosity verbosity
void qudaLoadKSLink(int precision, QudaFatLinkArgs_t fatlink_args, const double act_path_coeff[6], void *inlink, void *fatlink, void *longlink)
QudaExtLibType solver_ext_lib
QudaParity parity
Definition: covdev_test.cpp:53
void qudaGaugeForce(int precision, int num_loop_types, double milc_loop_coeff[3], double eb3, QudaMILCSiteArg_t *arg)
QudaPrecision prec
Definition: test_util.cpp:1615
enum QudaInverterType_s QudaInverterType
enum QudaMemoryType_s QudaMemoryType
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)
unsigned long long bytes
Definition: blas_quda.cu:43
enum QudaExtLibType_s QudaExtLibType
QudaLayout_t layout
void qudaDestroyGaugeField(void *gauge)