QUDA
v1.1.0
A library for QCD on GPUs
|
Functions | |
void | init () |
void | destroy () |
void | setParam (int kernel, int prec, int threads, int blocks) |
void | zero (ColorSpinorField &a) |
void | copy (ColorSpinorField &dst, const ColorSpinorField &src) |
void | ax (double a, ColorSpinorField &x) |
void | axpbyz (double a, ColorSpinorField &x, double b, ColorSpinorField &y, ColorSpinorField &z) |
void | xpy (ColorSpinorField &x, ColorSpinorField &y) |
void | mxpy (ColorSpinorField &x, ColorSpinorField &y) |
void | axpy (double a, ColorSpinorField &x, ColorSpinorField &y) |
void | axpby (double a, ColorSpinorField &x, double b, ColorSpinorField &y) |
void | xpay (ColorSpinorField &x, double a, ColorSpinorField &y) |
void | xpayz (ColorSpinorField &x, double a, ColorSpinorField &y, ColorSpinorField &z) |
void | axpyZpbx (double a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z, double b) |
void | axpyBzpcx (double a, ColorSpinorField &x, ColorSpinorField &y, double b, ColorSpinorField &z, double c) |
void | caxpby (const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y) |
void | caxpy (const Complex &a, ColorSpinorField &x, ColorSpinorField &y) |
void | caxpbypczw (const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y, const Complex &c, ColorSpinorField &z, ColorSpinorField &w) |
void | cxpaypbz (ColorSpinorField &, const Complex &b, ColorSpinorField &y, const Complex &c, ColorSpinorField &z) |
void | caxpbypzYmbw (const Complex &, ColorSpinorField &, const Complex &, ColorSpinorField &, ColorSpinorField &, ColorSpinorField &) |
void | caxpyBzpx (const Complex &, ColorSpinorField &, ColorSpinorField &, const Complex &, ColorSpinorField &) |
void | caxpyBxpz (const Complex &, ColorSpinorField &, ColorSpinorField &, const Complex &, ColorSpinorField &) |
void | cabxpyAx (double a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y) |
void | caxpyXmaz (const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
void | caxpyXmazMR (const double &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
void | tripleCGUpdate (double alpha, double beta, ColorSpinorField &q, ColorSpinorField &r, ColorSpinorField &x, ColorSpinorField &p) |
double | norm1 (const ColorSpinorField &b) |
double | norm2 (const ColorSpinorField &a) |
double | axpyReDot (double a, ColorSpinorField &x, ColorSpinorField &y) |
double | reDotProduct (ColorSpinorField &x, ColorSpinorField &y) |
double | axpbyzNorm (double a, ColorSpinorField &x, double b, ColorSpinorField &y, ColorSpinorField &z) |
double | axpyNorm (double a, ColorSpinorField &x, ColorSpinorField &y) |
double | xmyNorm (ColorSpinorField &x, ColorSpinorField &y) |
Complex | cDotProduct (ColorSpinorField &, ColorSpinorField &) |
double3 | cDotProductNormA (ColorSpinorField &a, ColorSpinorField &b) |
double3 | cDotProductNormB (ColorSpinorField &a, ColorSpinorField &b) |
Return (a,b) and ||b||^2 - implemented using cDotProductNormA. More... | |
double3 | caxpbypzYmbwcDotProductUYNormY (const Complex &a, ColorSpinorField &x, const Complex &b, ColorSpinorField &y, ColorSpinorField &z, ColorSpinorField &w, ColorSpinorField &u) |
double | caxpyNorm (const Complex &a, ColorSpinorField &x, ColorSpinorField &y) |
double | caxpyXmazNormX (const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
double | cabxpyzAxNorm (double a, const Complex &b, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
Complex | caxpyDotzy (const Complex &a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
Complex | axpyCGNorm (double a, ColorSpinorField &x, ColorSpinorField &y) |
double3 | HeavyQuarkResidualNorm (ColorSpinorField &x, ColorSpinorField &r) |
double3 | xpyHeavyQuarkResidualNorm (ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &r) |
double3 | tripleCGReduction (ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
double4 | quadrupleCGReduction (ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
double | quadrupleCG3InitNorm (double a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z, ColorSpinorField &w, ColorSpinorField &v) |
double | quadrupleCG3UpdateNorm (double a, double b, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z, ColorSpinorField &w, ColorSpinorField &v) |
void | axpy (const double *a, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &y) |
Compute the block "axpy" with over the set of ColorSpinorFields. E.g., it computes y = x * a + y The dimensions of a can be rectangular, e.g., the width of x and y need not be same. More... | |
void | axpy (const double *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "axpy" with a composite ColorSpinorField. E.g., it computes y = x * a + y. More... | |
void | axpy_U (const double *a, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &y) |
Compute the block "axpy_U" with over the set of ColorSpinorFields. E.g., it computes. More... | |
void | axpy_U (const double *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "axpy_U" with a composite ColorSpinorField. E.g., it computes. More... | |
void | axpy_L (const double *a, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &y) |
Compute the block "axpy_L" with over the set of ColorSpinorFields. E.g., it computes. More... | |
void | axpy_L (const double *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "axpy_U" with a composite ColorSpinorField. E.g., it computes. More... | |
void | caxpy (const Complex *a, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &y) |
Compute the block "caxpy" with over the set of ColorSpinorFields. E.g., it computes. More... | |
void | caxpy (const Complex *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "caxpy" with a composite ColorSpinorField. E.g., it computes. More... | |
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. More... | |
void | caxpy_U (const Complex *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "caxpy_U" with a composite ColorSpinorField. E.g., it computes. More... | |
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. More... | |
void | caxpy_L (const Complex *a, ColorSpinorField &x, ColorSpinorField &y) |
This is a wrapper for calling the block "caxpy_U" with a composite ColorSpinorField. E.g., it computes. More... | |
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. More... | |
void | caxpyz (const Complex *a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes. More... | |
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. More... | |
void | caxpyz_U (const Complex *a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes. More... | |
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. More... | |
void | caxpyz_L (const Complex *a, ColorSpinorField &x, ColorSpinorField &y, ColorSpinorField &z) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes. More... | |
void | axpyBzpcx (const double *a, std::vector< ColorSpinorField * > &x, std::vector< ColorSpinorField * > &y, const double *b, ColorSpinorField &z, const double *c) |
Compute the vectorized "axpyBzpcx" with over the set of ColorSpinorFields, where the third vector, z, is constant over the batch. E.g., it computes. More... | |
void | caxpyBxpz (const Complex *a_, std::vector< ColorSpinorField * > &x_, ColorSpinorField &y_, const Complex *b_, ColorSpinorField &z_) |
Compute the vectorized "caxpyBxpz" over the set of ColorSpinorFields, where the second and third vector, y and z, is constant over the batch. E.g., it computes. More... | |
void | reDotProduct (double *result, std::vector< ColorSpinorField * > &a, std::vector< ColorSpinorField * > &b) |
void | cDotProduct (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. More... | |
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. This routine is specifically for the case where the result matrix is guarantted to be Hermitian. Requires a.size()==b.size(). More... | |
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. This routine is specifically for the case where the result matrix is guarantted to be Hermitian. Uniquely defined for cases like (p, Ap) where the output is Hermitian, but there's an A-norm instead of an L2 norm. Requires a.size()==b.size(). More... | |
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, and copies b into c. More... | |
Variables | |
unsigned long long | flops |
unsigned long long | bytes |
void quda::blas::ax | ( | double | a, |
ColorSpinorField & | x | ||
) |
|
inline |
Definition at line 44 of file blas_quda.h.
void quda::blas::axpbyz | ( | double | a, |
ColorSpinorField & | x, | ||
double | b, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
double quda::blas::axpbyzNorm | ( | double | a, |
ColorSpinorField & | x, | ||
double | b, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::axpy | ( | const double * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "axpy" with a composite ColorSpinorField. E.g., it computes y = x * a + y.
a[in] | Matrix of real coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::axpy | ( | const double * | a, |
std::vector< ColorSpinorField * > & | x, | ||
std::vector< ColorSpinorField * > & | y | ||
) |
Compute the block "axpy" with over the set of ColorSpinorFields. E.g., it computes y = x * a + y The dimensions of a can be rectangular, e.g., the width of x and y need not be same.
a[in] | Matrix of real coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
|
inline |
Definition at line 43 of file blas_quda.h.
void quda::blas::axpy_L | ( | const double * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "axpy_U" with a composite ColorSpinorField. E.g., it computes.
y = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::axpy_L | ( | const double * | a, |
std::vector< ColorSpinorField * > & | x, | ||
std::vector< ColorSpinorField * > & | y | ||
) |
Compute the block "axpy_L" with over the set of ColorSpinorFields. E.g., it computes.
y = x * a + y
Where 'a' must be a square, lower triangular matrix.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
void quda::blas::axpy_U | ( | const double * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "axpy_U" with a composite ColorSpinorField. E.g., it computes.
y = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::axpy_U | ( | const double * | a, |
std::vector< ColorSpinorField * > & | x, | ||
std::vector< ColorSpinorField * > & | y | ||
) |
Compute the block "axpy_U" with over the set of ColorSpinorFields. E.g., it computes.
y = x * a + y
Where 'a' must be a square, upper triangular matrix.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
void quda::blas::axpyBzpcx | ( | const double * | a, |
std::vector< ColorSpinorField * > & | x, | ||
std::vector< ColorSpinorField * > & | y, | ||
const double * | b, | ||
ColorSpinorField & | z, | ||
const double * | c | ||
) |
Compute the vectorized "axpyBzpcx" with over the set of ColorSpinorFields, where the third vector, z, is constant over the batch. E.g., it computes.
y = a * x + y x = b * z + c * x
The dimensions of a, b, c are the same as the size of x and y, with a maximum size of 16.
a[in] | Array of coefficients |
b[in] | Array of coefficients |
c[in] | Array of coefficients |
x[in,out] | vector of ColorSpinorFields |
y[in,out] | vector of ColorSpinorFields |
z[in] | input ColorSpinorField |
void quda::blas::axpyBzpcx | ( | double | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
double | b, | ||
ColorSpinorField & | z, | ||
double | c | ||
) |
Complex quda::blas::axpyCGNorm | ( | double | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
|
inline |
Definition at line 78 of file blas_quda.h.
double quda::blas::axpyReDot | ( | double | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
void quda::blas::axpyZpbx | ( | double | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z, | ||
double | b | ||
) |
void quda::blas::cabxpyAx | ( | double | a, |
const Complex & | b, | ||
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
double quda::blas::cabxpyzAxNorm | ( | double | a, |
const Complex & | b, | ||
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::caxpby | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
const Complex & | b, | ||
ColorSpinorField & | y | ||
) |
void quda::blas::caxpbypczw | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
const Complex & | b, | ||
ColorSpinorField & | y, | ||
const Complex & | c, | ||
ColorSpinorField & | z, | ||
ColorSpinorField & | w | ||
) |
void quda::blas::caxpbypzYmbw | ( | const Complex & | , |
ColorSpinorField & | , | ||
const Complex & | , | ||
ColorSpinorField & | , | ||
ColorSpinorField & | , | ||
ColorSpinorField & | |||
) |
double3 quda::blas::caxpbypzYmbwcDotProductUYNormY | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
const Complex & | b, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z, | ||
ColorSpinorField & | w, | ||
ColorSpinorField & | u | ||
) |
void quda::blas::caxpy | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
void quda::blas::caxpy | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "caxpy" with a composite ColorSpinorField. E.g., it computes.
y = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::caxpy | ( | const Complex * | a, |
std::vector< ColorSpinorField * > & | x, | ||
std::vector< ColorSpinorField * > & | y | ||
) |
Compute the block "caxpy" with over the set of ColorSpinorFields. E.g., it computes.
y = x * a + y
The dimensions of a can be rectangular, e.g., the width of x and y need not be same.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
void quda::blas::caxpy_L | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "caxpy_U" with a composite ColorSpinorField. E.g., it computes.
y = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::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.
y = x * a + y
Where 'a' must be a square, lower triangular matrix.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
void quda::blas::caxpy_U | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
This is a wrapper for calling the block "caxpy_U" with a composite ColorSpinorField. E.g., it computes.
y = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in,out] | Computed output matrix |
void quda::blas::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.
y = x * a + y
Where 'a' must be a square, upper triangular matrix.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in,out] | vector of input/output ColorSpinorFields |
void quda::blas::caxpyBxpz | ( | const Complex & | , |
ColorSpinorField & | , | ||
ColorSpinorField & | , | ||
const Complex & | , | ||
ColorSpinorField & | |||
) |
void quda::blas::caxpyBxpz | ( | const Complex * | a_, |
std::vector< ColorSpinorField * > & | x_, | ||
ColorSpinorField & | y_, | ||
const Complex * | b_, | ||
ColorSpinorField & | z_ | ||
) |
Compute the vectorized "caxpyBxpz" over the set of ColorSpinorFields, where the second and third vector, y and z, is constant over the batch. E.g., it computes.
y = a * x + y z = b * x + z
The dimensions of a, b are the same as the size of x, with a maximum size of 16.
a[in] | Array of coefficients |
b[in] | Array of coefficients |
x[in] | vector of ColorSpinorFields |
y[in,out] | input ColorSpinorField |
z[in,out] | input ColorSpinorField |
void quda::blas::caxpyBzpx | ( | const Complex & | , |
ColorSpinorField & | , | ||
ColorSpinorField & | , | ||
const Complex & | , | ||
ColorSpinorField & | |||
) |
Complex quda::blas::caxpyDotzy | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
double quda::blas::caxpyNorm | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y | ||
) |
void quda::blas::caxpyXmaz | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::caxpyXmazMR | ( | const double & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
double quda::blas::caxpyXmazNormX | ( | const Complex & | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::caxpyz | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes.
z = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in] | Computed output matrix |
z[out] | vector of input/output ColorSpinorFields |
void quda::blas::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.
z = x * a + y
The dimensions of a can be rectangular, e.g., the width of x and y need not be same, though the maximum width for both is 16.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in] | vector of input ColorSpinorFields |
z[out] | vector of output ColorSpinorFields |
void quda::blas::caxpyz_L | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes.
z = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in] | Computed output matrix |
z[out] | vector of input/output ColorSpinorFields |
void quda::blas::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.
z = x * a + y
Where 'a' is assumed to be lower triangular
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in] | vector of input ColorSpinorFields |
z[out] | vector of output ColorSpinorFields |
void quda::blas::caxpyz_U | ( | const Complex * | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
This is a wrapper for calling the block "caxpyz" with a composite ColorSpinorField. E.g., it computes.
z = x * a + y
a[in] | Matrix of coefficients |
x[in] | Input matrix |
y[in] | Computed output matrix |
z[out] | vector of input/output ColorSpinorFields |
void quda::blas::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.
z = x * a + y
Where 'a' is assumed to be upper triangular.
a[in] | Matrix of coefficients |
x[in] | vector of input ColorSpinorFields |
y[in] | vector of input ColorSpinorFields |
z[out] | vector of output ColorSpinorFields |
Complex quda::blas::cDotProduct | ( | ColorSpinorField & | , |
ColorSpinorField & | |||
) |
void quda::blas::cDotProduct | ( | 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.
result[out] | Matrix of inner product result[i][j] = (a[j],b[i]) |
a[in] | set of input ColorSpinorFields |
b[in] | set of input ColorSpinorFields |
void quda::blas::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, and copies b into c.
result[out] | Matrix of inner product result[i][j] = (a[j],b[i]) |
a[in] | set of input ColorSpinorFields |
b[in] | set of input ColorSpinorFields |
c[out] | set of output ColorSpinorFields |
double3 quda::blas::cDotProductNormA | ( | ColorSpinorField & | a, |
ColorSpinorField & | b | ||
) |
|
inline |
Return (a,b) and ||b||^2 - implemented using cDotProductNormA.
Definition at line 87 of file blas_quda.h.
|
inline |
Definition at line 24 of file blas_quda.h.
void quda::blas::cxpaypbz | ( | ColorSpinorField & | , |
const Complex & | b, | ||
ColorSpinorField & | y, | ||
const Complex & | c, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::destroy | ( | ) |
void quda::blas::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. This routine is specifically for the case where the result matrix is guarantted to be Hermitian. Requires a.size()==b.size().
result[out] | Matrix of inner product result[i][j] = (a[j],b[i]) |
a[in] | set of input ColorSpinorFields |
b[in] | set of input ColorSpinorFields |
void quda::blas::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. This routine is specifically for the case where the result matrix is guarantted to be Hermitian. Uniquely defined for cases like (p, Ap) where the output is Hermitian, but there's an A-norm instead of an L2 norm. Requires a.size()==b.size().
result[out] | Matrix of inner product result[i][j] = (a[j],b[i]) |
a[in] | set of input ColorSpinorFields |
b[in] | set of input ColorSpinorFields |
double3 quda::blas::HeavyQuarkResidualNorm | ( | ColorSpinorField & | x, |
ColorSpinorField & | r | ||
) |
void quda::blas::init | ( | ) |
|
inline |
Definition at line 42 of file blas_quda.h.
double quda::blas::norm1 | ( | const ColorSpinorField & | b | ) |
double quda::blas::norm2 | ( | const ColorSpinorField & | a | ) |
double quda::blas::quadrupleCG3InitNorm | ( | double | a, |
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z, | ||
ColorSpinorField & | w, | ||
ColorSpinorField & | v | ||
) |
double quda::blas::quadrupleCG3UpdateNorm | ( | double | a, |
double | b, | ||
ColorSpinorField & | x, | ||
ColorSpinorField & | y, | ||
ColorSpinorField & | z, | ||
ColorSpinorField & | w, | ||
ColorSpinorField & | v | ||
) |
double4 quda::blas::quadrupleCGReduction | ( | ColorSpinorField & | x, |
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
double quda::blas::reDotProduct | ( | ColorSpinorField & | x, |
ColorSpinorField & | y | ||
) |
void quda::blas::reDotProduct | ( | double * | result, |
std::vector< ColorSpinorField * > & | a, | ||
std::vector< ColorSpinorField * > & | b | ||
) |
void quda::blas::setParam | ( | int | kernel, |
int | prec, | ||
int | threads, | ||
int | blocks | ||
) |
double3 quda::blas::tripleCGReduction | ( | ColorSpinorField & | x, |
ColorSpinorField & | y, | ||
ColorSpinorField & | z | ||
) |
void quda::blas::tripleCGUpdate | ( | double | alpha, |
double | beta, | ||
ColorSpinorField & | q, | ||
ColorSpinorField & | r, | ||
ColorSpinorField & | x, | ||
ColorSpinorField & | p | ||
) |
|
inline |
Definition at line 79 of file blas_quda.h.
|
inline |
Definition at line 45 of file blas_quda.h.
|
inline |
Definition at line 46 of file blas_quda.h.
|
inline |
Definition at line 41 of file blas_quda.h.
double3 quda::blas::xpyHeavyQuarkResidualNorm | ( | ColorSpinorField & | x, |
ColorSpinorField & | y, | ||
ColorSpinorField & | r | ||
) |
void quda::blas::zero | ( | ColorSpinorField & | a | ) |
|
extern |
|
extern |