QUDA  v1.1.0
A library for QCD on GPUs
Classes | Macros | Typedefs | Functions
quda.h File Reference

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

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

Go to the source code of this file.

Classes

struct  QudaGaugeParam_s
 
struct  QudaInvertParam_s
 
struct  QudaEigParam_s
 
struct  QudaMultigridParam_s
 
struct  QudaGaugeObservableParam_s
 
struct  QudaBLASParam_s
 

Macros

#define double_complex   double _Complex
 

Typedefs

typedef struct QudaGaugeParam_s QudaGaugeParam
 
typedef struct QudaInvertParam_s QudaInvertParam
 
typedef struct QudaEigParam_s QudaEigParam
 
typedef struct QudaMultigridParam_s QudaMultigridParam
 
typedef struct QudaGaugeObservableParam_s QudaGaugeObservableParam
 
typedef struct QudaBLASParam_s QudaBLASParam
 
typedef int(* QudaCommsMap) (const int *coords, void *fdata)
 

Functions

void setVerbosityQuda (QudaVerbosity verbosity, const char prefix[], FILE *outfile)
 
void qudaSetCommHandle (void *mycomm)
 
void initCommsGridQuda (int nDim, const int *dims, QudaCommsMap func, void *fdata)
 
void initQudaDevice (int device)
 
void initQudaMemory ()
 
void initQuda (int device)
 
void endQuda (void)
 
void updateR ()
 update the radius for halos. More...
 
QudaGaugeParam newQudaGaugeParam (void)
 
QudaInvertParam newQudaInvertParam (void)
 
QudaMultigridParam newQudaMultigridParam (void)
 
QudaEigParam newQudaEigParam (void)
 
QudaGaugeObservableParam newQudaGaugeObservableParam (void)
 
QudaBLASParam newQudaBLASParam (void)
 
void printQudaGaugeParam (QudaGaugeParam *param)
 
void printQudaInvertParam (QudaInvertParam *param)
 
void printQudaMultigridParam (QudaMultigridParam *param)
 
void printQudaEigParam (QudaEigParam *param)
 
void printQudaGaugeObservableParam (QudaGaugeObservableParam *param)
 
void printQudaBLASParam (QudaBLASParam *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 eigensolveQuda (void **h_evecs, double_complex *h_evals, QudaEigParam *param)
 
void invertQuda (void *h_x, void *h_b, QudaInvertParam *param)
 
void invertMultiSrcQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge, QudaGaugeParam *gauge_param)
 Perform the solve like @invertQuda but for multiple rhs by spliting the comm grid into sub-partitions: each sub-partition invert one or more rhs'. The QudaInvertParam object specifies how the solve should be performed on each sub-partition. Unlike @invertQuda, the interface also takes the host side gauge as input. The gauge pointer and gauge_param are used if for inv_param split_grid[0] * split_grid[1] * split_grid[2] * split_grid[3] is larger than 1, in which case gauge field is not required to be loaded beforehand; otherwise this interface would just work as @invertQuda, which requires gauge field to be loaded beforehand, and the gauge field pointer and gauge_param are not used. More...
 
void invertMultiSrcStaggeredQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, void *milc_fatlinks, void *milc_longlinks, QudaGaugeParam *gauge_param)
 Really the same with @invertMultiSrcQuda but for staggered-style fermions, by accepting pointers to fat links and long links. More...
 
void invertMultiSrcCloverQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, void *h_gauge, QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv)
 Really the same with @invertMultiSrcQuda but for clover-style fermions, by accepting pointers to direct and inverse clover field pointers. More...
 
void invertMultiShiftQuda (void **_hp_x, void *_hp_b, QudaInvertParam *param)
 
void * newMultigridQuda (QudaMultigridParam *param)
 
void destroyMultigridQuda (void *mg_instance)
 Free resources allocated by the multigrid solver. More...
 
void updateMultigridQuda (void *mg_instance, QudaMultigridParam *param)
 Updates the multigrid preconditioner for the new gauge / clover field. More...
 
void dumpMultigridQuda (void *mg_instance, QudaMultigridParam *param)
 Dump the null-space vectors to disk. More...
 
void dslashQuda (void *h_out, void *h_in, QudaInvertParam *inv_param, QudaParity parity)
 
void dslashMultiSrcQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge, QudaGaugeParam *gauge_param)
 Perform the solve like @dslashQuda but for multiple rhs by spliting the comm grid into sub-partitions: each sub-partition does one or more rhs'. The QudaInvertParam object specifies how the solve should be performed on each sub-partition. Unlike @invertQuda, the interface also takes the host side gauge as input - gauge field is not required to be loaded beforehand. More...
 
void dslashMultiSrcStaggeredQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *milc_fatlinks, void *milc_longlinks, QudaGaugeParam *gauge_param)
 Really the same with @dslashMultiSrcQuda but for staggered-style fermions, by accepting pointers to fat links and long links. More...
 
void dslashMultiSrcCloverQuda (void **_hp_x, void **_hp_b, QudaInvertParam *param, QudaParity parity, void *h_gauge, QudaGaugeParam *gauge_param, void *h_clover, void *h_clovinv)
 Really the same with @dslashMultiSrcQuda but for clover-style fermions, by accepting pointers to direct and inverse clover field pointers. More...
 
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 computeKSLinkQuda (void *fatlink, void *longlink, void *ulink, void *inlink, double *path_coeff, QudaGaugeParam *param)
 
void momResidentQuda (void *mom, QudaGaugeParam *param)
 
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)
 
void updateGaugeFieldQuda (void *gauge, void *momentum, double dt, int conj_mom, int exact, QudaGaugeParam *param)
 
void staggeredPhaseQuda (void *gauge_h, QudaGaugeParam *param)
 
void projectSU3Quda (void *gauge_h, double tol, QudaGaugeParam *param)
 
double momActionQuda (void *momentum, QudaGaugeParam *param)
 
void * createGaugeFieldQuda (void *gauge, int geometry, QudaGaugeParam *param)
 
void saveGaugeFieldQuda (void *outGauge, void *inGauge, QudaGaugeParam *param)
 
void destroyGaugeFieldQuda (void *gauge)
 
void createCloverQuda (QudaInvertParam *param)
 
void computeCloverForceQuda (void *mom, double dt, void **x, void **p, double *coeff, double kappa2, double ck, int nvector, double multiplicity, void *gauge, QudaGaugeParam *gauge_param, QudaInvertParam *inv_param)
 
void computeStaggeredForceQuda (void *mom, double dt, double delta, void *gauge, void **x, QudaGaugeParam *gauge_param, QudaInvertParam *invert_param)
 
void computeHISQForceQuda (void *momentum, double dt, const double level2_coeff[6], const double fat7_coeff[6], const void *const w_link, const void *const v_link, const void *const u_link, void **quark, int num, int num_naik, double **coeff, QudaGaugeParam *param)
 
void gaussGaugeQuda (unsigned long long seed, double sigma)
 Generate Gaussian distributed fields and store in the resident gauge field. We create a Gaussian-distributed su(n) field and exponentiate it, e.g., U = exp(sigma * H), where H is the distributed su(n) field and beta is the width of the distribution (beta = 0 results in a free field, and sigma = 1 has maximum disorder). More...
 
void plaqQuda (double plaq[3])
 
void copyExtendedResidentGaugeQuda (void *resident_gauge, QudaFieldLocation loc)
 
void performWuppertalnStep (void *h_out, void *h_in, QudaInvertParam *param, unsigned int n_steps, double alpha)
 
void performAPEnStep (unsigned int n_steps, double alpha, int meas_interval)
 
void performSTOUTnStep (unsigned int n_steps, double rho, int meas_interval)
 
void performOvrImpSTOUTnStep (unsigned int n_steps, double rho, double epsilon, int meas_interval)
 
void performWFlownStep (unsigned int n_steps, double step_size, int meas_interval, QudaWFlowType wflow_type)
 
void gaugeObservablesQuda (QudaGaugeObservableParam *param)
 Calculates a variety of gauge-field observables. If a smeared gauge field is presently loaded (in gaugeSmeared) the observables are computed on this, else the resident gauge field will be used. More...
 
void contractQuda (const void *x, const void *y, void *result, const QudaContractType cType, QudaInvertParam *param, const int *X)
 
