QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
quda::GaugeField Class Referenceabstract

#include <gauge_field.h>

Inheritance diagram for quda::GaugeField:
Inheritance graph
[legend]
Collaboration diagram for quda::GaugeField:
Collaboration graph
[legend]

Public Member Functions

 GaugeField (const GaugeFieldParam &param)
 
virtual ~GaugeField ()
 
virtual void exchangeGhost (QudaLinkDirection=QUDA_LINK_BACKWARDS)=0
 
virtual void injectGhost (QudaLinkDirection=QUDA_LINK_BACKWARDS)=0
 
size_t Length () const
 
int Ncolor () const
 
QudaReconstructType Reconstruct () const
 
QudaGaugeFieldOrder Order () const
 
double Anisotropy () const
 
double Tadpole () const
 
QudaTboundary TBoundary () const
 
QudaLinkType LinkType () const
 
QudaGaugeFixed GaugeFixed () const
 
QudaGaugeFieldOrder FieldOrder () const
 
QudaFieldGeometry Geometry () const
 
QudaStaggeredPhase StaggeredPhase () const
 
bool StaggeredPhaseApplied () const
 
void applyStaggeredPhase (QudaStaggeredPhase phase=QUDA_STAGGERED_PHASE_INVALID)
 
void removeStaggeredPhase ()
 
double iMu () const
 
const double & LinkMax () const
 
int Nface () const
 
virtual void exchangeExtendedGhost (const int *R, bool no_comms_fill=false)=0
 This does routine will populate the border / halo region of a gauge field that has been created using copyExtendedGauge. More...
 
virtual void exchangeExtendedGhost (const int *R, TimeProfile &profile, bool no_comms_fill=false)=0
 This does routine will populate the border / halo region of a gauge field that has been created using copyExtendedGauge. Overloaded variant that will start and stop a comms profile. More...
 
void checkField (const LatticeField &) const
 
bool isNative () const
 
size_t Bytes () const
 
size_t PhaseBytes () const
 
size_t PhaseOffset () const
 
virtual void * Gauge_p ()
 
virtual void * Even_p ()
 
virtual void * Odd_p ()
 
virtual const void * Gauge_p () const
 
virtual const void * Even_p () const
 
virtual const void * Odd_p () const
 
const void ** Ghost () const
 
void ** Ghost ()
 
size_t SiteOffset () const
 
size_t SiteSize () const
 
virtual void zero ()=0
 
virtual void copy (const GaugeField &src)=0
 
double norm1 (int dim=-1) const
 Compute the L1 norm of the field. More...
 
double norm2 (int dim=-1) const
 Compute the L2 norm squared of the field. More...
 
double abs_max (int dim=-1) const
 Compute the absolute maximum of the field (Linfinity norm) More...
 
double abs_min (int dim=-1) const
 Compute the absolute minimum of the field. More...
 
uint64_t checksum (bool mini=false) const
 
- Public Member Functions inherited from quda::LatticeField
 LatticeField (const LatticeFieldParam &param)
 
 LatticeField (const LatticeField &field)
 
virtual ~LatticeField ()
 
void allocateGhostBuffer (size_t ghost_bytes) const
 Allocate the static ghost buffers. More...
 
void createComms (bool no_comms_fill=false, bool bidir=true)
 
void destroyComms ()
 
void createIPCComms ()
 
bool ipcCopyComplete (int dir, int dim)
 
bool ipcRemoteCopyComplete (int dir, int dim)
 
const cudaEvent_t & getIPCCopyEvent (int dir, int dim) const
 
const cudaEvent_t & getIPCRemoteCopyEvent (int dir, int dim) const
 
int Ndim () const
 
const int * X () const
 
int Volume () const
 
int VolumeCB () const
 
const int * SurfaceCB () const
 
int SurfaceCB (const int i) const
 
int Stride () const
 
int Pad () const
 
const int * R () const
 
QudaGhostExchange GhostExchange () const
 
QudaPrecision Precision () const
 
