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.HMCType
HMC(
    U,
    integrator,
    trajectory,
    steps,
    friction = 0,
    numsmear = 0,
    ρ_stout = 0;
    hmc_logging = true,
    fermion_action = QuenchedFermionAction,
    heavy_flavours = 0,
    num_cv = 0,
    logdir = "",
)

Create an HMC object, that can be used as an update algorithm.

Arguments

  • U: The gauge field on which the update is performed.
  • integrator: The integrator used to evolve the field.
  • 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: Number of Stout smearing steps applied to the gauge action.
  • ρ_stout: Step length of the Stout smearing applied to the gauge action.
  • hmc_logging: If true, creates a logfile in logdir containing information

on the trajectories, unless logdir = ""

  • fermion_action: An AbstratFermionAction to initialize the appropriate fermion fields
  • heavy_flavours: The number of non-degenerate heavy flavours, again to initialize the

right number of fermion fields

  • num_cv: If bigger than 0, additional fields are initialized that are needed for Stout

force recursion when using a bias.

Supported Integrators

  • Leapfrog
  • OMF2
  • OMF2Slow
  • OMF4
  • OMF4Slow

Supported Fermion Actions

  • WilsonFermionAction
  • WilsonEOPreFermionAction
  • StaggeredFermionAction
  • StaggeredEOPreFermionAction
source
MetaQCD.Updates.MetropolisType
Metropolis(U::Gaugefield{B,T,A,GA}, eo, ϵ, numhits, target_acc, or_alg, numorelax) where {B,T,A,GA}

Create a Metropolis object.

Arguments

  • U::Gaugefield{B,T,A,GA}: Gauge field object.
  • 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.

source
MetaQCD.Updates.HeatbathType
Heatbath(U::Gaugefield{B,T,A,GA}, MAXIT, numheatbath, or_alg, numorelax) where {B,T,A,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.

source
MetaQCD.Updates.OverrelaxationType
Overrelaxation(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)
source