QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
copy_gauge_inc.cu
Go to the documentation of this file.
1 #include <gauge_field_order.h>
2 #include <copy_gauge_helper.cuh>
3 
4 namespace quda {
5 
6  template <typename FloatOut, typename FloatIn, int length, typename InOrder>
7  void copyGauge(const InOrder &inOrder, const GaugeField &out, const GaugeField &in,
8  QudaFieldLocation location, FloatOut *Out, FloatOut **outGhost, int type) {
9  if (out.isNative()) {
10  // this overrides the check that the texture maps to the gauge
11  // pointer - this is safe here since it only occurs when running
12  // the copier on the host when we will not be using texture
13  // reads
14  const bool override = true;
15  if (out.Reconstruct() == QUDA_RECONSTRUCT_NO) {
17  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
18  } else if (out.Reconstruct() == QUDA_RECONSTRUCT_12) {
19 #if QUDA_RECONSTRUCT & 2
21  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
22 #else
23  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-12", QUDA_RECONSTRUCT);
24 #endif
25  } else if (out.Reconstruct() == QUDA_RECONSTRUCT_8) {
26 #if QUDA_RECONSTRUCT & 1
28  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
29 #else
30  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-8", QUDA_RECONSTRUCT);
31 #endif
32 #ifdef GPU_STAGGERED_DIRAC
33  } else if (out.Reconstruct() == QUDA_RECONSTRUCT_13) {
34 #if QUDA_RECONSTRUCT & 2
36  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
37 #else
38  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-13", QUDA_RECONSTRUCT);
39 #endif
40  } else if (out.Reconstruct() == QUDA_RECONSTRUCT_9) {
41 #if QUDA_RECONSTRUCT & 1
44  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
45  } else if (out.StaggeredPhase() == QUDA_STAGGERED_PHASE_TIFR) {
46 #ifdef BUILD_TIFR_INTERFACE
48  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
49 #else
50  errorQuda("TIFR interface has not been built so TIFR phase type not enabled\n");
51 #endif
52  } else if (out.StaggeredPhase() == QUDA_STAGGERED_PHASE_NO) {
54  copyGauge<FloatOut, FloatIn, length>(G(out, Out, outGhost, override), inOrder, out, in, location, type);
55  } else {
56  errorQuda("Staggered phase type %d not supported", out.StaggeredPhase());
57  }
58 #else
59  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-9", QUDA_RECONSTRUCT);
60 #endif
61 #endif // GPU_STAGGERED_DIRAC
62  } else {
63  errorQuda("Reconstruction %d and order %d not supported", out.Reconstruct(), out.Order());
64  }
65  } else if (out.Order() == QUDA_QDP_GAUGE_ORDER) {
66 
67 #ifdef BUILD_QDP_INTERFACE
68  copyGauge<FloatOut,FloatIn,length>
69  (QDPOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
70 #else
71  errorQuda("QDP interface has not been built\n");
72 #endif
73 
74  } else if (out.Order() == QUDA_QDPJIT_GAUGE_ORDER) {
75 
76 #ifdef BUILD_QDPJIT_INTERFACE
77  copyGauge<FloatOut,FloatIn,length>
78  (QDPJITOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
79 #else
80  errorQuda("QDPJIT interface has not been built\n");
81 #endif
82 
83  } else if (out.Order() == QUDA_CPS_WILSON_GAUGE_ORDER) {
84 
85 #ifdef BUILD_CPS_INTERFACE
86  copyGauge<FloatOut,FloatIn,length>
87  (CPSOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
88 #else
89  errorQuda("CPS interface has not been built\n");
90 #endif
91 
92  } else if (out.Order() == QUDA_MILC_GAUGE_ORDER) {
93 
94 #ifdef BUILD_MILC_INTERFACE
95  copyGauge<FloatOut,FloatIn,length>
96  (MILCOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
97 #else
98  errorQuda("MILC interface has not been built\n");
99 #endif
100 
101  } else if (out.Order() == QUDA_MILC_SITE_GAUGE_ORDER) {
102 
103 #ifdef BUILD_MILC_INTERFACE
104  copyGauge<FloatOut,FloatIn,length>
105  (MILCSiteOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
106 #else
107  errorQuda("MILC interface has not been built\n");
108 #endif
109 
110  } else if (out.Order() == QUDA_BQCD_GAUGE_ORDER) {
111 
112 #ifdef BUILD_BQCD_INTERFACE
113  copyGauge<FloatOut,FloatIn,length>
114  (BQCDOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
115 #else
116  errorQuda("BQCD interface has not been built\n");
117 #endif
118 
119  } else if (out.Order() == QUDA_TIFR_GAUGE_ORDER) {
120 
121 #ifdef BUILD_TIFR_INTERFACE
122  copyGauge<FloatOut,FloatIn,length>
123  (TIFROrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
124 #else
125  errorQuda("TIFR interface has not been built\n");
126 #endif
127 
128  } else if (out.Order() == QUDA_TIFR_PADDED_GAUGE_ORDER) {
129 
130 #ifdef BUILD_TIFR_INTERFACE
131  copyGauge<FloatOut,FloatIn,length>
132  (TIFRPaddedOrder<FloatOut,length>(out, Out, outGhost), inOrder, out, in, location, type);
133 #else
134  errorQuda("TIFR interface has not been built\n");
135 #endif
136 
137  } else {
138  errorQuda("Gauge field %d order not supported", out.Order());
139  }
140 
141  }
142 
143  template <typename FloatOut, typename FloatIn, int length>
144  void copyGauge(GaugeField &out, const GaugeField &in, QudaFieldLocation location, FloatOut *Out, FloatIn *In,
145  FloatOut **outGhost, FloatIn **inGhost, int type)
146  {
147 
148  // reconstruction only supported on FloatN fields currently
149  if (in.isNative()) {
150  // this overrides the check that the texture maps to the gauge
151  // pointer - this is safe here since it only occurs when running
152  // the copier on the host when we will not be using texture
153  // reads
154  const bool override = true;
155  if (in.Reconstruct() == QUDA_RECONSTRUCT_NO) {
157  copyGauge<FloatOut, FloatIn, length>(G(in, In, inGhost, override), out, in, location, Out, outGhost, type);
158  } else if (in.Reconstruct() == QUDA_RECONSTRUCT_12) {
159 #if QUDA_RECONSTRUCT & 2
161  copyGauge<FloatOut,FloatIn,length> (G(in,In,inGhost,override), out, in, location, Out, outGhost, type);
162 #else
163  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-12", QUDA_RECONSTRUCT);
164 #endif
165  } else if (in.Reconstruct() == QUDA_RECONSTRUCT_8) {
166 #if QUDA_RECONSTRUCT & 1
168  copyGauge<FloatOut,FloatIn,length> (G(in,In,inGhost,override), out, in, location, Out, outGhost, type);
169 #else
170  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-8", QUDA_RECONSTRUCT);
171 #endif
172 #ifdef GPU_STAGGERED_DIRAC
173  } else if (in.Reconstruct() == QUDA_RECONSTRUCT_13) {
174 #if QUDA_RECONSTRUCT & 2
176  copyGauge<FloatOut,FloatIn,length> (G(in,In,inGhost,override), out, in, location, Out, outGhost, type);
177 #else
178  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-13", QUDA_RECONSTRUCT);
179 #endif
180  } else if (in.Reconstruct() == QUDA_RECONSTRUCT_9) {
181 #if QUDA_RECONSTRUCT & 1
184  copyGauge<FloatOut, FloatIn, length>(G(in, In, inGhost, override), out, in, location, Out, outGhost, type);
185  } else if (in.StaggeredPhase() == QUDA_STAGGERED_PHASE_NO) {
187  copyGauge<FloatOut, FloatIn, length>(G(in, In, inGhost, override), out, in, location, Out, outGhost, type);
188  } else {
189  errorQuda("Staggered phase type %d not supported", in.StaggeredPhase());
190  }
191 #else
192  errorQuda("QUDA_RECONSTRUCT=%d does not enable reconstruct-9", QUDA_RECONSTRUCT);
193 #endif
194 #endif // GPU_STAGGERED_DIRAC
195  } else {
196  errorQuda("Reconstruction %d and order %d not supported", in.Reconstruct(), in.Order());
197  }
198  } else if (in.Order() == QUDA_QDP_GAUGE_ORDER) {
199 
200 #ifdef BUILD_QDP_INTERFACE
201  copyGauge<FloatOut,FloatIn,length>(QDPOrder<FloatIn,length>(in, In, inGhost),
202  out, in, location, Out, outGhost, type);
203 #else
204  errorQuda("QDP interface has not been built\n");
205 #endif
206 
207  } else if (in.Order() == QUDA_QDPJIT_GAUGE_ORDER) {
208 
209 #ifdef BUILD_QDPJIT_INTERFACE
210  copyGauge<FloatOut,FloatIn,length>(QDPJITOrder<FloatIn,length>(in, In, inGhost),
211  out, in, location, Out, outGhost, type);
212 #else
213  errorQuda("QDPJIT interface has not been built\n");
214 #endif
215 
216  } else if (in.Order() == QUDA_CPS_WILSON_GAUGE_ORDER) {
217 
218 #ifdef BUILD_CPS_INTERFACE
219  copyGauge<FloatOut,FloatIn,length>(CPSOrder<FloatIn,length>(in, In, inGhost),
220  out, in, location, Out, outGhost, type);
221 #else
222  errorQuda("CPS interface has not been built\n");
223 #endif
224 
225  } else if (in.Order() == QUDA_MILC_GAUGE_ORDER) {
226 
227 #ifdef BUILD_MILC_INTERFACE
228  copyGauge<FloatOut,FloatIn,length>(MILCOrder<FloatIn,length>(in, In, inGhost),
229  out, in, location, Out, outGhost, type);
230 #else
231  errorQuda("MILC interface has not been built\n");
232 #endif
233 
234  } else if (in.Order() == QUDA_MILC_SITE_GAUGE_ORDER) {
235 
236 #ifdef BUILD_MILC_INTERFACE
237  copyGauge<FloatOut,FloatIn,length>(MILCSiteOrder<FloatIn,length>(in, In, inGhost),
238  out, in, location, Out, outGhost, type);
239 #else
240  errorQuda("MILC interface has not been built\n");
241 #endif
242 
243  } else if (in.Order() == QUDA_BQCD_GAUGE_ORDER) {
244 
245 #ifdef BUILD_BQCD_INTERFACE
246  copyGauge<FloatOut,FloatIn,length>(BQCDOrder<FloatIn,length>(in, In, inGhost),
247  out, in, location, Out, outGhost, type);
248 #else
249  errorQuda("BQCD interface has not been built\n");
250 #endif
251 
252  } else if (in.Order() == QUDA_TIFR_GAUGE_ORDER) {
253 
254 #ifdef BUILD_TIFR_INTERFACE
255  copyGauge<FloatOut,FloatIn,length>(TIFROrder<FloatIn,length>(in, In, inGhost),
256  out, in, location, Out, outGhost, type);
257 #else
258  errorQuda("TIFR interface has not been built\n");
259 #endif
260 
261  } else if (in.Order() == QUDA_TIFR_PADDED_GAUGE_ORDER) {
262 
263 #ifdef BUILD_TIFR_INTERFACE
264  copyGauge<FloatOut,FloatIn,length>(TIFRPaddedOrder<FloatIn,length>(in, In, inGhost),
265  out, in, location, Out, outGhost, type);
266 #else
267  errorQuda("TIFR interface has not been built\n");
268 #endif
269 
270  } else {
271  errorQuda("Gauge field order %d not supported", in.Order());
272  }
273  }
274 
275  void checkMomOrder(const GaugeField &u);
276 
277  template <typename FloatOut, typename FloatIn, int length, typename Out, typename In, typename Arg>
278  void copyMom(Arg &arg, const GaugeField &out, const GaugeField &in, QudaFieldLocation location) {
279  CopyGauge<FloatOut,FloatIn,length, Arg> momCopier(arg, out, in, location);
280  momCopier.apply(0);
281  }
282 
283  template <typename FloatOut, typename FloatIn>
284  void copyGauge(GaugeField &out, const GaugeField &in, QudaFieldLocation location, FloatOut *Out, FloatIn *In,
285  FloatOut **outGhost, FloatIn **inGhost, int type)
286  {
287 
288  if (in.Ncolor() != 3 && out.Ncolor() != 3) {
289  errorQuda("Unsupported number of colors; out.Nc=%d, in.Nc=%d", out.Ncolor(), in.Ncolor());
290  }
291 
292  if (out.Geometry() != in.Geometry()) {
293  errorQuda("Field geometries %d %d do not match", out.Geometry(), in.Geometry());
294  }
295 
297  // we are doing gauge field packing
298  copyGauge<FloatOut,FloatIn,18>(out, in, location, Out, In, outGhost, inGhost, type);
299  } else {
300  if (out.Geometry() != QUDA_VECTOR_GEOMETRY) errorQuda("Unsupported geometry %d", out.Geometry());
301 
302  checkMomOrder(in);
303  checkMomOrder(out);
304 
305  // this overrides the check that the texture maps to the gauge
306  // pointer - this is safe here since it only occurs when running
307  // the copier on the host when we will not be using texture
308  // reads
309  const bool override = true;
310 
311  // momentum only currently supported on MILC (10), TIFR (18) and Float2 (10) fields currently
312  if (out.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
313  if (in.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
314  typedef FloatNOrder<FloatOut,10,2,10> momOut;
315  typedef FloatNOrder<FloatIn,10,2,10> momIn;
316  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out, 0, override), momIn(in, In, 0, override), in);
317  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
318  } else if (in.Order() == QUDA_MILC_GAUGE_ORDER) {
319 #ifdef BUILD_MILC_INTERFACE
320  typedef FloatNOrder<FloatOut,10,2,10> momOut;
321  typedef MILCOrder<FloatIn,10> momIn;
322  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out, 0, override), momIn(in, In), in);
323  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
324 #else
325  errorQuda("MILC interface has not been built\n");
326 #endif
327  } else if (in.Order() == QUDA_MILC_SITE_GAUGE_ORDER) {
328 #ifdef BUILD_MILC_INTERFACE
329  typedef FloatNOrder<FloatOut,10,2,10> momOut;
330  typedef MILCSiteOrder<FloatIn,10> momIn;
331  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out, 0, override), momIn(in, In), in);
332  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
333 #else
334  errorQuda("MILC interface has not been built\n");
335 #endif
336  } else if (in.Order() == QUDA_TIFR_GAUGE_ORDER) {
337 #ifdef BUILD_TIFR_INTERFACE
338  typedef FloatNOrder<FloatOut,18,2,11> momOut;
339  typedef TIFROrder<FloatIn,18> momIn;
340  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out, 0, override), momIn(in, In), in);
341  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
342 #else
343  errorQuda("TIFR interface has not been built\n");
344 #endif
345  } else if (in.Order() == QUDA_TIFR_PADDED_GAUGE_ORDER) {
346 #ifdef BUILD_TIFR_INTERFACE
347  typedef FloatNOrder<FloatOut,18,2,11> momOut;
348  typedef TIFRPaddedOrder<FloatIn,18> momIn;
349  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out, 0, override), momIn(in, In), in);
350  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
351 #else
352  errorQuda("TIFR interface has not been built\n");
353 #endif
354  } else {
355  errorQuda("Gauge field orders %d not supported", in.Order());
356  }
357  } else if (out.Order() == QUDA_MILC_GAUGE_ORDER) {
358  typedef MILCOrder<FloatOut,10> momOut;
359 #ifdef BUILD_MILC_INTERFACE
360  if (in.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
361  typedef FloatNOrder<FloatIn,10,2,10> momIn;
362  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In, 0, override), in);
363  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
364  } else if (in.Order() == QUDA_MILC_GAUGE_ORDER) {
365  typedef MILCOrder<FloatIn,10> momIn;
366  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In), in);
367  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
368  } else {
369  errorQuda("Gauge field orders %d not supported", in.Order());
370  }
371 #else
372  errorQuda("MILC interface has not been built\n");
373 #endif
374  } else if (out.Order() == QUDA_MILC_SITE_GAUGE_ORDER) {
375  typedef MILCSiteOrder<FloatOut,10> momOut;
376 #ifdef BUILD_MILC_INTERFACE
377  if (in.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
378  typedef FloatNOrder<FloatIn,10,2,10> momIn;
379  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In, 0, override), in);
380  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
381  } else if (in.Order() == QUDA_MILC_SITE_GAUGE_ORDER) {
382  typedef MILCSiteOrder<FloatIn,10> momIn;
383  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In), in);
384  copyMom<FloatOut,FloatIn,10,momOut,momIn>(arg,out,in,location);
385  } else {
386  errorQuda("Gauge field orders %d not supported", in.Order());
387  }
388 #else
389  errorQuda("MILC interface has not been built\n");
390 #endif
391  } else if (out.Order() == QUDA_TIFR_GAUGE_ORDER) {
392  typedef TIFROrder<FloatOut,18> momOut;
393 #ifdef BUILD_TIFR_INTERFACE
394  if (in.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
395  // FIX ME - 11 is a misnomer to avoid confusion in template instantiation
396  typedef FloatNOrder<FloatIn,18,2,11> momIn;
397  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In, 0, override), in);
398  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
399  } else if (in.Order() == QUDA_TIFR_GAUGE_ORDER) {
400  typedef TIFROrder<FloatIn,18> momIn;
401  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In), in);
402  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
403  } else {
404  errorQuda("Gauge field orders %d not supported", in.Order());
405  }
406 #else
407  errorQuda("TIFR interface has not been built\n");
408 #endif
409  } else if (out.Order() == QUDA_TIFR_PADDED_GAUGE_ORDER) {
410  typedef TIFRPaddedOrder<FloatOut,18> momOut;
411 #ifdef BUILD_TIFR_INTERFACE
412  if (in.Order() == QUDA_FLOAT2_GAUGE_ORDER) {
413  // FIX ME - 11 is a misnomer to avoid confusion in template instantiation
414  typedef FloatNOrder<FloatIn,18,2,11> momIn;
415  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In, 0, override), in);
416  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
417  } else if (in.Order() == QUDA_TIFR_PADDED_GAUGE_ORDER) {
418  typedef TIFRPaddedOrder<FloatIn,18> momIn;
419  CopyGaugeArg<momOut,momIn> arg(momOut(out, Out), momIn(in, In), in);
420  copyMom<FloatOut,FloatIn,18,momOut,momIn>(arg,out,in,location);
421  } else {
422  errorQuda("Gauge field orders %d not supported", in.Order());
423  }
424 #else
425  errorQuda("TIFR interface has not been built\n");
426 #endif
427  } else {
428  errorQuda("Gauge field orders %d not supported", out.Order());
429  }
430  }
431  }
432 
433  // this is the function that is actually called, from here on down we instantiate all required templates
434  template <typename FloatOut>
435  void copyGenericGauge(GaugeField &out, const GaugeField &in, QudaFieldLocation location, void *Out, void *In,
436  void **ghostOut, void **ghostIn, int type)
437  {
438  if (in.Precision() == QUDA_DOUBLE_PRECISION) {
439  copyGauge(out, in, location, (FloatOut *)Out, (double *)In, (FloatOut **)ghostOut, (double **)ghostIn, type);
440  } else if (in.Precision() == QUDA_SINGLE_PRECISION) {
441 #if QUDA_PRECISION & 4
442  copyGauge(out, in, location, (FloatOut *)Out, (float *)In, (FloatOut **)ghostOut, (float **)ghostIn, type);
443 #else
444  errorQuda("QUDA_PRECISION=%d does not enable single precision", QUDA_PRECISION);
445 #endif
446  } else if (in.Precision() == QUDA_HALF_PRECISION) {
447 #if QUDA_PRECISION & 2
448  copyGauge(out, in, location, (FloatOut *)Out, (short *)In, (FloatOut **)ghostOut, (short **)ghostIn, type);
449 #else
450  errorQuda("QUDA_PRECISION=%d does not enable half precision", QUDA_PRECISION);
451 #endif
452  } else if (in.Precision() == QUDA_QUARTER_PRECISION) {
453 #if QUDA_PRECISION & 1
454  copyGauge(out, in, location, (FloatOut *)Out, (char *)In, (FloatOut **)ghostOut, (char **)ghostIn, type);
455 #else
456  errorQuda("QUDA_PRECISION=%d does not enable quarter precision", QUDA_PRECISION);
457 #endif
458  } else {
459  errorQuda("Unsupported precision %d", in.Precision());
460  }
461  }
462 
463 } // namespace quda
struct to define TIFR ordered gauge fields: [mu][parity][volumecb][col][row]
void apply(const cudaStream_t &stream)
void copyGenericGauge(GaugeField &out, const GaugeField &in, QudaFieldLocation location, void *Out=0, void *In=0, void **ghostOut=0, void **ghostIn=0, int type=0)
Definition: copy_gauge.cu:41
#define errorQuda(...)
Definition: util_quda.h:121
QudaLinkType LinkType() const
Definition: gauge_field.h:255
QudaFieldGeometry Geometry() const
Definition: gauge_field.h:258
void copyGauge(Arg &arg)
Definition: copy_gauge.cuh:32
struct to define gauge fields packed into an opaque MILC site struct:
int Ncolor() const
Definition: gauge_field.h:249
cpuColorSpinorField * in
Main header file for host and device accessors to GaugeFields.
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
void copyMom(Arg &arg, const GaugeField &out, const GaugeField &in, QudaFieldLocation location)
__host__ __device__ ValueType arg(const complex< ValueType > &z)
Returns the phase angle of z.
Accessor routine for CloverFields in native field order.
QudaReconstructType Reconstruct() const
Definition: gauge_field.h:250
QudaGaugeFieldOrder Order() const
Definition: gauge_field.h:251
QudaStaggeredPhase StaggeredPhase() const
Definition: gauge_field.h:259
void checkMomOrder(const GaugeField &u)
Definition: copy_gauge.cu:22
QudaPrecision Precision() const
bool isNative() const