-
Notifications
You must be signed in to change notification settings - Fork 5
Python API
def init(grid_size: List[int] = None, latt_size: List[int] = None, t_boundary: Literal[1, -1] = None, anisotropy: float = None, backend: Literal["numpy", "cupy", "torch"] = "cupy", **kwargs)
def readMPIFile( filename: str, dtype: str, offset: int, shape: Sequence[int], axes: Sequence[int]) -> numpy.ndarray
def writeMPIFile(filename: str, dtype: str, offset: int, shape: Sequence[int], axes: Sequence[int], buf: numpy.ndarray)
-
LatticeInfo(latt_size: List[int], t_boundary: Literal[1, -1] = 1, anisotropy: float = 1.0) -
mpi_comm: MPI.Comm -
mpi_size: int -
mpi_rank: int -
grid_size: List[int] -
grid_coord: List[int] -
global_size: List[int] -
global_volume: int -
size: List[int] -
volume: int -
size_cb2: List[int] -
volume_cb2: int -
t_boundary: Literat[-1, 1] -
anisotropy: float
-
latt_info: LatticeInfo -
data: numpy.ndarray | cupy.ndarray | torch.Tensor -
backend: "numpy" | "cupy" | "torch" -
location: "numpy" | "cupy" | "torch" -
def backup() -> numpy.ndarray | cupy.ndarray | torch.Tensor -
def copy() -> numpy.ndarray | cupy.ndarray | torch.Tensor -
def toDevice() -
def toHost() -
def getHost() -> numpy.ndarray
-
LatticeGauge(latt_info: LatticeInfo, value=None) -
def setAntiPeriodicT() -
def setAnisotropy() -
data_ptr: Pointer -
data_ptrs: Pointers -
def lexico() -> numpy.ndarray -
def ensurePureGauge() -
covDev(x: LatticeFerimon, covdev_mu: int) -> LatticeFermionApplies the covariant derivative on
xin the directioncovdev_mu.xshould beLatticeFermion. 0/1/2/3 represent +x/+y/+z/+t and 4/5/6/7 represent -x/-y/-z/-t. The covariant derivative is defined as$\psi'(x)=U_\mu(x)\psi(x+\hat{\mu})$ . -
def laplace(x: LatticeStaggeredFermion, laplace3D: int) -> LatticeStaggeredFermionApplies the Laplacian operator on
x, andlaplace3Dtakes 3 or 4 to apply Laplacian on spacial or all directions.xshould beLatticeStaggeredFermion. The Laplacian operator is defined as$\psi'(x)=\frac{1}{N_\mathrm{Lap}}\sum_\mu\psi(x)-\dfrac{1}{2}\left[U_\mu(x)\psi(x+\hat{\mu})+U_\mu^\dagger(x-\hat{\mu})\psi(x-\hat{\mu})\right]$ -
def staggeredPhase()Applies the staggered phase to the gauge field. The convention is controld by
LatticeGauge.pure_gauge.gauge_param.staggered_phase_type, which isQudaStaggeredPhase.QUDA_STAGGERED_PHASE_MILCby default. -
def projectSU3(tol: float)Projects the gauge field onto SU(3) matrix.
tolis the tolerance of how matrix deviates from SU(3).2e-15(which is 10x the epsilon of fp64) should be a good choice. -
def path(paths: List[List[int]]) -
def loop(loops: List[List[List[int]]], coeff: List[float])loopsis a list of length 4, which is[paths_x, paths_y, paths_z, paths_t]. Allpaths_*should have the same shape.paths_xis a list of any length. -
def loopTrace(loops: List[List[int]]) -> numpy.ndarrayloopsis similar topaths_xinLatticeGauge.path, but the function returns the traces of all loops. The traces is anumpy.ndarrayof complex number, and every element is defined as$\sum_x\mathrm{Tr}W(x)$ . -
def apeSmear(n_steps: int, alpha: float, dir_ignore: int)Applies the APE smearing to the gauge field.
alphais the smearing strength. -
def smearAPE(n_steps, factor, dir_ignore)Similar to
LatticeGauge.apeSmear()butfactormatchs Chroma. -
def stoutSmear(n_steps, rho, dir_ignore)Applies the stout smearing to the gauge field.
rhois the smearing strength. -
def smearSTOUT(n_steps, rho, dir_ignore)Similar to
LatticeGauge.stoutSmear(). -
def hypSmear(n_steps, alpha1, alpha2, alpha3, dir_ignore)Applies the stout smearing to the gauge field.
alpha1/alpha2/alpha3is the smearing strength on level 3/2/1. -
def smearHYP(n_steps, alpha1, alpha2, alpha3, dir_ignore)Similar to
LatticeGauge.hypSmear(). -
def wilsonFlow(n_steps, epsilon)Applies the Wilson flow to the gauge field. Returns the energy (all, spacial, temporal) for every step.
-
def wilsonFlowScale(max_steps, epsilon)Returns
$t_0$ and$w_0$ with up tomax_stepsWilson flow step. -
def symanzikFlow(n_steps, epsilon)Applies the Symanzik flow to the gauge field. Returns the energy (all, spacial, temporal) for every step.
-
def symanzikFlowScale(max_steps, epsilon)Returns
$t_0$ and$w_0$ with up tomax_stepsSymanzik flow step. -
def plaquette()Returns the plaquette (all, spatial, temporal) of the gauge field.
-
def polyakovLoop()Returns the Polyakov loop (real, image) of the gauge field.
-
def energy() -> List[int]Returns the energy (all, spacial, temporal) of the gauge field.
-
def qcharge() -> floatReturns the topological charge of the gauge field.
-
def qchargeDensity() -> numpy.ndarrayReturns the topological charge density with shape
(2, Lt, Lz, Ly, Lx // 2))`. -
def gauss(seed, sigma)Fills the gauge field with random SU(3) matrices.
sigma=1corresponds to the standard normal distribution. -
def fixingOVR(gauge_dir, Nsteps, verbose_interval, relax_boost, tolerance, reunit_interval, stopWtheta)Applies the gauge gixing to the gauge field and over relaxation is used to speed up the operation.
gauge_dir: {3, 4} 3 for Coulomb gauge fixing, 4 for Landau gauge fixing Nsteps: int maximum number of steps to perform gauge fixing verbose_interval: int print gauge fixing info when iteration count is a multiple of this relax_boost: float gauge fixing parameter of the overrelaxation method, most common value is 1.5 or 1.7. tolerance: float torelance value to stop the method, if this value is zero then the method stops when iteration reachs the maximum number of steps defined by Nsteps reunit_interval: int reunitarize gauge field when iteration count is a multiple of this stopWtheta: int 0 for MILC criterion and 1 to use the theta value -
def fixingFFT(gauge_dir, Nsteps, verbose_interval, alpha, autotune, tolerance, stopWtheta)Applies the gauge gixing to the gauge field and fast Fourier transform is used to speed up the operation.
gauge_dir: {3, 4} 3 for Coulomb gauge fixing, 4 for Landau gauge fixing Nsteps: int maximum number of steps to perform gauge fixing verbose_interval: int print gauge fixing info when iteration count is a multiple of this alpha: float gauge fixing parameter of the method, most common value is 0.08 autotune: int 1 to autotune the method, i.e., if the Fg inverts its tendency we decrease the alpha value tolerance: float torelance value to stop the method, if this value is zero then the method stops when iteration reachs the maximum number of steps defined by Nsteps stopWtheta: int 0 for MILC criterion and 1 to use the theta value