QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
dslash_clover_asym.cu
Go to the documentation of this file.
1 #include <cstdlib>
2 #include <cstdio>
3 #include <string>
4 #include <iostream>
5 
6 #include <color_spinor_field.h>
7 #include <clover_field.h>
8 
9 // these control the Wilson-type actions
10 #ifdef GPU_WILSON_DIRAC
11 //#define DIRECT_ACCESS_LINK
12 //#define DIRECT_ACCESS_WILSON_SPINOR
13 //#define DIRECT_ACCESS_WILSON_ACCUM
14 //#define DIRECT_ACCESS_WILSON_INTER
15 //#define DIRECT_ACCESS_WILSON_PACK_SPINOR
16 //#define DIRECT_ACCESS_CLOVER
17 #endif // GPU_WILSON_DIRAC
18 
19 #include <quda_internal.h>
20 #include <dslash_quda.h>
21 #include <sys/time.h>
22 #include <blas_quda.h>
23 #include <face_quda.h>
24 
25 #include <inline_ptx.h>
26 
27 namespace quda {
28 
29  namespace asym_clover {
30 
31 #undef GPU_STAGGERED_DIRAC
32 #include <dslash_constants.h>
33 #include <dslash_textures.h>
34 #include <dslash_index.cuh>
35 
36  // Enable shared memory dslash for Fermi architecture
37  //#define SHARED_WILSON_DSLASH
38  //#define SHARED_8_BYTE_WORD_SIZE // 8-byte shared memory access
39 
40 #ifdef GPU_CLOVER_DIRAC
41 #define DD_CLOVER 2
42 #include <wilson_dslash_def.h> // Wilson Dslash kernels (including clover)
43 #undef DD_CLOVER
44 #endif
45 
46 #ifndef DSLASH_SHARED_FLOATS_PER_THREAD
47 #define DSLASH_SHARED_FLOATS_PER_THREAD 0
48 #endif
49 
50 #include <dslash_quda.cuh>
51 
52  } // end namespace asym_clover
53 
54  // declare the dslash events
55 #include <dslash_events.cuh>
56 
57  using namespace asym_clover;
58 
59 #ifdef GPU_CLOVER_DIRAC
60  template <typename sFloat, typename gFloat, typename cFloat>
61  class AsymCloverDslashCuda : public SharedDslashCuda {
62 
63  private:
64  const gFloat *gauge0, *gauge1;
65  const cFloat *clover;
66  const float *cloverNorm;
67  const double a;
68 
69  protected:
70  unsigned int sharedBytesPerThread() const
71  {
72 #if (__COMPUTE_CAPABILITY__ >= 200)
73  if (dslashParam.kernel_type == INTERIOR_KERNEL) {
74  int reg_size = (typeid(sFloat)==typeid(double2) ? sizeof(double) : sizeof(float));
75  return DSLASH_SHARED_FLOATS_PER_THREAD * reg_size;
76  } else {
77  return 0;
78  }
79 #else
80  int reg_size = (typeid(sFloat)==typeid(double2) ? sizeof(double) : sizeof(float));
81  return DSLASH_SHARED_FLOATS_PER_THREAD * reg_size;
82 #endif
83  }
84 
85  public:
86  AsymCloverDslashCuda(cudaColorSpinorField *out, const gFloat *gauge0, const gFloat *gauge1,
87  const QudaReconstructType reconstruct, const cFloat *clover,
88  const float *cloverNorm, int cl_stride, const cudaColorSpinorField *in,
89  const cudaColorSpinorField *x, const double a, const int dagger)
90  : SharedDslashCuda(out, in, x, reconstruct, dagger), gauge0(gauge0), gauge1(gauge1), clover(clover),
91  cloverNorm(cloverNorm), a(a)
92  {
93  bindSpinorTex<sFloat>(in, out, x);
94  dslashParam.cl_stride = cl_stride;
95  if (!x) errorQuda("Asymmetric clover dslash only defined for Xpay");
96  }
97 
98  virtual ~AsymCloverDslashCuda() { unbindSpinorTex<sFloat>(in, out, x); }
99 
100  void apply(const cudaStream_t &stream)
101  {
102 #ifdef SHARED_WILSON_DSLASH
103  if (dslashParam.kernel_type == EXTERIOR_KERNEL_X)
104  errorQuda("Shared dslash does not yet support X-dimension partitioning");
105 #endif
106  TuneParam tp = tuneLaunch(*this, getTuning(), getVerbosity());
107  ASYM_DSLASH(asymCloverDslash, tp.grid, tp.block, tp.shared_bytes, stream, dslashParam,
108  (sFloat*)out->V(), (float*)out->Norm(), gauge0, gauge1, clover, cloverNorm,
109  (sFloat*)in->V(), (float*)in->Norm(), (sFloat*)x, (float*)x->Norm(), a);
110  }
111 
112  long long flops() const { return 1872ll * in->VolumeCB(); } // FIXME for multi-GPU
113  };
114 #endif // GPU_CLOVER_DIRAC
115 
116 #include <dslash_policy.cuh>
117 
119  const cudaColorSpinorField *in, const int parity, const int dagger,
120  const cudaColorSpinorField *x, const double &a, const int *commOverride,
121  TimeProfile &profile, const QudaDslashPolicy &dslashPolicy)
122  {
123  inSpinor = (cudaColorSpinorField*)in; // EVIL
124 
125 #ifdef GPU_CLOVER_DIRAC
126  int Npad = (in->Ncolor()*in->Nspin()*2)/in->FieldOrder(); // SPINOR_HOP in old code
127  for(int i=0;i<4;i++){
128  dslashParam.ghostDim[i] = commDimPartitioned(i); // determines whether to use regular or ghost indexing at boundary
129  dslashParam.ghostOffset[i] = Npad*(in->GhostOffset(i) + in->Stride());
130  dslashParam.ghostNormOffset[i] = in->GhostNormOffset(i) + in->Stride();
131  dslashParam.commDim[i] = (!commOverride[i]) ? 0 : commDimPartitioned(i); // switch off comms if override = 0
132  }
133 
134  void *cloverP, *cloverNormP;
135  QudaPrecision clover_prec = bindCloverTex(cloverInv, parity, &cloverP, &cloverNormP);
136 
137  void *gauge0, *gauge1;
138  bindGaugeTex(gauge, parity, &gauge0, &gauge1);
139 
140  if (in->Precision() != gauge.Precision())
141  errorQuda("Mixing gauge and spinor precision not supported");
142 
143  if (in->Precision() != clover_prec)
144  errorQuda("Mixing clover and spinor precision not supported");
145 
146  DslashCuda *dslash = 0;
147  size_t regSize = sizeof(float);
148 
149  if (in->Precision() == QUDA_DOUBLE_PRECISION) {
150 #if (__COMPUTE_CAPABILITY__ >= 130)
151  dslash = new AsymCloverDslashCuda<double2, double2, double2>
152  (out, (double2*)gauge0, (double2*)gauge1, gauge.Reconstruct(),
153  (double2*)cloverP, (float*)cloverNormP, cloverInv.stride, in, x, a, dagger);
154  regSize = sizeof(double);
155 #else
156  errorQuda("Double precision not supported on this GPU");
157 #endif
158  } else if (in->Precision() == QUDA_SINGLE_PRECISION) {
159  dslash = new AsymCloverDslashCuda<float4, float4, float4>
160  (out, (float4*)gauge0, (float4*)gauge1, gauge.Reconstruct(),
161  (float4*)cloverP, (float*)cloverNormP, cloverInv.stride, in, x, a, dagger);
162  } else if (in->Precision() == QUDA_HALF_PRECISION) {
163  dslash = new AsymCloverDslashCuda<short4, short4, short4>
164  (out, (short4*)gauge0, (short4*)gauge1, gauge.Reconstruct(),
165  (short4*)cloverP, (float*)cloverNormP, cloverInv.stride, in, x, a, dagger);
166  }
167 
168 #ifndef GPU_COMMS
169  DslashPolicyImp* dslashImp = DslashFactory::create(dslashPolicy);
170 #else
171  DslashPolicyImp* dslashImp = DslashFactory::create(QUDA_GPU_COMMS_DSLASH);
172 #endif
173  (*dslashImp)(*dslash, const_cast<cudaColorSpinorField*>(in), regSize, parity, dagger, in->Volume(), in->GhostFace(), profile);
174  delete dslashImp;
175 
176  delete dslash;
177  unbindGaugeTex(gauge);
178  unbindCloverTex(cloverInv);
179 
180  checkCudaError();
181 #else
182  errorQuda("Clover dslash has not been built");
183 #endif
184 
185  }
186 
187 }
void unbindGaugeTex(const cudaGaugeField &gauge)
enum QudaPrecision_s QudaPrecision
int commDimPartitioned(int dir)
QudaVerbosity getVerbosity()
Definition: util_quda.cpp:20
#define errorQuda(...)
Definition: util_quda.h:73
#define DSLASH_SHARED_FLOATS_PER_THREAD
cudaStream_t * stream
int GhostNormOffset(const int i) const
QudaDagType dagger
Definition: test_util.cpp:1558
QudaPrecision Precision() const
cpuColorSpinorField * in
enum QudaDslashPolicy_s QudaDslashPolicy
TuneParam & tuneLaunch(Tunable &tunable, QudaTune enabled, QudaVerbosity verbosity)
Definition: tune.cpp:271
QudaReconstructType Reconstruct() const
Definition: gauge_field.h:168
void unbindCloverTex(const FullClover clover)
void asymCloverDslashCuda(cudaColorSpinorField *out, const cudaGaugeField &gauge, const FullClover cloverInv, const cudaColorSpinorField *in, const int oddBit, const int daggerBit, const cudaColorSpinorField *x, const double &k, const int *commDim, TimeProfile &profile, const QudaDslashPolicy &dslashPolicy=QUDA_DSLASH2)
int x[4]
QudaFieldOrder FieldOrder() const
cpuColorSpinorField * out
enum QudaReconstructType_s QudaReconstructType
QudaPrecision Precision() const
void bindGaugeTex(const cudaGaugeField &gauge, const int oddBit, void **gauge0, void **gauge1)
QudaPrecision bindCloverTex(const FullClover clover, const int oddBit, void **cloverP, void **cloverNormP)
#define checkCudaError()
Definition: util_quda.h:110
QudaTune getTuning()
Definition: util_quda.cpp:32
int GhostOffset(const int i) const
const QudaParity parity
Definition: dslash_test.cpp:29
void * gauge[4]
Definition: su3_test.cpp:15
const int * GhostFace() const