QUDA  0.9.0
wilson_dslash_reference.cpp
Go to the documentation of this file.
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include <math.h>
4 
5 #include <util_quda.h>
6 
7 
8 #include <test_util.h>
9 #include <blas_reference.h>
11 
12 #include <gauge_field.h>
13 #include <color_spinor_field.h>
14 
15 #include <dslash_util.h>
16 #include <string.h>
17 
18 using namespace quda;
19 
20 static const double projector[8][4][4][2] = {
21  {
22  {{1,0}, {0,0}, {0,0}, {0,-1}},
23  {{0,0}, {1,0}, {0,-1}, {0,0}},
24  {{0,0}, {0,1}, {1,0}, {0,0}},
25  {{0,1}, {0,0}, {0,0}, {1,0}}
26  },
27  {
28  {{1,0}, {0,0}, {0,0}, {0,1}},
29  {{0,0}, {1,0}, {0,1}, {0,0}},
30  {{0,0}, {0,-1}, {1,0}, {0,0}},
31  {{0,-1}, {0,0}, {0,0}, {1,0}}
32  },
33  {
34  {{1,0}, {0,0}, {0,0}, {1,0}},
35  {{0,0}, {1,0}, {-1,0}, {0,0}},
36  {{0,0}, {-1,0}, {1,0}, {0,0}},
37  {{1,0}, {0,0}, {0,0}, {1,0}}
38  },
39  {
40  {{1,0}, {0,0}, {0,0}, {-1,0}},
41  {{0,0}, {1,0}, {1,0}, {0,0}},
42  {{0,0}, {1,0}, {1,0}, {0,0}},
43  {{-1,0}, {0,0}, {0,0}, {1,0}}
44  },
45  {
46  {{1,0}, {0,0}, {0,-1}, {0,0}},
47  {{0,0}, {1,0}, {0,0}, {0,1}},
48  {{0,1}, {0,0}, {1,0}, {0,0}},
49  {{0,0}, {0,-1}, {0,0}, {1,0}}
50  },
51  {
52  {{1,0}, {0,0}, {0,1}, {0,0}},
53  {{0,0}, {1,0}, {0,0}, {0,-1}},
54  {{0,-1}, {0,0}, {1,0}, {0,0}},
55  {{0,0}, {0,1}, {0,0}, {1,0}}
56  },
57  {
58  {{1,0}, {0,0}, {-1,0}, {0,0}},
59  {{0,0}, {1,0}, {0,0}, {-1,0}},
60  {{-1,0}, {0,0}, {1,0}, {0,0}},
61  {{0,0}, {-1,0}, {0,0}, {1,0}}
62  },
63  {
64  {{1,0}, {0,0}, {1,0}, {0,0}},
65  {{0,0}, {1,0}, {0,0}, {1,0}},
66  {{1,0}, {0,0}, {1,0}, {0,0}},
67  {{0,0}, {1,0}, {0,0}, {1,0}}
68  }
69 };
70 
71 
72 // todo pass projector
73 template <typename Float>
74 void multiplySpinorByDiracProjector(Float *res, int projIdx, Float *spinorIn) {
75  for (int i=0; i<4*3*2; i++) res[i] = 0.0;
76 
77  for (int s = 0; s < 4; s++) {
78  for (int t = 0; t < 4; t++) {
79  Float projRe = projector[projIdx][s][t][0];
80  Float projIm = projector[projIdx][s][t][1];
81 
82  for (int m = 0; m < 3; m++) {
83  Float spinorRe = spinorIn[t*(3*2) + m*(2) + 0];
84  Float spinorIm = spinorIn[t*(3*2) + m*(2) + 1];
85  res[s*(3*2) + m*(2) + 0] += projRe*spinorRe - projIm*spinorIm;
86  res[s*(3*2) + m*(2) + 1] += projRe*spinorIm + projIm*spinorRe;
87  }
88  }
89  }
90 }
91 
92 
93 //
94 // dslashReference()
95 //
96 // if oddBit is zero: calculate odd parity spinor elements (using even parity spinor)
97 // if oddBit is one: calculate even parity spinor elements
98 //
99 // if daggerBit is zero: perform ordinary dslash operator
100 // if daggerBit is one: perform hermitian conjugate of dslash
101 //
102 
103 #ifndef MULTI_GPU
104 
105 template <typename sFloat, typename gFloat>
106 void dslashReference(sFloat *res, gFloat **gaugeFull, sFloat *spinorField, int oddBit, int daggerBit) {
107  for (int i=0; i<Vh*mySpinorSiteSize; i++) res[i] = 0.0;
108 
109  gFloat *gaugeEven[4], *gaugeOdd[4];
110  for (int dir = 0; dir < 4; dir++) {
111  gaugeEven[dir] = gaugeFull[dir];
112  gaugeOdd[dir] = gaugeFull[dir]+Vh*gaugeSiteSize;
113  }
114 
115  for (int i = 0; i < Vh; i++) {
116  for (int dir = 0; dir < 8; dir++) {
117  gFloat *gauge = gaugeLink(i, dir, oddBit, gaugeEven, gaugeOdd, 1);
118  sFloat *spinor = spinorNeighbor(i, dir, oddBit, spinorField, 1);
119 
120  sFloat projectedSpinor[4*3*2], gaugedSpinor[4*3*2];
121  int projIdx = 2*(dir/2)+(dir+daggerBit)%2;
122  multiplySpinorByDiracProjector(projectedSpinor, projIdx, spinor);
123 
124  for (int s = 0; s < 4; s++) {
125  if (dir % 2 == 0) su3Mul(&gaugedSpinor[s*(3*2)], gauge, &projectedSpinor[s*(3*2)]);
126  else su3Tmul(&gaugedSpinor[s*(3*2)], gauge, &projectedSpinor[s*(3*2)]);
127  }
128 
129  sum(&res[i*(4*3*2)], &res[i*(4*3*2)], gaugedSpinor, 4*3*2);
130  }
131  }
132 }
133 
134 #else
135 
136 template <typename sFloat, typename gFloat>
137 void dslashReference(sFloat *res, gFloat **gaugeFull, gFloat **ghostGauge, sFloat *spinorField,
138  sFloat **fwdSpinor, sFloat **backSpinor, int oddBit, int daggerBit) {
139  for (int i=0; i<Vh*mySpinorSiteSize; i++) res[i] = 0.0;
140 
141  gFloat *gaugeEven[4], *gaugeOdd[4];
142  gFloat *ghostGaugeEven[4], *ghostGaugeOdd[4];
143  for (int dir = 0; dir < 4; dir++) {
144  gaugeEven[dir] = gaugeFull[dir];
145  gaugeOdd[dir] = gaugeFull[dir]+Vh*gaugeSiteSize;
146 
147  ghostGaugeEven[dir] = ghostGauge[dir];
148  ghostGaugeOdd[dir] = ghostGauge[dir] + (faceVolume[dir]/2)*gaugeSiteSize;
149  }
150 
151  for (int i = 0; i < Vh; i++) {
152 
153  for (int dir = 0; dir < 8; dir++) {
154  gFloat *gauge = gaugeLink_mg4dir(i, dir, oddBit, gaugeEven, gaugeOdd, ghostGaugeEven, ghostGaugeOdd, 1, 1);
155  sFloat *spinor = spinorNeighbor_mg4dir(i, dir, oddBit, spinorField, fwdSpinor, backSpinor, 1, 1);
156 
157  sFloat projectedSpinor[mySpinorSiteSize], gaugedSpinor[mySpinorSiteSize];
158  int projIdx = 2*(dir/2)+(dir+daggerBit)%2;
159  multiplySpinorByDiracProjector(projectedSpinor, projIdx, spinor);
160 
161  for (int s = 0; s < 4; s++) {
162  if (dir % 2 == 0) su3Mul(&gaugedSpinor[s*(3*2)], gauge, &projectedSpinor[s*(3*2)]);
163  else su3Tmul(&gaugedSpinor[s*(3*2)], gauge, &projectedSpinor[s*(3*2)]);
164  }
165 
166  sum(&res[i*(4*3*2)], &res[i*(4*3*2)], gaugedSpinor, 4*3*2);
167  }
168 
169  }
170 }
171 
172 #endif
173 
174 // this actually applies the preconditioned dslash, e.g., D_ee^{-1} D_eo or D_oo^{-1} D_oe
175 void wil_dslash(void *out, void **gauge, void *in, int oddBit, int daggerBit,
177 
178 #ifndef MULTI_GPU
179  if (precision == QUDA_DOUBLE_PRECISION)
180  dslashReference((double*)out, (double**)gauge, (double*)in, oddBit, daggerBit);
181  else
182  dslashReference((float*)out, (float**)gauge, (float*)in, oddBit, daggerBit);
183 #else
184 
185  GaugeFieldParam gauge_field_param(gauge, gauge_param);
186  gauge_field_param.ghostExchange = QUDA_GHOST_EXCHANGE_PAD;
187  cpuGaugeField cpu(gauge_field_param);
188  void **ghostGauge = (void**)cpu.Ghost();
189 
190  // Get spinor ghost fields
191  // First wrap the input spinor into a ColorSpinorField
193  csParam.v = in;
194  csParam.nColor = 3;
195  csParam.nSpin = 4;
196  csParam.nDim = 4;
197  for (int d=0; d<4; d++) csParam.x[d] = Z[d];
198  csParam.precision = precision;
199  csParam.pad = 0;
201  csParam.x[0] /= 2;
206 
207  cpuColorSpinorField inField(csParam);
208 
209  { // Now do the exchange
210  QudaParity otherParity = QUDA_INVALID_PARITY;
211  if (oddBit == QUDA_EVEN_PARITY) otherParity = QUDA_ODD_PARITY;
212  else if (oddBit == QUDA_ODD_PARITY) otherParity = QUDA_EVEN_PARITY;
213  else errorQuda("ERROR: full parity not supported in function %s", __FUNCTION__);
214  const int nFace = 1;
215 
216  inField.exchangeGhost(otherParity, nFace, daggerBit);
217  }
218  void** fwd_nbr_spinor = inField.fwdGhostFaceBuffer;
219  void** back_nbr_spinor = inField.backGhostFaceBuffer;
220 
221  if (precision == QUDA_DOUBLE_PRECISION) {
222  dslashReference((double*)out, (double**)gauge, (double**)ghostGauge, (double*)in,
223  (double**)fwd_nbr_spinor, (double**)back_nbr_spinor, oddBit, daggerBit);
224  } else{
225  dslashReference((float*)out, (float**)gauge, (float**)ghostGauge, (float*)in,
226  (float**)fwd_nbr_spinor, (float**)back_nbr_spinor, oddBit, daggerBit);
227  }
228 
229 #endif
230 
231 }
232 
233 // applies b*(1 + i*a*gamma_5)
234 template <typename sFloat>
235 void twistGamma5(sFloat *out, sFloat *in, const int dagger, const sFloat kappa, const sFloat mu,
236  const QudaTwistFlavorType flavor, const int V, QudaTwistGamma5Type twist) {
237 
238  sFloat a=0.0,b=0.0;
239  if (twist == QUDA_TWIST_GAMMA5_DIRECT) { // applying the twist
240  a = 2.0 * kappa * mu * flavor; // mu already includes the flavor
241  b = 1.0;
242  } else if (twist == QUDA_TWIST_GAMMA5_INVERSE) { // applying the inverse twist
243  a = -2.0 * kappa * mu * flavor;
244  b = 1.0 / (1.0 + a*a);
245  } else {
246  printf("Twist type %d not defined\n", twist);
247  exit(0);
248  }
249 
250  if (dagger) a *= -1.0;
251 
252  for(int i = 0; i < V; i++) {
253  sFloat tmp[24];
254  for(int s = 0; s < 4; s++)
255  for(int c = 0; c < 3; c++) {
256  sFloat a5 = ((s / 2) ? -1.0 : +1.0) * a;
257  tmp[s * 6 + c * 2 + 0] = b* (in[i * 24 + s * 6 + c * 2 + 0] - a5*in[i * 24 + s * 6 + c * 2 + 1]);
258  tmp[s * 6 + c * 2 + 1] = b* (in[i * 24 + s * 6 + c * 2 + 1] + a5*in[i * 24 + s * 6 + c * 2 + 0]);
259  }
260 
261  for (int j=0; j<24; j++) out[i*24+j] = tmp[j];
262  }
263 
264 }
265 
266 void twist_gamma5(void *out, void *in, int daggerBit, double kappa, double mu, QudaTwistFlavorType flavor,
267  int V, QudaTwistGamma5Type twist, QudaPrecision precision) {
268 
269  if (precision == QUDA_DOUBLE_PRECISION) {
270  twistGamma5((double*)out, (double*)in, daggerBit, kappa, mu, flavor, V, twist);
271  } else {
272  twistGamma5((float*)out, (float*)in, daggerBit, (float)kappa, (float)mu, flavor, V, twist);
273  }
274 }
275 
276 
277 void tm_dslash(void *res, void **gaugeFull, void *spinorField, double kappa, double mu,
278  QudaTwistFlavorType flavor, int oddBit, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision,
280 {
281 
282  if (daggerBit && (matpc_type == QUDA_MATPC_EVEN_EVEN || matpc_type == QUDA_MATPC_ODD_ODD)) twist_gamma5(spinorField, spinorField, daggerBit, kappa, mu,
283  flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
284 
285  wil_dslash(res, gaugeFull, spinorField, oddBit, daggerBit, precision, gauge_param);
286 
287  if (!daggerBit || (daggerBit && (matpc_type == QUDA_MATPC_EVEN_EVEN_ASYMMETRIC || matpc_type == QUDA_MATPC_ODD_ODD_ASYMMETRIC))) {
288  twist_gamma5(res, res, daggerBit, kappa, mu, flavor,
289  Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
290  } else {
291  twist_gamma5(spinorField, spinorField, daggerBit, kappa, mu, flavor,
292  Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
293  }
294 }
295 
296 void wil_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision,
298 
299  void *inEven = in;
300  void *inOdd = (char*)in + Vh*spinorSiteSize*precision;
301  void *outEven = out;
302  void *outOdd = (char*)out + Vh*spinorSiteSize*precision;
303 
304  wil_dslash(outOdd, gauge, inEven, 1, dagger_bit, precision, gauge_param);
305  wil_dslash(outEven, gauge, inOdd, 0, dagger_bit, precision, gauge_param);
306 
307  // lastly apply the kappa term
308  xpay(in, -kappa, out, V*spinorSiteSize, precision);
309 }
310 
311 void tm_mat(void *out, void **gauge, void *in, double kappa, double mu,
312  QudaTwistFlavorType flavor, int dagger_bit, QudaPrecision precision,
314 
315  void *inEven = in;
316  void *inOdd = (char*)in + Vh*spinorSiteSize*precision;
317  void *outEven = out;
318  void *outOdd = (char*)out + Vh*spinorSiteSize*precision;
319  void *tmp = malloc(V*spinorSiteSize*precision);
320 
321  wil_dslash(outOdd, gauge, inEven, 1, dagger_bit, precision, gauge_param);
322  wil_dslash(outEven, gauge, inOdd, 0, dagger_bit, precision, gauge_param);
323 
324  // apply the twist term to the full lattice
325  twist_gamma5(tmp, in, dagger_bit, kappa, mu, flavor, V, QUDA_TWIST_GAMMA5_DIRECT, precision);
326 
327  // combine
328  xpay(tmp, -kappa, (double*)out, V*spinorSiteSize, precision);
329 
330  free(tmp);
331 }
332 
333 // Apply the even-odd preconditioned Dirac operator
334 void wil_matpc(void *outEven, void **gauge, void *inEven, double kappa,
335  QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision,
337 
338  void *tmp = malloc(Vh*spinorSiteSize*precision);
339 
340  // FIXME: remove once reference clover is finished
341  // full dslash operator
343  wil_dslash(tmp, gauge, inEven, 1, daggerBit, precision, gauge_param);
344  wil_dslash(outEven, gauge, tmp, 0, daggerBit, precision, gauge_param);
345  } else {
346  wil_dslash(tmp, gauge, inEven, 0, daggerBit, precision, gauge_param);
347  wil_dslash(outEven, gauge, tmp, 1, daggerBit, precision, gauge_param);
348  }
349 
350  // lastly apply the kappa term
351  double kappa2 = -kappa*kappa;
352  xpay(inEven, kappa2, outEven, Vh*spinorSiteSize, precision);
353 
354  free(tmp);
355 }
356 
357 // Apply the even-odd preconditioned Dirac operator
358 void tm_matpc(void *outEven, void **gauge, void *inEven, double kappa, double mu, QudaTwistFlavorType flavor,
359  QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param) {
360 
361  void *tmp = malloc(Vh*spinorSiteSize*precision);
362 
364  wil_dslash(tmp, gauge, inEven, 1, daggerBit, precision, gauge_param);
365  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
366  wil_dslash(outEven, gauge, tmp, 0, daggerBit, precision, gauge_param);
367  twist_gamma5(tmp, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
369  wil_dslash(tmp, gauge, inEven, 0, daggerBit, precision, gauge_param);
370  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
371  wil_dslash(outEven, gauge, tmp, 1, daggerBit, precision, gauge_param);
372  twist_gamma5(tmp, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
373  } else if (!daggerBit) {
375  wil_dslash(tmp, gauge, inEven, 1, daggerBit, precision, gauge_param);
376  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
377  wil_dslash(outEven, gauge, tmp, 0, daggerBit, precision, gauge_param);
378  twist_gamma5(outEven, outEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
379  } else if (matpc_type == QUDA_MATPC_ODD_ODD) {
380  wil_dslash(tmp, gauge, inEven, 0, daggerBit, precision, gauge_param);
381  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
382  wil_dslash(outEven, gauge, tmp, 1, daggerBit, precision, gauge_param);
383  twist_gamma5(outEven, outEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
384  }
385  } else {
387  twist_gamma5(inEven, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
388  wil_dslash(tmp, gauge, inEven, 1, daggerBit, precision, gauge_param);
389  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
390  wil_dslash(outEven, gauge, tmp, 0, daggerBit, precision, gauge_param);
391  twist_gamma5(inEven, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
392  } else if (matpc_type == QUDA_MATPC_ODD_ODD) {
393  twist_gamma5(inEven, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
394  wil_dslash(tmp, gauge, inEven, 0, daggerBit, precision, gauge_param);
395  twist_gamma5(tmp, tmp, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
396  wil_dslash(outEven, gauge, tmp, 1, daggerBit, precision, gauge_param);
397  twist_gamma5(inEven, inEven, daggerBit, kappa, mu, flavor, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision); // undo
398  }
399  }
400  // lastly apply the kappa term
401  double kappa2 = -kappa*kappa;
403  xpay(inEven, kappa2, outEven, Vh*spinorSiteSize, precision);
404  } else {
405  xpay(tmp, kappa2, outEven, Vh*spinorSiteSize, precision);
406  }
407 
408  free(tmp);
409 }
410 
411 
412 //----- for non-degenerate dslash only----
413 template <typename sFloat>
414 void ndegTwistGamma5(sFloat *out1, sFloat *out2, sFloat *in1, sFloat *in2, const int dagger, const sFloat kappa, const sFloat mu,
415  const sFloat epsilon, const int V, QudaTwistGamma5Type twist) {
416 
417  sFloat a=0.0, b=0.0, d=0.0;
418  if (twist == QUDA_TWIST_GAMMA5_DIRECT) { // applying the twist
419  a = 2.0 * kappa * mu;
420  b = -2.0 * kappa * epsilon;
421  d = 1.0;
422  } else if (twist == QUDA_TWIST_GAMMA5_INVERSE) { // applying the inverse twist
423  a = -2.0 * kappa * mu;
424  b = 2.0 * kappa * epsilon;
425  d = 1.0 / (1.0 + a*a - b*b);
426  } else {
427  printf("Twist type %d not defined\n", twist);
428  exit(0);
429  }
430 
431  if (dagger) a *= -1.0;
432 
433  for(int i = 0; i < V; i++) {
434  sFloat tmp1[24];
435  sFloat tmp2[24];
436  for(int s = 0; s < 4; s++)
437  for(int c = 0; c < 3; c++) {
438  sFloat a5 = ((s / 2) ? -1.0 : +1.0) * a;
439  tmp1[s * 6 + c * 2 + 0] = d* (in1[i * 24 + s * 6 + c * 2 + 0] - a5*in1[i * 24 + s * 6 + c * 2 + 1] + b*in2[i * 24 + s * 6 + c * 2 + 0]);
440  tmp1[s * 6 + c * 2 + 1] = d* (in1[i * 24 + s * 6 + c * 2 + 1] + a5*in1[i * 24 + s * 6 + c * 2 + 0] + b*in2[i * 24 + s * 6 + c * 2 + 1]);
441  tmp2[s * 6 + c * 2 + 0] = d* (in2[i * 24 + s * 6 + c * 2 + 0] + a5*in2[i * 24 + s * 6 + c * 2 + 1] + b*in1[i * 24 + s * 6 + c * 2 + 0]);
442  tmp2[s * 6 + c * 2 + 1] = d* (in2[i * 24 + s * 6 + c * 2 + 1] - a5*in2[i * 24 + s * 6 + c * 2 + 0] + b*in1[i * 24 + s * 6 + c * 2 + 1]);
443  }
444  for (int j=0; j<24; j++) out1[i*24+j] = tmp1[j], out2[i*24+j] = tmp2[j];
445  }
446 
447 }
448 
449 void ndeg_twist_gamma5(void *outf1, void *outf2, void *inf1, void *inf2, const int dagger, const double kappa, const double mu,
450  const double epsilon, const int Vf, QudaTwistGamma5Type twist, QudaPrecision precision)
451 {
452  if (precision == QUDA_DOUBLE_PRECISION)
453  {
454  ndegTwistGamma5((double*)outf1, (double*)outf2, (double*)inf1, (double*)inf2, dagger, kappa, mu, epsilon, Vf, twist);
455  }
456  else //single precision dslash
457  {
458  ndegTwistGamma5((float*)outf1, (float*)outf2, (float*)inf1, (float*)inf2, dagger, (float)kappa, (float)mu, (float)epsilon, Vf, twist);
459  }
460 }
461 
462 void tm_ndeg_dslash(void *res1, void *res2, void **gauge, void *spinorField1, void *spinorField2, double kappa, double mu,
463  double epsilon, int oddBit, int daggerBit, QudaMatPCType matpc_type, QudaPrecision precision, QudaGaugeParam &gauge_param)
464 {
465  if (daggerBit && (matpc_type == QUDA_MATPC_EVEN_EVEN || matpc_type == QUDA_MATPC_ODD_ODD))
466  ndeg_twist_gamma5(spinorField1, spinorField2, spinorField1, spinorField2, daggerBit, kappa, -mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
467 
468  wil_dslash(res1, gauge, spinorField1, oddBit, daggerBit, precision, gauge_param);
469  wil_dslash(res2, gauge, spinorField2, oddBit, daggerBit, precision, gauge_param);
470 
472  ndeg_twist_gamma5(res1, res2, res1, res2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
473  }
474 }
475 
476 
477 void tm_ndeg_matpc(void *outEven1, void *outEven2, void **gauge, void *inEven1, void *inEven2, double kappa, double mu, double epsilon,
478  QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param) {
479 
480  void *tmp1 = malloc(Vh*spinorSiteSize*precision);
481  void *tmp2 = malloc(Vh*spinorSiteSize*precision);
482 
483  if (!daggerBit) {
485  wil_dslash(tmp1, gauge, inEven1, 1, daggerBit, precision, gauge_param);
486  wil_dslash(tmp2, gauge, inEven2, 1, daggerBit, precision, gauge_param);
487  ndeg_twist_gamma5(tmp1, tmp2, tmp1, tmp2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
488  wil_dslash(outEven1, gauge, tmp1, 0, daggerBit, precision, gauge_param);
489  wil_dslash(outEven2, gauge, tmp2, 0, daggerBit, precision, gauge_param);
490  ndeg_twist_gamma5(outEven1, outEven2, outEven1, outEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
491  } else if (matpc_type == QUDA_MATPC_ODD_ODD) {
492  wil_dslash(tmp1, gauge, inEven1, 0, daggerBit, precision, gauge_param);
493  wil_dslash(tmp2, gauge, inEven2, 0, daggerBit, precision, gauge_param);
494  ndeg_twist_gamma5(tmp1, tmp2, tmp1, tmp2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
495  wil_dslash(outEven1, gauge, tmp1, 1, daggerBit, precision, gauge_param);
496  wil_dslash(outEven2, gauge, tmp2, 1, daggerBit, precision, gauge_param);
497  ndeg_twist_gamma5(outEven1, outEven2, outEven1, outEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
498  }
499  } else {
501  ndeg_twist_gamma5(tmp1, tmp2, inEven1, inEven2, daggerBit, kappa, -mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
502  wil_dslash(outEven1, gauge, tmp1, 1, daggerBit, precision, gauge_param);
503  wil_dslash(outEven2, gauge, tmp2, 1, daggerBit, precision, gauge_param);
504  ndeg_twist_gamma5(tmp1, tmp2, outEven1, outEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
505  wil_dslash(outEven1, gauge, tmp1, 0, daggerBit, precision, gauge_param);
506  wil_dslash(outEven2, gauge, tmp2, 0, daggerBit, precision, gauge_param);
507  } else if (matpc_type == QUDA_MATPC_ODD_ODD) {
508  ndeg_twist_gamma5(tmp1, tmp2, inEven1, inEven2, daggerBit, kappa, -mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
509  wil_dslash(outEven1, gauge, tmp1, 0, daggerBit, precision, gauge_param);
510  wil_dslash(outEven2, gauge, tmp2, 0, daggerBit, precision, gauge_param);
511  ndeg_twist_gamma5(tmp1, tmp2, outEven1, outEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
512  wil_dslash(outEven1, gauge, tmp1, 1, daggerBit, precision, gauge_param);
513  wil_dslash(outEven2, gauge, tmp2, 1, daggerBit, precision, gauge_param);
514  }
515  }
516 
518  wil_dslash(tmp1, gauge, inEven1, 1, daggerBit, precision, gauge_param);
519  wil_dslash(tmp2, gauge, inEven2, 1, daggerBit, precision, gauge_param);
520  ndeg_twist_gamma5(tmp1, tmp2, tmp1, tmp2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
521  wil_dslash(outEven1, gauge, tmp1, 0, daggerBit, precision, gauge_param);
522  wil_dslash(outEven2, gauge, tmp2, 0, daggerBit, precision, gauge_param);
524  wil_dslash(tmp1, gauge, inEven1, 0, daggerBit, precision, gauge_param);
525  wil_dslash(tmp2, gauge, inEven2, 0, daggerBit, precision, gauge_param);
526  ndeg_twist_gamma5(tmp1, tmp2, tmp1, tmp2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_INVERSE, precision);
527  wil_dslash(outEven1, gauge, tmp1, 1, daggerBit, precision, gauge_param);
528  wil_dslash(outEven2, gauge, tmp2, 1, daggerBit, precision, gauge_param);
529  }
530 
531  // lastly apply the kappa term
532  double kappa2 = -kappa*kappa;
534  ndeg_twist_gamma5(inEven1, inEven2, inEven1, inEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
535  }
536 
537  xpay(inEven1, kappa2, outEven1, Vh*spinorSiteSize, precision);
538  xpay(inEven2, kappa2, outEven2, Vh*spinorSiteSize, precision);
539 
540  free(tmp1);
541  free(tmp2);
542 }
543 
544 
545 void tm_ndeg_mat(void *evenOut, void* oddOut, void **gauge, void *evenIn, void *oddIn, double kappa, double mu, double epsilon, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
546 {
547  //V-4d volume and Vh=V/2
548  void *inEven1 = evenIn;
549  void *inEven2 = (char*)evenIn + precision*Vh*spinorSiteSize;
550 
551  void *inOdd1 = oddIn;
552  void *inOdd2 = (char*)oddIn + precision*Vh*spinorSiteSize;
553 
554  void *outEven1 = evenOut;
555  void *outEven2 = (char*)evenOut + precision*Vh*spinorSiteSize;
556 
557  void *outOdd1 = oddOut;
558  void *outOdd2 = (char*)oddOut + precision*Vh*spinorSiteSize;
559 
560  void *tmpEven1 = malloc(Vh*spinorSiteSize*precision);
561  void *tmpEven2 = malloc(Vh*spinorSiteSize*precision);
562 
563  void *tmpOdd1 = malloc(Vh*spinorSiteSize*precision);
564  void *tmpOdd2 = malloc(Vh*spinorSiteSize*precision);
565 
566  // full dslash operator:
567  wil_dslash(outOdd1, gauge, inEven1, 1, daggerBit, precision, gauge_param);
568  wil_dslash(outOdd2, gauge, inEven2, 1, daggerBit, precision, gauge_param);
569 
570  wil_dslash(outEven1, gauge, inOdd1, 0, daggerBit, precision, gauge_param);
571  wil_dslash(outEven2, gauge, inOdd2, 0, daggerBit, precision, gauge_param);
572 
573  // apply the twist term
574  ndeg_twist_gamma5(tmpEven1, tmpEven2, inEven1, inEven2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
575  ndeg_twist_gamma5(tmpOdd1, tmpOdd2, inOdd1, inOdd2, daggerBit, kappa, mu, epsilon, Vh, QUDA_TWIST_GAMMA5_DIRECT, precision);
576  // combine
577  xpay(tmpOdd1, -kappa, outOdd1, Vh*spinorSiteSize, precision);
578  xpay(tmpOdd2, -kappa, outOdd2, Vh*spinorSiteSize, precision);
579 
580  xpay(tmpEven1, -kappa, outEven1, Vh*spinorSiteSize, precision);
581  xpay(tmpEven2, -kappa, outEven2, Vh*spinorSiteSize, precision);
582 
583  free(tmpOdd1);
584  free(tmpOdd2);
585  //
586  free(tmpEven1);
587  free(tmpEven2);
588 }
589 
590 //End of nondeg TM
QudaGhostExchange ghostExchange
Definition: lattice_field.h:60
void free(void *)
void xpay(ColorSpinorField &x, const double &a, ColorSpinorField &y)
Definition: blas_quda.cu:173
double mu
Definition: test_util.cpp:1643
enum QudaPrecision_s QudaPrecision
void tm_mat(void *out, void **gauge, void *in, double kappa, double mu, QudaTwistFlavorType flavor, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
void tm_dslash(void *res, void **gaugeFull, void *spinorField, double kappa, double mu, QudaTwistFlavorType flavor, int oddBit, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
#define errorQuda(...)
Definition: util_quda.h:90
cudaColorSpinorField * tmp
Definition: covdev_test.cpp:44
void tm_ndeg_mat(void *evenOut, void *oddOut, void **gauge, void *evenIn, void *oddIn, double kappa, double mu, double epsilon, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
QudaPrecision precision
Definition: lattice_field.h:54
QudaGaugeParam gauge_param
void tm_matpc(void *outEven, void **gauge, void *inEven, double kappa, double mu, QudaTwistFlavorType flavor, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
void twist_gamma5(void *out, void *in, int daggerBit, double kappa, double mu, QudaTwistFlavorType flavor, int V, QudaTwistGamma5Type twist, QudaPrecision precision)
#define spinorSiteSize
QudaSiteSubset siteSubset
Definition: lattice_field.h:55
void exit(int) __attribute__((noreturn))
static Float * spinorNeighbor(int i, int dir, int oddBit, Float *spinorField, int neighbor_distance)
Definition: dslash_util.h:127
#define b
int x[QUDA_MAX_DIM]
Definition: lattice_field.h:50
void * malloc(size_t __size) __attribute__((__warn_unused_result__)) __attribute__((alloc_size(1)))
int printf(const char *,...) __attribute__((__format__(__printf__
VOLATILE spinorFloat kappa
#define tmp2
Definition: tmc_core.h:16
__host__ __device__ void sum(double &a, double &b)
ColorSpinorParam csParam
Definition: pack_test.cpp:24
cpuColorSpinorField * in
#define mySpinorSiteSize
int V
Definition: test_util.cpp:28
enum QudaMatPCType_s QudaMatPCType
#define gaugeSiteSize
Definition: test_util.h:6
void multiplySpinorByDiracProjector(Float *res, int projIdx, Float *spinorIn)
static const double projector[8][4][4][2]
void ndegTwistGamma5(sFloat *out1, sFloat *out2, sFloat *in1, sFloat *in2, const int dagger, const sFloat kappa, const sFloat mu, const sFloat epsilon, const int V, QudaTwistGamma5Type twist)
int Z[4]
Definition: test_util.cpp:27
const void ** Ghost() const
Definition: gauge_field.h:254
#define tmp1
Definition: tmc_core.h:15
void tm_ndeg_matpc(void *outEven1, void *outEven2, void **gauge, void *inEven1, void *inEven2, double kappa, double mu, double epsilon, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
enum QudaParity_s QudaParity
QudaMatPCType matpc_type
Definition: test_util.cpp:1652
void tm_ndeg_dslash(void *res1, void *res2, void **gauge, void *spinorField1, void *spinorField2, double kappa, double mu, double epsilon, int oddBit, int daggerBit, QudaMatPCType matpc_type, QudaPrecision precision, QudaGaugeParam &gauge_param)
void ndeg_twist_gamma5(void *outf1, void *outf2, void *inf1, void *inf2, const int dagger, const double kappa, const double mu, const double epsilon, const int Vf, QudaTwistGamma5Type twist, QudaPrecision precision)
void dslashReference(sFloat *res, gFloat **gaugeFull, sFloat *spinorField, int oddBit, int daggerBit)
static Float * gaugeLink(int i, int dir, int oddBit, Float **gaugeEven, Float **gaugeOdd, int nbr_distance)
Definition: dslash_util.h:104
void wil_mat(void *out, void **gauge, void *in, double kappa, int dagger_bit, QudaPrecision precision, QudaGaugeParam &gauge_param)
cpuColorSpinorField * out
int Vh
Definition: test_util.cpp:29
int faceVolume[4]
Definition: test_util.cpp:32
const void * c
void wil_matpc(void *outEven, void **gauge, void *inEven, double kappa, QudaMatPCType matpc_type, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
static void su3Mul(sFloat *res, gFloat *mat, sFloat *vec)
Definition: dslash_util.h:80
enum QudaTwistGamma5Type_s QudaTwistGamma5Type
void wil_dslash(void *out, void **gauge, void *in, int oddBit, int daggerBit, QudaPrecision precision, QudaGaugeParam &gauge_param)
static void su3Tmul(sFloat *res, gFloat *mat, sFloat *vec)
Definition: dslash_util.h:85
static __inline__ size_t size_t d
#define a
cpuColorSpinorField * spinor
Definition: covdev_test.cpp:41
void twistGamma5(sFloat *out, sFloat *in, const int dagger, const sFloat kappa, const sFloat mu, const QudaTwistFlavorType flavor, const int V, QudaTwistGamma5Type twist)
enum QudaTwistFlavorType_s QudaTwistFlavorType