QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
gauge_alg_test.cpp
Go to the documentation of this file.
1 #include <stdio.h>
2 #include <stdlib.h>
3 #include <string.h>
4 
5 #include <quda.h>
6 #include <quda_internal.h>
7 #include <gauge_field.h>
8 
9 #include <comm_quda.h>
10 #include <test_util.h>
11 #include <gauge_tools.h>
12 
13 #include <pgauge_monte.h>
14 #include <random_quda.h>
15 #include <unitarization_links.h>
16 
17 #include <qio_field.h>
18 
19 #if defined(QMP_COMMS)
20 #include <qmp.h>
21 #elif defined(MPI_COMMS)
22 #include <mpi.h>
23 #endif
24 
25 #include <gtest/gtest.h>
26 
27 using namespace quda;
28 
29 extern int device;
30 extern int xdim;
31 extern int ydim;
32 extern int zdim;
33 extern int tdim;
34 extern int gridsize_from_cmdline[];
35 extern QudaPrecision prec;
39 extern double anisotropy;
40 extern char latfile[];
41 
44 
45 #define MAX(a,b) ((a)>(b)?(a):(b))
46 #define DABS(a) ((a)<(0.)?(-(a)):(a))
47 
51 
52 
54 
55  gauge_param.X[0] = xdim;
56  gauge_param.X[1] = ydim;
57  gauge_param.X[2] = zdim;
58  gauge_param.X[3] = tdim;
59 
60  gauge_param.anisotropy = anisotropy;
61  gauge_param.type = QUDA_WILSON_LINKS;
62  gauge_param.gauge_order = QUDA_QDP_GAUGE_ORDER;
63  gauge_param.t_boundary = QUDA_PERIODIC_T;
64 
65  gauge_param.cpu_prec = cpu_prec;
66 
67  gauge_param.cuda_prec = cuda_prec;
68  gauge_param.reconstruct = link_recon;
69 
70  gauge_param.cuda_prec_sloppy = cuda_prec_sloppy;
72 
73  gauge_param.gauge_fix = QUDA_GAUGE_FIXED_NO;
74 
75  gauge_param.ga_pad = 0;
76  // For multi-GPU, ga_pad must be large enough to store a time-slice
77 #ifdef MULTI_GPU
78  int x_face_size = gauge_param.X[1]*gauge_param.X[2]*gauge_param.X[3]/2;
79  int y_face_size = gauge_param.X[0]*gauge_param.X[2]*gauge_param.X[3]/2;
80  int z_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[3]/2;
81  int t_face_size = gauge_param.X[0]*gauge_param.X[1]*gauge_param.X[2]/2;
82  int pad_size =MAX(x_face_size, y_face_size);
83  pad_size = MAX(pad_size, z_face_size);
84  pad_size = MAX(pad_size, t_face_size);
85  gauge_param.ga_pad = pad_size;
86 #endif
87 }
88 
89 
90 class GaugeAlgTest : public ::testing::Test {
91  protected:
93  const double unitarize_eps = 1e-14;
94  const double max_error = 1e-10;
95  const int reunit_allow_svd = 1;
96  const int reunit_svd_only = 0;
97  const double svd_rel_error = 1e-6;
98  const double svd_abs_error = 1e-6;
99  setUnitarizeLinksConstants(unitarize_eps, max_error,
100  reunit_allow_svd, reunit_svd_only,
101  svd_rel_error, svd_abs_error);
102 
103  }
104 
107  return false;
108  }
109 
110  bool comparePlaquette(double3 a, double3 b){
111  double a0,a1,a2;
112  a0 = DABS(a.x - b.x);
113  a1=DABS(a.y - b.y);
114  a2=DABS(a.z - b.z);
115  double prec_val = 1.0e-5;
116  if(prec == QUDA_DOUBLE_PRECISION) prec_val = 1.0e-15;
117  if( (a0 < prec_val) && (a1 < prec_val) && (a2 < prec_val) ) return true;
118  return false;
119  }
120 
121  bool CheckDeterminant(double2 detu){
122  double prec_val = 5e-8;
123  if(prec == QUDA_DOUBLE_PRECISION) prec_val = 1.0e-15;
124  if(DABS(1.0 - detu.x) < prec_val && DABS(detu.y) < prec_val) return true;
125  return false;
126  }
127 
128 
129  void CallUnitarizeLinks(cudaGaugeField *cudaInGauge){
130  unitarizeLinks(*cudaInGauge, num_failures_dev);
131  cudaMemcpy(&num_failures, num_failures_dev, sizeof(int), cudaMemcpyDeviceToHost);
132  if(num_failures>0){
133  cudaFree(num_failures_dev);
134  errorQuda("Error in the unitarization\n");
135  exit(1);
136  }
137  cudaMemset(num_failures_dev, 0, sizeof(int));
138  }
139 
140  virtual void SetUp() {
142 
144 
145  //Setup Gauge container!!!!!!
146  param.cpu_prec = prec;
147  param.cpu_prec = prec;
148  param.cuda_prec = prec;
152 
155 
156  param.X[0] = xdim;
157  param.X[1] = ydim;
158  param.X[2] = zdim;
159  param.X[3] = tdim;
160  setDims(param.X);
161 
162  param.anisotropy = 1.0; //don't support anisotropy for now!!!!!!
165  param.ga_pad = 0;
166 
168  gParam.pad = 0;
171  gParam.link_type = param.type;
172  gParam.reconstruct = param.reconstruct;
174 
175 #ifdef MULTI_GPU
176  int y[4];
177  int R[4] = {0,0,0,0};
178  for(int dir=0; dir<4; ++dir) if(comm_dim_partitioned(dir)) R[dir] = 2;
179  for(int dir=0; dir<4; ++dir) y[dir] = param.X[dir] + 2 * R[dir];
180  int pad = 0;
181  GaugeFieldParam gParamEx(y, prec, link_recon,
183  gParamEx.create = QUDA_ZERO_FIELD_CREATE;
184  gParamEx.order = gParam.order;
185  gParamEx.siteSubset = QUDA_FULL_SITE_SUBSET;
186  gParamEx.t_boundary = gParam.t_boundary;
187  gParamEx.nFace = 1;
188  for(int dir=0; dir<4; ++dir) gParamEx.r[dir] = R[dir];
189  cudaInGauge = new cudaGaugeField(gParamEx);
190 #else
191  cudaInGauge = new cudaGaugeField(gParam);
192 #endif
193  // CURAND random generator initialization
194  randstates = new RNG(gParam, 1234);
195  randstates->Init();
196 
197  nsteps = 10;
198  nhbsteps = 4;
199  novrsteps = 4;
200  coldstart = false;
201  beta_value = 6.2;
202 
203  a0.Start(__func__, __FILE__, __LINE__);
204  a1.Start(__func__, __FILE__, __LINE__);
205 
206  cudaMalloc((void**)&num_failures_dev, sizeof(int));
207  cudaMemset(num_failures_dev, 0, sizeof(int));
208  if(num_failures_dev == NULL) errorQuda("cudaMalloc failed for dev_pointer\n");
209  if(link_recon != QUDA_RECONSTRUCT_8 && coldstart) InitGaugeField( *cudaInGauge);
210  else{
211  InitGaugeField( *cudaInGauge, *randstates );
212  }
213  // Reunitarization setup
214  SetReunitarizationConsts();
215  plaquette(*cudaInGauge);
216 
217  for(int step=1; step<=nsteps; ++step){
218  printfQuda("Step %d\n",step);
219  Monte( *cudaInGauge, *randstates, beta_value, nhbsteps, novrsteps);
220  //Reunitarize gauge links...
221  CallUnitarizeLinks(cudaInGauge);
222  plaquette(*cudaInGauge);
223  }
224  a1.Stop(__func__, __FILE__, __LINE__);
225 
226  printfQuda("Time Monte -> %.6f s\n", a1.Last());
227  plaq = plaquette(*cudaInGauge);
228  printfQuda("Plaq: %.16e , %.16e, %.16e\n", plaq.x, plaq.y, plaq.z);
229  }
230 
231  virtual void TearDown() {
232  detu = getLinkDeterminant(*cudaInGauge);
233  double2 tru = getLinkTrace(*cudaInGauge);
234  printfQuda("Det: %.16e:%.16e\n", detu.x, detu.y);
235  printfQuda("Tr: %.16e:%.16e\n", tru.x/3.0, tru.y/3.0);
236 
237 
238  delete cudaInGauge;
239  cudaFree(num_failures_dev);
240  //Release all temporary memory used for data exchange between GPUs in multi-GPU mode
242 
243  a0.Stop(__func__, __FILE__, __LINE__);
244  printfQuda("Time -> %.6f s\n", a0.Last());
245  randstates->Release();
246  delete randstates;
247  }
248 
249 
251 
252  Timer a0,a1;
253  double2 detu;// = getLinkDeterminant(*cudaInGauge);
254  double3 plaq;// = plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
256  int nsteps;
257  int nhbsteps;
259  bool coldstart;
260  double beta_value;
262 
263 };
264 
265 
266 TEST_F(GaugeAlgTest,Generation){
267  detu = getLinkDeterminant(*cudaInGauge);
268  plaq = plaquette(*cudaInGauge);
269  bool testgen = false;
270  //check plaquette value for beta = 6.2
271  if(plaq.x < 0.614 && plaq.x > 0.611 && plaq.y < 0.614 && plaq.y > 0.611) testgen = true;
272 
273  if(testgen){
274  ASSERT_TRUE(CheckDeterminant(detu));
275  }
276 }
277 
278 TEST_F(GaugeAlgTest,Landau_Overrelaxation){
279  const int reunit_interval = 10;
280  printfQuda("Landau gauge fixing with overrelaxation\n");
281  gaugefixingOVR(*cudaInGauge, 4, 100, 10, 1.5, 0, reunit_interval, 1);
282  ASSERT_TRUE(comparePlaquette(plaq, plaquette(*cudaInGauge)));
283 }
284 
285 TEST_F(GaugeAlgTest,Coulomb_Overrelaxation){
286  const int reunit_interval = 10;
287  printfQuda("Coulomb gauge fixing with overrelaxation\n");
288  gaugefixingOVR(*cudaInGauge, 3, 100, 10, 1.5, 0, reunit_interval, 1);
289  ASSERT_TRUE(comparePlaquette(plaq, plaquette(*cudaInGauge)));
290 }
291 
292 TEST_F(GaugeAlgTest,Landau_FFT){
293  if(!checkDimsPartitioned()){
294  printfQuda("Landau gauge fixing with steepest descent method with FFTs\n");
295  gaugefixingFFT(*cudaInGauge, 4, 100, 10, 0.08, 0, 0, 1);
296  ASSERT_TRUE(comparePlaquette(plaq, plaquette(*cudaInGauge)));
297  }
298 }
299 
300 TEST_F(GaugeAlgTest,Coulomb_FFT){
301  if(!checkDimsPartitioned()){
302  printfQuda("Coulomb gauge fixing with steepest descent method with FFTs\n");
303  gaugefixingFFT(*cudaInGauge, 3, 100, 10, 0.08, 0, 0, 1);
304  ASSERT_TRUE(comparePlaquette(plaq, plaquette(*cudaInGauge)));
305  }
306 }
307 
308 
309 int main(int argc, char **argv){
310  // initalize google test, includes command line options
311  ::testing::InitGoogleTest(&argc, argv);
312  // return code for google test
313  int test_rc = 0;
314  xdim=ydim=zdim=tdim=32;
315  int i;
316  for (i=1; i<argc; i++){
317  if(process_command_line_option(argc, argv, &i) == 0){
318  continue;
319  }
320 
321  fprintf(stderr, "ERROR: Invalid option:%s\n", argv[i]);
322  }
323 
324  initComms(argc, argv, gridsize_from_cmdline);
325 
326  initQuda(device);
327  test_rc = RUN_ALL_TESTS();
328  endQuda();
329 
330  finalizeComms();
331 
332  return test_rc;
333 }
static bool reunit_allow_svd
QudaTboundary t_boundary
Definition: gauge_field.h:20
QudaReconstructType reconstruct_sloppy
Definition: quda.h:53
double anisotropy
Definition: quda.h:38
double anisotropy
Definition: test_util.cpp:1650
QudaGhostExchange ghostExchange
Definition: lattice_field.h:76
void endQuda(void)
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:63
QudaPrecision & cuda_prec
bool comparePlaquette(double3 a, double3 b)
QudaGaugeFixed gauge_fix
Definition: quda.h:61
double2 getLinkDeterminant(cudaGaugeField &data)
Calculate the Determinant.
QudaLinkType type
Definition: quda.h:42
#define errorQuda(...)
Definition: util_quda.h:121
void setUnitarizeLinksConstants(double unitarize_eps, double max_error, bool allow_svd, bool svd_only, double svd_rel_error, double svd_abs_error)
int * num_failures_dev
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
virtual void SetUp()
static int R[4]
void finalizeComms()
Definition: test_util.cpp:128
QudaGaugeParam gauge_param
int tdim
Definition: test_util.cpp:1618
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
void gaugefixingOVR(cudaGaugeField &data, const int gauge_dir, const int Nsteps, const int verbose_interval, const double relax_boost, const double tolerance, const int reunit_interval, const int stopWtheta)
Gauge fixing with overrelaxation with support for single and multi GPU.
int ydim
Definition: test_util.cpp:1616
int num_failures
void setDims(int *)
Definition: test_util.cpp:151
virtual void TearDown()
int zdim
Definition: test_util.cpp:1617
QudaPrecision prec
Definition: test_util.cpp:1608
void cpuSetGaugeParam(QudaGaugeParam &gauge_param)
QudaGaugeParam param
Definition: pack_test.cpp:17
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1606
static double svd_rel_error
void initQuda(int device)
void unitarizeLinks(cudaGaugeField &outfield, const cudaGaugeField &infield, int *fails)
int device
Definition: test_util.cpp:1602
void PGaugeExchangeFree()
Release all allocated memory used to exchange data between nodes.
void CallUnitarizeLinks(quda::cudaGaugeField *cudaInGauge)
int xdim
Definition: test_util.cpp:1615
static bool reunit_svd_only
#define MAX(a, b)
QudaReconstructType link_recon
Definition: test_util.cpp:1605
Class declaration to initialize and hold CURAND RNG states.
Definition: random_quda.h:23
QudaPrecision prec_sloppy
Definition: test_util.cpp:1609
QudaGaugeFieldOrder order
Definition: gauge_field.h:17
void CallUnitarizeLinks(cudaGaugeField *cudaInGauge)
double2 getLinkTrace(cudaGaugeField &data)
Calculate the Trace.
QudaPrecision & cuda_prec_sloppy
QudaPrecision cuda_prec_sloppy
Definition: quda.h:52
QudaGaugeParam param
void InitGaugeField(cudaGaugeField &data)
Perform a cold start to the gauge field, identity SU(3) matrix, also fills the ghost links in multi-G...
static double svd_abs_error
bool checkDimsPartitioned()
QudaReconstructType reconstruct
Definition: quda.h:50
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
void SetReunitarizationConsts()
void gaugefixingFFT(cudaGaugeField &data, const int gauge_dir, const int Nsteps, const int verbose_interval, const double alpha, const int autotune, const double tolerance, const int stopWtheta)
Gauge fixing with Steepest descent method with FFTs with support for single GPU only.
QudaPrecision & cpu_prec
bool CheckDeterminant(double2 detu)
GaugeFieldParam gParam
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
static double unitarize_eps
TEST_F(GaugeAlgTest, Generation)
QudaLinkType link_type
Definition: gauge_field.h:19
int main(int argc, char **argv)
cudaGaugeField * cudaInGauge
#define printfQuda(...)
Definition: util_quda.h:115
QudaTboundary t_boundary
Definition: quda.h:45
char latfile[]
Definition: test_util.cpp:1623
QudaReconstructType reconstruct
Definition: gauge_field.h:16
QudaFieldCreate create
Definition: gauge_field.h:26
#define DABS(a)
void Monte(cudaGaugeField &data, RNG &rngstate, double Beta, int nhb, int nover)
Perform heatbath and overrelaxation. Performs nhb heatbath steps followed by nover overrelaxation ste...
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
void initComms(int argc, char **argv, int *const commDims)
Definition: test_util.cpp:88
void setVerbosity(QudaVerbosity verbosity)
Definition: util_quda.cpp:25
double3 plaquette(const GaugeField &U)
Compute the plaquette of the gauge field.
Definition: gauge_plaq.cu:65
QudaPrecision cpu_prec
Definition: quda.h:47
int comm_dim_partitioned(int dim)
QudaGaugeParam newQudaGaugeParam(void)