QUDA  0.9.0
gauge_field.cpp
Go to the documentation of this file.
1 #include <gauge_field.h>
2 #include <typeinfo>
3 #include <blas_quda.h>
4 
5 namespace quda {
6 
8  nColor(u.Ncolor()),
9  nFace(u.Nface()),
10  reconstruct(u.Reconstruct()),
11  order(u.Order()),
12  fixed(u.GaugeFixed()),
13  link_type(u.LinkType()),
14  t_boundary(u.TBoundary()),
15  anisotropy(u.Anisotropy()),
16  tadpole(u.Tadpole()),
17  scale(u.Scale()),
18  gauge(NULL),
19  create(QUDA_NULL_FIELD_CREATE),
20  geometry(u.Geometry()),
21  compute_fat_link_max(false),
22  staggeredPhaseType(u.StaggeredPhase()),
23  staggeredPhaseApplied(u.StaggeredPhaseApplied()),
24  i_mu(u.iMu()),
25  site_offset(u.SiteOffset()),
26  site_size(u.SiteSize())
27  { }
28 
29 
31  LatticeField(param), bytes(0), phase_offset(0), phase_bytes(0), nColor(param.nColor), nFace(param.nFace),
32  geometry(param.geometry), reconstruct(param.reconstruct),
33  nInternal(reconstruct != QUDA_RECONSTRUCT_NO ? reconstruct : nColor * nColor * 2),
34  order(param.order), fixed(param.fixed), link_type(param.link_type), t_boundary(param.t_boundary),
35  anisotropy(param.anisotropy), tadpole(param.tadpole), fat_link_max(0.0), scale(param.scale),
36  create(param.create),
37  staggeredPhaseType(param.staggeredPhaseType), staggeredPhaseApplied(param.staggeredPhaseApplied), i_mu(param.i_mu),
38  site_offset(param.site_offset), site_size(param.site_size)
39  {
40  if (link_type != QUDA_COARSE_LINKS && nColor != 3)
41  errorQuda("nColor must be 3, not %d for this link type", nColor);
42  if (nDim != 4)
43  errorQuda("Number of dimensions must be 4 not %d", nDim);
44  if (link_type != QUDA_WILSON_LINKS && anisotropy != 1.0)
45  errorQuda("Anisotropy only supported for Wilson links");
47  errorQuda("Temporal gauge fixing only supported for Wilson links");
48 
50  errorQuda("reconstruct %d only supported for staggered long links\n", reconstruct);
51 
53 
56  length = 2*stride*nInternal; // two comes from being full lattice
57  } else if (geometry == QUDA_VECTOR_GEOMETRY) {
59  length = 2*nDim*stride*nInternal; // two comes from being full lattice
60  } else if(geometry == QUDA_TENSOR_GEOMETRY){
62  length = 2*(nDim*(nDim-1)/2)*stride*nInternal; // two comes from being full lattice
63  } else if(geometry == QUDA_COARSE_GEOMETRY){
65  length = 2*2*nDim*stride*nInternal; //two comes from being full lattice
66  }
67 
69  // Need to adjust the phase alignment as well.
70  int half_phase_bytes = ((size_t)length/(2*reconstruct))*precision; // number of bytes needed to store phases for a single parity
71  int half_gauge_bytes = ((size_t)length/2)*precision - half_phase_bytes; // number of bytes needed to store the gauge field for a single parity excluding the phases
72  // Adjust the alignments for the gauge and phase separately
73  half_phase_bytes = ((half_phase_bytes + (512-1))/512)*512;
74  half_gauge_bytes = ((half_gauge_bytes + (512-1))/512)*512;
75 
76  phase_offset = half_gauge_bytes;
77  phase_bytes = half_phase_bytes*2;
78  bytes = (half_gauge_bytes + half_phase_bytes)*2;
79  } else {
81  if (isNative()) bytes = 2*ALIGNMENT_ADJUST(bytes/2);
82  }
84  }
85 
87 
88  }
89 
90  void GaugeField::createGhostZone(const int *R, bool no_comms_fill) const
91  {
92  if (typeid(*this) == typeid(cpuGaugeField)) return;
93 
94  // calculate size of ghost zone required
95  ghost_bytes = 0;
96  for (int i=0; i<nDim; i++) {
97  ghost_face_bytes[i] = 0;
98  if ( !(comm_dim_partitioned(i) || (no_comms_fill && R[i])) ) ghostFace[i] = 0;
99  else ghostFace[i] = surface[i] * R[i]; // includes the radius (unlike ColorSpinorField)
100 
101  ghostOffset[i][0] = (i == 0) ? 0 : ghostOffset[i-1][1] + ghostFace[i-1]*geometry*nInternal;
103 
105  ghost_bytes += 2*ghost_face_bytes[i]; // factor of two from direction
106  }
107 
109 
110  } // createGhostZone
111 
113  if (staggeredPhaseApplied) errorQuda("Staggered phases already applied");
114  applyGaugePhase(*this);
116  if (typeid(*this)==typeid(cudaGaugeField)) {
117  static_cast<cudaGaugeField&>(*this).exchangeGhost();
118  } else {
119  static_cast<cpuGaugeField&>(*this).exchangeGhost();
120  }
121  }
122  staggeredPhaseApplied = true;
123  }
124 
126  if (!staggeredPhaseApplied) errorQuda("No staggered phases to remove");
127  applyGaugePhase(*this);
129  if (typeid(*this)==typeid(cudaGaugeField)) {
130  static_cast<cudaGaugeField&>(*this).exchangeGhost();
131  } else {
132  static_cast<cpuGaugeField&>(*this).exchangeGhost();
133  }
134  }
135  staggeredPhaseApplied = false;
136  }
137 
138  bool GaugeField::isNative() const {
140  if (order == QUDA_FLOAT2_GAUGE_ORDER) return true;
141  } else if (precision == QUDA_SINGLE_PRECISION ||
144  if (order == QUDA_FLOAT2_GAUGE_ORDER) return true;
146  if (order == QUDA_FLOAT4_GAUGE_ORDER) return true;
148  if (order == QUDA_FLOAT4_GAUGE_ORDER) return true;
149  } else if (reconstruct == QUDA_RECONSTRUCT_10) {
150  if (order == QUDA_FLOAT2_GAUGE_ORDER) return true;
151  }
152  }
153  return false;
154  }
155 
156  void GaugeField::exchange(void **ghost_link, void **link_sendbuf, QudaDirection dir) const {
157  MsgHandle *mh_send[4];
158  MsgHandle *mh_recv[4];
159  size_t bytes[4];
160 
161  for (int i=0; i<nDimComms; i++) bytes[i] = 2*nFace*surfaceCB[i]*nInternal*precision;
162 
163  // in general (standard ghost exchange) we always do the exchange
164  // even if a dimension isn't partitioned. However, this breaks
165  // GaugeField::injectGhost(), so when transferring backwards we
166  // only exchange if a dimension is partitioned. FIXME: this
167  // should probably be cleaned up.
168  bool no_comms_fill = (dir == QUDA_BACKWARDS) ? false : true;
169 
170  void *send[4];
171  void *receive[4];
173  for (int i=0; i<nDimComms; i++) {
174  if (comm_dim_partitioned(i)) {
175  send[i] = link_sendbuf[i];
176  receive[i] = ghost_link[i];
177  } else {
178  if (no_comms_fill) memcpy(ghost_link[i], link_sendbuf[i], bytes[i]);
179  }
180  }
181  } else { // FIXME for CUDA field copy back to the CPU
182  for (int i=0; i<nDimComms; i++) {
183  if (comm_dim_partitioned(i)) {
184  send[i] = pool_pinned_malloc(bytes[i]);
185  receive[i] = pool_pinned_malloc(bytes[i]);
186  qudaMemcpy(send[i], link_sendbuf[i], bytes[i], cudaMemcpyDeviceToHost);
187  } else {
188  if (no_comms_fill) qudaMemcpy(ghost_link[i], link_sendbuf[i], bytes[i], cudaMemcpyDeviceToDevice);
189  }
190  }
191  }
192 
193  for (int i=0; i<nDimComms; i++) {
194  if (!comm_dim_partitioned(i)) continue;
195  if (dir == QUDA_FORWARDS) {
196  mh_send[i] = comm_declare_send_relative(send[i], i, +1, bytes[i]);
197  mh_recv[i] = comm_declare_receive_relative(receive[i], i, -1, bytes[i]);
198  } else if (dir == QUDA_BACKWARDS) {
199  mh_send[i] = comm_declare_send_relative(send[i], i, -1, bytes[i]);
200  mh_recv[i] = comm_declare_receive_relative(receive[i], i, +1, bytes[i]);
201  } else {
202  errorQuda("Unsuported dir=%d", dir);
203  }
204 
205  }
206 
207  for (int i=0; i<nDimComms; i++) {
208  if (!comm_dim_partitioned(i)) continue;
209  comm_start(mh_send[i]);
210  comm_start(mh_recv[i]);
211  }
212 
213  for (int i=0; i<nDimComms; i++) {
214  if (!comm_dim_partitioned(i)) continue;
215  comm_wait(mh_send[i]);
216  comm_wait(mh_recv[i]);
217  }
218 
220  for (int i=0; i<nDimComms; i++) {
221  if (!comm_dim_partitioned(i)) continue;
222  qudaMemcpy(ghost_link[i], receive[i], bytes[i], cudaMemcpyHostToDevice);
223  pool_pinned_free(send[i]);
224  pool_pinned_free(receive[i]);
225  }
226  }
227 
228  for (int i=0; i<nDimComms; i++) {
229  if (!comm_dim_partitioned(i)) continue;
230  comm_free(mh_send[i]);
231  comm_free(mh_recv[i]);
232  }
233 
234  }
235 
236  void GaugeField::checkField(const LatticeField &l) const {
238  try {
239  const GaugeField &g = dynamic_cast<const GaugeField&>(l);
240  if (g.link_type != link_type) errorQuda("link_type does not match %d %d", link_type, g.link_type);
241  if (g.nColor != nColor) errorQuda("nColor does not match %d %d", nColor, g.nColor);
242  if (g.nFace != nFace) errorQuda("nFace does not match %d %d", nFace, g.nFace);
243  if (g.fixed != fixed) errorQuda("fixed does not match %d %d", fixed, g.fixed);
244  if (g.t_boundary != t_boundary) errorQuda("t_boundary does not match %d %d", t_boundary, g.t_boundary);
245  if (g.anisotropy != anisotropy) errorQuda("anisotropy does not match %e %e", anisotropy, g.anisotropy);
246  if (g.tadpole != tadpole) errorQuda("tadpole does not match %e %e", tadpole, g.tadpole);
247  //if (a.scale != scale) errorQuda("scale does not match %e %e", scale, a.scale);
248  }
249  catch(std::bad_cast &e) {
250  errorQuda("Failed to cast reference to GaugeField");
251  }
252  }
253 
254  std::ostream& operator<<(std::ostream& output, const GaugeFieldParam& param) {
255  output << static_cast<const LatticeFieldParam &>(param);
256  output << "nColor = " << param.nColor << std::endl;
257  output << "nFace = " << param.nFace << std::endl;
258  output << "reconstruct = " << param.reconstruct << std::endl;
259  int nInternal = (param.reconstruct != QUDA_RECONSTRUCT_NO ?
260  param.reconstruct : param.nColor * param.nColor * 2);
261  output << "nInternal = " << nInternal << std::endl;
262  output << "order = " << param.order << std::endl;
263  output << "fixed = " << param.fixed << std::endl;
264  output << "link_type = " << param.link_type << std::endl;
265  output << "t_boundary = " << param.t_boundary << std::endl;
266  output << "anisotropy = " << param.anisotropy << std::endl;
267  output << "tadpole = " << param.tadpole << std::endl;
268  output << "scale = " << param.scale << std::endl;
269  output << "create = " << param.create << std::endl;
270  output << "geometry = " << param.geometry << std::endl;
271  output << "staggeredPhaseType = " << param.staggeredPhaseType << std::endl;
272  output << "staggeredPhaseApplied = " << param.staggeredPhaseApplied << std::endl;
273 
274  return output; // for multiple << operators.
275  }
276 
278  if (a.FieldOrder() == QUDA_QDP_GAUGE_ORDER || a.FieldOrder() == QUDA_QDPJIT_GAUGE_ORDER)
279  errorQuda("Not implemented for this order %d", a.FieldOrder());
280 
281  if (a.LinkType() == QUDA_COARSE_LINKS) errorQuda("Not implemented for coarse-link type");
282  if (a.Ncolor() != 3) errorQuda("Not implemented for Ncolor = %d", a.Ncolor());
283 
284  if (a.Precision() == QUDA_HALF_PRECISION)
285  errorQuda("Casting a GaugeField into ColorSpinorField not possible in half precision");
286 
287  ColorSpinorParam spinor_param;
288  spinor_param.nColor = (a.Geometry()*a.Reconstruct())/2;
289  spinor_param.nSpin = 1;
290  spinor_param.nDim = a.Ndim();
291  for (int d=0; d<a.Ndim(); d++) spinor_param.x[d] = a.X()[d];
292  spinor_param.precision = a.Precision();
293  spinor_param.pad = a.Pad();
294  spinor_param.siteSubset = QUDA_FULL_SITE_SUBSET;
295  spinor_param.siteOrder = QUDA_EVEN_ODD_SITE_ORDER;
296  spinor_param.fieldOrder = (a.Precision() == QUDA_DOUBLE_PRECISION || spinor_param.nSpin == 1) ?
298  spinor_param.gammaBasis = QUDA_UKQCD_GAMMA_BASIS;
299  spinor_param.create = QUDA_REFERENCE_FIELD_CREATE;
300  spinor_param.v = (void*)a.Gauge_p();
301  spinor_param.location = a.Location();
302  return spinor_param;
303  }
304 
305  // Return the L2 norm squared of the gauge field
306  double norm2(const GaugeField &a) {
308  double nrm2 = blas::norm2(*b);
309  delete b;
310  return nrm2;
311  }
312 
313  // Return the L1 norm of the gauge field
314  double norm1(const GaugeField &a) {
316  double nrm1 = blas::norm1(*b);
317  delete b;
318  return nrm1;
319  }
320 
321  // Return the L2 norm squared of the gauge field
322  void ax(const double &a, GaugeField &u) {
324  blas::ax(a, *b);
325  delete b;
326  }
327 
328  uint64_t GaugeField::checksum(bool mini) const {
329  return Checksum(*this, mini);
330  }
331 
332 
333 
334 } // namespace quda
#define qudaMemcpy(dst, src, count, kind)
Definition: quda_cuda_api.h:32
double anisotropy
Definition: quda.h:31
#define pool_pinned_free(ptr)
Definition: malloc_quda.h:116
virtual ~GaugeField()
Definition: gauge_field.cpp:86
#define errorQuda(...)
Definition: util_quda.h:90
double norm2(const ColorSpinorField &a)
Definition: reduce_quda.cu:241
QudaReconstructType reconstruct
Definition: gauge_field.h:135
bool staggeredPhaseApplied
Definition: gauge_field.h:161
uint64_t checksum(bool mini=false) const
static ColorSpinorField * Create(const ColorSpinorParam &param)
uint64_t Checksum(const GaugeField &u, bool mini=false)
Definition: checksum.cu:34
QudaPrecision precision
Definition: lattice_field.h:54
static int R[4]
void ax(const double &a, ColorSpinorField &x)
Definition: blas_quda.cu:209
double anisotropy
Definition: test_util.cpp:1644
void applyGaugePhase(GaugeField &u)
Definition: gauge_phase.cu:244
GaugeField(const GaugeFieldParam &param)
Definition: gauge_field.cpp:30
QudaSiteSubset siteSubset
Definition: lattice_field.h:55
std::ostream & operator<<(std::ostream &output, const CloverFieldParam &param)
double norm2(const CloverField &a, bool inverse=false)
QudaFieldGeometry geometry
Definition: gauge_field.h:133
QudaGaugeParam param
Definition: pack_test.cpp:17
#define b
void comm_free(MsgHandle *mh)
Definition: comm_mpi.cpp:252
int x[QUDA_MAX_DIM]
Definition: lattice_field.h:50
double norm1(const CloverField &u, bool inverse=false)
const int * R() const
enum QudaDirection_s QudaDirection
double scale
Definition: quda.h:33
void checkField(const LatticeField &) const
void ax(const double &a, GaugeField &u)
Scale the gauge field by the scalar a.
QudaFieldLocation location
#define comm_declare_send_relative(buffer, dim, dir, nbytes)
Definition: comm_quda.h:59
void checkField(const LatticeField &a) const
#define comm_declare_receive_relative(buffer, dim, dir, nbytes)
Definition: comm_quda.h:74
const int nColor
Definition: covdev_test.cpp:77
#define ALIGNMENT_ADJUST(n)
Definition: quda_internal.h:32
void comm_start(MsgHandle *mh)
Definition: comm_mpi.cpp:260
QudaGhostExchange ghostExchange
unsigned long long uint64_t
void exchange(void **recv, void **send, QudaDirection dir) const
Exchange the buffers across all dimensions in a given direction.
void exchangeGhost(QudaLinkDirection link_direction=QUDA_LINK_BACKWARDS)
Exchange the ghost and store store in the padded region.
QudaReconstructType reconstruct
Definition: quda.h:43
long unsigned int size_t
size_t ghost_face_bytes[QUDA_MAX_DIM]
void * memcpy(void *__dst, const void *__src, size_t __n)
QudaFieldLocation Location() const
int surface[QUDA_MAX_DIM]
#define pool_pinned_malloc(size)
Definition: malloc_quda.h:115
void exchangeGhost(QudaLinkDirection link_direction=QUDA_LINK_BACKWARDS)
Exchange the ghost and store store in the padded region.
void createGhostZone(const int *R, bool no_comms_fill) const
Definition: gauge_field.cpp:90
int ghostOffset[QUDA_MAX_DIM][2]
double norm1(const ColorSpinorField &b)
Definition: reduce_quda.cu:200
int Ncolor
Definition: blas_test.cu:46
QudaLinkType link_type
Definition: gauge_field.h:139
QudaTboundary t_boundary
Definition: gauge_field.h:140
GaugeFieldParam(void *const h_gauge=NULL)
Definition: gauge_field.h:50
void applyStaggeredPhase()
QudaTboundary t_boundary
Definition: quda.h:38
ColorSpinorParam colorSpinorParam(const CloverField &a, bool inverse)
int surfaceCB[QUDA_MAX_DIM]
void comm_wait(MsgHandle *mh)
Definition: comm_mpi.cpp:266
QudaGaugeFieldOrder order
Definition: gauge_field.h:137
static __inline__ size_t size_t d
bool isNative() const
QudaPrecision precision
#define a
QudaGaugeFixed fixed
Definition: gauge_field.h:138
float fat_link_max
unsigned long long bytes
Definition: blas_quda.cu:43
int comm_dim_partitioned(int dim)
int ghostFace[QUDA_MAX_DIM]
Definition: gauge_field.h:151
void removeStaggeredPhase()