QudaPrecision GhostPrecision () const
 
double Scale () const
 
void Scale (double scale_)
 Set the scale factor for a fixed-point field. More...
 
virtual QudaSiteSubset SiteSubset () const
 
virtual QudaMemoryType MemType () const
 
int Nvec () const
 
QudaFieldLocation Location () const
 
size_t GBytes () const
 
void checkField (const LatticeField &a) const
 
virtual void read (char *filename)
 
virtual void write (char *filename)
 
virtual void gather (int nFace, int dagger, int dir, cudaStream_t *stream_p=NULL)
 
virtual void commsStart (int nFace, int dir, int dagger=0, cudaStream_t *stream_p=NULL, bool gdr_send=false, bool gdr_recv=true)
 
virtual int commsQuery (int nFace, int dir, int dagger=0, cudaStream_t *stream_p=NULL, bool gdr_send=false, bool gdr_recv=true)
 
virtual void commsWait (int nFace, int dir, int dagger=0, cudaStream_t *stream_p=NULL, bool gdr_send=false, bool gdr_recv=true)
 
virtual void scatter (int nFace, int dagger, int dir)
 
const char * VolString () const
 
const char * AuxString () const
 
virtual void backup () const
 Backs up the LatticeField. More...
 
virtual void restore () const
 Restores the LatticeField. More...
 
- Public Member Functions inherited from quda::Object
 Object ()
 
virtual ~Object ()
 

Static Public Member Functions

static GaugeFieldCreate (const GaugeFieldParam &param)
 Create the gauge field, with meta data specified in the parameter struct. More...
 
- Static Public Member Functions inherited from quda::LatticeField
static void freeGhostBuffer (void)
 Free statically allocated ghost buffers. More...
 
static void destroyIPCComms ()
 
- Static Public Member Functions inherited from quda::Object
static void * operator new (std::size_t size)
 
static void operator delete (void *p)
 
static void * operator new[] (std::size_t size)
 
static void operator delete[] (void *p)
 

Protected Member Functions

void exchange (void **recv, void **send, QudaDirection dir) const
 Exchange the buffers across all dimensions in a given direction. More...
 
void createGhostZone (const int *R, bool no_comms_fill, bool bidir=true) const
 
void setTuningString ()
 Set the vol_string and aux_string for use in tuning. More...
 
- Protected Member Functions inherited from quda::LatticeField
void precisionCheck ()
 

Protected Attributes

size_t bytes
 
size_t phase_offset
 
size_t phase_bytes
 
size_t length
 
size_t real_length
 
int nColor
 
int nFace
 
QudaFieldGeometry geometry
 
QudaReconstructType reconstruct
 
int nInternal
 
QudaGaugeFieldOrder order
 
QudaGaugeFixed fixed
 
QudaLinkType link_type
 
QudaTboundary t_boundary
 
double anisotropy
 
double tadpole
 
double fat_link_max
 
QudaFieldCreate create
 
void * ghost [2 *QUDA_MAX_DIM]
 
int ghostFace [QUDA_MAX_DIM]
 
QudaStaggeredPhase staggeredPhaseType
 
bool staggeredPhaseApplied
 
double i_mu
 
size_t site_offset
 
size_t site_size
 
- Protected Attributes inherited from quda::LatticeField
int volume
 
int volumeCB
 
int stride
 
int pad
 
size_t total_bytes
 
int nDim
 
int x [QUDA_MAX_DIM]
 
int surface [QUDA_MAX_DIM]
 
int surfaceCB [QUDA_MAX_DIM]
 
int r [QUDA_MAX_DIM]
 
QudaPrecision precision
 
QudaPrecision ghost_precision
 
bool ghost_precision_reset
 
double scale
 
QudaSiteSubset siteSubset
 
QudaGhostExchange ghostExchange
 
int nDimComms
 
size_t ghost_bytes
 
size_t ghost_bytes_old
 
size_t ghost_face_bytes [QUDA_MAX_DIM]
 