int computeGaugeFixingOVRQuda (void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double relax_boost, const double tolerance, const unsigned int reunit_interval, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
 Gauge fixing with overrelaxation with support for single and multi GPU. More...
 
int computeGaugeFixingFFTQuda (void *gauge, const unsigned int gauge_dir, const unsigned int Nsteps, const unsigned int verbose_interval, const double alpha, const unsigned int autotune, const double tolerance, const unsigned int stopWtheta, QudaGaugeParam *param, double *timeinfo)
 Gauge fixing with Steepest descent method with FFTs with support for single GPU only. More...
 
void blasGEMMQuda (void *arrayA, void *arrayB, void *arrayC, QudaBoolean native, QudaBLASParam *param)
 Strided Batched GEMM. More...
 
void flushChronoQuda (int index)
 Flush the chronological history for the given index. More...
 
void openMagma ()
 
void closeMagma ()
 
void * newDeflationQuda (QudaEigParam *param)
 
void destroyDeflationQuda (void *df_instance)
 
void setMPICommHandleQuda (void *mycomm)
 

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.

Macro Definition Documentation

◆ double_complex

#define double_complex   double _Complex

Definition at line 18 of file quda.h.

Typedef Documentation

◆ QudaBLASParam

◆ QudaCommsMap

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 817 of file quda.h.

◆ QudaEigParam

typedef struct QudaEigParam_s QudaEigParam

◆ QudaGaugeObservableParam

◆ QudaGaugeParam

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

◆ QudaInvertParam

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

◆ QudaMultigridParam

Function Documentation

◆ blasGEMMQuda()

void blasGEMMQuda ( void *  arrayA,
void *  arrayB,
void *  arrayC,
QudaBoolean  native,
QudaBLASParam param 
)

Strided Batched GEMM.

Parameters
[in]arrayAThe array containing the A matrix data
[in]arrayBThe array containing the A matrix data
[in]arrayCThe array containing the A matrix data
[in]nativeboolean to use either the native or generic version
[in]paramThe data defining the problem execution.

Definition at line 12 of file blas_interface.cpp.

◆ closeMagma()

void closeMagma ( )

Definition at line 91 of file interface_quda.cpp.

◆ cloverQuda()

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

Definition at line 2289 of file interface_quda.cpp.

◆ computeCloverForceQuda()

void computeCloverForceQuda ( void *  mom,
double  dt,
void **  x,
void **  p,
double *  coeff,
double  kappa2,
double  ck,
int  nvector,
double  multiplicity,
void *  gauge,
QudaGaugeParam gauge_param,
QudaInvertParam inv_param 
)

Compute the clover force contributions in each dimension mu given the array of solution fields, and compute the resulting momentum field.

Parameters
momForce matrix
dtIntegrating step size
xArray of solution vectors
pArray of intermediate vectors
coeffArray of residues for each contribution (multiplied by stepsize)
kappa2-kappa*kappa parameter
ck-clover_coefficient * kappa / 8
nvecNumber of vectors
multiplicityNumber fermions this bilinear reresents
gaugeGauge Field
gauge_paramGauge field meta data
inv_paramDirac and solver meta data

Definition at line 4842 of file interface_quda.cpp.

◆ computeGaugeFixingFFTQuda()

int computeGaugeFixingFFTQuda ( void *  gauge,
const unsigned int  gauge_dir,
const unsigned int  Nsteps,
const unsigned int  verbose_interval,
const double  alpha,
const unsigned int  autotune,
const double  tolerance,
const unsigned int  stopWtheta,
QudaGaugeParam param,
double *  timeinfo 
)

Gauge fixing with Steepest descent method with FFTs with support for single GPU only.

Parameters
[in,out]gauge,gaugefield to be fixed
[in]gauge_dir,3for Coulomb gauge fixing, other for Landau gauge fixing
[in]Nsteps,maximumnumber of steps to perform gauge fixing
[in]verbose_interval,printgauge fixing info when iteration count is a multiple of this
[in]alpha,gaugefixing parameter of the method, most common value is 0.08
[in]autotune,1to autotune the method, i.e., if the Fg inverts its tendency we decrease the alpha value
[in]tolerance,torelancevalue to stop the method, if this value is zero then the method stops when iteration reachs the maximum number of steps defined by Nsteps
[in]stopWtheta,0for MILC criterium and 1 to use the theta value
[in]paramThe parameters of the external fields and the computation settings
[out]timeinfo

Definition at line 5921 of file interface_quda.cpp.

◆ computeGaugeFixingOVRQuda()

int computeGaugeFixingOVRQuda ( void *  gauge,
const unsigned int  gauge_dir,
const unsigned int  Nsteps,
const unsigned int  verbose_interval,
const double  relax_boost,
const double  tolerance,
const unsigned int  reunit_interval,
const unsigned int  stopWtheta,
QudaGaugeParam param,
double *  timeinfo 
)

Gauge fixing with overrelaxation with support for single and multi GPU.

Parameters
[in,out]gauge,gaugefield to be fixed
[in]gauge_dir,3for Coulomb gauge fixing, other for Landau gauge fixing
[in]Nsteps,maximumnumber of steps to perform gauge fixing
[in]verbose_interval,printgauge fixing info when iteration count is a multiple of this
[in]relax_boost,gaugefixing parameter of the overrelaxation method, most common value is 1.5 or 1.7.
[in]tolerance,torelancevalue to stop the method, if this value is zero then the method stops when iteration reachs the maximum number of steps defined by Nsteps
[in]reunit_interval,reunitarizegauge field when iteration count is a multiple of this
[in]stopWtheta,0for MILC criterium and 1 to use the theta value
[in]paramThe parameters of the external fields and the computation settings
[out]timeinfo

if (!param->use_resident_gauge) { // load fields onto the device

Definition at line 5846 of file interface_quda.cpp.

◆ computeGaugeForceQuda()

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 
)

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

Definition at line 4186 of file interface_quda.cpp.

◆ computeHISQForceQuda()

void computeHISQForceQuda ( void *  momentum,
double  dt,
const double  level2_coeff[6],
const double  fat7_coeff[6],
const void *const  w_link,
const void *const  v_link,
const void *const  u_link,
void **  quark,
int  num,
int  num_naik,
double **  coeff,
QudaGaugeParam param 
)

Compute the fermion force for the HISQ quark action and integrate the momentum.

Parameters
momentumThe momentum field we are integrating
dtThe stepsize used to integrate the momentum
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.
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.
quarkThe input fermion field.
numThe number of quark fields
num_naikThe number of naik contributions
coeffThe coefficient multiplying the fermion fields in the outer product
param.The field parameters.

Definition at line 4591 of file interface_quda.cpp.

◆ computeKSLinkQuda()

void computeKSLinkQuda ( void *  fatlink,
void *  longlink,
void *  ulink,
void *  inlink,
double *  path_coeff,
QudaGaugeParam param 
)

Definition at line 4071 of file interface_quda.cpp.

◆ computeStaggeredForceQuda()

void computeStaggeredForceQuda ( void *  mom,
double  dt,
double  delta,
void *  gauge,
void **  x,
QudaGaugeParam gauge_param,
QudaInvertParam invert_param 
)

Compute the naive staggered force. All fields must be in the same precision.

Parameters
momMomentum field
dtIntegrating step size
deltaAdditional scale factor when updating momentum (mom += delta * [force]_TA
gaugeGauge field (at present only supports resident gauge field)
xArray of single-parity solution vectors (at present only supports resident solutions)
gauge_paramGauge field meta data
invert_paramDirac and solver meta data

Definition at line 4436 of file interface_quda.cpp.

◆ contractQuda()

void contractQuda ( const void *  x,
const void *  y,
void *  result,
const QudaContractType  cType,
QudaInvertParam param,
const int *  X 
)

Public function to perform color contractions of the host spinors x and y.

Parameters
[in]xpointer to host data
[in]ypointer to host data
[out]resultpointer to the 16 spin projections per lattice site
[in]cTypeWhich type of contraction (open, degrand-rossi, etc)
[in]parammeta data for construction of ColorSpinorFields.
[in]Xspacetime data for construction of ColorSpinorFields.

Definition at line 5986 of file interface_quda.cpp.

◆ copyExtendedResidentGaugeQuda()

void copyExtendedResidentGaugeQuda ( void *  resident_gauge,
QudaFieldLocation  loc 
)

Performs a deep copy from the internal extendedGaugeResident field.

Parameters
Pointerto externalGaugeResident cudaGaugeField
Locationof gauge field

Definition at line 5600 of file interface_quda.cpp.

◆ createCloverQuda()

void createCloverQuda ( QudaInvertParam param)

Compute the clover field and its inverse from the resident gauge field.

Parameters
paramThe parameters of the clover field to create

Definition at line 4365 of file interface_quda.cpp.

◆ createGaugeFieldQuda()

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

Allocate a gauge (matrix) field on the device and optionally download a host gauge field.

Parameters
gaugeThe host gauge field (optional - if set to 0 then the gauge field zeroed)
geometryThe geometry of the matrix field to create (1 - scalar, 4 - vector, 6 - tensor)
paramThe parameters of the external field and the field to be created
Returns
Pointer to the gauge field (cast as a void*)

Definition at line 4394 of file interface_quda.cpp.

◆ destroyDeflationQuda()

void destroyDeflationQuda ( void *  df_instance)

Free resources allocated by the deflated solver

Definition at line 2830 of file interface_quda.cpp.

◆ destroyGaugeFieldQuda()

void destroyGaugeFieldQuda ( void *  gauge)

Reinterpret gauge as a pointer to cudaGaugeField and call destructor.

Parameters
gaugeGauge field to be freed

Definition at line 4430 of file interface_quda.cpp.

◆ destroyMultigridQuda()

void destroyMultigridQuda ( void *  mg_instance)

Free resources allocated by the multigrid solver.

Parameters
mg_instancePointer to instance of multigrid_solver
paramContains all metadata regarding host and device storage and solver parameters

Definition at line 2624 of file interface_quda.cpp.

◆ dslashMultiSrcCloverQuda()

void dslashMultiSrcCloverQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
QudaParity  parity,
void *  h_gauge,
QudaGaugeParam gauge_param,
void *  h_clover,
void *  h_clovinv 
)

Really the same with @dslashMultiSrcQuda but for clover-style fermions, by accepting pointers to direct and inverse clover field pointers.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
parityParity to apply dslash on
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field
h_cloverBase pointer to the direct clover field
h_clovinvBase pointer to the inverse clover field

Definition at line 3649 of file interface_quda.cpp.

◆ dslashMultiSrcQuda()

void dslashMultiSrcQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
QudaParity  parity,
void *  h_gauge,
QudaGaugeParam gauge_param 
)

