QUDA  0.9.0
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 
18 #include <gtest.h>
19 
20 using namespace quda;
21 
22 extern int device;
23 extern int xdim;
24 extern int ydim;
25 extern int zdim;
26 extern int tdim;
27 extern int gridsize_from_cmdline[];
29 extern QudaPrecision prec;
30 extern char latfile[];
31 
34 
35 #define MAX(a,b) ((a)>(b)?(a):(b))
36 #define DABS(a) ((a)<(0.)?(-(a)):(a))
37 
38 
39 
40 
41 
42 class GaugeAlgTest : public ::testing::Test {
43  protected:
45  const double unitarize_eps = 1e-14;
46  const double max_error = 1e-10;
47  const int reunit_allow_svd = 1;
48  const int reunit_svd_only = 0;
49  const double svd_rel_error = 1e-6;
50  const double svd_abs_error = 1e-6;
54 
55  }
56 
59  return false;
60  }
61 
62  bool comparePlaquette(double3 a, double3 b){
63  double a0,a1,a2;
64  a0 = DABS(a.x - b.x);
65  a1=DABS(a.y - b.y);
66  a2=DABS(a.z - b.z);
67  double prec_val = 1.0e-5;
68  if(prec == QUDA_DOUBLE_PRECISION) prec_val = 1.0e-15;
69  if( (a0 < prec_val) && (a1 < prec_val) && (a2 < prec_val) ) return true;
70  return false;
71  }
72 
73  bool CheckDeterminant(double2 detu){
74  double prec_val = 5e-8;
75  if(prec == QUDA_DOUBLE_PRECISION) prec_val = 1.0e-15;
76  if(DABS(1.0 - detu.x) < prec_val && DABS(detu.y) < prec_val) return true;
77  return false;
78  }
79 
80 
81  void CallUnitarizeLinks(cudaGaugeField *cudaInGauge){
82  unitarizeLinks(*cudaInGauge, num_failures_dev);
83  cudaMemcpy(&num_failures, num_failures_dev, sizeof(int), cudaMemcpyDeviceToHost);
84  if(num_failures>0){
85  cudaFree(num_failures_dev);
86  errorQuda("Error in the unitarization\n");
87  exit(1);
88  }
89  cudaMemset(num_failures_dev, 0, sizeof(int));
90  }
91 
92  virtual void SetUp() {
94 
96 
97  //Setup Gauge container!!!!!!
100  param.cuda_prec = prec;
104 
107 
108  param.X[0] = xdim;
109  param.X[1] = ydim;
110  param.X[2] = zdim;
111  param.X[3] = tdim;
112  setDims(param.X);
113 
114  param.anisotropy = 1.0; //don't support anisotropy for now!!!!!!
117  param.ga_pad = 0;
118 
120  gParam.pad = 0;
126 
127 
128 #ifdef MULTI_GPU
129  int y[4];
130  int R[4] = {0,0,0,0};
131  for(int dir=0; dir<4; ++dir) if(comm_dim_partitioned(dir)) R[dir] = 2;
132  for(int dir=0; dir<4; ++dir) y[dir] = param.X[dir] + 2 * R[dir];
133  int pad = 0;
134  GaugeFieldParam gParamEx(y, prec, link_recon,
136  gParamEx.create = QUDA_ZERO_FIELD_CREATE;
137  gParamEx.order = gParam.order;
138  gParamEx.siteSubset = QUDA_FULL_SITE_SUBSET;
139  gParamEx.t_boundary = gParam.t_boundary;
140  gParamEx.nFace = 1;
141  for(int dir=0; dir<4; ++dir) gParamEx.r[dir] = R[dir];
142  cudaInGauge = new cudaGaugeField(gParamEx);
143 #else
144  cudaInGauge = new cudaGaugeField(gParam);
145 #endif
146  int halfvolume = xdim*ydim*zdim*tdim >> 1;
147  printfQuda("xdim=%d\tydim=%d\tzdim=%d\ttdim=%d\trng_size=%d\n",xdim,ydim,zdim,tdim,halfvolume);
148  // CURAND random generator initialization
149  randstates = new RNG(halfvolume, 1234, param.X);
150  randstates->Init();
151 
152  nsteps = 10;
153  nhbsteps = 4;
154  novrsteps = 4;
155  coldstart = false;
156  beta_value = 6.2;
157 
158 
159 
160  a0.Start(__func__, __FILE__, __LINE__);
161  a1.Start(__func__, __FILE__, __LINE__);
162 
163  cudaMalloc((void**)&num_failures_dev, sizeof(int));
164  cudaMemset(num_failures_dev, 0, sizeof(int));
165  if(num_failures_dev == NULL) errorQuda("cudaMalloc failed for dev_pointer\n");
166  if(link_recon != QUDA_RECONSTRUCT_8 && coldstart) InitGaugeField( *cudaInGauge);
167  else{
168  InitGaugeField( *cudaInGauge, *randstates );
169  }
170  // Reunitarization setup
171  SetReunitarizationConsts();
172  plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
173 
174  for(int step=1; step<=nsteps; ++step){
175  printfQuda("Step %d\n",step);
176  Monte( *cudaInGauge, *randstates, beta_value, nhbsteps, novrsteps);
177  //Reunitarize gauge links...
178  CallUnitarizeLinks(cudaInGauge);
179  plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
180  }
181  a1.Stop(__func__, __FILE__, __LINE__);
182 
183  printfQuda("Time Monte -> %.6f s\n", a1.Last());
184  plaq = plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
185  printfQuda("Plaq: %.16e , %.16e, %.16e\n", plaq.x, plaq.y, plaq.z);
186 
187  }
188 
189  virtual void TearDown() {
190  detu = getLinkDeterminant(*cudaInGauge);
191  double2 tru = getLinkTrace(*cudaInGauge);
192  printfQuda("Det: %.16e:%.16e\n", detu.x, detu.y);
193  printfQuda("Tr: %.16e:%.16e\n", tru.x/3.0, tru.y/3.0);
194 
195 
196  delete cudaInGauge;
197  cudaFree(num_failures_dev);
198  //Release all temporary memory used for data exchange between GPUs in multi-GPU mode
200 
201  a0.Stop(__func__, __FILE__, __LINE__);
202  printfQuda("Time -> %.6f s\n", a0.Last());
203  randstates->Release();
204  delete randstates;
205  }
206 
207 
209 
210  Timer a0,a1;
211  double2 detu;// = getLinkDeterminant(*cudaInGauge);
212  double3 plaq;// = plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
214  int nsteps;
215  int nhbsteps;
217  bool coldstart;
218  double beta_value;
220 
221 };
222 
223 
224 TEST_F(GaugeAlgTest,Generation){
225  detu = getLinkDeterminant(*cudaInGauge);
226  plaq = plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION) ;
227  bool testgen = false;
228  //check plaquette value for beta = 6.2
229  if(plaq.x < 0.614 && plaq.x > 0.611 && plaq.y < 0.614 && plaq.y > 0.611) testgen = true;
230 
231  if(testgen){
232  ASSERT_TRUE(CheckDeterminant(detu));
233  }
234 }
235 
236 TEST_F(GaugeAlgTest,Landau_Overrelaxation){
237  const int reunit_interval = 10;
238  printfQuda("Landau gauge fixing with overrelaxation\n");
239  gaugefixingOVR(*cudaInGauge, 4, 100, 10, 1.5, 0, reunit_interval, 1);
240  ASSERT_TRUE(comparePlaquette(plaq, plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION)));
241 }
242 
243 TEST_F(GaugeAlgTest,Coulomb_Overrelaxation){
244  const int reunit_interval = 10;
245  printfQuda("Coulomb gauge fixing with overrelaxation\n");
246  gaugefixingOVR(*cudaInGauge, 3, 100, 10, 1.5, 0, reunit_interval, 1);
247  ASSERT_TRUE(comparePlaquette(plaq, plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION)));
248 }
249 
250 TEST_F(GaugeAlgTest,Landau_FFT){
251  if(!checkDimsPartitioned()){
252  printfQuda("Landau gauge fixing with steepest descent method with FFTs\n");
253  gaugefixingFFT(*cudaInGauge, 4, 100, 10, 0.08, 0, 0, 1);
254  ASSERT_TRUE(comparePlaquette(plaq, plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION)));
255  }
256 }
257 
258 TEST_F(GaugeAlgTest,Coulomb_FFT){
259  if(!checkDimsPartitioned()){
260  printfQuda("Coulomb gauge fixing with steepest descent method with FFTs\n");
261  gaugefixingFFT(*cudaInGauge, 3, 100, 10, 0.08, 0, 0, 1);
262  ASSERT_TRUE(comparePlaquette(plaq, plaquette( *cudaInGauge, QUDA_CUDA_FIELD_LOCATION)));
263  }
264 }
265 
266 
267 
268 
269 
270 
271 
272 int main(int argc, char **argv){
273  // initalize google test, includes command line options
274  ::testing::InitGoogleTest(&argc, argv);
275  // return code for google test
276  int test_rc = 0;
277  xdim=ydim=zdim=tdim=32;
278  int i;
279  for (i=1; i<argc; i++){
280  if(process_command_line_option(argc, argv, &i) == 0){
281  continue;
282  }
283 
284  fprintf(stderr, "ERROR: Invalid option:%s\n", argv[i]);
285  }
286 
287  initComms(argc, argv, gridsize_from_cmdline);
288 
289  initQuda(device);
290  test_rc = RUN_ALL_TESTS();
291  endQuda();
292 
293  finalizeComms();
294 
295  return test_rc;
296 }
QudaTboundary t_boundary
Definition: gauge_field.h:18
QudaReconstructType reconstruct_sloppy
Definition: quda.h:46
double anisotropy
Definition: quda.h:31
QudaGhostExchange ghostExchange
Definition: lattice_field.h:60
void endQuda(void)
double3 plaquette(const GaugeField &U, QudaFieldLocation location)
Definition: gauge_plaq.cu:138
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:53
bool comparePlaquette(double3 a, double3 b)
QudaGaugeFixed gauge_fix
Definition: quda.h:51
double2 getLinkDeterminant(cudaGaugeField &data)
Calculate the Determinant.
QudaLinkType type
Definition: quda.h:35
#define errorQuda(...)
Definition: util_quda.h:90
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:1795
virtual void SetUp()
void CallUnitarizeLinks(quda::cudaGaugeField *cudaInGauge)
static int R[4]
void finalizeComms()
Definition: test_util.cpp:107
int tdim
Definition: test_util.cpp:1623
QudaGaugeFieldOrder gauge_order
Definition: quda.h:36
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:1621
int num_failures
void exit(int) __attribute__((noreturn))
void setDims(int *)
Definition: test_util.cpp:130
virtual void TearDown()
int zdim
Definition: test_util.cpp:1622
QudaPrecision prec
Definition: test_util.cpp:1615
QudaGaugeParam param
Definition: pack_test.cpp:17
#define b
void initQuda(int device)
void unitarizeLinks(cudaGaugeField &outfield, const cudaGaugeField &infield, int *fails)
void PGaugeExchangeFree()
Release all allocated memory used to exchange data between nodes.
int xdim
Definition: test_util.cpp:1620
QudaReconstructType link_recon
Definition: test_util.cpp:1612
Class declaration to initialize and hold CURAND RNG states.
Definition: random_quda.h:23
QudaGaugeFieldOrder order
Definition: gauge_field.h:15
void CallUnitarizeLinks(cudaGaugeField *cudaInGauge)
double2 getLinkTrace(cudaGaugeField &data)
Calculate the Trace.
QudaPrecision cuda_prec_sloppy
Definition: quda.h:45
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...
bool checkDimsPartitioned()
QudaReconstructType reconstruct
Definition: quda.h:43
QudaPrecision cuda_prec
Definition: quda.h:42
int X[4]
Definition: quda.h:29
void SetReunitarizationConsts()
int fprintf(FILE *, const char *,...) __attribute__((__format__(__printf__
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.
bool CheckDeterminant(double2 detu)
GaugeFieldParam gParam
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
TEST_F(GaugeAlgTest, Generation)
QudaLinkType link_type
Definition: gauge_field.h:17
int main(int argc, char **argv)
cudaGaugeField * cudaInGauge
#define printfQuda(...)
Definition: util_quda.h:84
QudaTboundary t_boundary
Definition: quda.h:38
char latfile[]
Definition: test_util.cpp:1627
QudaReconstructType reconstruct
Definition: gauge_field.h:14
QudaFieldCreate create
Definition: gauge_field.h:25
#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:50
#define a
void initComms(int argc, char **argv, const int *commDims)
Definition: test_util.cpp:72
void setVerbosity(const QudaVerbosity verbosity)
Definition: util_quda.cpp:24
QudaPrecision cpu_prec
Definition: quda.h:40
int comm_dim_partitioned(int dim)
QudaGaugeParam newQudaGaugeParam(void)