int ghostOffset [QUDA_MAX_DIM][2]
 
int ghostNormOffset [QUDA_MAX_DIM][2]
 
void * my_face_h [2]
 
void * my_face_hd [2]
 
void * my_face_d [2]
 
void * my_face_dim_dir_h [2][QUDA_MAX_DIM][2]
 
void * my_face_dim_dir_hd [2][QUDA_MAX_DIM][2]
 
void * my_face_dim_dir_d [2][QUDA_MAX_DIM][2]
 
void * from_face_h [2]
 
void * from_face_hd [2]
 
void * from_face_d [2]
 
void * from_face_dim_dir_h [2][QUDA_MAX_DIM][2]
 
void * from_face_dim_dir_hd [2][QUDA_MAX_DIM][2]
 
void * from_face_dim_dir_d [2][QUDA_MAX_DIM][2]
 
MsgHandlemh_recv_fwd [2][QUDA_MAX_DIM]
 
MsgHandlemh_recv_back [2][QUDA_MAX_DIM]
 
MsgHandlemh_send_fwd [2][QUDA_MAX_DIM]
 
MsgHandlemh_send_back [2][QUDA_MAX_DIM]
 
MsgHandlemh_recv_rdma_fwd [2][QUDA_MAX_DIM]
 
MsgHandlemh_recv_rdma_back [2][QUDA_MAX_DIM]
 
MsgHandlemh_send_rdma_fwd [2][QUDA_MAX_DIM]
 
MsgHandlemh_send_rdma_back [2][QUDA_MAX_DIM]
 
bool initComms
 
char vol_string [TuneKey::volume_n]
 
char aux_string [TuneKey::aux_n]
 
QudaMemoryType mem_type
 
char * backup_h
 
char * backup_norm_h
 
bool backed_up
 

Additional Inherited Members

- Static Public Attributes inherited from quda::LatticeField
static int bufferIndex = 0
 
static bool ghost_field_reset = false
 
- Static Protected Attributes inherited from quda::LatticeField
static void * ghost_send_buffer_d [2] = {nullptr, nullptr}
 
static void * ghost_recv_buffer_d [2] = {nullptr, nullptr}
 
static void * ghost_pinned_send_buffer_h [2] = {nullptr, nullptr}
 
static void * ghost_pinned_recv_buffer_h [2] = {nullptr, nullptr}
 
static void * ghost_pinned_send_buffer_hd [2] = {nullptr, nullptr}
 
static void * ghost_pinned_recv_buffer_hd [2] = {nullptr, nullptr}
 
static void * ghost_remote_send_buffer_d [2][QUDA_MAX_DIM][2]
 
static size_t ghostFaceBytes = 0
 
static bool initGhostFaceBuffer = false
 
static MsgHandlemh_send_p2p_fwd [2][QUDA_MAX_DIM] { }
 
static MsgHandlemh_send_p2p_back [2][QUDA_MAX_DIM] { }
 
static MsgHandlemh_recv_p2p_fwd [2][QUDA_MAX_DIM] { }
 
static MsgHandlemh_recv_p2p_back [2][QUDA_MAX_DIM] { }
 
static int buffer_send_p2p_fwd [2][QUDA_MAX_DIM] { }
 
static int buffer_recv_p2p_fwd [2][QUDA_MAX_DIM] { }
 
static int buffer_send_p2p_back [2][QUDA_MAX_DIM] { }
 
static int buffer_recv_p2p_back [2][QUDA_MAX_DIM] { }
 
static cudaEvent_t ipcCopyEvent [2][2][QUDA_MAX_DIM]
 
static cudaEvent_t ipcRemoteCopyEvent [2][2][QUDA_MAX_DIM]
 
static bool initIPCComms = false
 

Detailed Description

Definition at line 164 of file gauge_field.h.

Constructor & Destructor Documentation

◆ GaugeField()

quda::GaugeField::GaugeField ( const GaugeFieldParam param)

◆ ~GaugeField()