Perform the solve like @dslashQuda but for multiple rhs by spliting the comm grid into sub-partitions: each sub-partition does one or more rhs'. The QudaInvertParam object specifies how the solve should be performed on each sub-partition. Unlike @invertQuda, the interface also takes the host side gauge as input - gauge field is not required to be loaded beforehand.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
parityParity to apply dslash on
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field

Definition at line 3634 of file interface_quda.cpp.

◆ dslashMultiSrcStaggeredQuda()

void dslashMultiSrcStaggeredQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
QudaParity  parity,
void *  milc_fatlinks,
void *  milc_longlinks,
QudaGaugeParam gauge_param 
)

Really the same with @dslashMultiSrcQuda but for staggered-style fermions, by accepting pointers to fat links and long links.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
parityParity to apply dslash on
milc_fatlinksBase pointer to host fat gauge field (regardless of dimensionality)
milc_longlinksBase pointer to host long gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field

Definition at line 3641 of file interface_quda.cpp.

◆ dslashQuda()

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 1934 of file interface_quda.cpp.

◆ dumpMultigridQuda()

void dumpMultigridQuda ( void *  mg_instance,
QudaMultigridParam param 
)

Dump the null-space vectors to disk.

Parameters
[in]mg_instancePointer to the instance of multigrid_solver
[in]paramContains all metadata regarding host and device storage and solver parameters (QudaMultigridParam::vec_outfile sets the output filename prefix).

