QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
7 #ifdef __cplusplus
8 extern "C" {
9 #endif
10 
11  typedef struct {
12  int max_iter;
13  QudaParity evenodd; // options are QUDA_EVEN_PARITY, QUDA_ODD_PARITY, QUDA_INVALID_PARITY
15  double boundary_phase[4];
17 
18 
19  typedef struct {
20  const int* latsize;
21  const int* machsize; // grid size
22  int device; // device number
23  } QudaLayout_t;
24 
25 
26  typedef struct {
29  } QudaInitArgs_t; // passed to the initialization struct
30 
31 
32  typedef struct {
37  double force_filter;
39 
40 
41  typedef struct {
42  int su3_source; // is the incoming gauge field su3?
43  int use_pinned_memory; // use page-locked memory in Quda?
45 
46 
47  void qudaInit(QudaInitArgs_t input);
48 
50 
51  void qudaFinalize();
52 
53 
54  void qudaHisqParamsInit(QudaHisqParams_t hisq_params);
55 
56 
57  void qudaLoadKSLink(int precision, QudaFatLinkArgs_t fatlink_args, const double act_path_coeff[6], void* inlink, void* fatlink, void* longlink);
58 
59 
60  void qudaLoadUnitarizedLink(int precision, QudaFatLinkArgs_t fatlink_args, const double path_coeff[6], void* inlink, void* fatlink, void* ulink);
61 
62 
63  void qudaInvert(int external_precision,
64  int quda_precision,
65  double mass,
66  QudaInvertArgs_t inv_args,
67  double target_resid,
68  double target_relresid,
69  const void* const milc_fatlink,
70  const void* const milc_longlink,
71  const double tadpole,
72  void* source,
73  void* solution,
74  double* const final_resid,
75  double* const final_rel_resid,
76  int* num_iters);
77 
78 
79  void qudaDDInvert(int external_precision,
80  int quda_precision,
81  double mass,
82  QudaInvertArgs_t inv_args,
83  double target_residual,
84  double target_fermilab_residual,
85  const int * const domain_overlap,
86  const void* const fatlink,
87  const void* const longlink,
88  void* source,
89  void* solution,
90  double* const final_residual,
91  double* const final_fermilab_residual,
92  int* num_iters);
93 
94 
96  int external_precision,
97  int precision,
98  int num_offsets,
99  double* const offset,
100  QudaInvertArgs_t inv_args,
101  const double* target_residual,
102  const double* target_relative_residual,
103  const void* const milc_fatlink,
104  const void* const milc_longlink,
105  const double tadpole,
106  void* source,
107  void** solutionArray,
108  double* const final_residual,
109  double* const final_relative_residual,
110  int* num_iters);
111 
112 
113  void qudaCloverInvert(int external_precision,
114  int quda_precision,
115  double kappa,
116  double clover_coeff,
117  QudaInvertArgs_t inv_args,
118  double target_residual,
119  double target_fermilab_residual,
120  const void* milc_link,
121  void* milc_clover,
122  void* milc_clover_inv,
123  void* source,
124  void* solution,
125  double* const final_residual,
126  double* const final_fermilab_residual,
127  int* num_iters
128  );
129 
130  void qudaLoadGaugeField(int external_precision,
131  int quda_precision,
132  QudaInvertArgs_t inv_args,
133  const void* milc_link) ;
134 
135  void qudaFreeGaugeField();
136 
137  void qudaLoadCloverField(int external_precision,
138  int quda_precision,
139  QudaInvertArgs_t inv_args,
140  void* milc_clover,
141  void* milc_clover_inv,
142  QudaSolutionType solution_type,
143  QudaSolveType solve_type,
144  double clover_coeff,
145  int compute_trlog,
146  double *trlog) ;
147 
148  void qudaFreeCloverField();
149 
150  void qudaCloverMultishiftInvert(int external_precision,
151  int quda_precision,
152  int num_offsets,
153  double* const offset,
154  double kappa,
155  double clover_coeff,
156  QudaInvertArgs_t inv_args,
157  const double* target_residual,
158  const void* milc_link,
159  void* milc_clover,
160  void* milc_clover_inv,
161  void* source,
162  void** solutionArray,
163  double* const final_residual,
164  int* num_iters
165  );
166 
167  void qudaCloverMultishiftMDInvert(int external_precision,
168  int quda_precision,
169  int num_offsets,
170  double* const offset,
171  double kappa,
172  double clover_coeff,
173  QudaInvertArgs_t inv_args,
174  const double* target_residual,
175  const void* milc_link,
176  void* milc_clover,
177  void* milc_clover_inv,
178  void* source,
179  void** psiEven,
180  void** psiOdd,
181  void** pEven,
182  void** pOdd,
183  double* const final_residual,
184  int* num_iters
185  );
186 
187  void qudaHisqForce(
188  int precision,
189  const double level2_coeff[6],
190  const double fat7_coeff[6],
191  const void* const staple_src[4],
192  const void* const one_link_src[4],
193  const void* const naik_src[4],
194  const void* const w_link,
195  const void* const v_link,
196  const void* const u_link,
197  void* const milc_momentum);
198 
199 
200  void qudaAsqtadForce(
201  int precision,
202  const double act_path_coeff[6],
203  const void* const one_link_src[4],
204  const void* const naik_src[4],
205  const void* const link,
206  void* const milc_momentum);
207 
208 
209  void qudaGaugeForce(int precision,
210  int num_loop_types,
211  double milc_loop_coeff[3],
212  double eb3,
213  void* milc_sitelink,
214  void* milc_momentum);
215 
216 
217  void qudaComputeOprod(int precision,
218  int num_terms,
219  double** coeff,
220  void** quark_field,
221  void* oprod[2]);
222 
223 
224  void qudaUpdateU(int precision,
225  double eps,
226  void* momentum,
227  void* link);
228 
229  void qudaCloverTrace(void* out, void* clover, int mu, int nu);
230 
231 
232  void qudaCloverDerivative(void* out, void* gauge, void* oprod,
233  int mu, int nu, double coeff, int precision, int parity, int conjugate);
234 
235 
236  void* qudaCreateExtendedGaugeField(void* gauge, int geometry, int precision);
237 
238  void* qudaCreateGaugeField(void* gauge, int geometry, int precision);
239 
240  void qudaSaveGaugeField(void* gauge, void* inGauge);
241 
242  void qudaDestroyGaugeField(void* gauge);
243 
244 
245 #ifdef __cplusplus
246 }
247 #endif
248 
249 
250 #endif // _QUDA_MILC_INTERFACE_H
void qudaGaugeForce(int precision, int num_loop_types, double milc_loop_coeff[3], double eb3, void *milc_sitelink, void *milc_momentum)
void qudaHisqParamsInit(QudaHisqParams_t hisq_params)
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 * qudaCreateExtendedGaugeField(void *gauge, int geometry, int precision)
void qudaComputeOprod(int precision, int num_terms, double **coeff, void **quark_field, void *oprod[2])
enum QudaSolveType_s QudaSolveType
__global__ void const RealA *const const RealA *const const RealA *const const RealB *const const RealB *const int int mu
void qudaInit(QudaInitArgs_t input)
void qudaLoadGaugeField(int external_precision, int quda_precision, QudaInvertArgs_t inv_args, const void *milc_link)
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 qudaSaveGaugeField(void *gauge, void *inGauge)
void qudaMultishiftInvert(int external_precision, int precision, int num_offsets, double *const offset, QudaInvertArgs_t inv_args, const double *target_residual, const double *target_relative_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_relative_residual, int *num_iters)
void qudaSetLayout(QudaLayout_t layout)
void * longlink[4]
void qudaCloverTrace(void *out, void *clover, int mu, int nu)
VOLATILE spinorFloat kappa
void qudaFreeCloverField()
void * qudaCreateGaugeField(void *gauge, int geometry, int precision)
void qudaUpdateU(int precision, double eps, void *momentum, void *link)
enum QudaSolutionType_s QudaSolutionType
const int * machsize
__constant__ double coeff
void qudaFreeGaugeField()
enum QudaParity_s QudaParity
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()
void qudaHisqForce(int precision, const double level2_coeff[6], const double fat7_coeff[6], const void *const staple_src[4], const void *const one_link_src[4], const void *const naik_src[4], const void *const w_link, const void *const v_link, const void *const u_link, void *const milc_momentum)
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 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)
cpuColorSpinorField * out
void qudaCloverMultishiftMDInvert(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 **psiEven, void **psiOdd, void **pEven, void **pOdd, double *const final_residual, int *num_iters)
Main header file for the QUDA library.
const int * latsize
void qudaInvert(int external_precision, int quda_precision, double mass, QudaInvertArgs_t inv_args, double target_resid, double target_relresid, 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)
QIO_Layout layout
Definition: gauge_qio.cpp:6
enum QudaVerbosity_s QudaVerbosity
QudaVerbosity verbosity
double mass
Definition: test_util.cpp:1569
void qudaLoadKSLink(int precision, QudaFatLinkArgs_t fatlink_args, const double act_path_coeff[6], void *inlink, void *fatlink, void *longlink)
const QudaParity parity
Definition: dslash_test.cpp:29
void * gauge[4]
Definition: su3_test.cpp:15
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)
QudaLayout_t layout
void qudaDestroyGaugeField(void *gauge)
void * fatlink[4]