quda::GaugeField::~GaugeField ( )
virtual

Definition at line 100 of file gauge_field.cpp.

Member Function Documentation

◆ abs_max()

double quda::GaugeField::abs_max ( int  dim = -1) const

Compute the absolute maximum of the field (Linfinity norm)

Parameters
[in]dimWhich dimension we are taking the norm of (dim=-1 mean all dimensions)
Returns
Absolute maximum value

Definition at line 84 of file max_gauge.cu.

References quda::ABS_MAX, errorQuda, QUDA_DOUBLE_PRECISION, QUDA_HALF_PRECISION, QUDA_QUARTER_PRECISION, QUDA_RECONSTRUCT_NO, and QUDA_SINGLE_PRECISION.

Referenced by quda::copyGaugeMG(), and quda::cpuGaugeField::cpuGaugeField().

Here is the caller graph for this function:

◆ abs_min()

double quda::GaugeField::abs_min ( int  dim = -1) const

Compute the absolute minimum of the field.

Parameters
[in]dimWhich dimension we are taking the norm of (dim=-1 mean all dimensions)
Returns
Absolute minimum value

Definition at line 97 of file max_gauge.cu.

References quda::ABS_MIN, errorQuda, QUDA_DOUBLE_PRECISION, QUDA_HALF_PRECISION, QUDA_QUARTER_PRECISION, QUDA_RECONSTRUCT_NO, and QUDA_SINGLE_PRECISION.

◆ Anisotropy()

double quda::GaugeField::Anisotropy ( ) const
inline

Definition at line 252 of file gauge_field.h.

References quda::GaugeFieldParam::anisotropy.

Referenced by cloverQuda(), quda::CoarseOp(), createExtendedGauge(), loadCloverQuda(), and quda::setDiracParam().

Here is the caller graph for this function:

◆ applyStaggeredPhase()

void quda::GaugeField::applyStaggeredPhase ( QudaStaggeredPhase  phase = QUDA_STAGGERED_PHASE_INVALID)

Apply the staggered phase factors to the gauge field.

Parameters
[in]phaseThe phase we will apply to the field. If this is QUDA_STAGGERED_PHASE_INVALID, the default value, then apply the phase set internal to the field.

Definition at line 139 of file gauge_field.cpp.

References quda::applyGaugePhase(), errorQuda, quda::cudaGaugeField::exchangeGhost(), quda::cpuGaugeField::exchangeGhost(), quda::LatticeField::ghostExchange, QUDA_GHOST_EXCHANGE_PAD, QUDA_STAGGERED_PHASE_INVALID, staggeredPhaseApplied, and staggeredPhaseType.

Referenced by apply_staggered_phase_quda_(), and staggeredPhaseQuda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Bytes()

size_t quda::GaugeField::Bytes ( ) const
inline

◆ checkField()

void quda::GaugeField::checkField ( const LatticeField l) const

Definition at line 265 of file gauge_field.cpp.

References anisotropy, quda::LatticeField::checkField(), errorQuda, fixed, link_type, nColor, nFace, t_boundary, and tadpole.

Referenced by quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), and quda::cudaGaugeField::saveCPUField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ checksum()

uint64_t quda::GaugeField::checksum ( bool  mini = false) const

Compute checksum of this gauge field: this uses a XOR-based checksum method

Parameters
[in]miniWhether to compute a mini checksum or global checksum. A mini checksum only computes the checksum over a subset of the lattice sites and is to be used for online comparisons, e.g., checking a field has changed with a global update algorithm.
Returns
checksum value

Definition at line 355 of file gauge_field.cpp.

References quda::Checksum().

Referenced by loadGaugeQuda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ copy()

virtual void quda::GaugeField::copy ( const GaugeField src)
pure virtual

Generic gauge field copy

Parameters
[in]srcSource from which we are copying

Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.

Referenced by quda::CoarseOp().

Here is the caller graph for this function:

◆ Create()

