QUDA  1.0.0
contract_test.cpp
Go to the documentation of this file.
1 #include <stdlib.h>
2 #include <stdio.h>
3 #include <time.h>
4 #include <math.h>
5 #include <string.h>
6 
7 #include <util_quda.h>
8 #include <test_util.h>
9 #include <dslash_util.h>
10 #include <contract_reference.h>
11 #include "misc.h"
12 
13 #if defined(QMP_COMMS)
14 #include <qmp.h>
15 #elif defined(MPI_COMMS)
16 #include <mpi.h>
17 #endif
18 
19 // google test
20 #include <gtest/gtest.h>
21 
22 // In a typical application, quda.h is the only QUDA header required.
23 #include <quda.h>
24 #include <color_spinor_field.h>
25 
26 extern int device;
27 extern int xdim;
28 extern int ydim;
29 extern int zdim;
30 extern int tdim;
31 extern int Lsdim;
32 extern int gridsize_from_cmdline[];
33 extern QudaPrecision prec;
37 
39 
40 // If you add a new contraction type, this must be updated++
41 const int NcontractType = 2;
42 
43 extern bool verify_results;
44 
45 extern void usage(char **);
46 
47 namespace quda
48 {
49  extern void setTransferGPU(bool);
50 }
51 
53 {
54  printfQuda("running the following test:\n");
55 
56  printfQuda("prec sloppy_prec S_dimension T_dimension Ls_dimension\n");
57  printfQuda("%s %s %d/%d/%d %d %d\n", get_prec_str(prec), get_prec_str(prec_sloppy),
58  xdim, ydim, zdim, tdim, Lsdim);
59 
60  printfQuda("Grid partition info: X Y Z T\n");
61  printfQuda(" %d %d %d %d\n", dimPartitioned(0), dimPartitioned(1), dimPartitioned(2),
62  dimPartitioned(3));
63  return;
64 }
65 
70 
72 {
73 
74  inv_param.Ls = 1;
75  inv_param.sp_pad = 0;
76  inv_param.cl_pad = 0;
77 
78  inv_param.cpu_prec = cpu_prec;
79  inv_param.cuda_prec = cuda_prec;
82 
84  inv_param.dirac_order = QUDA_DIRAC_ORDER;
85  // Quda performs contractions in Degrand-Rossi gamma basis,
86  // but the user may suppy vectors in any supported order.
88 
91 }
92 
93 const char *prec_str[] = {"single", "double"};
94 
95 // For googletest, names must be non-empty, unique, and may only contain ASCII
96 // alphanumeric characters or underscore.
97 const char *names[] = {"OpenSpin", "DegrandRossi"};
98 
99 int main(int argc, char **argv)
100 {
101 
102  // QUDA initialise
103  //-----------------------------------------------------------------------------
104  for (int i = 1; i < argc; i++) {
105  if (process_command_line_option(argc, argv, &i) == 0) { continue; }
106  printf("ERROR: Invalid option:%s\n", argv[i]);
107  usage(argv);
108  }
109 
110  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
111  initComms(argc, argv, gridsize_from_cmdline);
112  // call srand() with a rank-dependent seed
113  initRand();
115 
116  // initialize the QUDA library
117  initQuda(device);
118  int X[4] = {xdim, ydim, zdim, tdim};
119  setDims(X);
120  setSpinorSiteSize(24);
121  //-----------------------------------------------------------------------------
122 
123  // Start Google Test Suite
124  //-----------------------------------------------------------------------------
125  ::testing::InitGoogleTest(&argc, argv);
126 
128 
129  // Clear previous error state if it exists
130  cudaGetLastError();
131 
132  // Check for correctness
133  if (verify_results) {
134  ::testing::TestEventListeners &listeners = ::testing::UnitTest::GetInstance()->listeners();
135  if (comm_rank() != 0) { delete listeners.Release(listeners.default_result_printer()); }
136  int result = RUN_ALL_TESTS();
137  if (result) warningQuda("Google tests for QUDA contraction failed!");
138  }
139  //-----------------------------------------------------------------------------
140 
141  // finalize the QUDA library
142  endQuda();
143 
144  // finalize the communications layer
145  finalizeComms();
146 
147  return 0;
148 }
149 
150 // Functions used for Google testing
151 //-----------------------------------------------------------------------------
152 
153 // Performs the CPU GPU comparison with the given parameters
154 void test(int contractionType, int Prec)
155 {
156 
158  switch (Prec) {
159  case 0: testPrec = QUDA_SINGLE_PRECISION; break;
160  case 1: testPrec = QUDA_DOUBLE_PRECISION; break;
161  default: errorQuda("Undefined QUDA precision type %d\n", Prec);
162  }
163 
164  int X[4] = {xdim, ydim, zdim, tdim};
165 
167  setInvertParam(inv_param);
168  inv_param.cpu_prec = testPrec;
169  inv_param.cuda_prec = testPrec;
170  inv_param.cuda_prec_sloppy = testPrec;
171  inv_param.cuda_prec_precondition = testPrec;
172 
173  size_t sSize = (testPrec == QUDA_DOUBLE_PRECISION) ? sizeof(double) : sizeof(float);
174  void *spinorX = malloc(V * spinorSiteSize * sSize);
175  void *spinorY = malloc(V * spinorSiteSize * sSize);
176  void *d_result = malloc(2 * V * 16 * sSize);
177 
178  if (testPrec == QUDA_SINGLE_PRECISION) {
179  for (int i = 0; i < V * spinorSiteSize; i++) {
180  ((float *)spinorX)[i] = rand() / (float)RAND_MAX;
181  ((float *)spinorY)[i] = rand() / (float)RAND_MAX;
182  }
183  } else {
184  for (int i = 0; i < V * spinorSiteSize; i++) {
185  ((double *)spinorX)[i] = rand() / (double)RAND_MAX;
186  ((double *)spinorY)[i] = rand() / (double)RAND_MAX;
187  }
188  }
189 
190  // Host side spinor data and result passed to QUDA.
191  // QUDA will allocate GPU memory, transfer the data,
192  // perform the requested contraction, and return the
193  // result in the array 'result'
194  // We then compare the GPU result with a CPU refernce code
195 
197  switch (contractionType) {
198  case 0: cType = QUDA_CONTRACT_TYPE_OPEN; break;
199  case 1: cType = QUDA_CONTRACT_TYPE_DR; break;
200  default: errorQuda("Undefined contraction type %d\n", contractionType);
201  }
202 
203  // Perform GPU contraction.
204  contractQuda(spinorX, spinorY, d_result, cType, &inv_param, X);
205 
206  // Compare each site contraction from the host and device.
207  // It returns the number of faults it detects.
208  int faults = 0;
209  if (testPrec == QUDA_DOUBLE_PRECISION) {
210  faults = contraction_reference((double *)spinorX, (double *)spinorY, (double *)d_result, cType, X);
211  } else {
212  faults = contraction_reference((float *)spinorX, (float *)spinorY, (float *)d_result, cType, X);
213  }
214 
215  printfQuda("Contraction comparison for contraction type %s complete with %d/%d faults\n", get_contract_str(cType),
216  faults, V * 16 * 2);
217 
218  EXPECT_LE(faults, 0) << "CPU and GPU implementations do not agree";
219 
220  free(spinorX);
221  free(spinorY);
222  free(d_result);
223 }
224 
225 // The following tests gets each contraction type and precision using google testing framework
226 using ::testing::Bool;
227 using ::testing::Combine;
228 using ::testing::Range;
229 using ::testing::TestWithParam;
230 using ::testing::Values;
231 
232 class ContractionTest : public ::testing::TestWithParam<::testing::tuple<int, int>>
233 {
234 
235  protected:
236  ::testing::tuple<int, int> param;
237 
238  public:
239  virtual ~ContractionTest() {}
240  virtual void SetUp() { param = GetParam(); }
241 };
242 
243 // Sets up the Google test
245 {
246  int prec = ::testing::get<0>(GetParam());
247  int contractionType = ::testing::get<1>(GetParam());
248  test(contractionType, prec);
249 }
250 
251 // Helper function to construct the test name
252 std::string getContractName(testing::TestParamInfo<::testing::tuple<int, int>> param)
253 {
254  int prec = ::testing::get<0>(param.param);
255  int contractType = ::testing::get<1>(param.param);
256  std::string str(names[contractType]);
257  str += std::string("_");
258  str += std::string(prec_str[prec]);
259  return str; // names[contractType] + "_" + prec_str[prec];
260 }
261 
262 // Instantiate all test cases
263 INSTANTIATE_TEST_SUITE_P(QUDA, ContractionTest, Combine(Range(0, 2), Range(0, NcontractType)), getContractName);
int dimPartitioned(int dim)
Definition: test_util.cpp:1776
QudaDiracFieldOrder dirac_order
Definition: quda.h:219
int comm_rank(void)
Definition: comm_mpi.cpp:82
TEST_P(ContractionTest, verify)
void endQuda(void)
enum QudaPrecision_s QudaPrecision
const char * prec_str[]
#define errorQuda(...)
Definition: util_quda.h:121
QudaPrecision cuda_prec
Definition: quda.h:214
QudaPrecision cpu_prec
Definition: quda.h:213
QudaPrecision & cuda_prec
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
virtual void SetUp()
void setInvertParam(QudaInvertParam &inv_param)
int device
Definition: test_util.cpp:1602
void finalizeComms()
Definition: test_util.cpp:128
const char * get_contract_str(QudaContractType type)
Definition: misc.cpp:970
QudaPrecision prec_sloppy
Definition: test_util.cpp:1609
const char * names[]
const char * get_prec_str(QudaPrecision prec)
Definition: misc.cpp:701
#define spinorSiteSize
QudaPrecision prec_precondition
Definition: test_util.cpp:1611
QudaPrecision & cuda_prec_precondition
void setDims(int *)
Definition: test_util.cpp:151
QudaFieldLocation input_location
Definition: quda.h:99
QudaGaugeParam param
Definition: pack_test.cpp:17
QudaPrecision & cuda_prec_sloppy
int zdim
Definition: test_util.cpp:1617
void initQuda(int device)
QudaFieldLocation output_location
Definition: quda.h:100
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
QudaInvertParam inv_param
Definition: covdev_test.cpp:37
void test(int contractionType, int Prec)
void setTransferGPU(bool)
INSTANTIATE_TEST_SUITE_P(QUDA, ContractionTest, Combine(Range(0, 2), Range(0, NcontractType)), getContractName)
int contraction_reference(Float *spinorX, Float *spinorY, Float *d_result, QudaContractType cType, int X[])
virtual ~ContractionTest()
QudaPrecision cuda_prec_sloppy
Definition: quda.h:215
void setSpinorSiteSize(int n)
Definition: test_util.cpp:211
QudaInvertParam newQudaInvertParam(void)
QudaContractType contract_type
Definition: test_util.cpp:1772
#define warningQuda(...)
Definition: util_quda.h:133
QudaGammaBasis gamma_basis
Definition: quda.h:221
int X[4]
Definition: covdev_test.cpp:70
bool verify_results
Definition: test_util.cpp:1643
int tdim
Definition: test_util.cpp:1618
int Lsdim
Definition: test_util.cpp:1619
int V
Definition: test_util.cpp:27
int xdim
Definition: test_util.cpp:1615
QudaPrecision cuda_prec_precondition
Definition: quda.h:217
void usage(char **)
Definition: test_util.cpp:1783
QudaPrecision prec_null
Definition: test_util.cpp:1612
::testing::tuple< int, int > param
void contractQuda(const ColorSpinorField &x, const ColorSpinorField &y, void *result, QudaContractType cType)
Definition: contract.cu:107
Main header file for the QUDA library.
void initRand()
Definition: test_util.cpp:138
int main(int argc, char **argv)
void display_test_info()
#define printfQuda(...)
Definition: util_quda.h:115
const int NcontractType
int ydim
Definition: test_util.cpp:1616
enum QudaContractType_s QudaContractType
QudaPrecision & cpu_prec
QudaPrecision prec
Definition: test_util.cpp:1608
void initComms(int argc, char **argv, int *const commDims)
Definition: test_util.cpp:88
std::string getContractName(testing::TestParamInfo<::testing::tuple< int, int >> param)
QudaPreserveSource preserve_source
Definition: quda.h:211