QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 "misc.h"
9 #include "hisq_force_reference.h"
10 #include "ks_improved_force.h"
11 #include <sys/time.h>
12 #include <dslash_quda.h>
13 
14 using namespace quda;
15 extern void usage(char** argv);
18 
21 
24 
25 
27 
29 
30 
31 extern bool verify_results;
32 double accuracy = 1e-5;
33 int ODD_BIT = 1;
34 extern int device;
35 extern int xdim, ydim, zdim, tdim;
36 extern int gridsize_from_cmdline[];
37 
39 extern QudaPrecision prec;
44 
45 void setPrecision(QudaPrecision precision)
46 {
47  link_prec = precision;
48  return;
49 }
50 
51 
52 
53 
54 // Create a field of links that are not su3_matrices
55 void createNoisyLinkCPU(void** field, QudaPrecision prec, int seed)
56 {
57  createSiteLinkCPU(field, prec, 0);
58 
59  srand(seed);
60  for(int dir=0; dir<4; ++dir){
61  for(int i=0; i<V*18; ++i){
62  if(prec == QUDA_DOUBLE_PRECISION){
63  double* ptr = ((double**)field)[dir] + i;
64  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
65  }else if(prec == QUDA_SINGLE_PRECISION){
66  float* ptr = ((float**)field)[dir]+i;
67  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
68  }
69  }
70  }
71  return;
72 }
73 
74 
75 
76 // allocate memory
77 // set the layout, etc.
78 static void
80 {
82 
83  gaugeParam.X[0] = xdim;
84  gaugeParam.X[1] = ydim;
85  gaugeParam.X[2] = zdim;
86  gaugeParam.X[3] = tdim;
87 
88  setDims(gaugeParam.X);
89 
90  gaugeParam.cpu_prec = QUDA_DOUBLE_PRECISION;
91  gaugeParam.cuda_prec = link_prec;
92  gaugeParam.reconstruct = link_recon;
93  gaugeParam.gauge_order = QUDA_QDP_GAUGE_ORDER;
94  GaugeFieldParam gParam(0, gaugeParam);
98  gParam.anisotropy = 1;
99 
100  cpuFatLink = new cpuGaugeField(gParam);
101  cpuOprod = new cpuGaugeField(gParam);
102  cpuResult = new cpuGaugeField(gParam);
103  cpuReference = new cpuGaugeField(gParam);
104 
105  // create "gauge fields"
106  int seed=0;
107 #ifdef MULTI_GPU
108  seed += comm_rank();
109 #endif
110 
111  createNoisyLinkCPU((void**)cpuFatLink->Gauge_p(), gaugeParam.cpu_prec, seed);
112  createNoisyLinkCPU((void**)cpuOprod->Gauge_p(), gaugeParam.cpu_prec, seed+1);
113 
114  gParam.order = QUDA_FLOAT2_GAUGE_ORDER;
115  cudaFatLink = new cudaGaugeField(gParam);
116  cudaOprod = new cudaGaugeField(gParam);
117  cudaResult = new cudaGaugeField(gParam);
118  gParam.order = QUDA_QDP_GAUGE_ORDER;
119 
120  cudaFatLink->loadCPUField(*cpuFatLink);
121  cudaOprod->loadCPUField(*cpuOprod);
122 
123 
124  return;
125 }
126 
127 
128 static void
130 {
131  delete cpuFatLink;
132  delete cpuOprod;
133  delete cpuResult;
134 
135  delete cudaFatLink;
136  delete cudaOprod;
137  delete cudaResult;
138 
139  delete cpuReference;
140 
141  endQuda();
142  return;
143 }
144 
145 static void
147 {
148  hisq_force_init();
149 
150  double unitarize_eps = 1e-5;
151  const double hisq_force_filter = 5e-5;
152  const double max_det_error = 1e-12;
153  const bool allow_svd = true;
154  const bool svd_only = false;
155  const double svd_rel_err = 1e-8;
156  const double svd_abs_err = 1e-8;
157 
158  fermion_force::setUnitarizeForceConstants(unitarize_eps, hisq_force_filter, max_det_error, allow_svd, svd_only, svd_rel_err, svd_abs_err);
159 
160 
161 
162  int* num_failures_dev;
163  if(cudaMalloc(&num_failures_dev, sizeof(int)) != cudaSuccess){
164  errorQuda("cudaMalloc failed for num_failures_dev\n");
165  }
166  cudaMemset(num_failures_dev, 0, sizeof(int));
167 
168  printfQuda("Calling unitarizeForceCuda\n");
169  fermion_force::unitarizeForce(*cudaResult, *cudaOprod, *cudaFatLink, num_failures_dev);
170 
171 
172  if(verify_results){
173  printfQuda("Calling unitarizeForceCPU\n");
174  fermion_force::unitarizeForceCPU(*cpuResult, *cpuOprod, *cpuFatLink);
175  }
176 
177  cudaResult->saveCPUField(*cpuReference);
178 
179  printfQuda("Comparing CPU and GPU results\n");
180  for(int dir=0; dir<4; ++dir){
181  int res = compare_floats(((char**)cpuReference->Gauge_p())[dir], ((char**)cpuResult->Gauge_p())[dir], cpuReference->Volume()*gaugeSiteSize, accuracy, gaugeParam.cpu_prec);
182 #ifdef MULTI_GPU
183  comm_allreduce_int(&res);
184  res /= comm_size();
185 #endif
186  printfQuda("Dir:%d Test %s\n",dir,(1 == res) ? "PASSED" : "FAILED");
187  }
188 
189  hisq_force_end();
190 }
191 
192 
193 static void
195 {
196  printfQuda("running the following fermion force computation test:\n");
197 
198  printfQuda("link_precision link_reconstruct space_dim(x/y/z) T_dimension\n");
199  printfQuda("%s %s %d/%d/%d %d \n",
202  xdim, ydim, zdim, tdim);
203  return ;
204 
205 }
206 
207 int
208 main(int argc, char **argv)
209 {
210  int i;
211  for (i =1;i < argc; i++){
212 
213  if(process_command_line_option(argc, argv, &i) == 0){
214  continue;
215  }
216 
217  fprintf(stderr, "ERROR: Invalid option:%s\n", argv[i]);
218  usage(argv);
219  }
220 
221  initComms(argc, argv, gridsize_from_cmdline);
222 
224 
226 
227  hisq_force_test();
228 
229  finalizeComms();
230 
231  return EXIT_SUCCESS;
232 }
233 
static void hisq_force_test()
int comm_rank(void)
Definition: comm_mpi.cpp:82
QudaGhostExchange ghostExchange
Definition: lattice_field.h:76
static void hisq_force_init()
cudaGaugeField * cudaResult
void endQuda(void)
enum QudaPrecision_s QudaPrecision
void saveCPUField(cpuGaugeField &cpu) const
Upload from this field into a CPU field.
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)
Set the constant parameters for the force unitarization.
bool verify_results
Definition: test_util.cpp:1643
#define errorQuda(...)
Definition: util_quda.h:121
int main(int argc, char **argv)
int ydim
Definition: test_util.cpp:1616
int * num_failures_dev
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
static QudaGaugeParam gaugeParam
void unitarizeForce(cudaGaugeField &newForce, const cudaGaugeField &oldForce, const cudaGaugeField &gauge, int *unitarization_failed)
Unitarize the fermion force.
cpuGaugeField * cpuResult
void finalizeComms()
Definition: test_util.cpp:128
void setPrecision(QudaPrecision precision)
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
int compare_floats(void *a, void *b, int len, double epsilon, QudaPrecision precision)
Definition: test_util.cpp:434
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:701
void unitarizeForceCPU(cpuGaugeField &newForce, const cpuGaugeField &oldForce, const cpuGaugeField &gauge)
Unitarize the fermion force on CPU.
int zdim
Definition: test_util.cpp:1617
void createSiteLinkCPU(void **link, QudaPrecision precision, int phase)
Definition: test_util.cpp:1227
void setDims(int *)
Definition: test_util.cpp:151
QudaPrecision link_prec
cpuGaugeField * cpuReference
void loadCPUField(const cpuGaugeField &cpu)
Download into this field from a CPU field.
int comm_size(void)
Definition: comm_mpi.cpp:88
void initQuda(int device)
static void display_test_info()
cpuGaugeField * cpuFatLink
cpuGaugeField * cpuOprod
cudaGaugeField * cudaFatLink
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:768
QudaGaugeFieldOrder order
Definition: gauge_field.h:17
QudaReconstructType link_recon
Definition: test_util.cpp:1605
int Volume() const
double accuracy
QudaReconstructType reconstruct
Definition: quda.h:50
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
int device
Definition: test_util.cpp:1602
int V
Definition: test_util.cpp:27
cudaGaugeField * cudaOprod
static void hisq_force_end()
GaugeFieldParam gParam
QudaPrecision prec
Definition: test_util.cpp:1608
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
static double unitarize_eps
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
QudaLinkType link_type
Definition: gauge_field.h:19
void createNoisyLinkCPU(void **field, QudaPrecision prec, int seed)
int tdim
Definition: test_util.cpp:1618
#define printfQuda(...)
Definition: util_quda.h:115
void comm_allreduce_int(int *data)
Definition: comm_mpi.cpp:304
QudaFieldCreate create
Definition: gauge_field.h:26
QudaPrecision hw_prec
void usage(char **argv)
Definition: test_util.cpp:1783
void initComms(int argc, char **argv, int *const commDims)
Definition: test_util.cpp:88
QudaPrecision mom_prec
QudaPrecision cpu_hw_prec
int xdim
Definition: test_util.cpp:1615
QudaPrecision cpu_prec
Definition: quda.h:47
#define gaugeSiteSize
Definition: face_gauge.cpp:34