QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
misc.cpp
Go to the documentation of this file.
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include "quda.h"
4 #include <string.h>
5 #include "invert_quda.h"
6 #include "misc.h"
7  #include <assert.h>
8 #include "util_quda.h"
9 #include <test_util.h>
10 
11 
12 #define stSpinorSiteSize 6
13 template<typename Float>
15 {
16  printf("(%f,%f) (%f,%f) (%f,%f) \t",
17  spinor[0], spinor[1], spinor[2],
18  spinor[3], spinor[4], spinor[5]);
19 
20  printf("\n");
21  return;
22 }
23 
24 
25 
26 void display_spinor(void* spinor, int len, int precision)
27 {
28  int i;
29 
30  if (precision == QUDA_DOUBLE_PRECISION){
31  double* myspinor = (double*)spinor;
32  for (i = 0;i < len; i++){
34  }
35  }else if (precision == QUDA_SINGLE_PRECISION){
36  float* myspinor = (float*)spinor;
37  for (i=0;i < len ;i++){
39  }
40  }
41  return;
42 }
43 
44 
45 
46 template<typename Float>
48 {
49  int i, j;
50 
51  for (i = 0;i < 3; i++){
52  for(j=0;j < 3; j++){
53  printf("(%.10f,%.10f) \t", link[i*3*2 + j*2], link[i*3*2 + j*2 + 1]);
54  }
55  printf("\n");
56  }
57  printf("\n");
58  return;
59 }
60 
61 
62 
63 void display_link(void* link, int len, int precision)
64 {
65  int i;
66 
67  if (precision == QUDA_DOUBLE_PRECISION){
68  double* mylink = (double*)link;
69  for (i = 0;i < len; i++){
71  }
72  }else if (precision == QUDA_SINGLE_PRECISION){
73  float* mylink = (float*)link;
74  for (i=0;i < len ;i++){
76  }
77  }
78  return;
79 }
80 
81 
82 
83 template <typename Float>
84 void accumulateConjugateProduct(Float *a, Float *b, Float *c, int sign) {
85  a[0] += sign * (b[0]*c[0] - b[1]*c[1]);
86  a[1] -= sign * (b[0]*c[1] + b[1]*c[0]);
87 }
88 
89 
90 template<typename Float>
92 {
93  //printf("link sanity check is called\n");
94 
95  int ret =0;
96 
97  Float refc_buf[6];
98  Float* refc = &refc_buf[0];
99 
100  memset((void*)refc, 0, sizeof(refc_buf));
101 
102  Float* a = link;
103  Float* b = link + 6;
104  Float* c = link + 12;
105 
106  accumulateConjugateProduct(refc + 0*2, a + 1*2, b + 2*2, +1);
107  accumulateConjugateProduct(refc + 0*2, a + 2*2, b + 1*2, -1);
108  accumulateConjugateProduct(refc + 1*2, a + 2*2, b + 0*2, +1);
109  accumulateConjugateProduct(refc + 1*2, a + 0*2, b + 2*2, -1);
110  accumulateConjugateProduct(refc + 2*2, a + 0*2, b + 1*2, +1);
111  accumulateConjugateProduct(refc + 2*2, a + 1*2, b + 0*2, -1);
112 
113  int X1h=gaugeParam->X[0]/2;
114  int X1 =gaugeParam->X[0];
115  int X2 =gaugeParam->X[1];
116  int X3 =gaugeParam->X[2];
117  int X4 =gaugeParam->X[3];
118  double t_boundary = (gaugeParam->t_boundary ==QUDA_ANTI_PERIODIC_T)? -1.0:1.0;
119 
120  double u0 = gaugeParam->tadpole_coeff;
121  double coff= -u0*u0*24;
122  //coff = (dir < 6) ? coff : ( (ga_idx >= (X4-3)*X1h*X2*X3 )? t_boundary : 1);
123 
124  //float u0 = (dir < 6) ? gaugeParam->anisotropy : ( (ga_idx >= (X4-3)*X1h*X2*X3 )? t_boundary : 1);
125 
126 
127 #if 1
128 
129  {
130  int index = fullLatticeIndex(ga_idx, oddBit);
131  int i4 = index /(X3*X2*X1);
132  int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
133  int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
134  int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
135 
136  if (dir == 0) {
137  if (i4 % 2 == 1){
138  coff *= -1;
139  }
140  }
141 
142  if (dir == 2){
143  if ((i1+i4) % 2 == 1){
144  coff *= -1;
145  }
146  }
147  if (dir == 4){
148  if ( (i4+i1+i2) % 2 == 1){
149  coff *= -1;
150  }
151  }
152  if (dir == 6){
153  if (ga_idx >= (X4-3)*X1h*X2*X3 ){
154  coff *= -1;
155  }
156  }
157 
158  //printf("local ga_idx =%d, index=%d, i4,3,2,1 =%d %d %d %d\n", ga_idx, index, i4, i3, i2,i1);
159 
160  }
161 #endif
162 
163 
164  refc[0]*=coff; refc[1]*=coff; refc[2]*=coff; refc[3]*=coff; refc[4]*=coff; refc[5]*=coff;
165 
166 
167  double delta = 0.0001;
168  int i;
169  for (i =0;i < 6; i++){
170  double diff = refc[i] - c[i];
171  double absdiff = diff > 0? diff: (-diff);
172  if (absdiff > delta){
173  printf("ERROR: sanity check failed for link\n");
174  display_link_internal(link);
175  printf("refc = (%.10f,%.10f) (%.10f,%.10f) (%.10f,%.10f)\n",
176  refc[0], refc[1], refc[2], refc[3], refc[4], refc[5]);
177  printf("dir=%d, ga_idx=%d, coff=%f, t_boundary=%f\n",dir, ga_idx,coff, t_boundary);
178  printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
179  return -1;
180  }
181 
182  }
183 
184 
185  return ret;
186 }
187 
188 
189 template<typename Float>
191 {
192  int ret = 0;
193 
194  Float refc_buf[6];
195  Float* refc = &refc_buf[0];
196 
197  memset((void*)refc, 0, sizeof(refc_buf));
198 
199  Float* a = link;
200  Float* b = link + 6;
201  Float* c = link + 12;
202 
203  accumulateConjugateProduct(refc + 0*2, a + 1*2, b + 2*2, +1);
204  accumulateConjugateProduct(refc + 0*2, a + 2*2, b + 1*2, -1);
205  accumulateConjugateProduct(refc + 1*2, a + 2*2, b + 0*2, +1);
206  accumulateConjugateProduct(refc + 1*2, a + 0*2, b + 2*2, -1);
207  accumulateConjugateProduct(refc + 2*2, a + 0*2, b + 1*2, +1);
208  accumulateConjugateProduct(refc + 2*2, a + 1*2, b + 0*2, -1);
209 
210  int X1h=gaugeParam->X[0]/2;
211  int X1 =gaugeParam->X[0];
212  int X2 =gaugeParam->X[1];
213  int X3 =gaugeParam->X[2];
214  int X4 =gaugeParam->X[3];
215 
216  // only apply temporal boundary condition if I'm the last node in T
217 #ifdef MULTI_GPU
218  bool last_node_in_t = (commCoords(3) == commDim(3)-1);
219 #else
220  bool last_node_in_t = true;
221 #endif
222 
223 #if 1
224  double coeff= 1.0;
225 
226  {
227  int index = fullLatticeIndex(ga_idx, oddBit);
228  int i4 = index /(X3*X2*X1);
229  int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
230  int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
231  int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
232 
233  if (dir == XUP) {
234  if (i4 % 2 == 1){
235  coeff *= -1;
236  }
237  }
238 
239  if (dir == YUP){
240  if ((i1+i4) % 2 == 1){
241  coeff *= -1;
242  }
243  }
244  if (dir == ZUP){
245  if ( (i4+i1+i2) % 2 == 1){
246  coeff *= -1;
247  }
248  }
249  if (dir == TUP){
250  if (last_node_in_t && i4 == (X4-1) ){
251  coeff *= -1;
252  }
253  }
254  }
255 
256 
257  refc[0]*=coeff; refc[1]*=coeff; refc[2]*=coeff; refc[3]*=coeff; refc[4]*=coeff; refc[5]*=coeff;
258 #endif
259 
260 
261  double delta = 0.0001;
262  int i;
263  for (i =0;i < 6; i++){
264  double diff = refc[i] - c[i];
265  double absdiff = diff > 0? diff: (-diff);
266  if (absdiff > delta){
267  printf("ERROR: sanity check failed for site link\n");
268  display_link_internal(link);
269  printf("refc = (%.10f,%.10f) (%.10f,%.10f) (%.10f,%.10f)\n",
270  refc[0], refc[1], refc[2], refc[3], refc[4], refc[5]);
271  printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
272  return -1;
273  }
274 
275  }
276 
277 
278  return ret;
279 }
280 
281 
282 
283 
284 
285 
286 // a+=b
287 template <typename Float>
288 void complexAddTo(Float *a, Float *b) {
289  a[0] += b[0];
290  a[1] += b[1];
291 }
292 
293 // a = b*c
294 template <typename Float>
295 void complexProduct(Float *a, Float *b, Float *c) {
296  a[0] = b[0]*c[0] - b[1]*c[1];
297  a[1] = b[0]*c[1] + b[1]*c[0];
298 }
299 
300 // a = conj(b)*conj(c)
301 template <typename Float>
303  a[0] = b[0]*c[0] - b[1]*c[1];
304  a[1] = -b[0]*c[1] - b[1]*c[0];
305 }
306 
307 // a = conj(b)*c
308 template <typename Float>
309 void complexDotProduct(Float *a, Float *b, Float *c) {
310  a[0] = b[0]*c[0] + b[1]*c[1];
311  a[1] = b[0]*c[1] - b[1]*c[0];
312 }
313 
314 // a += b*c
315 template <typename Float>
317  a[0] += sign*(b[0]*c[0] - b[1]*c[1]);
318  a[1] += sign*(b[0]*c[1] + b[1]*c[0]);
319 }
320 
321 // a += conj(b)*c)
322 template <typename Float>
324  a[0] += b[0]*c[0] + b[1]*c[1];
325  a[1] += b[0]*c[1] - b[1]*c[0];
326 }
327 
328 
329 template<typename Float>
331 {
332  int ret =0;
333 
334  Float ref_link_buf[18];
335  Float* ref = & ref_link_buf[0];
336  memset(ref, 0, sizeof(ref_link_buf));
337 
338  ref[0] = atan2(link[1], link[0]);
339  ref[1] = atan2(link[13], link[12]);
340  for (int i=2; i<7; i++) {
341  ref[i] = link[i];
342  }
343 
344  int X1h=gaugeParam->X[0]/2;
345  int X2 =gaugeParam->X[1];
346  int X3 =gaugeParam->X[2];
347  int X4 =gaugeParam->X[3];
348  double t_boundary = (gaugeParam->t_boundary ==QUDA_ANTI_PERIODIC_T)? -1.0:1.0;
349 
350 
351  // First reconstruct first row
352  Float row_sum = 0.0;
353  row_sum += ref[2]*ref[2];
354  row_sum += ref[3]*ref[3];
355  row_sum += ref[4]*ref[4];
356  row_sum += ref[5]*ref[5];
357 
358 #define SMALL_NUM 1e-24
359  row_sum = (row_sum != 0)?row_sum: SMALL_NUM;
360 #if 1
361  Float u0= -gaugeParam->tadpole_coeff*gaugeParam->tadpole_coeff*24;
362  {
363  int X1h=gaugeParam->X[0]/2;
364  int X1 =gaugeParam->X[0];
365  int X2 =gaugeParam->X[1];
366  int X3 =gaugeParam->X[2];
367  int X4 =gaugeParam->X[3];
368 
369  int index = fullLatticeIndex(ga_idx, oddBit);
370  int i4 = index /(X3*X2*X1);
371  int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
372  int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
373  int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
374 
375  if (dir == 0) {
376  if (i4 % 2 == 1){
377  u0 *= -1;
378  }
379  }
380 
381  if (dir == 1){
382  if ((i1+i4) % 2 == 1){
383  u0 *= -1;
384  }
385  }
386  if (dir == 2){
387  if ( (i4+i1+i2) % 2 == 1){
388  u0 *= -1;
389  }
390  }
391  if (dir == 3){
392  if (ga_idx >= (X4-3)*X1h*X2*X3 ){
393  u0 *= -1;
394  }
395  }
396 
397  //printf("local ga_idx =%d, index=%d, i4,3,2,1 =%d %d %d %d\n", ga_idx, index, i4, i3, i2,i1);
398 
399  }
400 #endif
401 
402 
403  Float U00_mag = sqrt( (1.f/(u0*u0) - row_sum)>0? (1.f/(u0*u0)-row_sum):0);
404 
405  ref[14] = ref[0];
406  ref[15] = ref[1];
407 
408  ref[0] = U00_mag * cos(ref[14]);
409  ref[1] = U00_mag * sin(ref[14]);
410 
411  Float column_sum = 0.0;
412  for (int i=0; i<2; i++) column_sum += ref[i]*ref[i];
413  for (int i=6; i<8; i++) column_sum += ref[i]*ref[i];
414  Float U20_mag = sqrt( (1.f/(u0*u0) - column_sum) > 0? (1.f/(u0*u0)-column_sum) : 0);
415 
416  ref[12] = U20_mag * cos(ref[15]);
417  ref[13] = U20_mag * sin(ref[15]);
418 
419  // First column now restored
420 
421  // finally reconstruct last elements from SU(2) rotation
422  Float r_inv2 = 1.0/(u0*row_sum);
423 
424  // U11
425  Float A[2];
426  complexDotProduct(A, ref+0, ref+6);
427  complexConjugateProduct(ref+8, ref+12, ref+4);
428  accumulateComplexProduct(ref+8, A, ref+2, u0);
429  ref[8] *= -r_inv2;
430  ref[9] *= -r_inv2;
431 
432  // U12
433  complexConjugateProduct(ref+10, ref+12, ref+2);
434  accumulateComplexProduct(ref+10, A, ref+4, -u0);
435  ref[10] *= r_inv2;
436  ref[11] *= r_inv2;
437 
438  // U21
439  complexDotProduct(A, ref+0, ref+12);
440  complexConjugateProduct(ref+14, ref+6, ref+4);
441  accumulateComplexProduct(ref+14, A, ref+2, -u0);
442  ref[14] *= r_inv2;
443  ref[15] *= r_inv2;
444 
445  // U12
446  complexConjugateProduct(ref+16, ref+6, ref+2);
447  accumulateComplexProduct(ref+16, A, ref+4, u0);
448  ref[16] *= -r_inv2;
449  ref[17] *= -r_inv2;
450 
451  double delta = 0.0001;
452  int i;
453  for (i =0;i < 18; i++){
454 
455  double diff = ref[i] - link[i];
456  double absdiff = diff > 0? diff: (-diff);
457  if ( (ref[i] != ref[i]) || (absdiff > delta)){
458  printf("ERROR: sanity check failed for link\n");
459  display_link_internal(link);
460  printf("reconstructed link is\n");
462  printf("dir=%d, ga_idx=%d, u0=%f, t_boundary=%f\n",dir, ga_idx, u0, t_boundary);
463  printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
464  return -1;
465  }
466 
467  }
468 
469 
470  return ret;
471 }
472 
473 
474 //this len must be V
475 int
476 link_sanity_check(void* link, int len, int precision, int dir, QudaGaugeParam* gaugeParam)
477 {
478  int i;
479  int rc = 0;
480 
481  if (precision == QUDA_DOUBLE_PRECISION){
482  double* mylink = (double*)link;
483  //even
484  for (i = 0;i < len/2; i++){
485  rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
486  if (rc != 0){
487  printf("ERROR: even link sanity check failed, i=%d\n",i);
489  exit(1);
490  }
491  }
492 
493  mylink = mylink + gaugeSiteSize*len/2;
494  //odd
495  for (i = 0;i < len/2; i++){
496  rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 1);
497  if (rc != 0){
498  printf("ERROR: odd link sanity check failed, i=%d\n",i);
500  exit(1);
501  }
502  }
503 
504  }else if (precision == QUDA_SINGLE_PRECISION){
505  float* mylink = (float*)link;
506 
507  //even
508  for (i=0;i < len/2 ;i++){
509  rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
510  if (rc != 0){
511  printf("ERROR: even link sanity check 12 failed, i=%d\n",i);
512  exit(1);
513  }
514  /*
515  rc = link_sanity_check_internal_8(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
516  if (rc != 0){
517  printf("ERROR: even link sanity check 8 failed, i=%d\n",i);
518  exit(1);
519  }
520  */
521 
522  }
523  mylink = mylink + gaugeSiteSize*len/2;
524  //odd
525  for (i=0;i < len/2 ;i++){
526  rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 1);
527  if (rc != 0){
528  printf("ERROR: odd link sanity check 12 failed, i=%d\n", i);
529  exit(1);
530  }
531  /*
532  rc = link_sanity_check_internal_8(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
533  if (rc != 0){
534  printf("ERROR: even link sanity check 8 failed, i=%d\n",i);
535  exit(1);
536  }
537  */
538  }
539 
540  }
541 
542  return rc;
543 }
544 
545 
546 
547 //this len must be V
548 int
549 site_link_sanity_check(void* link, int len, int precision, QudaGaugeParam* gaugeParam)
550 {
551  int i;
552  int rc = 0;
553  int dir;
554 
555  if (precision == QUDA_DOUBLE_PRECISION){
556  double* mylink = (double*)link;
557  //even
558  for (i = 0;i < len/2; i++){
559  for(dir=XUP;dir <= TUP; dir++){
560  rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 0);
561  if (rc != 0){
562  printf("ERROR: even link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
564  exit(1);
565  }
566  }
567  }
568 
569  mylink = mylink + 4*gaugeSiteSize*len/2;
570  //odd
571  for (i = 0;i < len/2; i++){
572  for(dir=XUP;dir <= TUP; dir++){
573  rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 1);
574  if (rc != 0){
575  printf("ERROR: odd link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
577  exit(1);
578  }
579  }
580  }
581 
582  }else if (precision == QUDA_SINGLE_PRECISION){
583  float* mylink = (float*)link;
584 
585  //even
586  for (i=0;i < len/2 ;i++){
587  for(dir=XUP;dir <= TUP; dir++){
588  rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 0);
589  if (rc != 0){
590  printf("ERROR: even link sanity check 12 failed, i=%d, function %s\n",i, __FUNCTION__);
591  exit(1);
592  }
593  }
594  }
595  mylink = mylink + 4*gaugeSiteSize*len/2;
596  //odd
597  for (i=0;i < len/2 ;i++){
598  for(dir=XUP;dir <= TUP; dir++){
599  rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 1);
600  if (rc != 0){
601  printf("ERROR: odd link sanity check 12 failed, i=%d, function %s\n", i, __FUNCTION__);
602  exit(1);
603  }
604  }
605  }
606 
607  }
608 
609  return rc;
610 }
611 
612 
614 get_recon(char* s)
615 {
617 
618  if (strcmp(s, "8") == 0){
619  ret = QUDA_RECONSTRUCT_8;
620  }else if (strcmp(s, "9") == 0){
621  ret = QUDA_RECONSTRUCT_9;
622  }else if (strcmp(s, "12") == 0){
623  ret = QUDA_RECONSTRUCT_12;
624  }else if (strcmp(s, "13") == 0){
625  ret = QUDA_RECONSTRUCT_13;
626  }else if (strcmp(s, "18") == 0){
627  ret = QUDA_RECONSTRUCT_NO;
628  }else{
629  fprintf(stderr, "Error: invalid reconstruct type\n");
630  exit(1);
631  }
632 
633  return ret;
634 
635 
636 }
637 
639 get_prec(char* s)
640 {
642 
643  if (strcmp(s, "double") == 0){
644  ret = QUDA_DOUBLE_PRECISION;
645  }else if (strcmp(s, "single") == 0){
646  ret = QUDA_SINGLE_PRECISION;
647  }else if (strcmp(s, "half") == 0){
648  ret = QUDA_HALF_PRECISION;
649  }else{
650  fprintf(stderr, "Error: invalid precision type\n");
651  exit(1);
652  }
653 
654  return ret;
655 }
656 
657 const char*
659 {
660  const char* ret;
661 
662  switch( prec){
663 
665  ret= "double";
666  break;
668  ret= "single";
669  break;
670  case QUDA_HALF_PRECISION:
671  ret= "half";
672  break;
673  default:
674  ret = "unknown";
675  break;
676  }
677 
678 
679  return ret;
680 }
681 
682 
683 const char*
684 get_unitarization_str(bool svd_only)
685 {
686  const char* ret;
687 
688  if(svd_only){
689  ret = "SVD";
690  }else{
691  ret = "Cayley-Hamilton/SVD";
692  }
693  return ret;
694 }
695 
696 const char*
698 {
699  const char* ret;
700 
701  switch(order){
703  ret = "qdp";
704  break;
705 
707  ret = "milc";
708  break;
709 
711  ret = "cps_wilson";
712  break;
713 
714  default:
715  ret = "unknown";
716  break;
717  }
718 
719  return ret;
720 }
721 
722 
723 const char*
725 {
726  const char* ret;
727  switch(recon){
728  case QUDA_RECONSTRUCT_13:
729  ret="13";
730  break;
731  case QUDA_RECONSTRUCT_12:
732  ret= "12";
733  break;
734  case QUDA_RECONSTRUCT_9:
735  ret="9";
736  break;
737  case QUDA_RECONSTRUCT_8:
738  ret = "8";
739  break;
740  case QUDA_RECONSTRUCT_NO:
741  ret = "18";
742  break;
743  default:
744  ret="unknown";
745  break;
746  }
747 
748  return ret;
749 }
750 
751 const char*
753 {
754  const char* ret;
755  switch(t){
756  case 0:
757  ret = "even";
758  break;
759  case 1:
760  ret = "odd";
761  break;
762  case 2:
763  ret = "full";
764  break;
765  case 3:
766  ret = "mcg_even";
767  break;
768  case 4:
769  ret = "mcg_odd";
770  break;
771  case 5:
772  ret = "mcg_full";
773  break;
774  default:
775  ret = "unknown";
776  break;
777  }
778 
779  return ret;
780 }
781 
784 {
786 
787  if (strcmp(s, "wilson") == 0){
788  ret = QUDA_WILSON_DSLASH;
789  }else if (strcmp(s, "clover") == 0){
791  }else if (strcmp(s, "twisted_mass") == 0){
793  }else if (strcmp(s, "twisted_clover") == 0){
795  }else if (strcmp(s, "staggered") == 0){
796  ret = QUDA_STAGGERED_DSLASH;
797  }else if (strcmp(s, "asqtad") == 0){
798  ret = QUDA_ASQTAD_DSLASH;
799  }else if (strcmp(s, "domain_wall") == 0){
801  }else if (strcmp(s, "domain_wall_4d") == 0){
803  }else if (strcmp(s, "mobius") == 0){
805  }else{
806  fprintf(stderr, "Error: invalid dslash type\n");
807  exit(1);
808  }
809 
810  return ret;
811 }
812 
813 const char*
815 {
816  const char* ret;
817 
818  switch( type){
819  case QUDA_WILSON_DSLASH:
820  ret= "wilson";
821  break;
823  ret= "clover";
824  break;
826  ret= "twisted_mass";
827  break;
829  ret= "twisted_clover";
830  break;
832  ret = "staggered";
833  break;
834  case QUDA_ASQTAD_DSLASH:
835  ret = "asqtad";
836  break;
838  ret = "domain_wall";
839  break;
841  ret = "domain_wall_4d";
842  break;
844  ret = "mobius";
845  break;
846  default:
847  ret = "unknown";
848  break;
849  }
850 
851 
852  return ret;
853 
854 }
855 
858 {
860 
861  if (strcmp(s, "kappa") == 0){
863  }else if (strcmp(s, "mass") == 0){
865  }else if (strcmp(s, "asym_mass") == 0){
867  }else{
868  fprintf(stderr, "Error: invalid mass normalization\n");
869  exit(1);
870  }
871 
872  return ret;
873 }
874 
875 const char*
877 {
878  const char *s;
879 
880  switch (type) {
882  s = "kappa";
883  break;
885  s = "mass";
886  break;
888  s = "asym_mass";
889  break;
890  default:
891  fprintf(stderr, "Error: invalid mass normalization\n");
892  exit(1);
893  }
894 
895  return s;
896 }
897 
900 {
902 
903  if (strcmp(s, "even_even") == 0){
904  ret = QUDA_MATPC_EVEN_EVEN;
905  }else if (strcmp(s, "odd_odd") == 0){
906  ret = QUDA_MATPC_ODD_ODD;
907  }else if (strcmp(s, "even_even_asym") == 0){
909  }else if (strcmp(s, "odd_odd_asym") == 0){
911  }else{
912  fprintf(stderr, "Error: invalid matpc type\n");
913  exit(1);
914  }
915 
916  return ret;
917 }
918 
919 const char *
921 {
922  const char* ret;
923 
924  switch(type) {
926  ret = "even_even";
927  break;
928  case QUDA_MATPC_ODD_ODD:
929  ret = "odd_odd";
930  break;
932  ret = "even_even_asym";
933  break;
935  ret = "odd_odd_asym";
936  break;
937  default:
938  fprintf(stderr, "Error: invalid matpc type\n");
939  exit(1);
940  }
941 
942  return ret;
943 }
944 
947 {
949 
950  if (strcmp(s, "minus") == 0){
951  ret = QUDA_TWIST_MINUS;
952  }else if (strcmp(s, "plus") == 0){
953  ret = QUDA_TWIST_PLUS;
954  }else if (strcmp(s, "deg_doublet") == 0){
956  }else if (strcmp(s, "nondeg_doublet") == 0){
958  }else if (strcmp(s, "no") == 0){
959  ret = QUDA_TWIST_NO;
960  }else{
961  fprintf(stderr, "Error: invalid flavor type\n");
962  exit(1);
963  }
964 
965  return ret;
966 }
967 
968 const char*
970 {
971  const char* ret;
972 
973  switch(type) {
974  case QUDA_TWIST_MINUS:
975  ret = "minus";
976  break;
977  case QUDA_TWIST_PLUS:
978  ret = "plus";
979  break;
981  ret = "deg_doublet";
982  break;
984  ret = "nondeg_doublet";
985  break;
986  case QUDA_TWIST_NO:
987  ret = "no";
988  break;
989  default:
990  ret = "unknown";
991  break;
992  }
993 
994  return ret;
995 }
996 
999 {
1001 
1002  if (strcmp(s, "cg") == 0){
1003  ret = QUDA_CG_INVERTER;
1004  }else if (strcmp(s, "bicgstab") == 0){
1005  ret = QUDA_BICGSTAB_INVERTER;
1006  }else if (strcmp(s, "gcr") == 0){
1007  ret = QUDA_GCR_INVERTER;
1008  }else if (strcmp(s, "pcg") == 0){
1009  ret = QUDA_PCG_INVERTER;
1010  }else if (strcmp(s, "mpcg") == 0){
1011  ret = QUDA_MPCG_INVERTER;
1012  }else if (strcmp(s, "mpbicgstab") == 0){
1014  }else if (strcmp(s, "mr") == 0){
1015  ret = QUDA_MR_INVERTER;
1016  }else{
1017  fprintf(stderr, "Error: invalid solver type\n");
1018  exit(1);
1019  }
1020 
1021  return ret;
1022 }
1023 
1024 const char*
1026 {
1027  const char* ret;
1028 
1029  switch( type){
1030  case QUDA_CG_INVERTER:
1031  ret= "cg";
1032  break;
1034  ret= "bicgstab";
1035  break;
1036  case QUDA_GCR_INVERTER:
1037  ret= "gcr";
1038  break;
1039  default:
1040  ret = "unknown";
1041  break;
1042  }
1043 
1044 
1045  return ret;
1046 
1047 }
1048 
1049 const char*
1051 {
1052  static char vstr[32];
1053  int major_num = QUDA_VERSION_MAJOR;
1054  int minor_num = QUDA_VERSION_MINOR;
1055  int ext_num = QUDA_VERSION_SUBMINOR;
1056  sprintf(vstr, "%1d.%1d.%1d",
1057  major_num,
1058  minor_num,
1059  ext_num);
1060  return vstr;
1061 }
1062 
int commDim(int)
enum QudaMassNormalization_s QudaMassNormalization
#define SMALL_NUM
__constant__ int X1h
void complexDotProduct(Float *a, Float *b, Float *c)
Definition: misc.cpp:309
__constant__ int X2
enum QudaPrecision_s QudaPrecision
QudaTwistFlavorType get_flavor_type(char *s)
Definition: misc.cpp:946
QudaPrecision get_prec(char *s)
Definition: misc.cpp:639
void accumulateComplexDotProduct(Float *a, Float *b, Float *c)
Definition: misc.cpp:323
int site_link_sanity_check(void *link, int len, int precision, QudaGaugeParam *gaugeParam)
Definition: misc.cpp:549
#define stSpinorSiteSize
Definition: misc.cpp:12
QudaGaugeParam gaugeParam
__constant__ int X1
__host__ __device__ ValueType sqrt(ValueType x)
Definition: complex_quda.h:105
#define QUDA_VERSION_MINOR
Definition: quda_constants.h:2
#define gaugeSiteSize
const char * get_matpc_str(QudaMatPCType type)
Definition: misc.cpp:920
const char * get_gauge_order_str(QudaGaugeFieldOrder order)
Definition: misc.cpp:697
cpuColorSpinorField * spinor
Definition: dslash_test.cpp:40
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:658
QudaReconstructType get_recon(char *s)
Definition: misc.cpp:614
void accumulateConjugateProduct(Float *a, Float *b, Float *c, int sign)
Definition: misc.cpp:84
void accumulateComplexProduct(Float *a, Float *b, Float *c, Float sign)
Definition: misc.cpp:316
const char * get_flavor_str(QudaTwistFlavorType type)
Definition: misc.cpp:969
const char * get_test_type(int t)
Definition: misc.cpp:752
const char * get_solver_str(QudaInverterType type)
Definition: misc.cpp:1025
int fullLatticeIndex(int dim[4], int index, int oddBit)
Definition: test_util.cpp:400
void complexProduct(Float *a, Float *b, Float *c)
Definition: misc.cpp:295
const char * get_unitarization_str(bool svd_only)
Definition: misc.cpp:684
QudaInverterType get_solver_type(char *s)
Definition: misc.cpp:998
void display_link_internal(Float *link)
Definition: misc.cpp:47
__host__ __device__ ValueType sin(ValueType x)
Definition: complex_quda.h:40
void display_spinor(void *spinor, int len, int precision)
Definition: misc.cpp:26
FloatingPoint< float > Float
Definition: gtest.h:7350
const char * get_mass_normalization_str(QudaMassNormalization type)
Definition: misc.cpp:876
QudaDslashType get_dslash_type(char *s)
Definition: misc.cpp:783
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:724
__host__ __device__ ValueType atan2(ValueType x, ValueType y)
Definition: complex_quda.h:65
enum QudaMatPCType_s QudaMatPCType
int link_sanity_check(void *link, int len, int precision, int dir, QudaGaugeParam *gaugeParam)
Definition: misc.cpp:476
#define QUDA_VERSION_SUBMINOR
Definition: quda_constants.h:3
__constant__ double coeff
QudaMatPCType get_matpc_type(char *s)
Definition: misc.cpp:899
const char * get_dslash_str(QudaDslashType type)
Definition: misc.cpp:814
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
int link_sanity_check_internal_8(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
Definition: misc.cpp:330
int X[4]
Definition: quda.h:29
int commCoords(int)
void display_link(void *link, int len, int precision)
Definition: misc.cpp:63
double tadpole_coeff
Definition: quda.h:32
int link_sanity_check_internal_12(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
Definition: misc.cpp:91
void * memset(void *s, int c, size_t n)
cpuColorSpinorField * ref
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
__constant__ int X3
QudaTboundary t_boundary
Definition: quda.h:38
void complexConjugateProduct(Float *a, Float *b, Float *c)
Definition: misc.cpp:302
QudaMassNormalization get_mass_normalization_type(char *s)
Definition: misc.cpp:857
enum QudaDslashType_s QudaDslashType
void display_spinor_internal(Float *spinor)
Definition: misc.cpp:14
int site_link_sanity_check_internal_12(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
Definition: misc.cpp:190
__host__ __device__ ValueType cos(ValueType x)
Definition: complex_quda.h:35
__constant__ double t_boundary
VOLATILE spinorFloat * s
QudaPrecision prec
Definition: test_util.cpp:1551
const char * get_quda_ver_str()
Definition: misc.cpp:1050
#define QUDA_VERSION_MAJOR
Definition: quda_constants.h:1
void complexAddTo(Float *a, Float *b)
Definition: misc.cpp:288
enum QudaInverterType_s QudaInverterType
int oddBit
__constant__ int X4
enum QudaTwistFlavorType_s QudaTwistFlavorType