QUDA
1.0.0
|
#include <gauge_field.h>
Public Member Functions | |
GaugeField (const GaugeFieldParam ¶m) | |
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 |
![]() | |
LatticeField (const LatticeFieldParam ¶m) | |
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... | |
![]() | |
Object () | |
virtual | ~Object () |
Static Public Member Functions | |
static GaugeField * | Create (const GaugeFieldParam ¶m) |
Create the gauge field, with meta data specified in the parameter struct. More... | |
![]() | |
static void | freeGhostBuffer (void) |
Free statically allocated ghost buffers. More... | |
static void | destroyIPCComms () |
![]() | |
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... | |
![]() | |
void | precisionCheck () |
Additional Inherited Members | |
![]() | |
static int | bufferIndex = 0 |
static bool | ghost_field_reset = false |
![]() | |
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 MsgHandle * | mh_send_p2p_fwd [2][QUDA_MAX_DIM] { } |
static MsgHandle * | mh_send_p2p_back [2][QUDA_MAX_DIM] { } |
static MsgHandle * | mh_recv_p2p_fwd [2][QUDA_MAX_DIM] { } |
static MsgHandle * | mh_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 |
Definition at line 164 of file gauge_field.h.
quda::GaugeField::GaugeField | ( | const GaugeFieldParam & | param | ) |
Definition at line 30 of file gauge_field.cpp.
References ALIGNMENT_ADJUST, anisotropy, bytes, errorQuda, fixed, geometry, quda::LatticeField::ghost_precision, isNative(), length, link_type, nColor, quda::LatticeField::nDim, nInternal, phase_bytes, phase_offset, quda::LatticeField::precision, QUDA_COARSE_GEOMETRY, QUDA_COARSE_LINKS, QUDA_GAUGE_FIXED_YES, QUDA_RECONSTRUCT_13, QUDA_RECONSTRUCT_9, QUDA_SCALAR_GEOMETRY, QUDA_TENSOR_GEOMETRY, QUDA_VECTOR_GEOMETRY, QUDA_WILSON_LINKS, real_length, reconstruct, setTuningString(), quda::LatticeField::stride, quda::LatticeField::total_bytes, and quda::LatticeField::volume.
|
virtual |
Definition at line 100 of file gauge_field.cpp.
double quda::GaugeField::abs_max | ( | int | dim = -1 | ) | const |
Compute the absolute maximum of the field (Linfinity norm)
[in] | dim | Which dimension we are taking the norm of (dim=-1 mean all dimensions) |
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().
double quda::GaugeField::abs_min | ( | int | dim = -1 | ) | const |
Compute the absolute minimum of the field.
[in] | dim | Which dimension we are taking the norm of (dim=-1 mean all dimensions) |
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.
|
inline |
Definition at line 252 of file gauge_field.h.
References quda::GaugeFieldParam::anisotropy.
Referenced by cloverQuda(), quda::CoarseOp(), createExtendedGauge(), loadCloverQuda(), and quda::setDiracParam().
void quda::GaugeField::applyStaggeredPhase | ( | QudaStaggeredPhase | phase = QUDA_STAGGERED_PHASE_INVALID | ) |
Apply the staggered phase factors to the gauge field.
[in] | phase | The 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().
|
inline |
Definition at line 311 of file gauge_field.h.
References quda::blas::bytes.
Referenced by quda::GaugeGauss< Float, Arg >::bytes(), quda::DslashCoarsePolicyTune::bytes(), computeHISQForceQuda(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), and quda::cudaGaugeField::saveCPUField().
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().
uint64_t quda::GaugeField::checksum | ( | bool | mini = false | ) | const |
Compute checksum of this gauge field: this uses a XOR-based checksum method
[in] | mini | Whether 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. |
Definition at line 355 of file gauge_field.cpp.
References quda::Checksum().
Referenced by loadGaugeQuda().
|
pure virtual |
Generic gauge field copy
[in] | src | Source from which we are copying |
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Referenced by quda::CoarseOp().
|
static |
Create the gauge field, with meta data specified in the parameter struct.
param | Parameter struct specifying the 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().
|
protected |
Compute the required extended ghost zone sizes and offsets
[in] | R | Radius of the ghost zone |
[in] | no_comms_fill | If true we create a full halo regardless of partitioning |
[in] | bidir | Is 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().
|
inlinevirtual |
Reimplemented in quda::cudaGaugeField.
Definition at line 316 of file gauge_field.h.
References errorQuda.
|
inlinevirtual |
Reimplemented in quda::cudaGaugeField.
Definition at line 320 of file gauge_field.h.
References errorQuda.
|
protected |
Exchange the buffers across all dimensions in a given direction.
[out] | recv | Receive buffer |
[in] | send | Send buffer |
[in] | dir | Direction 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().
|
pure virtual |
This does routine will populate the border / halo region of a gauge field that has been created using copyExtendedGauge.
R | The thickness of the extended region in each dimension |
no_comms_fill | Do local exchange to fill out the extended region in non-partitioned dimensions |
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Referenced by quda::gaugeGauss().
|
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.
R | The thickness of the extended region in each dimension |
profile | TimeProfile intance which will record the time taken |
no_comms_fill | Do local exchange to fill out the extended region in non-partitioned dimensions |
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
|
pure virtual |
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Referenced by quda::applyGaugePhase(), and quda::gaugeGauss().
|
inline |
Definition at line 257 of file gauge_field.h.
References quda::GaugeFieldParam::order.
Referenced by quda::colorSpinorParam(), quda::CovDevArg< Float, nColor, reconstruct_ >::CovDevArg(), quda::LaplaceArg< Float, nColor, reconstruct_ >::LaplaceArg(), quda::norm(), quda::StaggeredArg< Float, nColor, reconstruct_u_, reconstruct_l_, improved_, phase_ >::StaggeredArg(), and quda::WilsonArg< Float, nColor, reconstruct_ >::WilsonArg().
|
inlinevirtual |
Reimplemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Definition at line 315 of file gauge_field.h.
References errorQuda.
Referenced by quda::gauge::Accessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, storeFloat, use_tex >::Accessor(), quda::gauge::Accessor< Float, nColor, QUDA_FLOAT2_GAUGE_ORDER, storeFloat, use_tex >::Accessor(), quda::colorSpinorParam(), quda::cudaGaugeField::copy(), quda::gauge::FloatNOrder< Float, length, N, reconLenParam, stag_phase, huge_alloc, ghostExchange_, use_inphase >::FloatNOrder(), quda::CalculateY< from_coarse, Float, fineSpin, fineColor, coarseSpin, coarseColor, Arg >::postTune(), quda::CalculateY< from_coarse, Float, fineSpin, fineColor, coarseSpin, coarseColor, Arg >::preTune(), quda::gauge::QDPJITOrder< Float, length >::QDPJITOrder(), and quda::gauge::QDPOrder< Float, length >::QDPOrder().
|
inlinevirtual |
Reimplemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Definition at line 319 of file gauge_field.h.
References errorQuda.
|
inline |
Definition at line 256 of file gauge_field.h.
References quda::GaugeFieldParam::fixed.
Referenced by quda::CoarseOp().
|
inline |
Definition at line 258 of file gauge_field.h.
References quda::GaugeFieldParam::geometry.
Referenced by quda::gauge::Accessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, storeFloat, use_tex >::Accessor(), quda::cloverDerivative(), quda::CoarseOp(), quda::colorSpinorParam(), quda::cudaGaugeField::copy(), quda::copyGauge(), quda::copyGaugeEx(), quda::copyGenericGauge(), createExtendedGauge(), quda::gauge::GhostAccessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_MILC_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_FLOAT2_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::cudaGaugeField::saveCPUField(), saveGaugeFieldQuda(), and quda::CopyGauge< FloatOut, FloatIn, length, Arg >::set_ghost().
|
inline |
Definition at line 323 of file gauge_field.h.
References errorQuda.
Referenced by quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), dslash_4_4d(), dw_dslash(), eigensolve_test(), quda::gauge::GhostAccessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_MILC_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), init(), invert_test(), quda::gauge::LegacyOrder< Float, length >::LegacyOrder(), quda::cudaGaugeField::saveCPUField(), and wil_dslash().
|
inline |
Definition at line 328 of file gauge_field.h.
References errorQuda.
|
inline |
Return the imaginary chemical potential applied to this field
Definition at line 278 of file gauge_field.h.
References quda::GaugeFieldParam::i_mu.
|
pure virtual |
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
bool quda::GaugeField::isNative | ( | ) | const |
This function returns true if the field is stored in an internal field order for the given precision.
Definition at line 167 of file gauge_field.cpp.
References order, quda::LatticeField::precision, QUDA_DOUBLE_PRECISION, QUDA_FLOAT2_GAUGE_ORDER, QUDA_FLOAT4_GAUGE_ORDER, QUDA_HALF_PRECISION, QUDA_QUARTER_PRECISION, QUDA_RECONSTRUCT_10, QUDA_RECONSTRUCT_12, QUDA_RECONSTRUCT_13, QUDA_RECONSTRUCT_8, QUDA_RECONSTRUCT_9, QUDA_RECONSTRUCT_NO, QUDA_SINGLE_PRECISION, and reconstruct.
Referenced by quda::APEStep(), quda::computeQCharge(), quda::computeQChargeDensity(), quda::cpuGaugeField::copy(), quda::copyGauge(), quda::copyGaugeEx(), quda::copyGaugeMG(), quda::CovDevArg< Float, nColor, reconstruct_ >::CovDevArg(), createGhostZone(), quda::cudaGaugeField::cudaGaugeField(), quda::cudaGaugeField::exchangeGhost(), quda::extractGhost(), quda::extractGhostEx(), quda::extractGhostMG(), GaugeField(), quda::gaugeGauss(), quda::cudaGaugeField::injectGhost(), quda::isUnitary(), quda::LaplaceArg< Float, nColor, reconstruct_ >::LaplaceArg(), quda::OvrImpSTOUTStep(), quda::StaggeredArg< Float, nColor, reconstruct_u_, reconstruct_l_, improved_, phase_ >::StaggeredArg(), quda::STOUTStep(), quda::WilsonArg< Float, nColor, reconstruct_ >::WilsonArg(), quda::cudaGaugeField::zeroPad(), and quda::cudaGaugeField::~cudaGaugeField().
|
inline |
Definition at line 248 of file gauge_field.h.
References length.
|
inline |
Definition at line 280 of file gauge_field.h.
References fat_link_max.
Referenced by quda::cudaGaugeField::copy(), and quda::cpuGaugeField::copy().
|
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().
|
inline |
Definition at line 249 of file gauge_field.h.
References quda::GaugeFieldParam::nColor.
Referenced by quda::Checksum(), quda::colorSpinorParam(), quda::copyGauge(), quda::copyGaugeEx(), quda::copyGaugeMG(), quda::copyGenericGauge(), quda::extractGaugeGhost(), quda::extractGhostMG(), quda::gaugeGauss(), quda::gauge::GhostAccessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_MILC_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::instantiate(), invert_multishift_quda_(), and quda::norm().
|
inline |
Definition at line 281 of file gauge_field.h.
References quda::exchangeExtendedGhost(), quda::GaugeFieldParam::nFace, and R.
Referenced by quda::CopyGaugeArg< OutOrder, InOrder >::CopyGaugeArg(), quda::copyGaugeEx(), quda::ExtractGhostArg< Float, nColor_, Order, nDim >::ExtractGhostArg(), quda::gauge::FloatNOrder< Float, length, N, reconLenParam, stag_phase, huge_alloc, ghostExchange_, use_inphase >::FloatNOrder(), quda::gauge::GhostAccessor< Float, nColor, QUDA_QDP_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_MILC_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), quda::gauge::GhostAccessor< Float, nColor, QUDA_FLOAT2_GAUGE_ORDER, native_ghost, storeFloat, use_tex >::GhostAccessor(), and quda::gauge::LegacyOrder< Float, length >::LegacyOrder().
double quda::GaugeField::norm1 | ( | int | dim = -1 | ) | const |
Compute the L1 norm of the field.
[in] | dim | Which dimension we are taking the norm of (dim=-1 mean all dimensions) |
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.
double quda::GaugeField::norm2 | ( | int | dim = -1 | ) | const |
Compute the L2 norm squared of the field.
[in] | dim | Which dimension we are taking the norm of (dim=-1 mean all dimensions) |
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.
|
inlinevirtual |
Reimplemented in quda::cudaGaugeField.
Definition at line 317 of file gauge_field.h.
References errorQuda.
|
inlinevirtual |
Reimplemented in quda::cudaGaugeField.
Definition at line 321 of file gauge_field.h.
References errorQuda.
|
inline |
Definition at line 251 of file gauge_field.h.
References quda::GaugeFieldParam::order.
Referenced by quda::APEStep(), quda::applyU(), quda::checkMomOrder(), quda::Checksum(), quda::computeCloverForce(), quda::computeCloverSigmaOprod(), quda::computeMomAction(), quda::computeQCharge(), quda::computeQChargeDensity(), quda::computeStaggeredOprod(), quda::cudaGaugeField::copy(), quda::copyGauge(), quda::copyGaugeEx(), quda::copyGaugeMG(), createExtendedGauge(), quda::extractGhost(), quda::extractGhostEx(), quda::extractGhostMG(), quda::forceRecord(), quda::gaugeGauss(), quda::isUnitary(), loadGaugeQuda(), quda::LatticeField::Nvec(), quda::OvrImpSTOUTStep(), quda::cudaGaugeField::saveCPUField(), quda::STOUTStep(), and quda::updateMomentum().
|
inline |
Definition at line 312 of file gauge_field.h.
|
inline |
Definition at line 313 of file gauge_field.h.
|
inline |
Definition at line 250 of file gauge_field.h.
References quda::GaugeFieldParam::reconstruct.
Referenced by quda::APEStep(), checkGauge(), quda::checkMomOrder(), quda::CoarseOp(), quda::colorSpinorParam(), quda::completeKSForce(), quda::computeCloverForce(), quda::computeKSLongLinkForce(), quda::computeMomAction(), quda::computeQCharge(), quda::computeQChargeDensity(), quda::cudaGaugeField::copy(), quda::copyGauge(), quda::copyGaugeEx(), quda::copyGaugeMG(), createCloverQuda(), createExtendedGauge(), quda::extractGhost(), quda::extractGhostEx(), quda::extractGhostMG(), quda::fatLongKSLink(), quda::gauge::FieldOrder< Float, nColor, nSpinCoarse, order, native_ghost, storeFloat, use_tex >::FieldOrder(), quda::forceRecord(), quda::gaugeGauss(), quda::instantiate(), quda::isUnitary(), loadGaugeQuda(), loadSloppyGaugeQuda(), quda::OvrImpSTOUTStep(), quda::projectSU3(), quda::cudaGaugeField::saveCPUField(), quda::STOUTStep(), quda::updateMomentum(), and quda::wuppertalStep().
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().
|
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().
|
inline |
Definition at line 337 of file gauge_field.h.
References quda::GaugeFieldParam::site_offset.
|
inline |
Definition at line 343 of file gauge_field.h.
References quda::copy(), quda::norm1(), quda::norm2(), quda::GaugeFieldParam::site_size, and quda::zero().
|
inline |
Definition at line 259 of file gauge_field.h.
References quda::GaugeFieldParam::staggeredPhaseType.
Referenced by computeStaggeredForceQuda(), quda::cudaGaugeField::copy(), quda::copyGauge(), quda::extractGhost(), and quda::StaggeredApply< Float, nColor, recon_u >::StaggeredApply().
|
inline |
Definition at line 260 of file gauge_field.h.
References QUDA_STAGGERED_PHASE_INVALID, and quda::GaugeFieldParam::staggeredPhaseApplied.
Referenced by computeStaggeredForceQuda(), quda::cudaGaugeField::copy(), quda::projectSU3(), and staggeredPhaseQuda().
|
inline |
Definition at line 253 of file gauge_field.h.
References quda::GaugeFieldParam::tadpole.
Referenced by createExtendedGauge().
|
inline |
Definition at line 254 of file gauge_field.h.
References quda::GaugeFieldParam::t_boundary.
Referenced by quda::CoarseOp(), and createExtendedGauge().
|
pure virtual |
Set all field elements to zero (virtual)
Implemented in quda::cpuGaugeField, and quda::cudaGaugeField.
Referenced by quda::calculateY().
|
protected |
Definition at line 183 of file gauge_field.h.
Referenced by checkField(), and GaugeField().
|
protected |
Definition at line 167 of file gauge_field.h.
Referenced by quda::cudaGaugeField::backup(), quda::cpuGaugeField::backup(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), exchange(), quda::cpuGaugeField::exchangeExtendedGhost(), GaugeField(), quda::cudaGaugeField::restore(), quda::cudaGaugeField::saveCPUField(), quda::cudaGaugeField::zero(), and quda::cudaGaugeField::zeroPad().
|
protected |
Definition at line 187 of file gauge_field.h.
Referenced by quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), quda::cudaGaugeField::setGauge(), quda::cpuGaugeField::setGauge(), quda::cpuGaugeField::~cpuGaugeField(), and quda::cudaGaugeField::~cudaGaugeField().
|
protected |
Definition at line 185 of file gauge_field.h.
Referenced by quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), and quda::cpuGaugeField::cpuGaugeField().
|
protected |
Definition at line 179 of file gauge_field.h.
Referenced by checkField(), and GaugeField().
|
protected |
Definition at line 174 of file gauge_field.h.
Referenced by quda::cpuGaugeField::backup(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), quda::cpuGaugeField::cpuGaugeField(), quda::create_gauge_buffer(), quda::create_ghost_buffer(), createGhostZone(), quda::cudaGaugeField::cudaGaugeField(), quda::cpuGaugeField::exchangeExtendedGhost(), quda::cudaGaugeField::exchangeGhost(), quda::cpuGaugeField::exchangeGhost(), quda::free_gauge_buffer(), quda::free_ghost_buffer(), GaugeField(), quda::cudaGaugeField::injectGhost(), quda::cpuGaugeField::injectGhost(), quda::cpuGaugeField::restore(), quda::cudaGaugeField::saveCPUField(), setTuningString(), quda::cudaGaugeField::zeroPad(), quda::cpuGaugeField::~cpuGaugeField(), and quda::cudaGaugeField::~cudaGaugeField().
|
mutableprotected |
Definition at line 189 of file gauge_field.h.
Referenced by quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), quda::cudaGaugeField::exchangeGhost(), quda::cpuGaugeField::exchangeGhost(), quda::cudaGaugeField::injectGhost(), quda::cpuGaugeField::injectGhost(), quda::cpuGaugeField::~cpuGaugeField(), and quda::cudaGaugeField::~cudaGaugeField().
|
mutableprotected |
Definition at line 191 of file gauge_field.h.
Referenced by createGhostZone().
|
protected |
Imaginary chemical potential
Definition at line 214 of file gauge_field.h.
|
protected |
Definition at line 170 of file gauge_field.h.
Referenced by GaugeField().
|
protected |
Definition at line 180 of file gauge_field.h.
Referenced by checkField(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), GaugeField(), and quda::cpuGaugeField::~cpuGaugeField().
|
protected |
Definition at line 172 of file gauge_field.h.
Referenced by checkField(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), GaugeField(), quda::cudaGaugeField::saveCPUField(), setTuningString(), and quda::cudaGaugeField::zeroPad().
|
protected |
Definition at line 173 of file gauge_field.h.
Referenced by checkField(), quda::cudaGaugeField::copy(), quda::cpuGaugeField::copy(), quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), exchange(), quda::cudaGaugeField::exchangeGhost(), quda::cpuGaugeField::exchangeGhost(), quda::cudaGaugeField::injectGhost(), quda::cpuGaugeField::injectGhost(), and quda::cudaGaugeField::saveCPUField().
|
protected |
Definition at line 177 of file gauge_field.h.
Referenced by quda::cpuGaugeField::cpuGaugeField(), createGhostZone(), quda::cudaGaugeField::cudaGaugeField(), exchange(), quda::cpuGaugeField::exchangeExtendedGhost(), quda::cpuGaugeField::exchangeGhost(), GaugeField(), quda::cpuGaugeField::injectGhost(), and quda::operator<<().
|
protected |
Definition at line 178 of file gauge_field.h.
Referenced by quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), quda::cpuGaugeField::exchangeExtendedGhost(), isNative(), quda::cudaGaugeField::zeroPad(), and quda::cpuGaugeField::~cpuGaugeField().
|
protected |
Definition at line 169 of file gauge_field.h.
Referenced by GaugeField(), and quda::cudaGaugeField::zeroPad().
|
protected |
Definition at line 168 of file gauge_field.h.
Referenced by quda::cudaGaugeField::cudaGaugeField(), and GaugeField().
|
protected |
Definition at line 171 of file gauge_field.h.
Referenced by GaugeField().
|
protected |
Definition at line 176 of file gauge_field.h.
Referenced by quda::cpuGaugeField::copy(), quda::cpuGaugeField::cpuGaugeField(), quda::cudaGaugeField::cudaGaugeField(), GaugeField(), isNative(), and quda::cudaGaugeField::zeroPad().
|
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.
|
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().
|
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().
|
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().
|
protected |
Definition at line 181 of file gauge_field.h.
Referenced by checkField().
|
protected |
Definition at line 184 of file gauge_field.h.
Referenced by checkField().