QUDA  0.9.0
blas_quda.h
Go to the documentation of this file.
1 #ifndef _QUDA_BLAS_H
2 #define _QUDA_BLAS_H
3 
4 #include <quda_internal.h>
5 #include <color_spinor_field.h>
6 
7 // ---------- blas_quda.cu ----------
8 
9 // these defitions are used to avoid calling
10 // std::complex<type>::real/imag which have C++11 ABI incompatibility
11 // issues with certain versions of GCC
12 
13 #define REAL(a) (*((double*)&a))
14 #define IMAG(a) (*((double*)&a+1))
15 
16 namespace quda {
17 
18  namespace blas {
19 
20  // creates and destroys reduction buffers
21  void init();
22  void end(void);
23 
24  void* getDeviceReduceBuffer();
26  void* getHostReduceBuffer();
27 
28  void setParam(int kernel, int prec, int threads, int blocks);
29 
30  extern unsigned long long flops;
31  extern unsigned long long bytes;
32 
33  double norm2(const ColorSpinorField &a);
34  double norm1(const ColorSpinorField &b);
35 
36  void zero(ColorSpinorField &a);
37  void copy(ColorSpinorField &dst, const ColorSpinorField &src);
38 
39  double axpyNorm(const double &a, ColorSpinorField &x, ColorSpinorField &y);
40  double axpyReDot(const double &a, ColorSpinorField &x, ColorSpinorField &y);
41 
44 
46 
47  void axpby(const double &a, ColorSpinorField &x, const double &b, ColorSpinorField &y);
48  void axpy(const double &a, ColorSpinorField &x, ColorSpinorField &y);
49  void ax(const double &a, ColorSpinorField &x);
51  void xpay(ColorSpinorField &x, const double &a, ColorSpinorField &y);
52  void xpayz(ColorSpinorField &x, const double &a, ColorSpinorField &y, ColorSpinorField &z);
54 
55  void axpyZpbx(const double &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z, const double &b);
56  void axpyBzpcx(const double &a, ColorSpinorField& x, ColorSpinorField& y, const double &b, ColorSpinorField& z, const double &c);
57 
58  void caxpby(const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y);
64 
67 
72 
73  void cabxpyAx(const double &a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y);
75  void caxpyXmaz(const Complex &a, ColorSpinorField &x,
77  void caxpyXmazMR(const Complex &a, ColorSpinorField &x,
79  double caxpyXmazNormX(const Complex &a, ColorSpinorField &x,
81  double cabxpyAxNorm(const double &a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y);
82 
83  void caxpbypz(const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &,
85  void caxpbypczpw(const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &,
92 
93  void tripleCGUpdate(const double &alpha, const double &beta, ColorSpinorField &q,
110  void caxpy(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y);
111 
123 
136  void caxpy_U(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y);
137 
149 
162  void caxpy_L(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y);
163 
175 
191  void caxpyz(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y, std::vector<ColorSpinorField*> &z);
192 
205 
219  void caxpyz_U(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y, std::vector<ColorSpinorField*> &z);
220 
233 
247  void caxpyz_L(const Complex *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y, std::vector<ColorSpinorField*> &z);
248 
261 
262 
281  void axpyBzpcx(const double *a, std::vector<ColorSpinorField*> &x, std::vector<ColorSpinorField*> &y,
282  const double *b, ColorSpinorField &z, const double *c);
283 
301  void caxpyBxpz(const Complex *a_, std::vector<ColorSpinorField*> &x_, ColorSpinorField &y_,
302  const Complex *b_, ColorSpinorField &z_);
303 
304  void reDotProduct(double* result, std::vector<ColorSpinorField*>& a, std::vector<ColorSpinorField*>& b);
305 
313  void cDotProduct(Complex* result, std::vector<ColorSpinorField*>& a, std::vector<ColorSpinorField*>& b);
314 
325  void hDotProduct(Complex* result, std::vector<ColorSpinorField*>& a, std::vector<ColorSpinorField*>& b);
326 
339  void hDotProduct_Anorm(Complex* result, std::vector<ColorSpinorField*>& a, std::vector<ColorSpinorField*>& b);
340 
341 
350  void cDotProductCopy(Complex* result, std::vector<ColorSpinorField*>& a, std::vector<ColorSpinorField*>& b, std::vector<ColorSpinorField*>& c);
351 
352  } // namespace blas
353 
354 } // namespace quda
355 
356 #endif // _QUDA_BLAS_H
void caxpyz(const Complex *a, std::vector< ColorSpinorField *> &x, std::vector< ColorSpinorField *> &y, std::vector< ColorSpinorField *> &z)
Compute the block "caxpyz" with over the set of ColorSpinorFields. E.g., it computes.
void caxpyz_U(const Complex *a, std::vector< ColorSpinorField *> &x, std::vector< ColorSpinorField *> &y, std::vector< ColorSpinorField *> &z)
Compute the block "caxpyz" with over the set of ColorSpinorFields. E.g., it computes.
void xpay(ColorSpinorField &x, const double &a, ColorSpinorField &y)
Definition: blas_quda.cu:173
void caxpyXmazMR(const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
Definition: blas_quda.cu:583
void * getHostReduceBuffer()
Definition: reduce_quda.cu:75
double3 cDotProductNormA(ColorSpinorField &a, ColorSpinorField &b)
Definition: reduce_quda.cu:572
double caxpyNorm(const Complex &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:402
void setParam(int kernel, int prec, int threads, int blocks)
const void * src
void end(void)
Definition: blas_quda.cu:70
double norm2(const ColorSpinorField &a)
Definition: reduce_quda.cu:241
void init()
Definition: blas_quda.cu:64
Complex cDotProduct(ColorSpinorField &, ColorSpinorField &)
Definition: reduce_quda.cu:500
void cDotProductCopy(Complex *result, std::vector< ColorSpinorField *> &a, std::vector< ColorSpinorField *> &b, std::vector< ColorSpinorField *> &c)
Computes the matrix of inner products between the vector set a and the vector set b...
std::complex< double > Complex
Definition: eig_variables.h:13
void xpayz(ColorSpinorField &x, const double &a, ColorSpinorField &y, ColorSpinorField &z)
Definition: blas_quda.cu:177
double3 xpyHeavyQuarkResidualNorm(ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &r)
Definition: reduce_quda.cu:742
double axpyNorm(const double &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:325
double reDotProduct(ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:277
void copy(ColorSpinorField &dst, const ColorSpinorField &src)
Definition: copy_quda.cu:263
void ax(const double &a, ColorSpinorField &x)
Definition: blas_quda.cu:209
void caxpy_U(const Complex *a, std::vector< ColorSpinorField *> &x, std::vector< ColorSpinorField *> &y)
Compute the block "caxpy_U" with over the set of ColorSpinorFields. E.g., it computes.
void * getMappedHostReduceBuffer()
Definition: reduce_quda.cu:74
double xmyNorm(ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:364
void caxpyBzpx(const Complex &, ColorSpinorField &, ColorSpinorField &, const Complex &, ColorSpinorField &)
Definition: blas_quda.cu:412
void caxpyBxpz(const Complex &, ColorSpinorField &, ColorSpinorField &, const Complex &, ColorSpinorField &)
Definition: blas_quda.cu:438
#define b
void cabxpyAx(const double &a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:484
double cabxpyAxNorm(const double &a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:449
void axpyZpbx(const double &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z, const double &b)
Definition: blas_quda.cu:384
int int int w
void caxpbypzYmbw(const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &, ColorSpinorField &, ColorSpinorField &)
Definition: blas_quda.cu:464
static __inline__ size_t p
Complex axpyCGNorm(const double &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:654
void tripleCGUpdate(const double &alpha, const double &beta, ColorSpinorField &q, ColorSpinorField &r, ColorSpinorField &x, ColorSpinorField &p)
Definition: blas_quda.cu:610
double4 quadrupleCGReduction(ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
double3 HeavyQuarkResidualNorm(ColorSpinorField &x, ColorSpinorField &r)
Definition: reduce_quda.cu:703
void hDotProduct_Anorm(Complex *result, std::vector< ColorSpinorField *> &a, std::vector< ColorSpinorField *> &b)
Computes the matrix of inner products between the vector set a and the vector set b...
double caxpyXmazNormX(const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
Definition: reduce_quda.cu:424
Complex caxpyDotzy(const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
Definition: reduce_quda.cu:544
void caxpy(const Complex &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:246
void zero(ColorSpinorField &a)
Definition: blas_quda.cu:45
void caxpy_L(const Complex *a, std::vector< ColorSpinorField *> &x, std::vector< ColorSpinorField *> &y)
Compute the block "caxpy_L" with over the set of ColorSpinorFields. E.g., it computes.
void caxpbypczpw(const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &, ColorSpinorField &)
Definition: blas_quda.cu:527
void axpy(const double &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:150
void axpby(const double &a, ColorSpinorField &x, const double &b, ColorSpinorField &y)
Definition: blas_quda.cu:106
double3 caxpbypzYmbwcDotProductUYNormY(const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y, ColorSpinorField &z, ColorSpinorField &w, ColorSpinorField &u)
Definition: reduce_quda.cu:619
void caxpyz_L(const Complex *a, std::vector< ColorSpinorField *> &x, std::vector< ColorSpinorField *> &y, std::vector< ColorSpinorField *> &z)
Compute the block "caxpyz" with over the set of ColorSpinorFields. E.g., it computes.
void caxpbypz(const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &, ColorSpinorField &)
Definition: blas_quda.cu:505
double norm1(const ColorSpinorField &b)
Definition: reduce_quda.cu:200
void hDotProduct(Complex *result, std::vector< ColorSpinorField *> &a, std::vector< ColorSpinorField *> &b)
Computes the matrix of inner products between the vector set a and the vector set b...
double axpyReDot(const double &a, ColorSpinorField &x, ColorSpinorField &y)
Definition: reduce_quda.cu:345
void axpyBzpcx(const double &a, ColorSpinorField &x, ColorSpinorField &y, const double &b, ColorSpinorField &z, const double &c)
Definition: blas_quda.cu:356
void caxpyXmaz(const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
Definition: blas_quda.cu:549
Complex xpaycDotzy(ColorSpinorField &x, const double &a, ColorSpinorField &y, ColorSpinorField &z)
Definition: reduce_quda.cu:521
void * getDeviceReduceBuffer()
Definition: reduce_quda.cu:73
unsigned long long flops
Definition: blas_quda.cu:42
void xpy(ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:128
void caxpby(const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y)
Definition: blas_quda.cu:292
const void * c
void mxpy(ColorSpinorField &x, ColorSpinorField &y)
Definition: blas_quda.cu:192
void cxpaypbz(ColorSpinorField &, const Complex &b, ColorSpinorField &y, const Complex &c, ColorSpinorField &z)
Definition: blas_quda.cu:335
double3 cDotProductNormB(ColorSpinorField &a, ColorSpinorField &b)
Definition: reduce_quda.cu:599
QudaPrecision prec
Definition: test_util.cpp:1615
#define a
double2 reDotProductNormA(ColorSpinorField &a, ColorSpinorField &b)
Definition: reduce_quda.cu:305
unsigned long long bytes
Definition: blas_quda.cu:43
double3 tripleCGReduction(ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z)
Definition: reduce_quda.cu:767