Definition at line 2733 of file interface_quda.cpp.

◆ eigensolveQuda()

void eigensolveQuda ( void **  h_evecs,
double_complex h_evals,
QudaEigParam param 
)

Perform the eigensolve. The problem matrix is defined by the invert param, the mode of solution is specified by the eig param. It is assumed that the gauge field has already been loaded via loadGaugeQuda().

Parameters
h_evecsArray of pointers to application eigenvectors
h_evalsHost side eigenvalues
paramContains all metadata regarding the type of solve.

◆ endQuda()

void endQuda ( void  )

Finalize the library.

Definition at line 1474 of file interface_quda.cpp.

◆ flushChronoQuda()

void flushChronoQuda ( int  index)

Flush the chronological history for the given index.

Parameters
[in]indexIndex for which we are flushing

Definition at line 1461 of file interface_quda.cpp.

◆ freeCloverQuda()

void freeCloverQuda ( void  )

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

Definition at line 1453 of file interface_quda.cpp.

◆ freeGaugeQuda()

void freeGaugeQuda ( void  )

Free QUDA's internal copy of the gauge field.

Definition at line 1190 of file interface_quda.cpp.

◆ gaugeObservablesQuda()

void gaugeObservablesQuda ( QudaGaugeObservableParam param)

Calculates a variety of gauge-field observables. If a smeared gauge field is presently loaded (in gaugeSmeared) the observables are computed on this, else the resident gauge field will be used.

Parameters
[in,out]paramParameter struct that defines which observables we are making and the resulting observables.

Definition at line 6042 of file interface_quda.cpp.

◆ gaussGaugeQuda()

void gaussGaugeQuda ( unsigned long long  seed,
double  sigma 
)

Generate Gaussian distributed fields and store in the resident gauge field. We create a Gaussian-distributed su(n) field and exponentiate it, e.g., U = exp(sigma * H), where H is the distributed su(n) field and beta is the width of the distribution (beta = 0 results in a free field, and sigma = 1 has maximum disorder).

Parameters
seedThe seed used for the RNG
sigmaWidth of Gaussian distrubution

Definition at line 5545 of file interface_quda.cpp.

◆ initCommsGridQuda()

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 371 of file interface_quda.cpp.

◆ initQuda()

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 536 of file interface_quda.cpp.

◆ initQudaDevice()

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 453 of file interface_quda.cpp.

◆ initQudaMemory()

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 503 of file interface_quda.cpp.

◆ invertMultiShiftQuda()

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_bSource 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.

For Wilson-type fermions, 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.

For Staggered-type fermions, the solution_type must be MATPC, and the solve type must be DIRECT_PC. This difference in convention is because preconditioned staggered operator is normal, unlike with Wilson-type fermions.

Definition at line 3668 of file interface_quda.cpp.

◆ invertMultiSrcCloverQuda()

void invertMultiSrcCloverQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
void *  h_gauge,
QudaGaugeParam gauge_param,
void *  h_clover,
void *  h_clovinv 
)

Really the same with @invertMultiSrcQuda but for clover-style fermions, by accepting pointers to direct and inverse clover field pointers.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field
h_cloverBase pointer to the direct clover field
h_clovinvBase pointer to the inverse clover field

Definition at line 3627 of file interface_quda.cpp.

◆ invertMultiSrcQuda()

void invertMultiSrcQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
void *  h_gauge,
QudaGaugeParam gauge_param 
)

Perform the solve like @invertQuda but for multiple rhs by spliting the comm grid into sub-partitions: each sub-partition invert one or more rhs'. The QudaInvertParam object specifies how the solve should be performed on each sub-partition. Unlike @invertQuda, the interface also takes the host side gauge as input. The gauge pointer and gauge_param are used if for inv_param split_grid[0] * split_grid[1] * split_grid[2] * split_grid[3] is larger than 1, in which case gauge field is not required to be loaded beforehand; otherwise this interface would just work as @invertQuda, which requires gauge field to be loaded beforehand, and the gauge field pointer and gauge_param are not used.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
h_gaugeBase pointer to host gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field

Definition at line 3614 of file interface_quda.cpp.

◆ invertMultiSrcStaggeredQuda()

void invertMultiSrcStaggeredQuda ( void **  _hp_x,
void **  _hp_b,
QudaInvertParam param,
void *  milc_fatlinks,
void *  milc_longlinks,
QudaGaugeParam gauge_param 
)

Really the same with @invertMultiSrcQuda but for staggered-style fermions, by accepting pointers to fat links and long links.

Parameters
_hp_xArray of solution spinor fields
_hp_bArray of source spinor fields
paramContains all metadata regarding host and device storage and solver parameters
milc_fatlinksBase pointer to host fat gauge field (regardless of dimensionality)
milc_longlinksBase pointer to host long gauge field (regardless of dimensionality)
gauge_paramContains all metadata regarding host and device storage for gauge field

Definition at line 3620 of file interface_quda.cpp.

◆ invertQuda()

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 2837 of file interface_quda.cpp.

◆ lanczosQuda()

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

◆ loadCloverQuda()

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 849 of file interface_quda.cpp.

◆ loadGaugeQuda()

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 553 of file interface_quda.cpp.

◆ MatDagMatQuda()

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

Apply M^{\dag}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 2099 of file interface_quda.cpp.

◆ MatQuda()

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 2029 of file interface_quda.cpp.

