QUDA  v0.7.0
A library for QCD on GPUs
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Typedefs | Functions
quda.h File Reference

Main header file for the QUDA library. More...

#include <enum_quda.h>
#include <stdio.h>
#include <quda_constants.h>

Go to the source code of this file.

Classes

struct  QudaGaugeParam_s
 
struct  QudaInvertParam_s
 
struct  QudaEigParam_s
 

Typedefs

typedef struct QudaGaugeParam_s QudaGaugeParam
 
typedef struct QudaInvertParam_s QudaInvertParam
 
typedef struct QudaEigParam_s QudaEigParam
 
typedef int(* QudaCommsMap )(const int *coords, void *fdata)
 

Functions

void setVerbosityQuda (QudaVerbosity verbosity, const char prefix[], FILE *outfile)
 
void initCommsGridQuda (int nDim, const int *dims, QudaCommsMap func, void *fdata)
 
void initQudaDevice (int device)
 
void initQudaMemory ()
 
void initQuda (int device)
 
void endQuda (void)
 
QudaGaugeParam newQudaGaugeParam (void)
 
QudaInvertParam newQudaInvertParam (void)
 
QudaEigParam newQudaEigParam (void)
 
void printQudaGaugeParam (QudaGaugeParam *param)
 
void printQudaInvertParam (QudaInvertParam *param)
 
void printQudaEigParam (QudaEigParam *param)
 
void loadGaugeQuda (void *h_gauge, QudaGaugeParam *param)
 
void freeGaugeQuda (void)
 
void saveGaugeQuda (void *h_gauge, QudaGaugeParam *param)
 
void loadCloverQuda (void *h_clover, void *h_clovinv, QudaInvertParam *inv_param)
 
void freeCloverQuda (void)
 
void lanczosQuda (int k0, int m, void *hp_Apsi, void *hp_r, void *hp_V, void *hp_alpha, void *hp_beta, QudaEigParam *eig_param)
 
void invertQuda (void *h_x, void *h_b, QudaInvertParam *param)
 
void invertMultiShiftMDQuda (void **_hp_xe, void **_hp_xo, void **_hp_ye, void **_hp_yo, void *_hp_b, QudaInvertParam *param)
 
void invertMultiShiftQuda (void **_hp_x, void *_hp_b, QudaInvertParam *param)
 
void incrementalEigQuda (void *_h_x, void *_h_b, QudaInvertParam *param, void *_h_u, double *inv_eigenvals, int last_rhs)
 
void dslashQuda (void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
 
void dslashQuda_4dpc (void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
 
void dslashQuda_mdwf (void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity, int test_type)
 
void cloverQuda (void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity *parity, int inverse)
 
void MatQuda (void *h_out, void *h_in, QudaInvertParam *inv_param)
 
void MatDagMatQuda (void *h_out, void *h_in, QudaInvertParam *inv_param)
 
void set_dim (int *)
 
void pack_ghost (void **cpuLink, void **cpuGhost, int nFace, QudaPrecision precision)
 
void setFatLinkPadding (QudaComputeFatMethod method, QudaGaugeParam *param)
 
void computeKSLinkQuda (void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param, QudaComputeFatMethod method)
 
int computeGaugeForceQuda (void *mom, void *sitelink, int ***input_path_buf, int *path_length, double *loop_coeff, int num_paths, int max_length, double dt, QudaGaugeParam *qudaGaugeParam, double *timeinfo)
 
void updateGaugeFieldQuda (void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
 
void * createExtendedGaugeField (void *gauge, int geometry, QudaGaugeParam *param)
 
void * createGaugeField (void *gauge, int geometry, QudaGaugeParam *param)
 
void saveGaugeField (void *outGauge, void *inGauge, QudaGaugeParam *param)
 
void extendGaugeField (void *outGauge, void *inGauge)
 
void destroyQudaGaugeField (void *gauge)
 
void createCloverQuda (QudaInvertParam *param)
 
void computeCloverTraceQuda (void *out, void *clover, int mu, int nu, int dim[4])
 
void computeCloverDerivativeQuda (void *out, void *gauge, void *oprod, int mu, int nu, double coeff, QudaParity parity, QudaGaugeParam *param, int conjugate)
 
void computeStaggeredOprodQuda (void **oprod, void **quark, int num, double **coeff, QudaGaugeParam *param)
 
void computeStaggeredForceQuda (void *mom, void *quark, double *coeff)
 
void computeAsqtadForceQuda (void *const momentum, long long *flops, const double act_path_coeff[6], const void *const one_link_src[4], const void *const naik_src[4], const void *const link, const QudaGaugeParam *param)
 
void computeHISQForceQuda (void *momentum, long long *flops, const double level2_coeff[6], const double fat7_coeff[6], const void *const staple_src[4], const void *const one_link_src[4], const void *const naik_src[4], const void *const w_link, const void *const v_link, const void *const u_link, const QudaGaugeParam *param)
 
void computeHISQForceCompleteQuda (void *momentum, const double level2_coeff[6], const double fat7_coeff[6], void **quark_array, int num_terms, double **quark_coeff, const void *const w_link, const void *const v_link, const void *const u_link, const QudaGaugeParam *param)
 
void performAPEnStep (unsigned int nSteps, double alpha)
 
double plaqCuda ()
 
void openMagma ()
 
void closeMagma ()
 

Detailed Description

Main header file for the QUDA library.

Note to QUDA developers: When adding new members to QudaGaugeParam and QudaInvertParam, be sure to update lib/check_params.h as well as the Fortran interface in lib/quda_fortran.F90.

Definition in file quda.h.

Typedef Documentation

typedef int(* QudaCommsMap)(const int *coords, void *fdata)

initCommsGridQuda() takes an optional "rank_from_coords" argument that should be a pointer to a user-defined function with this prototype.

Parameters
coordsNode coordinates
fdataAny auxiliary data needed by the function
Returns
MPI rank or QMP node ID cooresponding to the node coordinates
See Also
initCommsGridQuda

Definition at line 320 of file quda.h.

typedef struct QudaEigParam_s QudaEigParam

Parameters having to do with the gauge field or the interpretation of the gauge field by various Dirac operators

Parameters relating to the solver and the choice of Dirac operator.

Function Documentation

void closeMagma ( )

Definition at line 98 of file interface_quda.cpp.

void cloverQuda ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param,
QudaParity parity,
int  inverse 
)

Apply the clover operator or its inverse.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage
parityThe source and destination parity of the field
inverseWhether to apply the inverse of the clover term
void computeAsqtadForceQuda ( void *const  momentum,
long long *  flops,
const double  act_path_coeff[6],
const void *const  one_link_src[4],
const void *const  naik_src[4],
const void *const  link,
const QudaGaugeParam param 
)

Compute the fermion force for the asqtad quark action.

Parameters
momentumThe momentum contribution from the quark action.
act_path_coeffThe coefficients that define the asqtad action.
one_link_srcThe quark field outer product corresponding to the one-link term in the action.
naik_srcThe quark field outer product corresponding to the naik term in the action.
linkThe gauge field.
paramThe field parameters.

Definition at line 4235 of file interface_quda.cpp.

void computeCloverDerivativeQuda ( void *  out,
void *  gauge,
void *  oprod,
int  mu,
int  nu,
double  coeff,
QudaParity  parity,
QudaGaugeParam param,
int  conjugate 
)

Definition at line 4064 of file interface_quda.cpp.

void computeCloverTraceQuda ( void *  out,
void *  clover,
int  mu,
int  nu,
int  dim[4] 
)

Definition at line 4041 of file interface_quda.cpp.

int computeGaugeForceQuda ( void *  mom,
void *  sitelink,
int ***  input_path_buf,
int *  path_length,
double *  loop_coeff,
int  num_paths,
int  max_length,
double  dt,
QudaGaugeParam qudaGaugeParam,
double *  timeinfo 
)

Compute the gauge force and update the mometum field

Parameters
momThe momentum field to be updated
sitelinkThe gauge field from which we compute the force
input_path_buf[dim][num_paths][path_length]
path_lengthOne less that the number of links in a loop (e.g., 3 for a staple)
loop_coeffCoefficients of the different loops in the Symanzik action
num_pathsHow many contributions from path_length different "staples"
max_lengthThe maximum number of non-zero of links in any path in the action
dtThe integration step size (for MILC this is dt*beta/3)
paramThe parameters of the external fields and the computation settings
timeinfo

Definition at line 3608 of file interface_quda.cpp.

void computeHISQForceCompleteQuda ( void *  momentum,
const double  level2_coeff[6],
const double  fat7_coeff[6],
void **  quark_array,
int  num_terms,
double **  quark_coeff,
const void *const  w_link,
const void *const  v_link,
const void *const  u_link,
const QudaGaugeParam param 
)

Definition at line 4408 of file interface_quda.cpp.

void computeHISQForceQuda ( void *  momentum,
long long *  flops,
const double  level2_coeff[6],
const double  fat7_coeff[6],
const void *const  staple_src[4],
const void *const  one_link_src[4],
const void *const  naik_src[4],
const void *const  w_link,
const void *const  v_link,
const void *const  u_link,
const QudaGaugeParam param 
)

Compute the fermion force for the HISQ quark action.

Parameters
momentumThe momentum contribution from the quark action.
level2_coeffThe coefficients for the second level of smearing in the quark action.
fat7_coeffThe coefficients for the first level of smearing (fat7) in the quark action.
staple_srcQuark outer-product for the staple.
one_link_srcQuark outer-product for the one-link term in the action.
naik_srcQuark outer-product for the three-hop term in the action.
w_linkUnitarized link variables obtained by applying fat7 smearing and unitarization to the original links.
v_linkFat7 link variables.
u_linkSU(3) think link variables.
param.The field parameters.

Definition at line 4448 of file interface_quda.cpp.

void computeKSLinkQuda ( void *  fatlink,
void *  longlink,
void *  ulink,
void *  inlink,
double *  path_coeff,
QudaGaugeParam param,
QudaComputeFatMethod  method 
)
void computeStaggeredForceQuda ( void *  mom,
void *  quark,
double *  coeff 
)
void computeStaggeredOprodQuda ( void **  oprod,
void **  quark,
int  num,
double **  coeff,
QudaGaugeParam param 
)

Compute the quark-field outer product needed for gauge generation

Parameters
oprodThe outer product to be computed.
quarkThe input fermion field.
displacementThe fermion-field displacement in the outer product.
coeffThe coefficient multiplying the fermion fields in the outer product
paramThe parameters of the outer-product field.

Definition at line 4705 of file interface_quda.cpp.

void createCloverQuda ( QudaInvertParam param)

Definition at line 3776 of file interface_quda.cpp.

void* createExtendedGaugeField ( void *  gauge,
int  geometry,
QudaGaugeParam param 
)

Take a gauge field on the host, extend it and load it onto the device. Return a pointer to the extended gauge field.

Definition at line 3939 of file interface_quda.cpp.

void* createGaugeField ( void *  gauge,
int  geometry,
QudaGaugeParam param 
)

Definition at line 3891 of file interface_quda.cpp.

void destroyQudaGaugeField ( void *  gauge)

Reinterpret gauge as a pointer to cudaGaugeField and call destructor.

Definition at line 4035 of file interface_quda.cpp.

void dslashQuda ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param,
QudaParity  parity 
)

Apply the Dslash operator (D_{eo} or D_{oe}).

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage
parityThe destination parity of the field

Definition at line 1329 of file interface_quda.cpp.

void dslashQuda_4dpc ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param,
QudaParity  parity,
int  test_type 
)

