QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 sFloatOut, typename sFloatIn, int Nc, typename InOrder>
10  void copyGaugeMG(const InOrder &inOrder, GaugeField &out, const GaugeField &in,
11  QudaFieldLocation location, sFloatOut *Out, sFloatOut **outGhost, int type) {
12 
13  typedef typename mapper<sFloatOut>::type FloatOut;
14  typedef typename mapper<sFloatIn>::type FloatIn;
15  constexpr int length = 2*Nc*Nc;
16 
17  if (out.Reconstruct() != QUDA_RECONSTRUCT_NO)
18  errorQuda("Reconstruct type %d not supported", out.Reconstruct());
19 
20 #ifdef FINE_GRAINED_ACCESS
21  if (out.Precision() == QUDA_HALF_PRECISION) {
22  if (in.Precision() == QUDA_HALF_PRECISION) {
23  out.Scale(in.Scale());
24  } else {
25  InOrder in_(const_cast<GaugeField&>(in));
26  out.Scale( in.abs_max() );
27  }
28  }
29 #endif
30 
31  if (out.isNative()) {
32 
33 #ifdef FINE_GRAINED_ACCESS
34  if (outGhost) {
36  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out, in, location, type);
37  } else {
39  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out, in, location, type);
40  }
41 #else
43  copyGauge<FloatOut,FloatIn,length>
44  (G(out,Out,outGhost), inOrder, out, in, location, type);
45 #endif
46 
47  } else if (out.Order() == QUDA_QDP_GAUGE_ORDER) {
48 
49 #ifdef FINE_GRAINED_ACCESS
51  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out, in, location, type);
52 #else
53  typedef typename QDPOrder<FloatOut,length> G;
54  copyGauge<FloatOut,FloatIn,length>(G(out, Out, outGhost), inOrder, out, in, location, type);
55 #endif
56 
57  } else if (out.Order() == QUDA_MILC_GAUGE_ORDER) {
58 
59 #ifdef FINE_GRAINED_ACCESS
61  copyGauge<FloatOut,FloatIn,length>(G(out,(void*)Out,(void**)outGhost), inOrder, out, in, location, type);
62 #else
63  typedef typename MILCOrder<FloatOut,length> G;
64  copyGauge<FloatOut,FloatIn,length>(G(out, Out, outGhost), inOrder, out, in, location, type);
65 #endif
66 
67  } else {
68  errorQuda("Gauge field %d order not supported", out.Order());
69  }
70 
71  }
72 
73  template <typename sFloatOut, typename sFloatIn, int Nc>
75  sFloatOut *Out, sFloatIn *In, sFloatOut **outGhost, sFloatIn **inGhost, int type) {
76 
77  typedef typename mapper<sFloatOut>::type FloatOut;
78  typedef typename mapper<sFloatIn>::type FloatIn;
79 #ifndef FINE_GRAINED_ACCESS
80  constexpr int length = 2*Nc*Nc;
81 #endif
82 
83  if (in.Reconstruct() != QUDA_RECONSTRUCT_NO)
84  errorQuda("Reconstruct type %d not supported", in.Reconstruct());
85 
86  if (in.isNative()) {
87 #ifdef FINE_GRAINED_ACCESS
88  if (inGhost) {
90  copyGaugeMG<sFloatOut,sFloatIn,Nc> (G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
91  } else {
93  copyGaugeMG<sFloatOut,sFloatIn,Nc> (G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
94  }
95 #else
97  copyGaugeMG<FloatOut,FloatIn,Nc> (G(in, In,inGhost), out, in, location, Out, outGhost, type);
98 #endif
99  } else if (in.Order() == QUDA_QDP_GAUGE_ORDER) {
100 
101 #ifdef FINE_GRAINED_ACCESS
103  copyGaugeMG<sFloatOut,sFloatIn,Nc>(G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
104 #else
105  typedef typename QDPOrder<FloatIn,length> G;
106  copyGaugeMG<FloatOut,FloatIn,Nc>(G(in, In, inGhost), out, in, location, Out, outGhost, type);
107 #endif
108 
109  } else if (in.Order() == QUDA_MILC_GAUGE_ORDER) {
110 
111 #ifdef FINE_GRAINED_ACCESS
113  copyGaugeMG<sFloatOut,sFloatIn,Nc>(G(const_cast<GaugeField&>(in),(void*)In,(void**)inGhost), out, in, location, Out, outGhost, type);
114 #else
115  typedef typename MILCOrder<FloatIn,length> G;
116  copyGaugeMG<FloatOut,FloatIn,Nc>(G(in, In, inGhost), out, in, location, Out, outGhost, type);
117 #endif
118 
119  } else {
120  errorQuda("Gauge field %d order not supported", in.Order());
121  }
122 
123  }
124 
125  template <typename FloatOut, typename FloatIn>
126  void copyGaugeMG(GaugeField &out, const GaugeField &in, QudaFieldLocation location, FloatOut *Out,
127  FloatIn *In, FloatOut **outGhost, FloatIn **inGhost, int type) {
128 
129  switch(in.Ncolor()) {
130 #ifdef GPU_MULTIGRID
131  case 8: copyGaugeMG<FloatOut,FloatIn, 8>(out, in, location, Out, In, outGhost, inGhost, type); break;
132  case 12: copyGaugeMG<FloatOut,FloatIn,12>(out, in, location, Out, In, outGhost, inGhost, type); break;
133  case 16: copyGaugeMG<FloatOut,FloatIn,16>(out, in, location, Out, In, outGhost, inGhost, type); break;
134  case 24: copyGaugeMG<FloatOut,FloatIn,24>(out, in, location, Out, In, outGhost, inGhost, type); break;
135  case 32: copyGaugeMG<FloatOut,FloatIn,32>(out, in, location, Out, In, outGhost, inGhost, type); break;
136  case 40: copyGaugeMG<FloatOut,FloatIn,40>(out, in, location, Out, In, outGhost, inGhost, type); break;
137  case 48: copyGaugeMG<FloatOut,FloatIn,48>(out, in, location, Out, In, outGhost, inGhost, type); break;
138  case 56: copyGaugeMG<FloatOut,FloatIn,56>(out, in, location, Out, In, outGhost, inGhost, type); break;
139  case 64: copyGaugeMG<FloatOut,FloatIn,64>(out, in, location, Out, In, outGhost, inGhost, type); break;
140 #endif // GPU_MULTIGRID
141  default: errorQuda("Unsupported number of colors; out.Nc=%d, in.Nc=%d", out.Ncolor(), in.Ncolor());
142  }
143  }
144 
145  // this is the function that is actually called, from here on down we instantiate all required templates
147  void *Out, void *In, void **ghostOut, void **ghostIn, int type) {
148 
149 #ifndef FINE_GRAINED_ACCESS
151  errorQuda("Precision format not supported");
152 #endif
153 
154  if (out.Precision() == QUDA_DOUBLE_PRECISION) {
155 #ifdef GPU_MULTIGRID_DOUBLE
156  if (in.Precision() == QUDA_DOUBLE_PRECISION) {
157  copyGaugeMG(out, in, location, (double*)Out, (double*)In, (double**)ghostOut, (double**)ghostIn, type);
158  } else if (in.Precision() == QUDA_SINGLE_PRECISION) {
159  copyGaugeMG(out, in, location, (double*)Out, (float*)In, (double**)ghostOut, (float**)ghostIn, type);
160  } else if (in.Precision() == QUDA_HALF_PRECISION) {
161  copyGaugeMG(out, in, location, (double*)Out, (short*)In, (double**)ghostOut, (short**)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) {
169  if (in.Precision() == QUDA_DOUBLE_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 if (in.Precision() == QUDA_HALF_PRECISION) {
178  copyGaugeMG(out, in, location, (float*)Out, (short*)In, (float**)ghostOut, (short**)ghostIn, type);
179  } else {
180  errorQuda("Precision %d not supported", in.Precision());
181  }
182  } else if (out.Precision() == QUDA_HALF_PRECISION) {
183  if (in.Precision() == QUDA_DOUBLE_PRECISION) {
184 #ifdef GPU_MULTIGRID_DOUBLE
185  copyGaugeMG(out, in, location, (short*)Out, (double*)In, (short**)ghostOut, (double**)ghostIn, type);
186 #else
187  errorQuda("Double precision multigrid has not been enabled");
188 #endif
189  } else if (in.Precision() == QUDA_SINGLE_PRECISION) {
190  copyGaugeMG(out, in, location, (short*)Out, (float*)In, (short**)ghostOut, (float**)ghostIn, type);
191  } else if (in.Precision() == QUDA_HALF_PRECISION) {
192  copyGaugeMG(out, in, location, (short*)Out, (short*)In, (short**)ghostOut, (short**)ghostIn, type);
193  } else {
194  errorQuda("Precision %d not supported", in.Precision());
195  }
196  } else {
197  errorQuda("Precision %d not supported", out.Precision());
198  }
199  }
200 
201 
202 
203 } // namespace quda
#define errorQuda(...)
Definition: util_quda.h:121
double Scale() const
double abs_max(int dim=-1) const
Compute the absolute maximum of the field (Linfinity norm)
Definition: max_gauge.cu:84
int length[]
int Ncolor() const
Definition: gauge_field.h:249
void copyGaugeMG(const InOrder &inOrder, GaugeField &out, const GaugeField &in, QudaFieldLocation location, sFloatOut *Out, sFloatOut **outGhost, int type)
cpuColorSpinorField * in
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)
QudaReconstructType Reconstruct() const
Definition: gauge_field.h:250
QudaGaugeFieldOrder Order() const
Definition: gauge_field.h:251
QudaPrecision Precision() const
bool isNative() const