◆ momActionQuda()

double momActionQuda ( void *  momentum,
QudaGaugeParam param 
)

Evaluate the momentum contribution to the Hybrid Monte Carlo action.

Parameters
momentumThe momentum field
paramThe parameters of the external fields and the computation settings
Returns
momentum action

Definition at line 5242 of file interface_quda.cpp.

◆ momResidentQuda()

void momResidentQuda ( void *  mom,
QudaGaugeParam param 
)

Either downloads and sets the resident momentum field, or uploads and returns the resident momentum field

Parameters
[in,out]momThe external momentum field
[in]paramThe parameters of the external field

Definition at line 4310 of file interface_quda.cpp.

◆ newDeflationQuda()

void* newDeflationQuda ( QudaEigParam param)

Create deflation solver resources.

Definition at line 2816 of file interface_quda.cpp.

◆ newMultigridQuda()

void* newMultigridQuda ( QudaMultigridParam param)

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

Parameters
paramContains all metadata regarding host and device storage and solver parameters

Definition at line 2607 of file interface_quda.cpp.

◆ newQudaBLASParam()

QudaBLASParam newQudaBLASParam ( void  )

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

QudaBLASParam blas_param = newQudaBLASParam();

◆ newQudaEigParam()

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();

◆ newQudaGaugeObservableParam()

QudaGaugeObservableParam newQudaGaugeObservableParam ( void  )

A new QudaGaugeObservableParam 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 obs_param = newQudaGaugeObservableParam();

◆ newQudaGaugeParam()

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();

◆ newQudaInvertParam()

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();

◆ newQudaMultigridParam()

QudaMultigridParam newQudaMultigridParam ( void  )

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

QudaMultigridParam mg_param = newQudaMultigridParam();

◆ openMagma()

void openMagma ( )

Open/Close MAGMA library

Definition at line 80 of file interface_quda.cpp.

◆ pack_ghost()

void pack_ghost ( void **  cpuLink,
void **  cpuGhost,
int  nFace,
QudaPrecision  precision 
)

◆ performAPEnStep()

void performAPEnStep ( unsigned int  n_steps,
double  alpha,
int  meas_interval 
)

Performs APE smearing on gaugePrecise and stores it in gaugeSmeared

Parameters
n_stepsNumber of steps to apply.
alphaAlpha coefficient for APE smearing.
meas_intervalMeasure the Q charge every Nth step

Definition at line 5691 of file interface_quda.cpp.

◆ performOvrImpSTOUTnStep()

void performOvrImpSTOUTnStep ( unsigned int  n_steps,
double  rho,
double  epsilon,
int  meas_interval 
)

Performs Over Imroved STOUT smearing on gaugePrecise and stores it in gaugeSmeared

Parameters
n_stepsNumber of steps to apply.
rhoRho coefficient for STOUT smearing.
epsilonEpsilon coefficient for Over Improved STOUT smearing.
meas_intervalMeasure the Q charge every Nth step

Definition at line 5759 of file interface_quda.cpp.

◆ performSTOUTnStep()

void performSTOUTnStep ( unsigned int  n_steps,
double  rho,
int  meas_interval 
)

Performs STOUT smearing on gaugePrecise and stores it in gaugeSmeared

Parameters
n_stepsNumber of steps to apply.
rhoRho coefficient for STOUT smearing.
meas_intervalMeasure the Q charge every Nth step

Definition at line 5725 of file interface_quda.cpp.

◆ performWFlownStep()

void performWFlownStep ( unsigned int  n_steps,
double  step_size,
int  meas_interval,
QudaWFlowType  wflow_type 
)

Performs Wilson Flow on gaugePrecise and stores it in gaugeSmeared

Parameters
n_stepsNumber of steps to apply.
step_sizeSize of Wilson Flow step
meas_intervalMeasure the Q charge and field energy every Nth step
wflow_type1x1 Wilson or 2x1 Symanzik flow type

Definition at line 5793 of file interface_quda.cpp.

◆ performWuppertalnStep()

void performWuppertalnStep ( void *  h_out,
void *  h_in,
QudaInvertParam param,
unsigned int  n_steps,
double  alpha 
)

Performs Wuppertal smearing on a given spinor using the gauge field gaugeSmeared, if it exist, or gaugePrecise if no smeared field is present.