GaugeField * quda::GaugeField::Create ( const GaugeFieldParam param)
static

Create the gauge field, with meta data specified in the parameter struct.

Parameters
paramParameter struct specifying the gauge field
Returns
Pointer to allcoated gauge field

Definition at line 359 of file gauge_field.cpp.

References errorQuda, quda::GaugeFieldParam::location, QUDA_CPU_FIELD_LOCATION, and QUDA_CUDA_FIELD_LOCATION.

Referenced by quda::CoarseCoarseOp(), quda::CoarseOp(), and computeGaugeForceQuda().

Here is the caller graph for this function:

◆ createGhostZone()

void quda::GaugeField::createGhostZone ( const int *  R,
bool  no_comms_fill,
bool  bidir = true 
) const
protected

Compute the required extended ghost zone sizes and offsets

Parameters
[in]RRadius of the ghost zone
[in]no_comms_fillIf true we create a full halo regardless of partitioning
[in]bidirIs this a bi-directional exchange - if not then we alias the fowards and backwards offsetss

Definition at line 112 of file gauge_field.cpp.

References ALIGNMENT_ADJUST, comm_dim_partitioned(), geometry, quda::LatticeField::ghost_bytes, quda::LatticeField::ghost_bytes_old, quda::LatticeField::ghost_face_bytes, quda::LatticeField::ghost_precision, ghostFace, quda::LatticeField::ghostOffset, isNative(), quda::LatticeField::nDim, nInternal, QUDA_COARSE_GEOMETRY, QUDA_SCALAR_GEOMETRY, QUDA_VECTOR_GEOMETRY, and quda::LatticeField::surface.

Referenced by quda::cudaGaugeField::allocateGhostBuffer().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Even_p() [1/2]

virtual void* quda::GaugeField::Even_p ( )
inlinevirtual

Reimplemented in quda::cudaGaugeField.

Definition at line 316 of file gauge_field.h.

References errorQuda.

◆ Even_p() [2/2]

virtual const void* quda::GaugeField::Even_p ( ) const
inlinevirtual

Reimplemented in quda::cudaGaugeField.

Definition at line 320 of file gauge_field.h.

References errorQuda.

◆ exchange()

void quda::GaugeField::exchange ( void **  recv,
void **  send,
QudaDirection  dir 
) const
protected

Exchange the buffers across all dimensions in a given direction.

Parameters
[out]recvReceive buffer
[in]sendSend buffer
[in]dirDirection in which we are sending (forwards OR backwards only)

Definition at line 185 of file gauge_field.cpp.

References bytes, comm_declare_receive_relative, comm_declare_send_relative, comm_dim_partitioned(), comm_free(), comm_start(), comm_wait(), errorQuda, quda::LatticeField::Location(), quda::LatticeField::nDimComms, nFace, nInternal, pool_pinned_free, pool_pinned_malloc, quda::LatticeField::precision, QUDA_BACKWARDS, QUDA_CPU_FIELD_LOCATION, QUDA_CUDA_FIELD_LOCATION, QUDA_FORWARDS, qudaMemcpy, and quda::LatticeField::surfaceCB.

Referenced by quda::cpuGaugeField::exchangeGhost(), and quda::cpuGaugeField::injectGhost().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ exchangeExtendedGhost() [1/2]

virtual void quda::GaugeField::exchangeExtendedGhost ( const int *  R,
bool  no_comms_fill = false 
)
pure virtual

This does routine will populate the border / halo region of a gauge field that has been created using copyExtendedGauge.

Parameters
RThe thickness of the extended region in each dimension
no_comms_fillDo local exchange to fill out the extended region in non-partitioned dimensions

Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.

Referenced by quda::gaugeGauss().

Here is the caller graph for this function:

◆ exchangeExtendedGhost() [2/2]

virtual void quda::GaugeField::exchangeExtendedGhost ( const int *  R,
TimeProfile profile,
bool  no_comms_fill = false 
)
pure virtual

