QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
hisq_unitarize_force_test.cpp
Go to the documentation of this file.
1 #include <cstdio>
2 #include <cstdlib>
3 #include <cstring>
4 
5 #include <quda.h>
6 #include "test_util.h"
7 #include "gauge_field.h"
8 #include "fat_force_quda.h"
9 #include "misc.h"
10 #include "hisq_force_reference.h"
11 #include "ks_improved_force.h"
12 #include "hw_quda.h"
13 #include <sys/time.h>
14 #include <dslash_quda.h>
15 
16 using namespace quda;
17 extern void usage(char** argv);
20 
23 
26 
27 
29 
31 
32 
33 extern bool verify_results;
34 double accuracy = 1e-5;
35 int ODD_BIT = 1;
36 extern int device;
37 extern int xdim, ydim, zdim, tdim;
38 extern int gridsize_from_cmdline[];
39 
41 extern QudaPrecision prec;
46 
47 void setPrecision(QudaPrecision precision)
48 {
49  link_prec = precision;
50  return;
51 }
52 
53 
54 
55 
56 // Create a field of links that are not su3_matrices
57 void createNoisyLinkCPU(void** field, QudaPrecision prec, int seed)
58 {
59  createSiteLinkCPU(field, prec, 0);
60 
61  srand(seed);
62  for(int dir=0; dir<4; ++dir){
63  for(int i=0; i<V*18; ++i){
64  if(prec == QUDA_DOUBLE_PRECISION){
65  double* ptr = ((double**)field)[dir] + i;
66  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
67  }else if(prec == QUDA_SINGLE_PRECISION){
68  float* ptr = ((float**)field)[dir]+i;
69  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
70  }
71  }
72  }
73  return;
74 }
75 
76 
77 
78 // allocate memory
79 // set the layout, etc.
80 static void
81 hisq_force_init()
82 {
84 
85  gaugeParam.X[0] = xdim;
86  gaugeParam.X[1] = ydim;
87  gaugeParam.X[2] = zdim;
88  gaugeParam.X[3] = tdim;
89 
91 
100  gParam.anisotropy = 1;
101 
106 
107  // create "gauge fields"
108  int seed=0;
109 #ifdef MULTI_GPU
110  seed += comm_rank();
111 #endif
112 
115 
121 
124 
125 
126  return;
127 }
128 
129 
130 static void
131 hisq_force_end()
132 {
133  delete cpuFatLink;
134  delete cpuOprod;
135  delete cpuResult;
136 
137  delete cudaFatLink;
138  delete cudaOprod;
139  delete cudaResult;
140 
141  delete cpuReference;
142 
143  endQuda();
144  return;
145 }
146 
147 static void
148 hisq_force_test()
149 {
150  hisq_force_init();
151 
152  double unitarize_eps = 1e-5;
153  const double hisq_force_filter = 5e-5;
154  const double max_det_error = 1e-12;
155  const bool allow_svd = true;
156  const bool svd_only = false;
157  const double svd_rel_err = 1e-8;
158  const double svd_abs_err = 1e-8;
159 
160  fermion_force::setUnitarizeForceConstants(unitarize_eps, hisq_force_filter, max_det_error, allow_svd, svd_only, svd_rel_err, svd_abs_err);
161 
162 
163 
164  int* num_failures_dev;
165  if(cudaMalloc(&num_failures_dev, sizeof(int)) != cudaSuccess){
166  errorQuda("cudaMalloc failed for num_failures_dev\n");
167  }
168  cudaMemset(num_failures_dev, 0, sizeof(int));
169 
170  printfQuda("Calling unitarizeForceCuda\n");
172 
173 
174  if(verify_results){
175  printfQuda("Calling unitarizeForceCPU\n");
177  }
178 
180 
181  printfQuda("Comparing CPU and GPU results\n");
182  for(int dir=0; dir<4; ++dir){
183  int res = compare_floats(((char**)cpuReference->Gauge_p())[dir], ((char**)cpuResult->Gauge_p())[dir], cpuReference->Volume()*gaugeSiteSize, accuracy, gaugeParam.cpu_prec);
184 #ifdef MULTI_GPU
185  comm_allreduce_int(&res);
186  res /= comm_size();
187 #endif
188  printfQuda("Dir:%d Test %s\n",dir,(1 == res) ? "PASSED" : "FAILED");
189  }
190 
191  hisq_force_end();
192 }
193 
194 
195 static void
197 {
198  printfQuda("running the following fermion force computation test:\n");
199 
200  printfQuda("link_precision link_reconstruct space_dim(x/y/z) T_dimension\n");
201  printfQuda("%s %s %d/%d/%d %d \n",
204  xdim, ydim, zdim, tdim);
205  return ;
206 
207 }
208 
209 int
210 main(int argc, char **argv)
211 {
212  int i;
213  for (i =1;i < argc; i++){
214 
215  if(process_command_line_option(argc, argv, &i) == 0){
216  continue;
217  }
218 
219  fprintf(stderr, "ERROR: Invalid option:%s\n", argv[i]);
220  usage(argv);
221  }
222 
223  initComms(argc, argv, gridsize_from_cmdline);
224 
226 
228 
229  hisq_force_test();
230 
231  finalizeComms();
232 
233  return EXIT_SUCCESS;
234 }
235 
int comm_rank(void)
Definition: comm_mpi.cpp:80
cudaGaugeField * cudaResult
void endQuda(void)
enum QudaPrecision_s QudaPrecision
int V
Definition: test_util.cpp:29
void setUnitarizeForceConstants(double unitarize_eps, double hisq_force_filter, double max_det_error, bool allow_svd, bool svd_only, double svd_rel_error, double svd_abs_error)
void display_test_info()
Definition: blas_test.cu:56
void unitarizeForceCuda(cudaGaugeField &cudaOldForce, cudaGaugeField &cudaGauge, cudaGaugeField *cudaNewForce, int *unitarization_failed, long long *flops=NULL)
bool verify_results
Definition: test_util.cpp:1568
void saveCPUField(cpuGaugeField &, const QudaFieldLocation &) const
#define errorQuda(...)
Definition: util_quda.h:73
void unitarizeForceCPU(cpuGaugeField &cpuOldForce, cpuGaugeField &cpuGauge, cpuGaugeField *cpuNewForce)
int main(int argc, char **argv)
QudaGaugeParam gaugeParam
int ydim
Definition: test_util.cpp:1554
void setDims(int *)
Definition: test_util.cpp:88
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:1635
#define gaugeSiteSize
cpuGaugeField * cpuResult
void finalizeComms()
Definition: test_util.cpp:65
void setPrecision(QudaPrecision precision)
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
int compare_floats(void *a, void *b, int len, double epsilon, QudaPrecision precision)
Definition: test_util.cpp:395
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:658
QudaPrecision cpu_prec
Definition: dslash_test.cpp:34
int zdim
Definition: test_util.cpp:1555
void createSiteLinkCPU(void **link, QudaPrecision precision, int phase)
Definition: test_util.cpp:1166
QudaPrecision link_prec
cpuGaugeField * cpuReference
int comm_size(void)
Definition: comm_mpi.cpp:86
void initQuda(int device)
cpuGaugeField * cpuFatLink
cpuGaugeField * cpuOprod
cudaGaugeField * cudaFatLink
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:724
QudaGaugeFieldOrder order
Definition: gauge_field.h:15
QudaReconstructType link_recon
Definition: test_util.cpp:1549
void loadCPUField(const cpuGaugeField &, const QudaFieldLocation &)
double accuracy
QudaReconstructType reconstruct
Definition: quda.h:43
QudaPrecision cuda_prec
Definition: quda.h:42
int X[4]
Definition: quda.h:29
int device
Definition: test_util.cpp:1546
cudaGaugeField * cudaOprod
QudaGhostExchange ghostExchange
Definition: gauge_field.h:40
GaugeFieldParam gParam
QudaPrecision prec
Definition: test_util.cpp:1551
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
int gridsize_from_cmdline[]
Definition: test_util.cpp:1559
QudaLinkType link_type
Definition: gauge_field.h:17
void createNoisyLinkCPU(void **field, QudaPrecision prec, int seed)
int tdim
Definition: test_util.cpp:1556
#define printfQuda(...)
Definition: util_quda.h:67
void comm_allreduce_int(int *data)
Definition: comm_mpi.cpp:225
QudaFieldCreate create
Definition: gauge_field.h:26
QudaPrecision hw_prec
void usage(char **argv)
Definition: test_util.cpp:1584
QudaPrecision mom_prec
void initComms(int argc, char **argv, const int *commDims)
Definition: test_util.cpp:48
QudaPrecision cpu_hw_prec
int xdim
Definition: test_util.cpp:1553
QudaPrecision cpu_prec
Definition: quda.h:40