Parameters
h_outResult spinor field
h_inInput spinor field
paramContains all metadata regarding host and device storage and operator which will be applied to the spinor
n_stepsNumber of steps to apply.
alphaAlpha coefficient for Wuppertal smearing.

Definition at line 5616 of file interface_quda.cpp.

◆ plaqQuda()

void plaqQuda ( double  plaq[3])

Computes the total, spatial and temporal plaquette averages of the loaded gauge configuration.

Parameters
Arrayfor storing the averages (total, spatial, temporal)

Definition at line 5578 of file interface_quda.cpp.

◆ printQudaBLASParam()

void printQudaBLASParam ( QudaBLASParam param)

Print the members of QudaBLASParam.

Parameters
paramThe QudaBLASParam whose elements we are to print.

Definition at line 981 of file check_params.h.

◆ printQudaEigParam()

void printQudaEigParam ( QudaEigParam param)

Print the members of QudaEigParam.

Parameters
paramThe QudaEigParam whose elements we are to print.

Definition at line 158 of file check_params.h.

◆ printQudaGaugeObservableParam()

void printQudaGaugeObservableParam ( QudaGaugeObservableParam param)

Print the members of QudaGaugeObservableParam.

Parameters
paramThe QudaGaugeObservableParam whose elements we are to print.

Definition at line 943 of file check_params.h.

◆ printQudaGaugeParam()

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.

◆ printQudaInvertParam()

void printQudaInvertParam ( QudaInvertParam param)

Print the members of QudaInvertParam.

Parameters
paramThe QudaInvertParam whose elements we are to print.

Whether to use a pipelined solver

< Number of offsets in the multi-shift solver

< Number of offsets in the multi-shift solver

< width of domain overlaps

< Number of sources per sub-partitions

Definition at line 342 of file check_params.h.

◆ printQudaMultigridParam()

void printQudaMultigridParam ( QudaMultigridParam param)

Print the members of QudaMultigridParam.

Parameters
paramThe QudaMultigridParam whose elements we are to print.

Definition at line 689 of file check_params.h.

◆ projectSU3Quda()

void projectSU3Quda ( void *  gauge_h,
double  tol,
QudaGaugeParam param 
)

Project the input field on the SU(3) group. If the target tolerance is not met, this routine will give a runtime error.

Parameters
gauge_hThe gauge field to be updated
tolThe tolerance to which we iterate
paramThe parameters of the gauge field

Definition at line 5126 of file interface_quda.cpp.

◆ qudaSetCommHandle()

void qudaSetCommHandle ( void *  mycomm)
Parameters
mycommUser provided MPI communicator in place of MPI_COMM_WORLD

◆ saveGaugeFieldQuda()

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

Copy the QUDA gauge (matrix) field on the device to the CPU

Parameters
outGaugePointer to the host gauge field
inGaugePointer to the device gauge field (QUDA device field)
paramThe parameters of the host and device fields

Definition at line 4417 of file interface_quda.cpp.

◆ saveGaugeQuda()

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 808 of file interface_quda.cpp.

◆ set_dim()

void set_dim ( int *  )

◆ setMPICommHandleQuda()

void setMPICommHandleQuda ( void *  mycomm)

Definition at line 361 of file interface_quda.cpp.

◆ setVerbosityQuda()

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 316 of file interface_quda.cpp.

◆ staggeredPhaseQuda()

void staggeredPhaseQuda ( void *  gauge_h,
QudaGaugeParam param 
)

Apply the staggered phase factors to the gauge field. If the imaginary chemical potential is non-zero then the phase factor exp(imu/T) will be applied to the links in the temporal direction.

Parameters
gauge_hThe gauge field
paramThe parameters of the gauge field

Definition at line 5187 of file interface_quda.cpp.

◆ updateGaugeFieldQuda()

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 5026 of file interface_quda.cpp.

◆ updateMultigridQuda()

void updateMultigridQuda ( void *  mg_instance,
QudaMultigridParam param 
)

Updates the multigrid preconditioner for the new gauge / clover field.

Parameters
mg_instancePointer to instance of multigrid_solver
paramContains all metadata regarding host and device storage and solver parameters, of note contains a flag specifying whether to do a full update or a thin update.

Definition at line 2628 of file interface_quda.cpp.

◆ updateR()

void updateR ( )

update the radius for halos.

This should only be needed for automated testing when different partitioning is applied within a single run.

Definition at line 531 of file interface_quda.cpp.