QUDA v0.4.0
A library for QCD on GPUs
quda/tests/misc.cpp
Go to the documentation of this file.
00001 #include <stdio.h>
00002 #include <stdlib.h>
00003 #include "quda.h"
00004 #include <string.h>
00005 #include "invert_quda.h"
00006 #include "misc.h"
00007  #include <assert.h>
00008 #include "util_quda.h"
00009 #include <test_util.h>
00010 
00011 
00012 extern int verbose;
00013 
00014 #define stSpinorSiteSize 6
00015 template<typename Float>
00016 void display_spinor_internal(Float* spinor)
00017 {
00018     printf("(%f,%f) (%f,%f) (%f,%f) \t", 
00019            spinor[0], spinor[1], spinor[2], 
00020            spinor[3], spinor[4], spinor[5]);
00021     
00022     printf("\n");
00023     return;
00024 }
00025 
00026 
00027 
00028 void display_spinor(void* spinor, int len, int precision)
00029 {    
00030     int i;
00031     
00032     if (precision == QUDA_DOUBLE_PRECISION){
00033         double* myspinor = (double*)spinor;
00034         for (i = 0;i < len; i++){
00035             display_spinor_internal(myspinor + stSpinorSiteSize*i);
00036         }
00037     }else if (precision == QUDA_SINGLE_PRECISION){
00038         float* myspinor = (float*)spinor;
00039         for (i=0;i < len ;i++){
00040             display_spinor_internal(myspinor + stSpinorSiteSize*i);
00041         }
00042     }
00043     return;
00044 }
00045 
00046 
00047 
00048 template<typename Float>
00049 void display_link_internal(Float* link)
00050 {
00051     int i, j;
00052     
00053     for (i = 0;i < 3; i++){
00054         for(j=0;j < 3; j++){
00055             printf("(%.10f,%.10f) \t", link[i*3*2 + j*2], link[i*3*2 + j*2 + 1]);
00056         }
00057         printf("\n");
00058     }
00059     printf("\n");
00060     return;
00061 }
00062 
00063 
00064 
00065 void display_link(void* link, int len, int precision)
00066 {    
00067     int i;
00068     
00069     if (precision == QUDA_DOUBLE_PRECISION){
00070         double* mylink = (double*)link;
00071         for (i = 0;i < len; i++){
00072             display_link_internal(mylink + gaugeSiteSize*i);
00073         }
00074     }else if (precision == QUDA_SINGLE_PRECISION){
00075         float* mylink = (float*)link;
00076         for (i=0;i < len ;i++){
00077             display_link_internal(mylink + gaugeSiteSize*i);
00078         }
00079     }
00080     return;
00081 }
00082 
00083 
00084 
00085 template <typename Float>
00086 void accumulateConjugateProduct(Float *a, Float *b, Float *c, int sign) {
00087   a[0] += sign * (b[0]*c[0] - b[1]*c[1]);
00088   a[1] -= sign * (b[0]*c[1] + b[1]*c[0]);
00089 }
00090 
00091 
00092 template<typename Float>
00093 int link_sanity_check_internal_12(Float* link, int dir, int ga_idx, QudaGaugeParam* gaugeParam, int oddBit)
00094 {
00095     //printf("link sanity check is called\n");
00096     
00097     int ret =0;
00098     
00099     Float refc_buf[6];
00100     Float* refc = &refc_buf[0];
00101 
00102     memset((void*)refc, 0, sizeof(refc_buf));
00103 
00104     Float* a = link;
00105     Float* b = link + 6;
00106     Float* c = link + 12;
00107     
00108     accumulateConjugateProduct(refc + 0*2, a + 1*2, b + 2*2, +1);
00109     accumulateConjugateProduct(refc + 0*2, a + 2*2, b + 1*2, -1);
00110     accumulateConjugateProduct(refc + 1*2, a + 2*2, b + 0*2, +1);
00111     accumulateConjugateProduct(refc + 1*2, a + 0*2, b + 2*2, -1);
00112     accumulateConjugateProduct(refc + 2*2, a + 0*2, b + 1*2, +1);
00113     accumulateConjugateProduct(refc + 2*2, a + 1*2, b + 0*2, -1);
00114     
00115     int X1h=gaugeParam->X[0]/2;
00116     int X1 =gaugeParam->X[0];    
00117     int X2 =gaugeParam->X[1];
00118     int X3 =gaugeParam->X[2];
00119     int X4 =gaugeParam->X[3];
00120     double t_boundary = (gaugeParam->t_boundary ==QUDA_ANTI_PERIODIC_T)? -1.0:1.0;
00121 
00122    double u0 = gaugeParam->tadpole_coeff;
00123    double coff= -u0*u0*24;
00124    //coff = (dir < 6) ? coff : ( (ga_idx >= (X4-3)*X1h*X2*X3 )? t_boundary : 1); 
00125 
00126    //float u0 = (dir < 6) ? gaugeParam->anisotropy : ( (ga_idx >= (X4-3)*X1h*X2*X3 )? t_boundary : 1); 
00127 
00128   
00129 #if 1
00130    
00131    {
00132        int index = fullLatticeIndex(ga_idx, oddBit);
00133        int i4 = index /(X3*X2*X1);
00134        int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
00135        int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
00136        int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
00137        
00138        if (dir == 0) {
00139            if (i4 % 2 == 1){
00140                coff *= -1;
00141            }
00142        }
00143 
00144        if (dir == 2){
00145            if ((i1+i4) % 2 == 1){
00146                coff *= -1;
00147            }
00148        }
00149        if (dir == 4){
00150            if ( (i4+i1+i2) % 2 == 1){
00151                coff *= -1;
00152            }
00153        }
00154        if (dir == 6){
00155            if (ga_idx >= (X4-3)*X1h*X2*X3 ){
00156                coff *= -1;
00157            } 
00158        }
00159 
00160        //printf("local ga_idx =%d, index=%d, i4,3,2,1 =%d %d %d %d\n", ga_idx, index, i4, i3, i2,i1);
00161  
00162    }
00163 #endif
00164  
00165 
00166    refc[0]*=coff; refc[1]*=coff; refc[2]*=coff; refc[3]*=coff; refc[4]*=coff; refc[5]*=coff;
00167    
00168     
00169     double delta = 0.0001;
00170     int i;
00171     for (i =0;i < 6; i++){
00172         double diff =  refc[i] -  c[i];
00173         double absdiff = diff > 0? diff: (-diff);
00174         if (absdiff  > delta){
00175             printf("ERROR: sanity check failed for link\n");
00176             display_link_internal(link);
00177             printf("refc = (%.10f,%.10f) (%.10f,%.10f) (%.10f,%.10f)\n", 
00178                    refc[0], refc[1], refc[2], refc[3], refc[4], refc[5]);
00179             printf("dir=%d, ga_idx=%d, coff=%f, t_boundary=%f\n",dir, ga_idx,coff, t_boundary);
00180             printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
00181             return -1;
00182         }
00183         
00184     }
00185     
00186 
00187     return ret;
00188 }
00189 
00190 
00191 template<typename Float>
00192 int site_link_sanity_check_internal_12(Float* link, int dir, int ga_idx, QudaGaugeParam* gaugeParam, int oddBit)
00193 {
00194     
00195     int ret =0;
00196     
00197     Float refc_buf[6];
00198     Float* refc = &refc_buf[0];
00199 
00200     memset((void*)refc, 0, sizeof(refc_buf));
00201 
00202     Float* a = link;
00203     Float* b = link + 6;
00204     Float* c = link + 12;
00205     
00206     accumulateConjugateProduct(refc + 0*2, a + 1*2, b + 2*2, +1);
00207     accumulateConjugateProduct(refc + 0*2, a + 2*2, b + 1*2, -1);
00208     accumulateConjugateProduct(refc + 1*2, a + 2*2, b + 0*2, +1);
00209     accumulateConjugateProduct(refc + 1*2, a + 0*2, b + 2*2, -1);
00210     accumulateConjugateProduct(refc + 2*2, a + 0*2, b + 1*2, +1);
00211     accumulateConjugateProduct(refc + 2*2, a + 1*2, b + 0*2, -1);
00212 
00213 
00214     int X1h=gaugeParam->X[0]/2;
00215     int X1 =gaugeParam->X[0];    
00216     int X2 =gaugeParam->X[1];
00217     int X3 =gaugeParam->X[2];
00218     int X4 =gaugeParam->X[3];
00219 
00220 #if 1        
00221     double coeff= 1.0;
00222    
00223    {
00224        int index = fullLatticeIndex(ga_idx, oddBit);
00225        int i4 = index /(X3*X2*X1);
00226        int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
00227        int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
00228        int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
00229        
00230        if (dir == XUP) {
00231            if (i4 % 2 == 1){
00232                coeff *= -1;
00233            }
00234        }
00235 
00236        if (dir == YUP){
00237            if ((i1+i4) % 2 == 1){
00238                coeff *= -1;
00239            }
00240        }
00241        if (dir == ZUP){
00242            if ( (i4+i1+i2) % 2 == 1){
00243                coeff *= -1;
00244            }
00245        }
00246        if (dir == TUP){
00247          if ((commCoords(3) == commDim(3) -1) && i4 == (X4-1) ){
00248            coeff *= -1;
00249          } 
00250        }       
00251    }
00252  
00253    
00254    refc[0]*=coeff; refc[1]*=coeff; refc[2]*=coeff; refc[3]*=coeff; refc[4]*=coeff; refc[5]*=coeff;
00255 #endif
00256    
00257     
00258     double delta = 0.0001;
00259     int i;
00260     for (i =0;i < 6; i++){
00261         double diff =  refc[i] -  c[i];
00262         double absdiff = diff > 0? diff: (-diff);
00263         if (absdiff  > delta){
00264             printf("ERROR: sanity check failed for site link\n");
00265             display_link_internal(link);
00266             printf("refc = (%.10f,%.10f) (%.10f,%.10f) (%.10f,%.10f)\n", 
00267                    refc[0], refc[1], refc[2], refc[3], refc[4], refc[5]);
00268             printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
00269             return -1;
00270         }
00271         
00272     }
00273     
00274 
00275     return ret;
00276 }
00277 
00278 
00279 
00280 
00281 
00282 
00283 // a+=b
00284 template <typename Float>
00285 void complexAddTo(Float *a, Float *b) {
00286   a[0] += b[0];
00287   a[1] += b[1];
00288 }
00289 
00290 // a = b*c
00291 template <typename Float>
00292 void complexProduct(Float *a, Float *b, Float *c) {
00293     a[0] = b[0]*c[0] - b[1]*c[1];
00294     a[1] = b[0]*c[1] + b[1]*c[0];
00295 }
00296 
00297 // a = conj(b)*conj(c)
00298 template <typename Float>
00299 void complexConjugateProduct(Float *a, Float *b, Float *c) {
00300     a[0] = b[0]*c[0] - b[1]*c[1];
00301     a[1] = -b[0]*c[1] - b[1]*c[0];
00302 }
00303 
00304 // a = conj(b)*c
00305 template <typename Float>
00306 void complexDotProduct(Float *a, Float *b, Float *c) {
00307     a[0] = b[0]*c[0] + b[1]*c[1];
00308     a[1] = b[0]*c[1] - b[1]*c[0];
00309 }
00310 
00311 // a += b*c
00312 template <typename Float>
00313 void accumulateComplexProduct(Float *a, Float *b, Float *c, Float sign) {
00314   a[0] += sign*(b[0]*c[0] - b[1]*c[1]);
00315   a[1] += sign*(b[0]*c[1] + b[1]*c[0]);
00316 }
00317 
00318 // a += conj(b)*c)
00319 template <typename Float>
00320 void accumulateComplexDotProduct(Float *a, Float *b, Float *c) {
00321     a[0] += b[0]*c[0] + b[1]*c[1];
00322     a[1] += b[0]*c[1] - b[1]*c[0];
00323 }
00324 
00325 
00326 template<typename Float>
00327 int link_sanity_check_internal_8(Float* link, int dir, int ga_idx, QudaGaugeParam* gaugeParam, int oddBit)
00328 {
00329     int ret =0;
00330     
00331     Float ref_link_buf[18];
00332     Float* ref = & ref_link_buf[0];    
00333     memset(ref, 0, sizeof(ref_link_buf));
00334     
00335     ref[0] = atan2(link[1], link[0]);
00336     ref[1] = atan2(link[13], link[12]);
00337     for (int i=2; i<7; i++) {
00338         ref[i] = link[i];
00339     }
00340 
00341     int X1h=gaugeParam->X[0]/2;
00342     int X2 =gaugeParam->X[1];
00343     int X3 =gaugeParam->X[2];
00344     int X4 =gaugeParam->X[3];
00345     double t_boundary = (gaugeParam->t_boundary ==QUDA_ANTI_PERIODIC_T)? -1.0:1.0;
00346     
00347     
00348     // First reconstruct first row
00349     Float row_sum = 0.0;
00350     row_sum += ref[2]*ref[2];
00351     row_sum += ref[3]*ref[3];
00352     row_sum += ref[4]*ref[4];
00353     row_sum += ref[5]*ref[5];
00354 
00355 #define SMALL_NUM 1e-24
00356     row_sum = (row_sum != 0)?row_sum: SMALL_NUM;
00357 #if 1
00358     Float u0= -gaugeParam->tadpole_coeff*gaugeParam->tadpole_coeff*24;
00359     {
00360         int X1h=gaugeParam->X[0]/2;
00361         int X1 =gaugeParam->X[0];
00362         int X2 =gaugeParam->X[1];
00363         int X3 =gaugeParam->X[2];
00364         int X4 =gaugeParam->X[3];
00365       
00366         int index = fullLatticeIndex(ga_idx, oddBit);
00367         int i4 = index /(X3*X2*X1);
00368         int i3 = (index - i4*(X3*X2*X1))/(X2*X1);
00369         int i2 = (index - i4*(X3*X2*X1) - i3*(X2*X1))/X1;
00370         int i1 = index - i4*(X3*X2*X1) - i3*(X2*X1) - i2*X1;
00371       
00372         if (dir == 0) {
00373             if (i4 % 2 == 1){
00374                 u0 *= -1;
00375             }
00376         }
00377       
00378         if (dir == 1){
00379             if ((i1+i4) % 2 == 1){
00380                 u0 *= -1;
00381             }
00382         }
00383         if (dir == 2){
00384             if ( (i4+i1+i2) % 2 == 1){
00385                 u0 *= -1;
00386             }
00387         }
00388         if (dir == 3){
00389             if (ga_idx >= (X4-3)*X1h*X2*X3 ){
00390                 u0 *= -1;
00391             }
00392         }
00393        
00394         //printf("local ga_idx =%d, index=%d, i4,3,2,1 =%d %d %d %d\n", ga_idx, index, i4, i3, i2,i1);
00395        
00396     }
00397 #endif
00398     
00399 
00400     Float U00_mag = sqrt( (1.f/(u0*u0) - row_sum)>0? (1.f/(u0*u0)-row_sum):0);
00401   
00402     ref[14] = ref[0];
00403     ref[15] = ref[1];
00404 
00405     ref[0] = U00_mag * cos(ref[14]);
00406     ref[1] = U00_mag * sin(ref[14]);
00407 
00408     Float column_sum = 0.0;
00409     for (int i=0; i<2; i++) column_sum += ref[i]*ref[i];
00410     for (int i=6; i<8; i++) column_sum += ref[i]*ref[i];
00411     Float U20_mag = sqrt( (1.f/(u0*u0) - column_sum) > 0? (1.f/(u0*u0)-column_sum) : 0);
00412 
00413     ref[12] = U20_mag * cos(ref[15]);
00414     ref[13] = U20_mag * sin(ref[15]);
00415 
00416     // First column now restored
00417 
00418     // finally reconstruct last elements from SU(2) rotation
00419     Float r_inv2 = 1.0/(u0*row_sum);
00420 
00421     // U11
00422     Float A[2];
00423     complexDotProduct(A, ref+0, ref+6);
00424     complexConjugateProduct(ref+8, ref+12, ref+4);
00425     accumulateComplexProduct(ref+8, A, ref+2, u0);
00426     ref[8] *= -r_inv2;
00427     ref[9] *= -r_inv2;
00428 
00429     // U12
00430     complexConjugateProduct(ref+10, ref+12, ref+2);
00431     accumulateComplexProduct(ref+10, A, ref+4, -u0);
00432     ref[10] *= r_inv2;
00433     ref[11] *= r_inv2;
00434 
00435     // U21
00436     complexDotProduct(A, ref+0, ref+12);
00437     complexConjugateProduct(ref+14, ref+6, ref+4);
00438     accumulateComplexProduct(ref+14, A, ref+2, -u0);
00439     ref[14] *= r_inv2;
00440     ref[15] *= r_inv2;
00441 
00442     // U12
00443     complexConjugateProduct(ref+16, ref+6, ref+2);
00444     accumulateComplexProduct(ref+16, A, ref+4, u0);
00445     ref[16] *= -r_inv2;
00446     ref[17] *= -r_inv2;
00447 
00448     double delta = 0.0001;
00449     int i;
00450     for (i =0;i < 18; i++){
00451 
00452         double diff =  ref[i] -  link[i];
00453         double absdiff = diff > 0? diff: (-diff);
00454         if ( (ref[i] !=  ref[i]) || (absdiff  > delta)){
00455             printf("ERROR: sanity check failed for link\n");
00456             display_link_internal(link);
00457             printf("reconstructed link is\n");
00458             display_link_internal(ref);
00459             printf("dir=%d, ga_idx=%d, u0=%f, t_boundary=%f\n",dir, ga_idx, u0, t_boundary);
00460             printf("X=%d %d %d %d, X1h=%d\n", gaugeParam->X[0], X2, X3, X4, X1h);
00461             return -1;
00462         }
00463         
00464     }
00465     
00466 
00467     return ret;
00468 }
00469 
00470 
00471 //this len must be V
00472 int
00473 link_sanity_check(void* link, int len, int precision, int dir, QudaGaugeParam* gaugeParam)
00474 {
00475     int i;
00476     int rc = 0;
00477     
00478     if (precision == QUDA_DOUBLE_PRECISION){
00479         double* mylink = (double*)link;
00480         //even
00481         for (i = 0;i < len/2; i++){
00482             rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
00483             if (rc != 0){
00484                 printf("ERROR: even link sanity check failed, i=%d\n",i);
00485                 display_link_internal(mylink+gaugeSiteSize*i);
00486                 exit(1);
00487             }
00488         }
00489         
00490         mylink = mylink + gaugeSiteSize*len/2;
00491         //odd
00492         for (i = 0;i < len/2; i++){
00493             rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 1);
00494             if (rc != 0){
00495                 printf("ERROR: odd link sanity check failed, i=%d\n",i);
00496                 display_link_internal(mylink+gaugeSiteSize*i);
00497                 exit(1);
00498             }
00499         }       
00500         
00501     }else if (precision == QUDA_SINGLE_PRECISION){
00502         float* mylink = (float*)link;
00503 
00504         //even
00505         for (i=0;i < len/2 ;i++){
00506             rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
00507             if (rc != 0){
00508                 printf("ERROR: even link sanity check 12 failed, i=%d\n",i);
00509                 exit(1);
00510             }
00511             /*
00512             rc = link_sanity_check_internal_8(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
00513             if (rc != 0){
00514                 printf("ERROR: even link sanity check 8 failed, i=%d\n",i);
00515                 exit(1);
00516             }
00517             */
00518             
00519         }
00520         mylink = mylink + gaugeSiteSize*len/2;
00521         //odd
00522         for (i=0;i < len/2 ;i++){
00523             rc = link_sanity_check_internal_12(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 1);
00524             if (rc != 0){
00525                 printf("ERROR: odd link sanity check 12 failed, i=%d\n", i);
00526                 exit(1);
00527             }   
00528             /*
00529             rc = link_sanity_check_internal_8(mylink + gaugeSiteSize*i, dir, i, gaugeParam, 0);
00530             if (rc != 0){
00531                 printf("ERROR: even link sanity check 8 failed, i=%d\n",i);
00532                 exit(1);
00533             }
00534             */
00535         }       
00536 
00537     }
00538     
00539     return rc;
00540 }
00541 
00542 
00543 
00544 //this len must be V
00545 int
00546 site_link_sanity_check(void* link, int len, int precision, QudaGaugeParam* gaugeParam)
00547 {
00548     int i;
00549     int rc = 0;
00550     int dir;
00551     
00552     if (precision == QUDA_DOUBLE_PRECISION){
00553         double* mylink = (double*)link;
00554         //even  
00555         for (i = 0;i < len/2; i++){
00556             for(dir=XUP;dir <= TUP; dir++){
00557                 rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 0);
00558                 if (rc != 0){
00559                     printf("ERROR: even link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
00560                     display_link_internal(mylink+gaugeSiteSize*i);
00561                     exit(1);
00562                 }
00563             }
00564         }
00565         
00566         mylink = mylink + 4*gaugeSiteSize*len/2;
00567         //odd
00568         for (i = 0;i < len/2; i++){
00569             for(dir=XUP;dir <= TUP; dir++){         
00570                 rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 1);
00571                 if (rc != 0){
00572                     printf("ERROR: odd link sanity check failed, i=%d, function %s\n",i, __FUNCTION__);
00573                     display_link_internal(mylink+gaugeSiteSize*i);
00574                     exit(1);
00575                 }
00576             }
00577         }       
00578         
00579     }else if (precision == QUDA_SINGLE_PRECISION){
00580         float* mylink = (float*)link;
00581 
00582         //even
00583         for (i=0;i < len/2 ;i++){
00584             for(dir=XUP;dir <= TUP; dir++){
00585                 rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 0);
00586                 if (rc != 0){
00587                     printf("ERROR: even link sanity check 12 failed, i=%d, function %s\n",i, __FUNCTION__);
00588                     exit(1);
00589                 }
00590             }
00591         }
00592         mylink = mylink + 4*gaugeSiteSize*len/2;
00593         //odd
00594         for (i=0;i < len/2 ;i++){
00595             for(dir=XUP;dir <= TUP; dir++){
00596                 rc = site_link_sanity_check_internal_12(mylink + gaugeSiteSize*(4*i+dir), dir, i, gaugeParam, 1);
00597                 if (rc != 0){
00598                     printf("ERROR: odd link sanity check 12 failed, i=%d, function %s\n", i, __FUNCTION__);
00599                     exit(1);
00600                 }       
00601             }
00602         }       
00603 
00604     }
00605     
00606     return rc;
00607 }
00608 
00609 
00610 QudaReconstructType
00611 get_recon(char* s)
00612 {
00613     QudaReconstructType  ret;
00614     
00615     if (strcmp(s, "8") == 0){
00616         ret =  QUDA_RECONSTRUCT_8;
00617     }else if (strcmp(s, "12") == 0){
00618         ret =  QUDA_RECONSTRUCT_12;
00619     }else if (strcmp(s, "18") == 0){
00620         ret =  QUDA_RECONSTRUCT_NO;
00621     }else{
00622         fprintf(stderr, "Error: invalid reconstruct type\n");
00623         exit(1);
00624     }
00625     
00626     return ret;
00627     
00628     
00629 }
00630 
00631 QudaPrecision
00632 get_prec(char* s)
00633 {
00634     QudaPrecision ret = QUDA_DOUBLE_PRECISION;
00635     
00636     if (strcmp(s, "double") == 0){
00637         ret = QUDA_DOUBLE_PRECISION;
00638     }else if (strcmp(s, "single") == 0){
00639         ret = QUDA_SINGLE_PRECISION;
00640     }else if (strcmp(s, "half") == 0){
00641         ret = QUDA_HALF_PRECISION;
00642     }else{
00643         fprintf(stderr, "Error: invalid precision type\n");     
00644         exit(1);
00645     }
00646     
00647     return ret;
00648 }
00649 
00650 const char* 
00651 get_prec_str(QudaPrecision prec)
00652 {
00653     const char* ret;
00654     
00655     switch( prec){
00656         
00657     case QUDA_DOUBLE_PRECISION:
00658         ret=  "double";
00659         break;
00660     case QUDA_SINGLE_PRECISION:
00661         ret= "single";
00662         break;
00663     case QUDA_HALF_PRECISION:
00664         ret= "half";
00665         break;
00666     default:
00667         ret = "unknown";        
00668         break;
00669     }
00670     
00671     
00672     return ret;
00673 }
00674 
00675 
00676 const char* 
00677 get_unitarization_str(bool svd_only)
00678 {
00679   const char* ret;
00680  
00681   if(svd_only){
00682     ret = "SVD";
00683   }else{
00684     ret = "Cayley-Hamilton/SVD";
00685   }
00686   return ret;
00687 }
00688 
00689 const char* 
00690 get_gauge_order_str(QudaGaugeFieldOrder order)
00691 {
00692   const char* ret;
00693 
00694   switch(order){
00695     case QUDA_QDP_GAUGE_ORDER:
00696         ret = "qdp";
00697         break;
00698 
00699     case QUDA_MILC_GAUGE_ORDER:
00700         ret = "milc";
00701         break;
00702 
00703     case QUDA_CPS_WILSON_GAUGE_ORDER:
00704         ret = "cps_wilson";
00705         break;
00706 
00707     default:
00708         ret = "unknown";
00709         break;
00710   }     
00711 
00712   return ret;
00713 }
00714 
00715 
00716 const char* 
00717 get_recon_str(QudaReconstructType recon)
00718 {
00719     const char* ret;
00720     switch(recon){
00721     case QUDA_RECONSTRUCT_12:
00722         ret= "12";
00723         break;
00724     case QUDA_RECONSTRUCT_8:
00725         ret = "8";
00726         break;
00727     case QUDA_RECONSTRUCT_NO:
00728         ret = "18";
00729         break;
00730     default:
00731         ret="unknown";
00732         break;
00733     }
00734     
00735     return ret;
00736 }
00737 
00738 const char*
00739 get_test_type(int t)
00740 {
00741     const char* ret;
00742     switch(t){
00743     case 0:
00744         ret = "even";
00745         break;
00746     case 1:
00747         ret = "odd";
00748         break;
00749     case 2:
00750         ret = "full";
00751         break;
00752     case 3:
00753         ret = "mcg_even";
00754         break;  
00755     case 4:
00756         ret = "mcg_odd";
00757         break;  
00758     case 5:
00759         ret = "mcg_full";
00760         break;  
00761     default:
00762         ret = "unknown";
00763         break;
00764     }
00765     
00766     return ret;
00767 }
00768 
00769 QudaDslashType
00770 get_dslash_type(char* s)
00771 {
00772   QudaDslashType ret =  QUDA_INVALID_DSLASH;
00773   
00774   if (strcmp(s, "wilson") == 0){
00775     ret = QUDA_WILSON_DSLASH;
00776   }else if (strcmp(s, "clover") == 0){
00777     ret = QUDA_CLOVER_WILSON_DSLASH;
00778   }else if (strcmp(s, "twisted_mass") == 0){
00779     ret = QUDA_TWISTED_MASS_DSLASH;
00780   }else if (strcmp(s, "asqtad") == 0){
00781     ret =  QUDA_ASQTAD_DSLASH;
00782   }else if (strcmp(s, "domain_wall") == 0){
00783     ret =  QUDA_DOMAIN_WALL_DSLASH;
00784   }else{
00785     fprintf(stderr, "Error: invalid dslash type\n");    
00786     exit(1);
00787   }
00788   
00789   return ret;
00790 }
00791 
00792 const char* 
00793 get_dslash_type_str(QudaDslashType type)
00794 {
00795   const char* ret;
00796   
00797   switch( type){        
00798   case QUDA_WILSON_DSLASH:
00799     ret=  "wilson";
00800     break;
00801   case QUDA_CLOVER_WILSON_DSLASH:
00802     ret= "clover";
00803     break;
00804   case QUDA_TWISTED_MASS_DSLASH:
00805     ret= "twisted_mass";
00806     break;
00807   case QUDA_ASQTAD_DSLASH:
00808     ret = "asqtad";
00809     break;
00810   case QUDA_DOMAIN_WALL_DSLASH:
00811     ret = "domain_wall";
00812       break;
00813   default:
00814     ret = "unknown";    
00815     break;
00816   }
00817   
00818   
00819   return ret;
00820     
00821 }
00822 
00823 const char* 
00824 get_quda_ver_str()
00825 {
00826   static char vstr[32];
00827   int major_num = QUDA_VERSION_MAJOR;
00828   int minor_num = QUDA_VERSION_MINOR;
00829   int ext_num = QUDA_VERSION_SUBMINOR;
00830   sprintf(vstr, "%1d.%1d.%1d", 
00831           major_num,
00832           minor_num,
00833           ext_num);
00834   return vstr;
00835 }
00836 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines