QUDA  0.9.0
copy_gauge_mg.cu
Go to the documentation of this file.
1 #include <gauge_field_order.h>
2 
3 #define FINE_GRAINED_ACCESS
4 
5 #include <copy_gauge_helper.cuh>
6 
7 namespace quda {
8 
9  template <typename FloatOut, typename FloatIn, int length, typename InOrder>
10  void copyGaugeMG(const InOrder &inOrder, GaugeField &out, const GaugeField &in,
11  QudaFieldLocation location, FloatOut *Out, FloatOut **outGhost, int type) {
12  if (out.Reconstruct() != QUDA_RECONSTRUCT_NO)
13  errorQuda("Reconstruct type %d not supported", out.Reconstruct());
14 
15  int faceVolumeCB[QUDA_MAX_DIM];
16  for (int i=0; i<4; i++) faceVolumeCB[i] = out.SurfaceCB(i) * out.Nface();
17  if (out.isNative()) {
18 
19 #ifdef FINE_GRAINED_ACCESS
20  if (outGhost) {
22  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out.Volume(), faceVolumeCB,
23  out.Ndim(), out.Geometry(), out, in, location, type);
24  } else {
26  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out.Volume(), faceVolumeCB,
27  out.Ndim(), out.Geometry(), out, in, location, type);
28  }
29 #else
31  copyGauge<FloatOut,FloatIn,length>
32  (G(out,Out,outGhost), inOrder, out.Volume(), faceVolumeCB,
33  out.Ndim(), out.Geometry(), out, in, location, type);
34 #endif
35 
36  } else if (out.Order() == QUDA_QDP_GAUGE_ORDER) {
37 
38 #ifdef FINE_GRAINED_ACCESS
40  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out.Volume(),
41  faceVolumeCB, out.Ndim(), out.Geometry(), out, in, location, type);
42 #else
43  typedef typename QDPOrder<FloatOut,length> G;
44  copyGauge<FloatOut,FloatIn,length>(G(out, Out, outGhost), inOrder, out.Volume(),
45  faceVolumeCB, out.Ndim(), out.Geometry(), out, in, location, type);
46 #endif
47 
48  } else if (out.Order() == QUDA_MILC_GAUGE_ORDER) {
49 
50 #ifdef FINE_GRAINED_ACCESS
52  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out.Volume(),
53  faceVolumeCB, out.Ndim(), out.Geometry(), out, in, location, type);
54 #else
55  typedef typename MILCOrder<FloatOut,length> G;
56  copyGauge<FloatOut,FloatIn,length>(G(out, Out, outGhost), inOrder, out.Volume(),
57  faceVolumeCB, out.Ndim(), out.Geometry(), out, in, location, type);
58 #endif
59 
60  } else {
61  errorQuda("Gauge field %d order not supported", out.Order());
62  }
63 
64  }
65 
66  template <typename FloatOut, typename FloatIn, int length>
68  FloatOut *Out, FloatIn *In, FloatOut **outGhost, FloatIn **inGhost, int type) {
69 
70  if (in.Reconstruct() != QUDA_RECONSTRUCT_NO)
71  errorQuda("Reconstruct type %d not supported", in.Reconstruct());
72 
73  // reconstruction only supported on FloatN fields currently
74  if (in.isNative()) {
75 #ifdef FINE_GRAINED_ACCESS
76  if (inGhost) {
78  copyGaugeMG<FloatOut,FloatIn,length> (G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
79  } else {
81  copyGaugeMG<FloatOut,FloatIn,length> (G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
82  }
83 #else
85  copyGaugeMG<FloatOut,FloatIn,length> (G(in, In,inGhost), out, in, location, Out, outGhost, type);
86 #endif
87  } else if (in.Order() == QUDA_QDP_GAUGE_ORDER) {
88 
89 #ifdef FINE_GRAINED_ACCESS
91  copyGaugeMG<FloatOut,FloatIn,length>(G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
92 #else
93  typedef typename QDPOrder<FloatIn,length> G;
94  copyGaugeMG<FloatOut,FloatIn,length>(G(in, In, inGhost), out, in, location, Out, outGhost, type);
95 #endif
96 
97  } else if (in.Order() == QUDA_MILC_GAUGE_ORDER) {
98 
99 #ifdef FINE_GRAINED_ACCESS
101  copyGaugeMG<FloatOut,FloatIn,length>(G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
102 #else
103  typedef typename MILCOrder<FloatIn,length> G;
104  copyGaugeMG<FloatOut,FloatIn,length>(G(in, In, inGhost), out, in, location, Out, outGhost, type);
105 #endif
106 
107  } else {
108  errorQuda("Gauge field %d order not supported", in.Order());
109  }
110 
111  }
112 
113  template <typename FloatOut, typename FloatIn>
114  void copyGaugeMG(GaugeField &out, const GaugeField &in, QudaFieldLocation location, FloatOut *Out,
115  FloatIn *In, FloatOut **outGhost, FloatIn **inGhost, int type) {
116 
117 #ifdef GPU_MULTIGRID
118  if (in.Ncolor() == 4) {
119  const int Nc = 4;
120  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
121  } else if (in.Ncolor() == 8) {
122  const int Nc = 8;
123  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
124  } else if (in.Ncolor() == 16) {
125  const int Nc = 16;
126  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
127  } else if (in.Ncolor() == 24) {
128  const int Nc = 24;
129  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
130  } else if (in.Ncolor() == 32) {
131  const int Nc = 32;
132  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
133  } else if (in.Ncolor() == 40) {
134  const int Nc = 40;
135  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
136  } else if (in.Ncolor() == 48) {
137  const int Nc = 48;
138  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
139  } else if (in.Ncolor() == 56) {
140  const int Nc = 56;
141  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
142  } else if (in.Ncolor() == 64) {
143  const int Nc = 64;
144  copyGaugeMG<FloatOut,FloatIn,2*Nc*Nc>(out, in, location, Out, In, outGhost, inGhost, type);
145  } else
146 #endif // GPU_MULTIGRID
147  {
148  errorQuda("Unsupported number of colors; out.Nc=%d, in.Nc=%d", out.Ncolor(), in.Ncolor());
149  }
150  }
151 
152  // this is the function that is actually called, from here on down we instantiate all required templates
154  void *Out, void *In, void **ghostOut, void **ghostIn, int type) {
155 
157 #ifdef GPU_MULTIGRID_DOUBLE
159  copyGaugeMG(out, in, location, (double*)Out, (double*)In, (double**)ghostOut, (double**)ghostIn, type);
160  } else if (in.Precision() == QUDA_SINGLE_PRECISION) {
161  copyGaugeMG(out, in, location, (double*)Out, (float*)In, (double**)ghostOut, (float**)ghostIn, type);
162  } else {
163  errorQuda("Precision %d not supported", in.Precision());
164  }
165 #else
166  errorQuda("Double precision multigrid has not been enabled");
167 #endif
168  } else if (out.Precision() == QUDA_SINGLE_PRECISION) {
170 #ifdef GPU_MULTIGRID_DOUBLE
171  copyGaugeMG(out, in, location, (float*)Out, (double*)In, (float**)ghostOut, (double**)ghostIn, type);
172 #else
173  errorQuda("Double precision multigrid has not been enabled");
174 #endif
175  } else if (in.Precision() == QUDA_SINGLE_PRECISION) {
176  copyGaugeMG(out, in, location, (float*)Out, (float*)In, (float**)ghostOut, (float**)ghostIn, type);
177  } else {
178  errorQuda("Precision %d not supported", in.Precision());
179  }
180  } else {
181  errorQuda("Precision %d not supported", out.Precision());
182  }
183  }
184 
185 
186 
187 } // namespace quda
#define errorQuda(...)
Definition: util_quda.h:90
const int * SurfaceCB() const
cpuColorSpinorField * in
void copyGaugeMG(const InOrder &inOrder, GaugeField &out, const GaugeField &in, QudaFieldLocation location, FloatOut *Out, FloatOut **outGhost, int type)
Main header file for host and device accessors to GaugeFields.
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
void copyGenericGaugeMG(GaugeField &out, const GaugeField &in, QudaFieldLocation location, void *Out, void *In, void **ghostOut, void **ghostIn, int type)
#define QUDA_MAX_DIM
Maximum number of dimensions supported by QUDA. In practice, no routines make use of more than 5...
QudaPrecision Precision() const