10 template <
typename Float,
int nColor, QudaReconstructType reconstruct_,
bool twist_ = false>
14 static constexpr
bool twist = twist_;
24 WilsonArg<Float,
nColor, reconstruct_>(out, in, U, a, x, parity, dagger, comm_override),
27 b(dagger ? -0.5 * b : 0.5 * b)
37 template <
typename Float,
int nDim,
int nColor,
int nParity,
bool dagger, KernelType kernel_type,
typename Arg>
48 int x_cb = getCoords<nDim, QUDA_4D_PC, kernel_type>(coord,
arg, idx,
parity, thread_dim);
54 applyWilson<Float, nDim, nColor, nParity, dagger, kernel_type>(
55 out,
arg, coord, x_cb, 0,
parity, idx, thread_dim, active);
58 Vector
x = arg.x(x_cb, my_spinor_parity);
64 for (
int chirality = 0; chirality < 2; chirality++) {
65 constexpr
int n =
nColor * Arg::nSpin / 2;
67 HalfVector x_chi = x.chiral_project(chirality);
68 HalfVector Ax_chi = A * x_chi;
70 const complex<real>
b(0.0, chirality == 0 ? static_cast<real>(arg.b) : -static_cast<real>(arg.b));
73 tmp += Ax_chi.chiral_reconstruct(chirality);
78 out = tmp + arg.a *
out;
80 Vector
x = arg.out(x_cb, my_spinor_parity);
81 out = x + arg.a *
out;
88 template <
typename Float,
int nDim,
int nColor,
int nParity,
bool dagger,
bool xpay, KernelType kernel_type,
typename Arg>
95 for (
int x_cb = 0; x_cb < arg.threads; x_cb++) {
96 wilsonClover<Float, nDim, nColor, nParity, dagger, kernel_type>(
arg, x_cb,
parity);
102 template <
typename Float,
int nDim,
int nColor,
int nParity,
bool dagger,
bool xpay, KernelType kernel_type,
typename Arg>
105 int x_cb = blockIdx.x * blockDim.x + threadIdx.x;
106 if (x_cb >= arg.threads)
return;
109 int parity =
nParity == 2 ? blockDim.z * blockIdx.z + threadIdx.z : arg.parity;
112 case 0: wilsonClover<Float, nDim, nColor, nParity, dagger, kernel_type>(
arg, x_cb, 0);
break;
113 case 1: wilsonClover<Float, nDim, nColor, nParity, dagger, kernel_type>(
arg, x_cb, 1);
break;
WilsonCloverArg(ColorSpinorField &out, const ColorSpinorField &in, const GaugeField &U, const CloverField &A, double a, double b, const ColorSpinorField &x, int parity, bool dagger, const int *comm_override)
clover_mapper< Float, length, true >::type C
cudaColorSpinorField * tmp
void wilsonCloverCPU(Arg arg)
Main header file for host and device accessors to CloverFields.
__global__ void wilsonCloverGPU(Arg arg)
Parameter structure for driving the Wilson operator.
__device__ __host__ void wilsonClover(Arg &arg, int idx, int parity)
Apply the Wilson-clover dslash out(x) = M*in = A(x)*x(x) + D * in(x-mu) Note this routine only exists...
Specialized container for Hermitian matrices (e.g., used for wrapping clover matrices) ...
static constexpr int length
mapper< Float >::type real
static constexpr bool twist
__host__ __device__ ValueType arg(const complex< ValueType > &z)
Returns the phase angle of z.
static constexpr int nSpin