QUDA  0.9.0
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 namespace quda {
8 
10  bool direct; // whether to create the direct clover
11  bool inverse; // whether to create the inverse clover
12  void *clover;
13  void *norm;
14  void *cloverInv;
15  void *invNorm;
16  double csw;
17  bool twisted; // whether to create twisted mass clover
18  double mu2;
19  double rho;
20 
24  this->precision = precision;
27  }
28 
30  direct(true), inverse(true), clover(nullptr), norm(nullptr),
31  cloverInv(nullptr), invNorm(nullptr), twisted(false), mu2(0.0), rho(0.0) { }
32 
38 
39  CloverFieldParam(const CloverField &field);
40  };
41 
42  std::ostream& operator<<(std::ostream& output, const CloverFieldParam& param);
43 
44  class CloverField : public LatticeField {
45 
46  protected:
47  size_t bytes; // bytes allocated per clover full field
48  size_t norm_bytes; // sizeof each norm full field
49  int length;
51  int nColor;
52  int nSpin;
53 
54  void *clover;
55  void *norm;
56  void *cloverInv;
57  void *invNorm;
58 
59  double csw;
60  bool twisted;
61  double mu2;
62  double rho;
63 
66 
67  mutable double trlog[2];
68 
69  public:
71  virtual ~CloverField();
72 
73  void* V(bool inverse=false) { return inverse ? cloverInv : clover; }
74  void* Norm(bool inverse=false) { return inverse ? invNorm : norm; }
75  const void* V(bool inverse=false) const { return inverse ? cloverInv : clover; }
76  const void* Norm(bool inverse=false) const { return inverse ? invNorm : norm; }
77 
82  bool isNative() const;
83 
87  double* TrLog() const { return trlog; }
88 
92  QudaCloverFieldOrder Order() const { return order; }
93 
97  size_t Bytes() const { return bytes; }
98 
102  size_t NormBytes() const { return norm_bytes; }
103 
107  bool Csw() const { return csw; }
108 
112  bool Twisted() const { return twisted; }
113 
117  double Mu2() const { return mu2; }
118 
123  double Rho() const { return rho; }
124 
129  void setRho(double rho);
130  };
131 
132  class cudaCloverField : public CloverField {
133 
134  private:
135  void *even, *odd;
136  void *evenNorm, *oddNorm;
137 
138  void *evenInv, *oddInv;
140 
141  // computes the clover field given the input gauge field
142  void compute(const cudaGaugeField &gauge);
143 
144 #ifdef USE_TEXTURE_OBJECTS
145  cudaTextureObject_t tex;
146  cudaTextureObject_t normTex;
147  cudaTextureObject_t invTex;
148  cudaTextureObject_t invNormTex;
149  cudaTextureObject_t evenTex;
150  cudaTextureObject_t evenNormTex;
151  cudaTextureObject_t oddTex;
152  cudaTextureObject_t oddNormTex;
153  cudaTextureObject_t evenInvTex;
154  cudaTextureObject_t evenInvNormTex;
155  cudaTextureObject_t oddInvTex;
156  cudaTextureObject_t oddInvNormTex;
157  void createTexObject(cudaTextureObject_t &tex, cudaTextureObject_t &texNorm, void *field, void *norm, bool full);
158  void destroyTexObject();
159 #endif
160 
161  public:
162  // create a cudaCloverField from a CloverFieldParam
164 
165  virtual ~cudaCloverField();
166 
167 #ifdef USE_TEXTURE_OBJECTS
168  const cudaTextureObject_t& Tex() const { return tex; }
169  const cudaTextureObject_t& NormTex() const { return normTex; }
170  const cudaTextureObject_t& InvTex() const { return invTex; }
171  const cudaTextureObject_t& InvNormTex() const { return invNormTex; }
172  const cudaTextureObject_t& EvenTex() const { return evenTex; }
173  const cudaTextureObject_t& EvenNormTex() const { return evenNormTex; }
174  const cudaTextureObject_t& OddTex() const { return oddTex; }
175  const cudaTextureObject_t& OddNormTex() const { return oddNormTex; }
176  const cudaTextureObject_t& EvenInvTex() const { return evenInvTex; }
177  const cudaTextureObject_t& EvenInvNormTex() const { return evenInvNormTex; }
178  const cudaTextureObject_t& OddInvTex() const { return oddInvTex; }
179  const cudaTextureObject_t& OddInvNormTex() const { return oddInvNormTex; }
180 #endif
181 
187  void copy(const CloverField &src, bool inverse=true);
188 
193  void loadCPUField(const cpuCloverField &cpu);
194 
195 
200  void saveCPUField(cpuCloverField &cpu) const;
201 
202  friend class DiracClover;
203  friend class DiracCloverPC;
204  friend struct FullClover;
205  };
206 
207  // this is a place holder for a future host-side clover object
208  class cpuCloverField : public CloverField {
209 
210  private:
211 
212  public:
214  virtual ~cpuCloverField();
215  };
216 
223  double norm1(const CloverField &u, bool inverse=false);
224 
231  double norm2(const CloverField &a, bool inverse=false);
232 
233 
234  // lightweight struct used to send pointers to cuda driver code
235  struct FullClover {
236  void *even;
237  void *odd;
238  void *evenNorm;
239  void *oddNorm;
241  size_t bytes; // sizeof each clover field (per parity)
242  size_t norm_bytes; // sizeof each norm field (per parity)
243  int stride; // stride (volume + pad)
244  double rho; // rho additive factor
245 
246 #ifdef USE_TEXTURE_OBJECTS
247  const cudaTextureObject_t &evenTex;
248  const cudaTextureObject_t &evenNormTex;
249  const cudaTextureObject_t &oddTex;
250  const cudaTextureObject_t &oddNormTex;
251  const cudaTextureObject_t& EvenTex() const { return evenTex; }
252  const cudaTextureObject_t& EvenNormTex() const { return evenNormTex; }
253  const cudaTextureObject_t& OddTex() const { return oddTex; }
254  const cudaTextureObject_t& OddNormTex() const { return oddNormTex; }
255 #endif
256 
257  FullClover(const cudaCloverField &clover, bool inverse=false) :
260 #ifdef USE_TEXTURE_OBJECTS
261  , evenTex(inverse ? clover.evenInvTex : clover.evenTex)
262  , evenNormTex(inverse ? clover.evenInvNormTex : clover.evenNormTex)
263  , oddTex(inverse ? clover.oddInvTex : clover.oddTex)
264  , oddNormTex(inverse ? clover.oddInvNormTex : clover.oddNormTex)
265 #endif
266  {
267  if (inverse) {
268  even = clover.evenInv;
269  evenNorm = clover.evenInvNorm;
270  odd = clover.oddInv;
271  oddNorm = clover.oddInvNorm;
272  } else {
273  even = clover.even;
274  evenNorm = clover.evenNorm;
275  odd = clover.odd;
276  oddNorm = clover.oddNorm;
277  }
278  }
279  };
280 
281 
282  // driver for computing the clover field from the gauge field
283  void computeClover(CloverField &clover, const GaugeField &gauge, double coeff, QudaFieldLocation location);
284 
285 
299  void copyGenericClover(CloverField &out, const CloverField &in, bool inverse,
300  QudaFieldLocation location, void *Out=0, void *In=0, void *outNorm=0, void *inNorm=0);
301 
302 
303 
312  void cloverInvert(CloverField &clover, bool computeTraceLog, QudaFieldLocation location);
313 
320  void cloverRho(CloverField &clover, double rho);
321 
339  void computeCloverForce(GaugeField& force, const GaugeField& U,
340  std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &p,
341  std::vector<double> &coeff);
352  void computeCloverSigmaOprod(GaugeField& oprod,
353  std::vector<ColorSpinorField*> &x,
354  std::vector<ColorSpinorField*> &p,
355  std::vector< std::vector<double> > &coeff);
364  void computeCloverSigmaTrace(GaugeField &output, const CloverField &clover, double coeff);
365 
377  void cloverDerivative(cudaGaugeField &force, cudaGaugeField& gauge, cudaGaugeField& oprod, double coeff, QudaParity parity);
378 
379 } // namespace quda
380 
381 #endif // _CLOVER_QUDA_H
virtual ~CloverField()
QudaCloverFieldOrder order
Definition: clover_field.h:21
void setRho(double rho)
Bakes in the rho factor into the clover field, (for real diagonal additive Hasenbusch), e.g., A + rho.
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...
double * TrLog() const
Definition: clover_field.h:87
cudaCloverField(const CloverFieldParam &param)
enum QudaPrecision_s QudaPrecision
void * V(bool inverse=false)
Definition: clover_field.h:73
const void * src
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.
QudaFieldCreate create
Definition: clover_field.h:22
size_t Bytes() const
Definition: clover_field.h:97
void saveCPUField(cpuCloverField &cpu) const
FullClover(const cudaCloverField &clover, bool inverse=false)
Definition: clover_field.h:257
void loadCPUField(const cpuCloverField &cpu)
QudaPrecision precision
Definition: lattice_field.h:54
QudaCloverFieldOrder order
Definition: clover_field.h:64
QudaCloverFieldOrder Order() const
Definition: clover_field.h:92
CloverField(const CloverFieldParam &param)
QudaFieldCreate create
Definition: clover_field.h:65
std::ostream & operator<<(std::ostream &output, const CloverFieldParam &param)
void computeCloverSigmaTrace(GaugeField &output, const CloverField &clover, double coeff)
Compute the matrix tensor field necessary for the force calculation from the clover trace action...
double norm2(const CloverField &a, bool inverse=false)
QudaGaugeParam param
Definition: pack_test.cpp:17
double norm1(const CloverField &u, bool inverse=false)
enum QudaCloverFieldOrder_s QudaCloverFieldOrder
void cloverInvert(CloverField &clover, bool computeTraceLog, QudaFieldLocation location)
This function compute the Cholesky decomposition of each clover matrix and stores the clover inverse ...
QudaPrecision precision
Definition: clover_field.h:240
bool Csw() const
Definition: clover_field.h:107
void compute(const cudaGaugeField &gauge)
cpuColorSpinorField * in
void setPrecision(QudaPrecision precision)
Definition: clover_field.h:23
static __inline__ size_t p
enum QudaParity_s QudaParity
cpuCloverField(const CloverFieldParam &param)
bool isNative() const
void copy(const CloverField &src, bool inverse=true)
Copy into this CloverField from the generic CloverField src.
enum QudaFieldLocation_s QudaFieldLocation
cpuColorSpinorField * out
bool twisted
Clover coefficient.
Definition: clover_field.h:17
void cloverRho(CloverField &clover, double rho)
This function adds a real scalar onto the clover diagonal (only to the direct field not the inverse) ...
size_t NormBytes() const
Definition: clover_field.h:102
void * Norm(bool inverse=false)
Definition: clover_field.h:74
double Mu2() const
Definition: clover_field.h:117
const void * V(bool inverse=false) const
Definition: clover_field.h:75
static __inline__ dim3 dim3 void size_t cudaStream_t int enum cudaTextureReadMode readMode static __inline__ const struct texture< T, dim, readMode > & tex
double Rho() const
Definition: clover_field.h:123
enum QudaFieldCreate_s QudaFieldCreate
const void * Norm(bool inverse=false) const
Definition: clover_field.h:76
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...
void computeClover(CloverField &clover, const GaugeField &gauge, double coeff, QudaFieldLocation location)
Definition: clover_quda.cu:204
CloverFieldParam(const CloverFieldParam &param)
Definition: clover_field.h:33
QudaParity parity
Definition: covdev_test.cpp:53
#define a
bool Twisted() const
Definition: clover_field.h:112
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...
Definition: copy_clover.cu:175