QUDA  v1.1.0
A library for QCD on GPUs
copy_clover_offset.cu
Go to the documentation of this file.
1 #include <clover_field_order.h>
2 #include <instantiate.h>
3 #include <kernels/copy_field_offset.cuh>
4 
5 namespace quda
6 {
7 
8  using namespace clover;
9 
10  template <typename Float> struct CopyCloverOffset {
11  CopyCloverOffset(CloverField &out, const CloverField &in, CommKey offset, bool inverse)
12  {
13  constexpr int length = 72;
14  using Field = CloverField;
15  using real = typename mapper<Float>::type;
16  using Element = typename mapper<Float>::type;
17 
18  if (in.isNative()) {
19  using C = typename clover_mapper<Float>::type;
20  C out_accessor(out, inverse);
21  C in_accessor(in, inverse);
22  CopyFieldOffsetArg<Field, Element, C> arg(out_accessor, out, in_accessor, in, offset);
23  CopyFieldOffset<decltype(arg)> copier(arg, in);
24  } else if (in.Order() == QUDA_PACKED_CLOVER_ORDER) {
25  using C = QDPOrder<Float, length>;
26  C out_accessor(out, inverse);
27  C in_accessor(in, inverse);
28  CopyFieldOffsetArg<Field, Element, C> arg(out_accessor, out, in_accessor, in, offset);
29  CopyFieldOffset<decltype(arg)> copier(arg, in);
30  } else if (in.Order() == QUDA_QDPJIT_CLOVER_ORDER) {
31 #ifdef BUILD_QDPJIT_INTERFACE
32  using C = QDPJITOrder<Float, length>;
33  C out_accessor(out, inverse);
34  C in_accessor(in, inverse);
35  CopyFieldOffsetArg<Field, Element, C> arg(out_accessor, out, in_accessor, in, offset);
36  CopyFieldOffset<decltype(arg)> copier(arg, in);
37 #else
38  errorQuda("QDPJIT interface has not been built\n");
39 #endif
40  } else if (in.Order() == QUDA_BQCD_CLOVER_ORDER) {
41 #ifdef BUILD_BQCD_INTERFACE
42  using C = BQCDOrder<Float, length>;
43  C out_accessor(out, inverse);
44  C in_accessor(in, inverse);
45  CopyFieldOffsetArg<Field, Element, C> arg(out_accessor, out, in_accessor, in, offset);
46  CopyFieldOffset<decltype(arg)> copier(arg, in);
47 #else
48  errorQuda("BQCD interface has not been built\n");
49 #endif
50  } else {
51  errorQuda("Clover field %d order not supported", in.Order());
52  }
53  }
54  };
55 
56  void copyFieldOffset(CloverField &out, const CloverField &in, CommKey offset, QudaPCType pc_type)
57  {
58 #ifdef GPU_CLOVER_DIRAC
59  if (out.Precision() < QUDA_SINGLE_PRECISION && out.Order() > 4) {
60  errorQuda("Fixed-point precision not supported for order %d", out.Order());
61  }
62 
63  if (in.Precision() < QUDA_SINGLE_PRECISION && in.Order() > 4) {
64  errorQuda("Fixed-point precision not supported for order %d", in.Order());
65  }
66 
67  if (out.Precision() != in.Precision()) {
68  errorQuda("Precision mismatch: %d (out) vs %d (in)", out.Precision(), in.Precision());
69  }
70  if (out.Order() != in.Order()) { errorQuda("Order mismatch: %d (out) vs %d (in)", out.Order(), in.Order()); }
71 
72  if (pc_type != QUDA_4D_PC) { errorQuda("Gauge field copy must use 4d even-odd preconditioning."); }
73 
74  if (in.V(true)) { instantiate<CopyCloverOffset>(out, in, offset, true); }
75  if (in.V(false)) { instantiate<CopyCloverOffset>(out, in, offset, false); }
76 #else
77  errorQuda("Clover has not been built");
78 #endif
79  }
80 
81 } // namespace quda