Updating a Gaugefield
To update a Gaugefield simply use the update! function that takes 2 positional arguments and however many keyword arguments specific to the update algorithm. The first positional is the actual update algorithm update_alg and the second is the Gaugefield U.
U = Gaugefield(...)
random_gauges!(U)
MAXIT = 100
numHB = 1
or_alg = MeatQCD.Updates.Subgroups
numOR = 4
update_alg = MeatQCD.Updates.Heatbath(U, MAXIT, numHB, or_alg, numOR)
update!(update_alg, U; ...) Supported Update Algorithms
MetaQCD.Updates.HMC — TypeHMC(
U,
hmc_levels,
trajectory,
friction=0.0,
numsmear_gauge=0,
numsmear_fermion=0,
rho_stout_gauge=0.0,
rho_stout_fermion=0.0;
rafriction=0.0,
hmc_logging=true,
fermion_action="quenched",
numfermions=0,
numcv=0,
logdir="",
instance=MPI_INSTANCE[],
)Create an HMC object, that can be used as an update algorithm.
Arguments
U: The gauge field on which the update is performed.levels: A vector ofDicts that define the levels of the hmc integration scheme.
To see what parameters are needed see the file "./hmc_levels".
trajectory: The length of the HMC trajectory.steps: The number of integrator steps within the trajectory.friction: Friction factor in the GHMC algorithm. Has to be in the range [0, 1].numsmear_gauge: Number of Stout smearing steps applied to the gauge action.numsmear_fermion: Number of Stout smearing steps applied to the fermion action.rho_stout_gauge: Step length of the Stout smearing applied to the gauge action.rho_stout_fermion: Step length of the Stout smearing applied to the fermion action.rafriction: Friction parameter for the repell-attract HMC.hmc_logging: If true, creates a logfile inlogdircontaining information
on the trajectories, unless logdir = ""
fermion_action: An String that identifies the fermion action type to initialize the appropriate fermion fieldsnumfermions: The number of non-degenerate heavy flavours, again to initialize the
right number of fermion fields
numcv: If bigger than 0, additional fields are initialized that are needed for Stout
force recursion when using a bias.
logdir: Directory that hmc data should be written into.instance: Integer identifier of current instance (for parallel tempering and multiple walkers).
Supported Fermion Actions
quenchedstaggeredstaggered_eowilsonwilson_eo
MetaQCD.Updates.Metropolis — TypeMetropolis(U::Gaugefield, eo, ϵ, numhits, target_acc, or_alg, numorelax) where {B,T,A,GA}Create a Metropolis object.
Arguments
U::Gaugefield: Gauge field.eo: Even-odd preconditioning.ϵ: Step size for the update.numhits: Number of Metropolis hits.target_acc: Target acceptance rate.or_alg: Overrelaxation algorithm.numorelax: Number of overrelaxation sweeps.
Returns
A Metropolis object with the specified parameters. The gauge action GA of the field U determines the iterator used. For the plaquette or Wilson action it uses a Checkerboard iterator and for rectangular actions it partitions the lattice into four sublattices.
MetaQCD.Updates.Heatbath — TypeHeatbath(U::Gaugefield, MAXIT, numheatbath, or_alg, numorelax) where {B,T,GA}Create a Heatbath` object.
Arguments
U: The gauge field on which the update is performed.MAXIT: The maximum iteration count in the Heatbath update.numheatbath: The number of Heatbath sweeps.or_alg: The overrelaxation algorithm used.numorelax: The number of overrelaxation sweeps.
Returns
A Heatbath object with the specified parameters. The gauge action GA of the field U determines the iterator used. For the plaquette or Wilson action it uses a Checkerboard iterator and for rectangular actions it partitions the lattice into four sublattices.
MetaQCD.Updates.Overrelaxation — TypeOverrelaxation(algorithm)Create an Overrelaxation object, that can be used within a Metropolis or Heatbath update step.
Supported Algorithms
"subgroups": Cabibbo-Marinari SU(2) subgroup embedding scheme"kenney_laub": Kenney-Laub projection onto SU(3)