This does routine will populate the border / halo region of a gauge field that has been created using copyExtendedGauge. Overloaded variant that will start and stop a comms profile.

Parameters
RThe thickness of the extended region in each dimension
profileTimeProfile intance which will record the time taken
no_comms_fillDo local exchange to fill out the extended region in non-partitioned dimensions

Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.

◆ exchangeGhost()

virtual void quda::GaugeField::exchangeGhost ( QudaLinkDirection  = QUDA_LINK_BACKWARDS)
pure virtual

Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.

Referenced by quda::applyGaugePhase(), and quda::gaugeGauss().

Here is the caller graph for this function:

◆ FieldOrder()

QudaGaugeFieldOrder quda::GaugeField::FieldOrder ( ) const
inline

◆ Gauge_p() [1/2]

virtual void* quda::GaugeField::Gauge_p ( )
inlinevirtual

◆ Gauge_p() [2/2]

virtual const void* quda::GaugeField::Gauge_p ( ) const
inlinevirtual

Reimplemented in quda::cpuGaugeField, and quda::cudaGaugeField.

Definition at line 319 of file gauge_field.h.

References errorQuda.

◆ GaugeFixed()

QudaGaugeFixed quda::GaugeField::GaugeFixed ( ) const
inline

Definition at line 256 of file gauge_field.h.

References quda::GaugeFieldParam::fixed.

Referenced by quda::CoarseOp().

Here is the caller graph for this function:

◆ Geometry()

QudaFieldGeometry quda::GaugeField::Geometry ( ) const
inline

◆ Ghost() [1/2]

const void** quda::GaugeField::Ghost ( ) const
inline

◆ Ghost() [2/2]

void** quda::GaugeField::Ghost ( )
inline

Definition at line 328 of file gauge_field.h.

References errorQuda.

◆ iMu()

double quda::GaugeField::iMu ( ) const
inline

Return the imaginary chemical potential applied to this field

Definition at line 278 of file gauge_field.h.

References quda::GaugeFieldParam::i_mu.

◆ injectGhost()

virtual void quda::GaugeField::injectGhost ( QudaLinkDirection  = QUDA_LINK_BACKWARDS)
pure virtual

◆ isNative()

bool quda::GaugeField::isNative ( ) const

◆ Length()

size_t quda::GaugeField::Length ( ) const
inline

Definition at line 248 of file gauge_field.h.

References length.

◆ LinkMax()

const double& quda::GaugeField::LinkMax ( ) const
inline

Definition at line 280 of file gauge_field.h.

References fat_link_max.

Referenced by quda::cudaGaugeField::copy(), and quda::cpuGaugeField::copy().

Here is the caller graph for this function:

◆ LinkType()

QudaLinkType quda::GaugeField::LinkType ( ) const
inline

Definition at line 255 of file gauge_field.h.

References quda::GaugeFieldParam::link_type.

Referenced by quda::CoarseOp(), quda::colorSpinorParam(), quda::copyGauge(), quda::copyGaugeEx(), quda::extractGhostMG(), and quda::gaugeGauss().

Here is the caller graph for this function:

◆ Ncolor()

int quda::GaugeField::Ncolor ( ) const
inline

◆ Nface()

int quda::GaugeField::Nface ( ) const
inline

◆ norm1()

double quda::GaugeField::norm1 ( int  dim = -1) const

Compute the L1 norm of the field.

Parameters
[in]dimWhich dimension we are taking the norm of (dim=-1 mean all dimensions)
Returns
L1 norm

Definition at line 58 of file max_gauge.cu.

References errorQuda, quda::NORM1, QUDA_DOUBLE_PRECISION, QUDA_HALF_PRECISION, QUDA_QUARTER_PRECISION, QUDA_RECONSTRUCT_NO, and QUDA_SINGLE_PRECISION.

◆ norm2()

double quda::GaugeField::norm2 ( int  dim = -1) const

Compute the L2 norm squared of the field.

Parameters
[in]dimWhich dimension we are taking the norm of (dim=-1 mean all dimensions)
Returns
L2 norm squared

Definition at line 71 of file max_gauge.cu.

References errorQuda, quda::NORM2, QUDA_DOUBLE_PRECISION, QUDA_HALF_PRECISION, QUDA_QUARTER_PRECISION, QUDA_RECONSTRUCT_NO, and QUDA_SINGLE_PRECISION.

◆ Odd_p() [1/2]

virtual void* quda::GaugeField::Odd_p ( )
inlinevirtual

Reimplemented in quda::cudaGaugeField.

Definition at line 317 of file gauge_field.h.

References errorQuda.

◆ Odd_p() [2/2]

virtual const void* quda::GaugeField::Odd_p ( ) const
inlinevirtual

Reimplemented in quda::cudaGaugeField.

Definition at line 321 of file gauge_field.h.

References errorQuda.

◆ Order()

QudaGaugeFieldOrder quda::GaugeField::Order ( ) const
inline

◆ PhaseBytes()

size_t quda::GaugeField::PhaseBytes ( ) const
inline

Definition at line 312 of file gauge_field.h.

◆ PhaseOffset()

size_t quda::GaugeField::PhaseOffset ( ) const
inline

Definition at line 313 of file gauge_field.h.

◆ Reconstruct()

QudaReconstructType quda::GaugeField::Reconstruct ( ) const
inline

◆ removeStaggeredPhase()

void quda::GaugeField::removeStaggeredPhase ( )

Remove the staggered phase factors from the gauge field.

Definition at line 154 of file gauge_field.cpp.

References quda::applyGaugePhase(), errorQuda, quda::cudaGaugeField::exchangeGhost(), quda::cpuGaugeField::exchangeGhost(), quda::LatticeField::ghostExchange, QUDA_GHOST_EXCHANGE_PAD, and staggeredPhaseApplied.

Referenced by remove_staggered_phase_quda_(), and staggeredPhaseQuda().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setTuningString()

void quda::GaugeField::setTuningString ( )
protectedvirtual

Set the vol_string and aux_string for use in tuning.

Reimplemented from quda::LatticeField.

Definition at line 104 of file gauge_field.cpp.

References quda::TuneKey::aux_n, quda::LatticeField::aux_string, errorQuda, geometry, nColor, quda::LatticeField::precision, quda::LatticeField::setTuningString(), quda::LatticeField::stride, and quda::LatticeField::volume.

Referenced by GaugeField().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SiteOffset()

size_t quda::GaugeField::SiteOffset ( ) const
inline
Returns
The offset into the struct to the start of the gauge field (only for order = QUDA_MILC_SITE_GAUGE_ORDER)

Definition at line 337 of file gauge_field.h.

References quda::GaugeFieldParam::site_offset.

◆ SiteSize()

size_t quda::GaugeField::SiteSize ( ) const
inline
Returns
The size of the struct into which the gauge field is packed (only for order = QUDA_MILC_SITE_GAUGE_ORDER)

Definition at line 343 of file gauge_field.h.

References quda::copy(), quda::norm1(), quda::norm2(), quda::GaugeFieldParam::site_size, and quda::zero().

Here is the call graph for this function:

◆ StaggeredPhase()

QudaStaggeredPhase quda::GaugeField::StaggeredPhase ( ) const
inline

◆ StaggeredPhaseApplied()

bool quda::GaugeField::StaggeredPhaseApplied ( ) const
inline

◆ Tadpole()

double quda::GaugeField::Tadpole ( ) const
inline

Definition at line 253 of file gauge_field.h.

References quda::GaugeFieldParam::tadpole.

Referenced by createExtendedGauge().

Here is the caller graph for this function:

◆ TBoundary()

QudaTboundary quda::GaugeField::TBoundary ( ) const
inline

Definition at line 254 of file gauge_field.h.

References quda::GaugeFieldParam::t_boundary.

Referenced by quda::CoarseOp(), and createExtendedGauge().

Here is the caller graph for this function:

◆ zero()

