QUDA  v1.1.0
A library for QCD on GPUs
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 "host_utils.h"
7 #include <command_line_params.h>
8 #include "gauge_field.h"
9 #include "misc.h"
10 #include "hisq_force_reference.h"
11 #include "ks_improved_force.h"
12 #include <sys/time.h>
13 #include <dslash_quda.h>
14 
15 using namespace quda;
16 
19 
22 
25 
27 
28 static QudaGaugeParam gaugeParam;
29 
30 // extern bool verify_results;
31 double accuracy = 1e-5;
32 int ODD_BIT = 1;
33 
35 
36 void setPrecision(QudaPrecision precision) { link_prec = precision; }
37 
38 // Create a field of links that are not su3_matrices
39 void createNoisyLinkCPU(void** field, QudaPrecision prec, int seed)
40 {
41  createSiteLinkCPU(field, prec, 0);
42 
43  srand(seed);
44  for(int dir=0; dir<4; ++dir){
45  for(int i=0; i<V*18; ++i){
47  double* ptr = ((double**)field)[dir] + i;
48  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
49  }else if(prec == QUDA_SINGLE_PRECISION){
50  float* ptr = ((float**)field)[dir]+i;
51  *ptr += (rand() - RAND_MAX/2.0)/(20.0*RAND_MAX);
52  }
53  }
54  }
55 }
56 
57 // allocate memory
58 // set the layout, etc.
59 static void hisq_force_init()
60 {
62 
63  gaugeParam.X[0] = xdim;
64  gaugeParam.X[1] = ydim;
65  gaugeParam.X[2] = zdim;
66  gaugeParam.X[3] = tdim;
67 
68  setDims(gaugeParam.X);
69 
70  gaugeParam.cpu_prec = QUDA_DOUBLE_PRECISION;
71  gaugeParam.cuda_prec = link_prec;
72  gaugeParam.reconstruct = link_recon;
73  gaugeParam.gauge_order = QUDA_QDP_GAUGE_ORDER;
74  GaugeFieldParam gParam(0, gaugeParam);
78  gParam.anisotropy = 1;
79 
84 
85  // create "gauge fields"
86  int seed=0;
87 #ifdef MULTI_GPU
88  seed += comm_rank();
89 #endif
90 
91  createNoisyLinkCPU((void**)cpuFatLink->Gauge_p(), gaugeParam.cpu_prec, seed);
92  createNoisyLinkCPU((void**)cpuOprod->Gauge_p(), gaugeParam.cpu_prec, seed+1);
93 
94  gParam.setPrecision(gaugeParam.cuda_prec, true);
95 
99 
101 
104 }
105 
106 static void hisq_force_end()
107 {
108  delete cpuFatLink;
109  delete cpuOprod;
110  delete cpuResult;
111 
112  delete cudaFatLink;
113  delete cudaOprod;
114  delete cudaResult;
115 
116  delete cpuReference;
117 
118  endQuda();
119 }
120 
121 static void hisq_force_test()
122 {
123  hisq_force_init();
124 
125  double unitarize_eps = 1e-5;
126  const double hisq_force_filter = 5e-5;
127  const double max_det_error = 1e-12;
128  const bool allow_svd = true;
129  const bool svd_only = false;
130  const double svd_rel_err = 1e-8;
131  const double svd_abs_err = 1e-8;
132 
133  fermion_force::setUnitarizeForceConstants(unitarize_eps, hisq_force_filter, max_det_error, allow_svd, svd_only, svd_rel_err, svd_abs_err);
134 
135  int *num_failures_dev = (int *)device_malloc(sizeof(int));
136  qudaMemset(num_failures_dev, 0, sizeof(int));
137 
138  printfQuda("Calling unitarizeForce\n");
140 
141  device_free(num_failures_dev);
142 
143  if (verify_results) {
144  printfQuda("Calling unitarizeForceCPU\n");
146  }
147 
149 
150  printfQuda("Comparing CPU and GPU results\n");
151  for(int dir=0; dir<4; ++dir){
152  int res = compare_floats(((char **)cpuReference->Gauge_p())[dir], ((char **)cpuResult->Gauge_p())[dir],
154 #ifdef MULTI_GPU
155  comm_allreduce_int(&res);
156  res /= comm_size();
157 #endif
158  printfQuda("Dir:%d Test %s\n",dir,(1 == res) ? "PASSED" : "FAILED");
159  }
160 
161  hisq_force_end();
162 }
163 
164 static void display_test_info()
165 {
166  printfQuda("running the following fermion force computation test:\n");
167 
168  printfQuda("link_precision link_reconstruct space_dim(x/y/z) T_dimension\n");
169  printfQuda("%s %s %d/%d/%d %d \n",
172  xdim, ydim, zdim, tdim);
173 }
174 
175 int main(int argc, char **argv)
176 {
177  auto app = make_app();
178  // app->get_formatter()->column_width(40);
179  // add_eigen_option_group(app);
180  // add_deflation_option_group(app);
181  // add_multigrid_option_group(app);
182  try {
183  app->parse(argc, argv);
184  } catch (const CLI::ParseError &e) {
185  return app->exit(e);
186  }
187 
188  initComms(argc, argv, gridsize_from_cmdline);
189 
191 
193 
194  hisq_force_test();
195 
196  finalizeComms();
197 
198  return EXIT_SUCCESS;
199 }
200 
void display_test_info()
size_t Volume() const
void loadCPUField(const cpuGaugeField &cpu)
Download into this field from a CPU field.
void saveCPUField(cpuGaugeField &cpu) const
Upload from this field into a CPU field.
int comm_rank(void)
int comm_size(void)
void comm_allreduce_int(int *data)
std::shared_ptr< QUDAApp > make_app(std::string app_description, std::string app_name)
QudaReconstructType link_recon
int device_ordinal
int & ydim
bool verify_results
int & zdim
QudaPrecision prec
int & tdim
int & xdim
std::array< int, 4 > gridsize_from_cmdline
int V
Definition: host_utils.cpp:37
void setDims(int *)
Definition: host_utils.cpp:315
enum QudaPrecision_s QudaPrecision
@ QUDA_GHOST_EXCHANGE_NO
Definition: enum_quda.h:508
@ QUDA_DOUBLE_PRECISION
Definition: enum_quda.h:65
@ QUDA_SINGLE_PRECISION
Definition: enum_quda.h:64
@ QUDA_QDP_GAUGE_ORDER
Definition: enum_quda.h:44
@ QUDA_ZERO_FIELD_CREATE
Definition: enum_quda.h:361
@ QUDA_GENERAL_LINKS
Definition: enum_quda.h:25
#define gauge_site_size
Definition: face_gauge.cpp:34
GaugeFieldParam gParam
QudaPrecision link_prec
int main(int argc, char **argv)
cpuGaugeField * cpuOprod
cudaGaugeField * cudaResult
cpuGaugeField * cpuResult
double accuracy
cudaGaugeField * cudaFatLink
cpuGaugeField * cpuFatLink
cudaGaugeField * cudaOprod
cpuGaugeField * cpuReference
void setPrecision(QudaPrecision precision)
void createNoisyLinkCPU(void **field, QudaPrecision prec, int seed)
int compare_floats(void *a, void *b, int len, double epsilon, QudaPrecision precision)
Definition: host_utils.cpp:889
void initComms(int argc, char **argv, std::array< int, 4 > &commDims)
Definition: host_utils.cpp:255
void finalizeComms()
Definition: host_utils.cpp:292
void createSiteLinkCPU(void **link, QudaPrecision precision, int phase)
#define device_malloc(size)
Definition: malloc_quda.h:102
#define device_free(ptr)
Definition: malloc_quda.h:110
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:26
const char * get_recon_str(QudaReconstructType recon)
Definition: misc.cpp:68
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.
void unitarizeForceCPU(GaugeField &newForce, const GaugeField &oldForce, const GaugeField &gauge)
Unitarize the fermion force on CPU.
void unitarizeForce(GaugeField &newForce, const GaugeField &oldForce, const GaugeField &gauge, int *unitarization_failed)
Unitarize the fermion force.
Main header file for the QUDA library.
void initQuda(int device)
void endQuda(void)
#define qudaMemset(ptr, value, count)
Definition: quda_api.h:218
QudaReconstructType reconstruct
Definition: quda.h:49
QudaGaugeFieldOrder gauge_order
Definition: quda.h:42
QudaPrecision cuda_prec
Definition: quda.h:48
int X[4]
Definition: quda.h:35
QudaPrecision cpu_prec
Definition: quda.h:46
QudaGaugeFieldOrder order
Definition: gauge_field.h:51
void setPrecision(QudaPrecision precision, bool force_native=false)
Helper function for setting the precision and corresponding field order for QUDA internal fields.
Definition: gauge_field.h:173
QudaLinkType link_type
Definition: gauge_field.h:53
QudaFieldCreate create
Definition: gauge_field.h:60
QudaGhostExchange ghostExchange
Definition: lattice_field.h:77
#define printfQuda(...)
Definition: util_quda.h:114