QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
dslash_twisted_mass.cuh
Go to the documentation of this file.
1 #pragma once
2 
4 
5 namespace quda
6 {
7 
8  template <typename Float, int nColor, QudaReconstructType reconstruct_>
9  struct TwistedMassArg : WilsonArg<Float, nColor, reconstruct_> {
10  typedef typename mapper<Float>::type real;
11  real a;
12  real b;
14  TwistedMassArg(ColorSpinorField &out, const ColorSpinorField &in, const GaugeField &U, double a, double b,
15  const ColorSpinorField &x, int parity, bool dagger, const int *comm_override) :
16  WilsonArg<Float, nColor, reconstruct_>(out, in, U, a, x, parity, dagger, comm_override),
17  a(a),
18  b(dagger ? -b : b) // if dagger flip the twist
19  {
20  }
21  };
22 
28  template <typename Float, int nDim, int nColor, int nParity, bool dagger, KernelType kernel_type, typename Arg>
29  __device__ __host__ inline void twistedMass(Arg &arg, int idx, int parity)
30  {
31  typedef typename mapper<Float>::type real;
33  typedef ColorSpinor<real, nColor, 2> HalfVector;
34 
35  bool active
36  = kernel_type == EXTERIOR_KERNEL_ALL ? false : true; // is thread active (non-trival for fused kernel only)
37  int thread_dim; // which dimension is thread working on (fused kernel only)
38  int coord[nDim];
39  int x_cb = getCoords<nDim, QUDA_4D_PC, kernel_type>(coord, arg, idx, parity, thread_dim);
40 
41  const int my_spinor_parity = nParity == 2 ? parity : 0;
42  Vector out;
43 
44  // defined in dslash_wilson.cuh
45  applyWilson<Float, nDim, nColor, nParity, dagger, kernel_type>(
46  out, arg, coord, x_cb, 0, parity, idx, thread_dim, active);
47 
49  Vector x = arg.x(x_cb, my_spinor_parity);
50  x += arg.b * x.igamma(4);
51  out = x + arg.a * out;
52  } else if (active) {
53  Vector x = arg.out(x_cb, my_spinor_parity);
54  out = x + arg.a * out;
55  }
56 
57  if (kernel_type != EXTERIOR_KERNEL_ALL || active) arg.out(x_cb, my_spinor_parity) = out;
58  }
59 
60  // CPU kernel for applying the twisted-mass operator to a vector
61  template <typename Float, int nDim, int nColor, int nParity, bool dagger, KernelType kernel_type, typename Arg>
63  {
64  for (int parity = 0; parity < nParity; parity++) {
65  // for full fields then set parity from loop else use arg setting
66  parity = nParity == 2 ? parity : arg.parity;
67 
68  for (int x_cb = 0; x_cb < arg.threads; x_cb++) { // 4-d volume
69  twistedMass<Float, nDim, nColor, nParity, dagger, kernel_type>(arg, x_cb, parity);
70  } // 4-d volumeCB
71  } // parity
72  }
73 
74  // GPU Kernel for applying the twisted-mass operator to a vector
75  template <typename Float, int nDim, int nColor, int nParity, bool dagger, bool xpay, KernelType kernel_type, typename Arg>
76  __global__ void twistedMassGPU(Arg arg)
77  {
78  int x_cb = blockIdx.x * blockDim.x + threadIdx.x;
79  if (x_cb >= arg.threads) return;
80 
81  // for full fields set parity from z thread index else use arg setting
82  int parity = nParity == 2 ? blockDim.z * blockIdx.z + threadIdx.z : arg.parity;
83 
84  switch (parity) {
85  case 0: twistedMass<Float, nDim, nColor, nParity, dagger, kernel_type>(arg, x_cb, 0); break;
86  case 1: twistedMass<Float, nDim, nColor, nParity, dagger, kernel_type>(arg, x_cb, 1); break;
87  }
88  }
89 
90 } // namespace quda
KernelType kernel_type
__global__ void twistedMassGPU(Arg arg)
__device__ __host__ void twistedMass(Arg &arg, int idx, int parity)
Apply the twisted-mass dslash out(x) = M*in = a * D * in + (1 + i*b*gamma_5)*x Note this routine only...
TwistedMassArg(ColorSpinorField &out, const ColorSpinorField &in, const GaugeField &U, double a, double b, const ColorSpinorField &x, int parity, bool dagger, const int *comm_override)
Parameter structure for driving the Wilson operator.
const int nColor
Definition: covdev_test.cpp:75
mapper< Float >::type real
void twistedMassCPU(Arg arg)
__host__ __device__ ValueType arg(const complex< ValueType > &z)
Returns the phase angle of z.
VectorXcd Vector