virtual void quda::GaugeField::zero ( )
pure virtual

Set all field elements to zero (virtual)

Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.

Referenced by quda::calculateY().

Here is the caller graph for this function:

Member Data Documentation

◆ anisotropy

double quda::GaugeField::anisotropy
protected

Definition at line 183 of file gauge_field.h.

Referenced by checkField(), and GaugeField().

◆ bytes

size_t quda::GaugeField::bytes
protected

◆ create

QudaFieldCreate quda::GaugeField::create
protected

◆ fat_link_max

double quda::GaugeField::fat_link_max
protected

◆ fixed

QudaGaugeFixed quda::GaugeField::fixed
protected

Definition at line 179 of file gauge_field.h.

Referenced by checkField(), and GaugeField().

◆ geometry

QudaFieldGeometry quda::GaugeField::geometry
protected

◆ ghost

void* quda::GaugeField::ghost[2 *QUDA_MAX_DIM]
mutableprotected

◆ ghostFace

int quda::GaugeField::ghostFace[QUDA_MAX_DIM]
mutableprotected

Definition at line 191 of file gauge_field.h.

Referenced by createGhostZone().

◆ i_mu

double quda::GaugeField::i_mu
protected

Imaginary chemical potential

Definition at line 214 of file gauge_field.h.

◆ length

size_t quda::GaugeField::length
protected

Definition at line 170 of file gauge_field.h.

Referenced by GaugeField().

◆ link_type

QudaLinkType quda::GaugeField::link_type
protected

◆ nColor

int quda::GaugeField::nColor
protected

◆ nFace

int quda::GaugeField::nFace
protected

◆ nInternal

int quda::GaugeField::nInternal
protected

◆ order

QudaGaugeFieldOrder quda::GaugeField::order
protected

◆ phase_bytes

size_t quda::GaugeField::phase_bytes
protected

Definition at line 169 of file gauge_field.h.

Referenced by GaugeField(), and quda::cudaGaugeField::zeroPad().

◆ phase_offset

size_t quda::GaugeField::phase_offset
protected

Definition at line 168 of file gauge_field.h.

Referenced by quda::cudaGaugeField::cudaGaugeField(), and GaugeField().

◆ real_length

size_t quda::GaugeField::real_length
protected

Definition at line 171 of file gauge_field.h.

Referenced by GaugeField().

◆ reconstruct

QudaReconstructType quda::GaugeField::reconstruct
protected

◆ site_offset

size_t quda::GaugeField::site_offset
protected

Offset into MILC site struct to the desired matrix field (only if gauge_order=MILC_SITE_GAUGE_ORDER)

Definition at line 219 of file gauge_field.h.

◆ site_size

size_t quda::GaugeField::site_size
protected

Size of MILC site struct (only if gauge_order=MILC_SITE_GAUGE_ORDER)

Definition at line 224 of file gauge_field.h.

Referenced by quda::cpuGaugeField::cpuGaugeField().

◆ staggeredPhaseApplied

bool quda::GaugeField::staggeredPhaseApplied
protected

Whether the staggered phase factor has been applied

Definition at line 201 of file gauge_field.h.

Referenced by applyStaggeredPhase(), quda::cudaGaugeField::copy(), removeStaggeredPhase(), and quda::cudaGaugeField::saveCPUField().

◆ staggeredPhaseType

QudaStaggeredPhase quda::GaugeField::staggeredPhaseType
protected

The staggered phase convention to use

Definition at line 196 of file gauge_field.h.

Referenced by applyStaggeredPhase(), quda::cudaGaugeField::copy(), and quda::cudaGaugeField::saveCPUField().

◆ t_boundary

QudaTboundary quda::GaugeField::t_boundary
protected

Definition at line 181 of file gauge_field.h.

Referenced by checkField().

◆ tadpole

double quda::GaugeField::tadpole
protected

Definition at line 184 of file gauge_field.h.

Referenced by checkField().


The documentation for this class was generated from the following files: