12 #define stSpinorSiteSize 6
13 template<
typename Float>
16 printf(
"(%f,%f) (%f,%f) (%f,%f) \t",
17 spinor[0], spinor[1], spinor[2],
18 spinor[3], spinor[4], spinor[5]);
31 double* myspinor = (
double*)spinor;
32 for (i = 0;i < len; i++){
36 float* myspinor = (
float*)spinor;
37 for (i=0;i < len ;i++){
46 template<
typename Float>
51 for (i = 0;i < 3; i++){
53 printf(
"(%.10f,%.10f) \t", link[i*3*2 + j*2], link[i*3*2 + j*2 + 1]);
68 double* mylink = (
double*)link;
69 for (i = 0;i < len; i++){
73 float* mylink = (
float*)link;
74 for (i=0;i < len ;i++){
83 template <
typename Float>
85 a[0] += sign * (b[0]*c[0] - b[1]*c[1]);
86 a[1] -= sign * (b[0]*c[1] + b[1]*c[0]);
90 template<
typename Float>
98 Float* refc = &refc_buf[0];
100 memset((
void*)refc, 0,
sizeof(refc_buf));
104 Float* c = link + 12;
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];
121 double coff= -u0*u0*24;
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;
143 if ((i1+i4) % 2 == 1){
148 if ( (i4+i1+i2) % 2 == 1){
153 if (ga_idx >= (X4-3)*X1h*X2*
X3 ){
164 refc[0]*=coff; refc[1]*=coff; refc[2]*=coff; refc[3]*=coff; refc[4]*=coff; refc[5]*=coff;
167 double delta = 0.0001;
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");
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);
189 template<
typename Float>
195 Float* refc = &refc_buf[0];
197 memset((
void*)refc, 0,
sizeof(refc_buf));
201 Float* c = link + 12;
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];
220 bool last_node_in_t =
true;
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;
240 if ((i1+i4) % 2 == 1){
245 if ( (i4+i1+i2) % 2 == 1){
250 if (last_node_in_t && i4 == (X4-1) ){
261 double delta = 0.0001;
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");
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);
287 template <
typename Float>
294 template <
typename Float>
296 a[0] = b[0]*c[0] - b[1]*c[1];
297 a[1] = b[0]*c[1] + b[1]*c[0];
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];
308 template <
typename Float>
310 a[0] = b[0]*c[0] + b[1]*c[1];
311 a[1] = b[0]*c[1] - b[1]*c[0];
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]);
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];
329 template<
typename Float>
334 Float ref_link_buf[18];
336 memset(ref, 0,
sizeof(ref_link_buf));
338 ref[0] =
atan2(link[1], link[0]);
339 ref[1] =
atan2(link[13], link[12]);
340 for (
int i=2; i<7; i++) {
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];
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];
358 #define SMALL_NUM 1e-24
359 row_sum = (row_sum != 0)?row_sum:
SMALL_NUM;
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];
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;
382 if ((i1+i4) % 2 == 1){
387 if ( (i4+i1+i2) % 2 == 1){
392 if (ga_idx >= (X4-3)*X1h*X2*
X3 ){
403 Float U00_mag =
sqrt( (1.f/(u0*u0) - row_sum)>0? (1.f/(u0*u0)-row_sum):0);
408 ref[0] = U00_mag *
cos(ref[14]);
409 ref[1] = U00_mag *
sin(ref[14]);
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);
416 ref[12] = U20_mag *
cos(ref[15]);
417 ref[13] = U20_mag *
sin(ref[15]);
422 Float r_inv2 = 1.0/(u0*row_sum);
451 double delta = 0.0001;
453 for (i =0;i < 18; i++){
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");
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);
482 double* mylink = (
double*)link;
484 for (i = 0;i < len/2; i++){
487 printf(
"ERROR: even link sanity check failed, i=%d\n",i);
495 for (i = 0;i < len/2; i++){
498 printf(
"ERROR: odd link sanity check failed, i=%d\n",i);
505 float* mylink = (
float*)link;
508 for (i=0;i < len/2 ;i++){
511 printf(
"ERROR: even link sanity check 12 failed, i=%d\n",i);
525 for (i=0;i < len/2 ;i++){
528 printf(
"ERROR: odd link sanity check 12 failed, i=%d\n", i);
556 double* mylink = (
double*)link;
558 for (i = 0;i < len/2; i++){
559 for(dir=
XUP;dir <=
TUP; dir++){
562 printf(
"ERROR: even link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
571 for (i = 0;i < len/2; i++){
572 for(dir=
XUP;dir <=
TUP; dir++){
575 printf(
"ERROR: odd link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
583 float* mylink = (
float*)link;
586 for (i=0;i < len/2 ;i++){
587 for(dir=
XUP;dir <=
TUP; dir++){
590 printf(
"ERROR: even link sanity check 12 failed, i=%d, function %s\n",i, __FUNCTION__);
597 for (i=0;i < len/2 ;i++){
598 for(dir=
XUP;dir <=
TUP; dir++){
601 printf(
"ERROR: odd link sanity check 12 failed, i=%d, function %s\n", i, __FUNCTION__);
618 if (strcmp(s,
"8") == 0){
620 }
else if (strcmp(s,
"9") == 0){
622 }
else if (strcmp(s,
"12") == 0){
624 }
else if (strcmp(s,
"13") == 0){
626 }
else if (strcmp(s,
"18") == 0){
629 fprintf(stderr,
"Error: invalid reconstruct type\n");
643 if (strcmp(s,
"double") == 0){
645 }
else if (strcmp(s,
"single") == 0){
647 }
else if (strcmp(s,
"half") == 0){
650 fprintf(stderr,
"Error: invalid precision type\n");
691 ret =
"Cayley-Hamilton/SVD";
787 if (strcmp(s,
"wilson") == 0){
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){
797 }
else if (strcmp(s,
"asqtad") == 0){
799 }
else if (strcmp(s,
"domain_wall") == 0){
801 }
else if (strcmp(s,
"domain_wall_4d") == 0){
803 }
else if (strcmp(s,
"mobius") == 0){
806 fprintf(stderr,
"Error: invalid dslash type\n");
829 ret=
"twisted_clover";
841 ret =
"domain_wall_4d";
861 if (strcmp(s,
"kappa") == 0){
863 }
else if (strcmp(s,
"mass") == 0){
865 }
else if (strcmp(s,
"asym_mass") == 0){
868 fprintf(stderr,
"Error: invalid mass normalization\n");
891 fprintf(stderr,
"Error: invalid mass normalization\n");
903 if (strcmp(s,
"even_even") == 0){
905 }
else if (strcmp(s,
"odd_odd") == 0){
907 }
else if (strcmp(s,
"even_even_asym") == 0){
909 }
else if (strcmp(s,
"odd_odd_asym") == 0){
912 fprintf(stderr,
"Error: invalid matpc type\n");
932 ret =
"even_even_asym";
935 ret =
"odd_odd_asym";
938 fprintf(stderr,
"Error: invalid matpc type\n");
950 if (strcmp(s,
"minus") == 0){
952 }
else if (strcmp(s,
"plus") == 0){
954 }
else if (strcmp(s,
"deg_doublet") == 0){
956 }
else if (strcmp(s,
"nondeg_doublet") == 0){
958 }
else if (strcmp(s,
"no") == 0){
961 fprintf(stderr,
"Error: invalid flavor type\n");
984 ret =
"nondeg_doublet";
1002 if (strcmp(s,
"cg") == 0){
1004 }
else if (strcmp(s,
"bicgstab") == 0){
1006 }
else if (strcmp(s,
"gcr") == 0){
1008 }
else if (strcmp(s,
"pcg") == 0){
1010 }
else if (strcmp(s,
"mpcg") == 0){
1012 }
else if (strcmp(s,
"mpbicgstab") == 0){
1014 }
else if (strcmp(s,
"mr") == 0){
1017 fprintf(stderr,
"Error: invalid solver type\n");
1052 static char vstr[32];
1056 sprintf(vstr,
"%1d.%1d.%1d",
enum QudaMassNormalization_s QudaMassNormalization
void complexDotProduct(Float *a, Float *b, Float *c)
enum QudaPrecision_s QudaPrecision
QudaTwistFlavorType get_flavor_type(char *s)
QudaPrecision get_prec(char *s)
void accumulateComplexDotProduct(Float *a, Float *b, Float *c)
int site_link_sanity_check(void *link, int len, int precision, QudaGaugeParam *gaugeParam)
QudaGaugeParam gaugeParam
__host__ __device__ ValueType sqrt(ValueType x)
#define QUDA_VERSION_MINOR
const char * get_matpc_str(QudaMatPCType type)
const char * get_gauge_order_str(QudaGaugeFieldOrder order)
cpuColorSpinorField * spinor
const char * get_prec_str(QudaPrecision prec)
QudaReconstructType get_recon(char *s)
void accumulateConjugateProduct(Float *a, Float *b, Float *c, int sign)
void accumulateComplexProduct(Float *a, Float *b, Float *c, Float sign)
const char * get_flavor_str(QudaTwistFlavorType type)
const char * get_test_type(int t)
const char * get_solver_str(QudaInverterType type)
int fullLatticeIndex(int dim[4], int index, int oddBit)
void complexProduct(Float *a, Float *b, Float *c)
const char * get_unitarization_str(bool svd_only)
QudaInverterType get_solver_type(char *s)
void display_link_internal(Float *link)
__host__ __device__ ValueType sin(ValueType x)
void display_spinor(void *spinor, int len, int precision)
FloatingPoint< float > Float
const char * get_mass_normalization_str(QudaMassNormalization type)
QudaDslashType get_dslash_type(char *s)
const char * get_recon_str(QudaReconstructType recon)
__host__ __device__ ValueType atan2(ValueType x, ValueType y)
enum QudaMatPCType_s QudaMatPCType
int link_sanity_check(void *link, int len, int precision, int dir, QudaGaugeParam *gaugeParam)
#define QUDA_VERSION_SUBMINOR
__constant__ double coeff
QudaMatPCType get_matpc_type(char *s)
const char * get_dslash_str(QudaDslashType type)
enum QudaGaugeFieldOrder_s QudaGaugeFieldOrder
int link_sanity_check_internal_8(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
void display_link(void *link, int len, int precision)
int link_sanity_check_internal_12(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
void * memset(void *s, int c, size_t n)
cpuColorSpinorField * ref
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
void complexConjugateProduct(Float *a, Float *b, Float *c)
QudaMassNormalization get_mass_normalization_type(char *s)
enum QudaDslashType_s QudaDslashType
void display_spinor_internal(Float *spinor)
int site_link_sanity_check_internal_12(Float *link, int dir, int ga_idx, QudaGaugeParam *gaugeParam, int oddBit)
__host__ __device__ ValueType cos(ValueType x)
__constant__ double t_boundary
const char * get_quda_ver_str()
#define QUDA_VERSION_MAJOR
void complexAddTo(Float *a, Float *b)
enum QudaInverterType_s QudaInverterType
enum QudaTwistFlavorType_s QudaTwistFlavorType