Apply the Dslash operator (D_{eo} or D_{oe}) for 4D EO preconditioned DWF.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage
parityThe destination parity of the field
test_typeChoose a type of dslash operators

Definition at line 1411 of file interface_quda.cpp.

void dslashQuda_mdwf ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param,
QudaParity  parity,
int  test_type 
)

Apply the Dslash operator (D_{eo} or D_{oe}) for Mobius DWF.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage
parityThe destination parity of the field
test_typeChoose a type of dslash operators

Definition at line 1487 of file interface_quda.cpp.

void endQuda ( void  )

Finalize the library.

Definition at line 1018 of file interface_quda.cpp.

void extendGaugeField ( void *  outGauge,
void *  inGauge 
)

Definition at line 4021 of file interface_quda.cpp.

void freeCloverQuda ( void  )

Free QUDA's internal copy of the clover term and/or clover inverse.

Definition at line 996 of file interface_quda.cpp.

void freeGaugeQuda ( void  )

Free QUDA's internal copy of the gauge field.

Definition at line 929 of file interface_quda.cpp.

void incrementalEigQuda ( void *  _h_x,
void *  _h_b,
QudaInvertParam param,
void *  _h_u,
double *  inv_eigenvals,
int  last_rhs 
)

Deflated solvers interface (e.g., based on invremental deflation space constructors, like incremental eigCG).

Parameters
_h_xOutnput: array of solution spinor fields (typically O(10))
_h_bInput: array of source spinor fields (typically O(10))
_h_uInput/Output: array of Ritz spinor fields (typically O(100))
_h_hInput/Output: complex projection mutirx (typically O(100))
paramContains all metadata regarding host and device storage and solver parameters

half precision Dirac field (for the initCG)

Definition at line 3071 of file interface_quda.cpp.

void initCommsGridQuda ( int  nDim,
const int *  dims,
QudaCommsMap  func,
void *  fdata 
)

