QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
plaq_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 "misc.h"
11 
12 #include <qio_field.h>
13 
14 #if defined(QMP_COMMS)
15 #include <qmp.h>
16 #elif defined(MPI_COMMS)
17 #include <mpi.h>
18 #endif
19 
20 // In a typical application, quda.h is the only QUDA header required.
21 #include <quda.h>
22 
23 extern bool tune;
24 extern int device;
25 extern int xdim;
26 extern int ydim;
27 extern int zdim;
28 extern int tdim;
29 extern int gridsize_from_cmdline[];
32 extern QudaPrecision prec;
34 extern double anisotropy;
35 extern double epsilon;
36 
37 extern bool verify_results;
38 
39 extern char latfile[];
40 extern double gaussian_sigma;
41 
43 
44 #define MAX(a, b) ((a) > (b) ? (a) : (b))
45 
49 
51 {
52  gauge_param.X[0] = xdim;
53  gauge_param.X[1] = ydim;
54  gauge_param.X[2] = zdim;
55  gauge_param.X[3] = tdim;
56 
57  gauge_param.anisotropy = anisotropy;
58  gauge_param.type = QUDA_WILSON_LINKS;
59  gauge_param.gauge_order = QUDA_QDP_GAUGE_ORDER;
60  gauge_param.t_boundary = QUDA_PERIODIC_T;
61 
62  gauge_param.cpu_prec = cpu_prec;
63 
64  gauge_param.cuda_prec = cuda_prec;
65  gauge_param.reconstruct = link_recon;
66 
67  gauge_param.cuda_prec_sloppy = cuda_prec_sloppy;
69 
70  gauge_param.gauge_fix = QUDA_GAUGE_FIXED_NO;
71 
72 #ifdef MULTI_GPU
73  int x_face_size = gauge_param.X[1] * gauge_param.X[2] * gauge_param.X[3] / 2;
74  int y_face_size = gauge_param.X[0] * gauge_param.X[2] * gauge_param.X[3] / 2;
75  int z_face_size = gauge_param.X[0] * gauge_param.X[1] * gauge_param.X[3] / 2;
76  int t_face_size = gauge_param.X[0] * gauge_param.X[1] * gauge_param.X[2] / 2;
77  int pad_size = MAX(x_face_size, y_face_size);
78  pad_size = MAX(pad_size, z_face_size);
79  pad_size = MAX(pad_size, t_face_size);
80  gauge_param.ga_pad = pad_size;
81 #else
82  gauge_param.ga_pad = 0;
83 #endif
84 }
85 
86 extern void usage(char **);
87 
88 void plaq_test(int argc, char **argv)
89 {
90 
91  for (int i = 1; i < argc; i++) {
92  if (process_command_line_option(argc, argv, &i) == 0) { continue; }
93  printf("ERROR: Invalid option:%s\n", argv[i]);
94  usage(argv);
95  }
96 
97  // initialize QMP/MPI, QUDA comms grid and RNG (test_util.cpp)
98  initComms(argc, argv, gridsize_from_cmdline);
99 
103 
104  setGaugeParam(gauge_param);
105  setDims(gauge_param.X);
106  size_t gSize = gauge_param.cpu_prec;
107 
108  void *gauge[4];
109 
110  for (int dir = 0; dir < 4; dir++) { gauge[dir] = malloc(V * gaugeSiteSize * gSize); }
111 
112  initQuda(device);
113 
115 
116  // call srand() with a rank-dependent seed
117  initRand();
118 
119  bool load_gauge = strcmp(latfile, "");
120  // load in the command line supplied gauge field
121  if (load_gauge) {
122  read_gauge_field(latfile, gauge, gauge_param.cpu_prec, gauge_param.X, argc, argv);
123  construct_gauge_field(gauge, 2, gauge_param.cpu_prec, &gauge_param);
124  }
125 
126  loadGaugeQuda(gauge, &gauge_param);
127 
128  if (!load_gauge) gaussGaugeQuda(1234, gaussian_sigma);
129 
130  double plaq[3];
131  plaqQuda(plaq);
132  printfQuda("Computed plaquette gauge precise is %16.15e (spatial = %16.15e, temporal = %16.15e)\n", plaq[0], plaq[1],
133  plaq[2]);
134 
135  freeGaugeQuda();
136  endQuda();
137 
138  // release memory
139  for (int dir = 0; dir < 4; dir++) { free(gauge[dir]); }
140 
141  finalizeComms();
142 }
143 
144 int main(int argc, char **argv)
145 {
146 
147  plaq_test(argc, argv);
148 
149  return 0;
150 }
double anisotropy
Definition: test_util.cpp:1650
static size_t gSize
QudaReconstructType reconstruct_sloppy
Definition: quda.h:53
double anisotropy
Definition: quda.h:38
void endQuda(void)
void construct_gauge_field(void **gauge, int type, QudaPrecision precision, QudaGaugeParam *param)
Definition: test_util.cpp:1047
enum QudaPrecision_s QudaPrecision
int ga_pad
Definition: quda.h:63
void usage(char **)
Definition: test_util.cpp:1783
QudaGaugeFixed gauge_fix
Definition: quda.h:61
QudaLinkType type
Definition: quda.h:42
QudaVerbosity verbosity
Definition: test_util.cpp:1614
void setGaugeParam(QudaGaugeParam &gauge_param)
Definition: plaq_test.cpp:50
void loadGaugeQuda(void *h_gauge, QudaGaugeParam *param)
char latfile[]
Definition: test_util.cpp:1623
int process_command_line_option(int argc, char **argv, int *idx)
Definition: test_util.cpp:2019
int ydim
Definition: test_util.cpp:1616
bool verify_results
Definition: test_util.cpp:1643
void plaqQuda(double plaq[3])
void finalizeComms()
Definition: test_util.cpp:128
QudaGaugeParam gauge_param
QudaGaugeFieldOrder gauge_order
Definition: quda.h:43
int main(int argc, char **argv)
Definition: plaq_test.cpp:144
void setDims(int *)
Definition: test_util.cpp:151
void freeGaugeQuda(void)
QudaReconstructType link_recon_sloppy
Definition: test_util.cpp:1606
void initQuda(int device)
void gaussGaugeQuda(unsigned long long seed, double sigma)
Generate Gaussian distributed fields and store in the resident gauge field. We create a Gaussian-dist...
int zdim
Definition: test_util.cpp:1617
QudaReconstructType link_recon
Definition: test_util.cpp:1605
int gridsize_from_cmdline[]
Definition: test_util.cpp:49
QudaPrecision cuda_prec_sloppy
Definition: quda.h:52
#define MAX(a, b)
Definition: plaq_test.cpp:44
QudaReconstructType reconstruct
Definition: quda.h:50
QudaPrecision cuda_prec
Definition: quda.h:49
int X[4]
Definition: quda.h:36
int V
Definition: test_util.cpp:27
int device
Definition: test_util.cpp:1602
void plaq_test(int argc, char **argv)
Definition: plaq_test.cpp:88
bool tune
enum QudaReconstructType_s QudaReconstructType
Main header file for the QUDA library.
void initRand()
Definition: test_util.cpp:138
QudaPrecision prec_sloppy
Definition: test_util.cpp:1609
QudaPrecision & cuda_prec_sloppy
Definition: plaq_test.cpp:48
#define printfQuda(...)
Definition: util_quda.h:115
QudaTboundary t_boundary
Definition: quda.h:45
int tdim
Definition: test_util.cpp:1618
QudaPrecision & cuda_prec
Definition: plaq_test.cpp:47
int xdim
Definition: test_util.cpp:1615
enum QudaVerbosity_s QudaVerbosity
double gaussian_sigma
Definition: test_util.cpp:1625
double epsilon
Definition: test_util.cpp:1649
QudaPrecision cpu_prec
Definition: plaq_test.cpp:46
void initComms(int argc, char **argv, int *const commDims)
Definition: test_util.cpp:88
void read_gauge_field(const char *filename, void *gauge[], QudaPrecision prec, const int *X, int argc, char *argv[])
Definition: qio_field.h:14
void setVerbosity(QudaVerbosity verbosity)
Definition: util_quda.cpp:25
QudaPrecision cpu_prec
Definition: quda.h:47
#define gaugeSiteSize
Definition: face_gauge.cpp:34
QudaGaugeParam newQudaGaugeParam(void)
QudaPrecision prec
Definition: test_util.cpp:1608