QUDA  v1.1.0
A library for QCD on GPUs
clover_field.h
Go to the documentation of this file.
1 #ifndef _CLOVER_QUDA_H
2 #define _CLOVER_QUDA_H
3 
4 #include <quda_internal.h>
5 #include <lattice_field.h>
6 
7 #include <comm_key.h>
8 
9 namespace quda {
10 
11  namespace clover
12  {
13 
14  inline bool isNative(QudaCloverFieldOrder order, QudaPrecision precision)
15  {
16  if (precision == QUDA_DOUBLE_PRECISION) {
17  if (order == QUDA_FLOAT2_CLOVER_ORDER) return true;
18  } else if (precision == QUDA_SINGLE_PRECISION || precision == QUDA_HALF_PRECISION
19  || precision == QUDA_QUARTER_PRECISION) {
20  if (order == QUDA_FLOAT4_CLOVER_ORDER) return true;
21  }
22  return false;
23  }
24 
25  } // namespace clover
26 
27  // Prefetch type
28  enum class CloverPrefetchType {
29  BOTH_CLOVER_PREFETCH_TYPE, // clover and inverse
30  CLOVER_CLOVER_PREFETCH_TYPE, // clover only
31  INVERSE_CLOVER_PREFETCH_TYPE, // inverse clover only
33  };
34 
36  bool direct; // whether to create the direct clover
37  bool inverse; // whether to create the inverse clover
38  void *clover;
39  void *norm;
40  void *cloverInv;
41  void *invNorm;
42  double csw;
43  double coeff;
44  bool twisted; // whether to create twisted mass clover
45  double mu2;
46  double rho;
47 
50 
52 
59  void setPrecision(QudaPrecision precision, bool force_native = false)
60  {
61  // is the current status in native field order?
62  bool native = force_native ? true : clover::isNative(order, this->precision);
63  this->precision = precision;
64  this->ghost_precision = precision;
65 
66  if (native) {
68  }
69  }
70 
73  direct(true),
74  inverse(true),
75  clover(nullptr),
76  norm(nullptr),
77  cloverInv(nullptr),
78  invNorm(nullptr),
79  twisted(false),
80  mu2(0.0),
81  rho(0.0),
83  {
84  }
85 
91  norm(param.norm),
95  mu2(param.mu2),
96  rho(param.rho),
98  {
99  }
100 
101  CloverFieldParam(const CloverField &field);
102  };
103 
104  std::ostream& operator<<(std::ostream& output, const CloverFieldParam& param);
105 
106  class CloverField : public LatticeField {
107 
108  protected:
109  size_t bytes; // bytes allocated per clover full field
110  size_t norm_bytes; // sizeof each norm full field
111  size_t length;
112  size_t real_length;
113  int nColor;
114  int nSpin;
115 
116  void *clover;
117  void *norm;
118  void *cloverInv;
119  void *invNorm;
120 
121  double csw;
122  double coeff;
123  bool twisted;
124  double mu2;
125  double rho;
126 
129 
130  mutable double trlog[2];
131 
132  public:
134  virtual ~CloverField();
135 
136  static CloverField *Create(const CloverFieldParam &param);
137 
138  void* V(bool inverse=false) { return inverse ? cloverInv : clover; }
139  void* Norm(bool inverse=false) { return inverse ? invNorm : norm; }
140  const void* V(bool inverse=false) const { return inverse ? cloverInv : clover; }
141  const void* Norm(bool inverse=false) const { return inverse ? invNorm : norm; }
142 
147 
152  bool isNative() const { return clover::isNative(order, precision); }
153 
157  double* TrLog() const { return trlog; }
158 
162  QudaCloverFieldOrder Order() const { return order; }
163 
167  size_t Bytes() const { return bytes; }
168 
172  size_t NormBytes() const { return norm_bytes; }
173 
177  size_t TotalBytes() const
178  {
179  int direct = V(false) ? 1 : 0;
180  int inverse = V(true) ? 1 : 0;
181  return (direct + inverse) * (bytes + norm_bytes);
182  }
183 
187  int Ncolor() const { return nColor; }
188 
192  int Nspin() const { return nSpin; }
193 
197  double Csw() const { return csw; }
198 
202  double Coeff() const { return coeff; }
203 
207  bool Twisted() const { return twisted; }
208 
212  double Mu2() const { return mu2; }
213 
218  double Rho() const { return rho; }
219 
224  void setRho(double rho);
225 
230  double norm1(bool inverse = false) const;
231 
236  double norm2(bool inverse = false) const;
237 
242  double abs_max(bool inverse = false) const;
243 
248  double abs_min(bool inverse = false) const;
249 
250  virtual int full_dim(int d) const { return x[d]; }
251  };
252 
253  class cudaCloverField : public CloverField {
254 
255  private:
256  void *even, *odd;
257  void *evenNorm, *oddNorm;
258 
259  void *evenInv, *oddInv;
260  void *evenInvNorm, *oddInvNorm;
261 
262  // computes the clover field given the input gauge field
263  void compute(const cudaGaugeField &gauge);
264 
265  public:
266  // create a cudaCloverField from a CloverFieldParam
268 
269  virtual ~cudaCloverField();
270 
276  void copy(const CloverField &src, bool inverse=true);
277 
282  void loadCPUField(const cpuCloverField &cpu);
283 
284 
289  void saveCPUField(cpuCloverField &cpu) const;
290 
298  void prefetch(QudaFieldLocation mem_space, qudaStream_t stream = 0) const;
299 
311 
316  virtual void copy_to_buffer(void *buffer) const;
317 
322  virtual void copy_from_buffer(void *buffer);
323 
324  friend class DiracClover;
325  friend class DiracCloverPC;
326  friend class DiracTwistedClover;
327  friend class DiracTwistedCloverPC;
328  friend struct FullClover;
329  };
330 
331  // this is a place holder for a future host-side clover object
332  class cpuCloverField : public CloverField {
333 
334  private:
335 
336  public:
338  virtual ~cpuCloverField();
339 
344  virtual void copy_to_buffer(void *buffer) const;
345 
350  virtual void copy_from_buffer(void *buffer);
351  };
352 
359  double norm1(const CloverField &u, bool inverse=false);
360 
367  double norm2(const CloverField &a, bool inverse=false);
368 
369 
370  // lightweight struct used to send pointers to cuda driver code
371  struct FullClover {
372  void *even;
373  void *odd;
374  void *evenNorm;
375  void *oddNorm;
377  size_t bytes; // sizeof each clover field (per parity)
378  size_t norm_bytes; // sizeof each norm field (per parity)
379  int stride; // stride (volume + pad)
380  double rho; // rho additive factor
381 
382  FullClover(const cudaCloverField &clover, bool inverse = false) :
383  precision(clover.precision),
384  bytes(clover.bytes),
385  norm_bytes(clover.norm_bytes),
386  stride(clover.stride),
387  rho(clover.rho)
388  {
389  if (inverse) {
390  even = clover.evenInv;
391  evenNorm = clover.evenInvNorm;
392  odd = clover.oddInv;
393  oddNorm = clover.oddInvNorm;
394  } else {
395  even = clover.even;
396  evenNorm = clover.evenNorm;
397  odd = clover.odd;
398  oddNorm = clover.oddNorm;
399  }
400  }
401  };
402 
410  void computeClover(CloverField &clover, const GaugeField &fmunu, double coeff);
411 
425  void copyGenericClover(CloverField &out, const CloverField &in, bool inverse,
426  QudaFieldLocation location, void *Out=0, void *In=0, void *outNorm=0, void *inNorm=0);
427 
428 
429 
437  void cloverInvert(CloverField &clover, bool computeTraceLog);
438 
445  void cloverRho(CloverField &clover, double rho);
446 
464  void computeCloverForce(GaugeField& force, const GaugeField& U,
465  std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &p,
466  std::vector<double> &coeff);
478  std::vector<ColorSpinorField*> &x,
479  std::vector<ColorSpinorField*> &p,
480  std::vector< std::vector<double> > &coeff);
489  void computeCloverSigmaTrace(GaugeField &output, const CloverField &clover, double coeff);
490 
502  void cloverDerivative(cudaGaugeField &force, cudaGaugeField& gauge, cudaGaugeField& oprod, double coeff, QudaParity parity);
503 
512  void copyFieldOffset(CloverField &out, const CloverField &in, CommKey offset, QudaPCType pc_type);
513 
518  constexpr bool dynamic_clover_inverse()
519  {
520 #ifdef DYNAMIC_CLOVER
521  return true;
522 #else
523  return false;
524 #endif
525  }
526 
527 } // namespace quda
528 
529 #endif // _CLOVER_QUDA_H
double Rho() const
Definition: clover_field.h:218
void * Norm(bool inverse=false)
Definition: clover_field.h:139
QudaCloverFieldOrder order
Definition: clover_field.h:127
double abs_min(bool inverse=false) const
Compute the absolute minimum of the field.
bool Twisted() const
Definition: clover_field.h:207
double norm1(bool inverse=false) const
Compute the L1 norm of the field.
double * TrLog() const
Definition: clover_field.h:157
QudaCloverFieldOrder Order() const
Definition: clover_field.h:162
double Coeff() const
Definition: clover_field.h:202
static CloverField * Create(const CloverFieldParam &param)
size_t NormBytes() const
Definition: clover_field.h:172
int Nspin() const
Definition: clover_field.h:192
void setRho(double rho)
Bakes in the rho factor into the clover field, (for real diagonal additive Hasenbusch),...
void * V(bool inverse=false)
Definition: clover_field.h:138
bool isNative() const
Definition: clover_field.h:152
QudaFieldCreate create
Definition: clover_field.h:128
virtual int full_dim(int d) const
Definition: clover_field.h:250
const void * V(bool inverse=false) const
Definition: clover_field.h:140
const void * Norm(bool inverse=false) const
Definition: clover_field.h:141
CloverField(const CloverFieldParam &param)
int Ncolor() const
Definition: clover_field.h:187
double Mu2() const
Definition: clover_field.h:212
double abs_max(bool inverse=false) const
Compute the absolute maximum of the field (Linfinity norm)
double norm2(bool inverse=false) const
Compute the L2 norm squared of the field.
virtual ~CloverField()
double Csw() const
Definition: clover_field.h:197
size_t Bytes() const
Definition: clover_field.h:167
size_t TotalBytes() const
Definition: clover_field.h:177
int x[QUDA_MAX_DIM]
QudaPrecision precision
virtual void copy_to_buffer(void *buffer) const
Copy all contents of the field to a host buffer.
cpuCloverField(const CloverFieldParam &param)
virtual void copy_from_buffer(void *buffer)
Copy all contents of the field from a host buffer to this field.
void copy(const CloverField &src, bool inverse=true)
Copy into this CloverField from the generic CloverField src.
void loadCPUField(const cpuCloverField &cpu)
cudaCloverField(const CloverFieldParam &param)
virtual void copy_from_buffer(void *buffer)
Copy all contents of the field from a host buffer to this field.
void saveCPUField(cpuCloverField &cpu) const
virtual void copy_to_buffer(void *buffer) const
Copy all contents of the field to a host buffer.
void prefetch(QudaFieldLocation mem_space, qudaStream_t stream=0) const
If managed memory and prefetch is enabled, prefetch the clover, the norm field (as appropriate),...
QudaParity parity
Definition: covdev_test.cpp:40
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
@ QUDA_FLOAT4_CLOVER_ORDER
Definition: enum_quda.h:257
@ QUDA_FLOAT2_CLOVER_ORDER
Definition: enum_quda.h:256
enum QudaPrecision_s QudaPrecision
@ QUDA_INVALID_FIELD_LOCATION
Definition: enum_quda.h:327
enum QudaPCType_s QudaPCType
@ QUDA_INVALID_PARITY
Definition: enum_quda.h:284
#define QUDA_INVALID_ENUM
Definition: enum_quda.h:4
enum QudaFieldLocation_s QudaFieldLocation
enum QudaFieldCreate_s QudaFieldCreate
@ QUDA_DOUBLE_PRECISION
Definition: enum_quda.h:65
@ QUDA_SINGLE_PRECISION
Definition: enum_quda.h:64
@ QUDA_QUARTER_PRECISION
Definition: enum_quda.h:62
@ QUDA_HALF_PRECISION
Definition: enum_quda.h:63
enum QudaParity_s QudaParity
bool isNative(QudaCloverFieldOrder order, QudaPrecision precision)
Definition: clover_field.h:14
double norm2(const CloverField &a, bool inverse=false)
void computeCloverSigmaTrace(GaugeField &output, const CloverField &clover, double coeff)
Compute the matrix tensor field necessary for the force calculation from the clover trace action....
__device__ __host__ Matrix< T, 3 > inverse(const Matrix< T, 3 > &u)
Definition: quda_matrix.h:605
void copyGenericClover(CloverField &out, const CloverField &in, bool inverse, QudaFieldLocation location, void *Out=0, void *In=0, void *outNorm=0, void *inNorm=0)
This generic function is used for copying the clover field where in the input and output can be in an...
void computeCloverSigmaOprod(GaugeField &oprod, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &p, std::vector< std::vector< double > > &coeff)
Compute the outer product from the solver solution fields arising from the diagonal term of the fermi...
qudaStream_t * stream
double norm1(const CloverField &u, bool inverse=false)
void computeClover(CloverField &clover, const GaugeField &fmunu, double coeff)
Driver for computing the clover field from the field strength tensor.
CloverPrefetchType
Definition: clover_field.h:28
void cloverDerivative(cudaGaugeField &force, cudaGaugeField &gauge, cudaGaugeField &oprod, double coeff, QudaParity parity)
Compute the derivative of the clover matrix in the direction mu,nu and compute the resulting force gi...
void cloverInvert(CloverField &clover, bool computeTraceLog)
This function compute the Cholesky decomposition of each clover matrix and stores the clover inverse ...
constexpr bool dynamic_clover_inverse()
Helper function that returns whether we have enabled dyanmic clover inversion or not.
Definition: clover_field.h:518
void cloverRho(CloverField &clover, double rho)
This function adds a real scalar onto the clover diagonal (only to the direct field not the inverse)
void copyFieldOffset(CloverField &out, const CloverField &in, CommKey offset, QudaPCType pc_type)
This function is used for copying from a source clover field to a destination clover field with an of...
void computeCloverForce(GaugeField &force, const GaugeField &U, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &p, std::vector< double > &coeff)
Compute the force contribution from the solver solution fields.
std::ostream & operator<<(std::ostream &output, const CloverFieldParam &param)
QudaGaugeParam param
Definition: pack_test.cpp:18
cudaStream_t qudaStream_t
Definition: quda_api.h:9
CloverFieldParam(const CloverFieldParam &param)
Definition: clover_field.h:86
bool twisted
Overall clover coefficient.
Definition: clover_field.h:44
double coeff
C_sw clover coefficient.
Definition: clover_field.h:43
QudaCloverFieldOrder order
Definition: clover_field.h:48
QudaFieldLocation location
Definition: clover_field.h:51
QudaFieldCreate create
Definition: clover_field.h:49
void setPrecision(QudaPrecision precision, bool force_native=false)
Helper function for setting the precision and corresponding field order for QUDA internal fields.
Definition: clover_field.h:59
QudaPrecision precision
Definition: clover_field.h:376
FullClover(const cudaCloverField &clover, bool inverse=false)
Definition: clover_field.h:382
QudaPrecision precision
Definition: lattice_field.h:52
QudaPrecision ghost_precision
Definition: lattice_field.h:55