Declare the grid mapping ("logical topology" in QMP parlance) used for communications in a multi-GPU grid. This function should be called prior to initQuda(). The only case in which it's optional is when QMP is used for communication and the logical topology has already been declared by the application.

Parameters
nDimNumber of grid dimensions. "4" is the only supported value currently.
dimsArray of grid dimensions. dims[0]*dims[1]*dims[2]*dims[3] must equal the total number of MPI ranks or QMP nodes.
funcPointer to a user-supplied function that maps coordinates in the communication grid to MPI ranks (or QMP node IDs). If the pointer is NULL, the default mapping depends on whether QMP or MPI is being used for communication. With QMP, the existing logical topology is used if it's been declared. With MPI or as a fallback with QMP, the default ordering is lexicographical with the fourth ("t") index varying fastest.
fdataPointer to any data required by "func" (may be NULL)
See Also
QudaCommsMap

Definition at line 254 of file interface_quda.cpp.

void initQuda ( int  device)

Initialize the library. This function is actually a wrapper around calls to initQudaDevice() and initQudaMemory().

Parameters
deviceCUDA device number to use. In a multi-GPU build, this parameter may either be set explicitly on a per-process basis or set to -1 to enable a default allocation of devices to processes.

Definition at line 428 of file interface_quda.cpp.

void initQudaDevice ( int  device)

Initialize the library. This is a low-level interface that is called by initQuda. Calling initQudaDevice requires that the user also call initQudaMemory before using QUDA.

Parameters
deviceCUDA device number to use. In a multi-GPU build, this parameter may either be set explicitly on a per-process basis or set to -1 to enable a default allocation of devices to processes.

Definition at line 320 of file interface_quda.cpp.

void initQudaMemory ( )

Initialize the library persistant memory allocations (both host and device). This is a low-level interface that is called by initQuda. Calling initQudaMemory requires that the user has previously called initQudaDevice.

Definition at line 398 of file interface_quda.cpp.

void invertMultiShiftMDQuda ( void **  _hp_xe,
void **  _hp_xo,
void **  _hp_ye,
void **  _hp_yo,
void *  _hp_b,
QudaInvertParam param 
)

Solve for multiple shifts (e.g., masses). This is a special variant of the multi-shift solver where the additional vectors required for force computation are also returned.

Parameters
_hp_xeArray of solution spinor fields
_hp_xoArray of fields with A_oo^{-1} D_oe * x
_hp_yeArray of fields with M_ee * x
_hp_yoArray of fields with A_oo^{-1} D_oe * M_ee * x
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters

Definition at line 2767 of file interface_quda.cpp.

void invertMultiShiftQuda ( void **  _hp_x,
void *  _hp_b,
QudaInvertParam param 
)

Solve for multiple shifts (e.g., masses).

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters

Generic version of the multi-shift solver. Should work for most fermions. Note that offset[0] is not folded into the mass parameter.

At present, the solution_type must be MATDAG_MAT or MATPCDAG_MATPC, and solve_type must be NORMOP or NORMOP_PC. The solution and solve preconditioning have to match.

Definition at line 2473 of file interface_quda.cpp.

void invertQuda ( void *  h_x,
void *  h_b,
QudaInvertParam param 
)

Perform the solve, according to the parameters set in param. It is assumed that the gauge field has already been loaded via loadGaugeQuda().

Parameters
h_xSolution spinor field
h_bSource spinor field
paramContains all metadata regarding host and device storage and solver parameters

Definition at line 1977 of file interface_quda.cpp.

void lanczosQuda ( int  k0,
int  m,
void *  hp_Apsi,
void *  hp_r,
void *  hp_V,
void *  hp_alpha,
void *  hp_beta,
QudaEigParam eig_param 
)

Perform the solve, according to the parameters set in param. It is assumed that the gauge field has already been loaded via loadGaugeQuda().

Parameters
h_xSolution spinor field
h_bSource spinor field
paramContains all metadata regarding host and device storage and solver parameters

Definition at line 1824 of file interface_quda.cpp.

void loadCloverQuda ( void *  h_clover,
void *  h_clovinv,
QudaInvertParam inv_param 
)

Load the clover term and/or the clover inverse from the host. Either h_clover or h_clovinv may be set to NULL.

Parameters
h_cloverBase pointer to host clover field
h_cloverinvBase pointer to host clover inverse field
inv_paramContains all metadata regarding host and device storage

Definition at line 659 of file interface_quda.cpp.

void loadGaugeQuda ( void *  h_gauge,
QudaGaugeParam param 
)

