QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
gauge_plaq.cu
Go to the documentation of this file.
1 #include <tune_quda.h>
2 #include <gauge_field.h>
3 #include <jitify_helper.cuh>
4 #include <kernels/gauge_plaq.cuh>
5 
6 namespace quda {
7 
8  template<typename Float, typename Gauge>
10 
12  const GaugeField &meta;
13 
14  private:
15  bool tuneGridDim() const { return true; }
16 
17  public:
19  : TunableLocalParity(), arg(arg), meta(meta) {
20 #ifdef JITIFY
21  create_jitify_program("kernels/gauge_plaq.cuh");
22 #endif
23  strcpy(aux,compile_type_str(meta));
24  }
25 
26  ~GaugePlaq () { }
27 
28  void apply(const cudaStream_t &stream){
29  if (meta.Location() == QUDA_CUDA_FIELD_LOCATION){
30  for (int i=0; i<2; i++) ((double*)arg.result_h)[i] = 0.0;
31  TuneParam tp = tuneLaunch(*this, getTuning(), getVerbosity());
32 #ifdef JITIFY
33  using namespace jitify::reflection;
34  jitify_error = program->kernel("quda::computePlaq")
35  .instantiate((int)tp.block.x,Type<Float>(),Type<Gauge>())
36  .configure(tp.grid,tp.block,tp.shared_bytes,stream).launch(arg);
37 #else
38  LAUNCH_KERNEL_LOCAL_PARITY(computePlaq, tp, stream, arg, Float, Gauge);
39 #endif
40  } else {
41  errorQuda("CPU not supported yet\n");
42  }
43  }
44 
45  TuneKey tuneKey() const { return TuneKey(meta.VolString(), typeid(*this).name(), aux); }
46  long long flops() const { return 6ll*2*arg.threads*(3*198+3); }
47  long long bytes() const { return 6ll*2*arg.threads*4*arg.dataOr.Bytes(); }
48  };
49 
50  template<typename Float, typename Gauge>
51  void plaquette(const Gauge dataOr, const GaugeField& data, double2 &plq, QudaFieldLocation location) {
52  GaugePlaqArg<Gauge> arg(dataOr, data);
53  GaugePlaq<Float,Gauge> gaugePlaq(arg, data);
54  gaugePlaq.apply(0);
56  comm_allreduce_array((double*)arg.result_h, 2);
57  for (int i=0; i<2; i++) ((double*)&plq)[i] = ((double*)arg.result_h)[i] / (9.*2*arg.threads*comm_size());
58  }
59 
60  template<typename Float>
61  void plaquette(const GaugeField& data, double2 &plq, QudaFieldLocation location) {
62  INSTANTIATE_RECONSTRUCT(plaquette<Float>, data, plq, location);
63  }
64 
65  double3 plaquette(const GaugeField &data)
66  {
67  double2 plq;
68  QudaFieldLocation location = data.Location();
69  INSTANTIATE_PRECISION(plaquette, data, plq, location);
70  double3 plaq = make_double3(0.5*(plq.x + plq.y), plq.x, plq.y);
71  return plaq;
72  }
73 
74 } // namespace quda
bool tuneGridDim() const
Definition: gauge_plaq.cu:15
#define LAUNCH_KERNEL_LOCAL_PARITY(kernel, tp, stream, arg,...)
QudaVerbosity getVerbosity()
Definition: util_quda.cpp:21
#define errorQuda(...)
Definition: util_quda.h:121
Helper file when using jitify run-time compilation. This file should be included in source code...
void apply(const cudaStream_t &stream)
Definition: gauge_plaq.cu:28
cudaStream_t * stream
long long bytes() const
Definition: gauge_plaq.cu:47
void comm_allreduce_array(double *data, size_t size)
Definition: comm_mpi.cpp:272
const char * VolString() const
const char * compile_type_str(const LatticeField &meta, QudaFieldLocation location_=QUDA_INVALID_FIELD_LOCATION)
Helper function for setting auxilary string.
const GaugeField & meta
Definition: gauge_plaq.cu:12
GaugePlaq(GaugePlaqArg< Gauge > &arg, const GaugeField &meta)
Definition: gauge_plaq.cu:18
int comm_size(void)
Definition: comm_mpi.cpp:88
#define qudaDeviceSynchronize()
#define INSTANTIATE_PRECISION(func, lat,...)
__global__ void computePlaq(GaugePlaqArg< Gauge > arg)
Definition: gauge_plaq.cuh:49
TuneParam & tuneLaunch(Tunable &tunable, QudaTune enabled, QudaVerbosity verbosity)
Definition: tune.cpp:643
CUresult jitify_error
Definition: tune_quda.h:276
QudaFieldLocation Location() const
long long flops() const
Definition: gauge_plaq.cu:46
TuneKey tuneKey() const
Definition: gauge_plaq.cu:45
GaugePlaqArg< Gauge > arg
Definition: gauge_plaq.cu:11
enum QudaFieldLocation_s QudaFieldLocation
#define INSTANTIATE_RECONSTRUCT(func, g,...)
QudaTune getTuning()
Query whether autotuning is enabled or not. Default is enabled but can be overridden by setting QUDA_...
Definition: util_quda.cpp:52
char aux[TuneKey::aux_n]
Definition: tune_quda.h:265
double3 plaquette(const GaugeField &U)
Compute the plaquette of the gauge field.
Definition: gauge_plaq.cu:65