Load the gauge field from the host.

Parameters
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
paramContains all metadata regarding host and device storage

Definition at line 452 of file interface_quda.cpp.

void MatDagMatQuda ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param 
)

Apply M^{}M, possibly even/odd preconditioned.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage

Definition at line 1643 of file interface_quda.cpp.

void MatQuda ( void *  h_out,
void *  h_in,
QudaInvertParam inv_param 
)

Apply the full Dslash matrix, possibly even/odd preconditioned.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage

Definition at line 1567 of file interface_quda.cpp.

QudaEigParam newQudaEigParam ( void  )

A new QudaEigParam should always be initialized immediately after it's defined (and prior to explicitly setting its members) using this function. Typical usage is as follows:

QudaEigParam eig_param = newQudaEigParam();

QudaGaugeParam newQudaGaugeParam ( void  )

A new QudaGaugeParam should always be initialized immediately after it's defined (and prior to explicitly setting its members) using this function. Typical usage is as follows:

QudaGaugeParam gauge_param = newQudaGaugeParam();

QudaInvertParam newQudaInvertParam ( void  )

A new QudaInvertParam should always be initialized immediately after it's defined (and prior to explicitly setting its members) using this function. Typical usage is as follows:

QudaInvertParam invert_param = newQudaInvertParam();

void openMagma ( )

Open/Close MAGMA library

Definition at line 87 of file interface_quda.cpp.

void pack_ghost ( void **  cpuLink,
void **  cpuGhost,
int  nFace,
QudaPrecision  precision 
)
void performAPEnStep ( unsigned int  nSteps,
double  alpha 
)

Definition at line 5210 of file interface_quda.cpp.

double plaqCuda ( )

Definition at line 5176 of file interface_quda.cpp.

void printQudaEigParam ( QudaEigParam param)

Print the members of QudaEigParam.

Parameters
paramThe QudaEigParam whose elements we are to print.

Definition at line 126 of file check_params.h.

void printQudaGaugeParam ( QudaGaugeParam param)

Print the members of QudaGaugeParam.

Parameters
paramThe QudaGaugeParam whose elements we are to print.

Definition at line 40 of file check_params.h.

void printQudaInvertParam ( QudaInvertParam param)

Print the members of QudaGaugeParam.

Parameters
paramThe QudaGaugeParam whose elements we are to print.

Whether to use a pipelined solver

< Number of offsets in the multi-shift solver

< width of domain overlaps

Definition at line 162 of file check_params.h.

void saveGaugeField ( void *  outGauge,
void *  inGauge,
QudaGaugeParam param 
)

Definition at line 3921 of file interface_quda.cpp.

void saveGaugeQuda ( void *  h_gauge,
QudaGaugeParam param 
)

Save the gauge field to the host.

Parameters
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
paramContains all metadata regarding host and device storage

Definition at line 623 of file interface_quda.cpp.

void set_dim ( int *  )
void setFatLinkPadding ( QudaComputeFatMethod  method,
QudaGaugeParam param 
)
void setVerbosityQuda ( QudaVerbosity  verbosity,
const char  prefix[],
FILE *  outfile 
)

Set parameters related to status reporting.

In typical usage, this function will be called once (or not at all) just before the call to initQuda(), but it's valid to call it any number of times at any point during execution. Prior to the first time it's called, the parameters take default values as indicated below.

Parameters
verbosityDefault verbosity, ranging from QUDA_SILENT to QUDA_DEBUG_VERBOSE. Within a solver, this parameter is overridden by the "verbosity" member of QudaInvertParam. The default value is QUDA_SUMMARIZE.
prefixString to prepend to all messages from QUDA. This defaults to the empty string (""), but you may wish to specify something like "QUDA: " to distinguish QUDA's output from that of your application.
outfileFile pointer (such as stdout, stderr, or a handle returned by fopen()) where messages should be printed. The default is stdout.

Definition at line 214 of file interface_quda.cpp.

void updateGaugeFieldQuda ( void *  gauge,
void *  momentum,
double  dt,
int  conj_mom,
int  exact,
QudaGaugeParam param 
)

Evolve the gauge field by step size dt, using the momentum field I.e., Evalulate U(t+dt) = e(dt pi) U(t)

Parameters
gaugeThe gauge field to be updated
momentumThe momentum field
dtThe integration step size step
conj_momWhether to conjugate the momentum matrix
exactWhether to use an exact exponential or Taylor expand
paramThe parameters of the external fields and the computation settings

Definition at line 4936 of file interface_quda.cpp.