Gradient flow
The gradient flow equations can be integrated in two different ways:
- Using a fixed step-size integrator. In this approach one fixes the step size $\epsilon$ and the links are evolved from $V_\mu(t)$ to $V_\mu(t +\epsilon)$ using some integration scheme.
- Using an adaptive step-size integrator. In this approach one fixes the tolerance $h$ and the links are evolved for a time $t_{\rm end}$ (i.e. from $V_\mu(t)$ to $V_\mu(t +t_{\rm end})$) with the condition that the maximum error while advancing is not larger than $h$.
In general adaptive step size integrators are much more efficient, but one loses the possibility to measure flow quantities at the intermediate times $\epsilon, 2\epsilon, 3\epsilon,...$. Adaptive step size integrators are ideal for finite size scaling studies, while a mix of both integrators is the most efficient approach in scale setting applications.
Integration schemes
LatticeGPU.YM.FlowIntr
— Type struct FlowIntr{N,T}
Structure containing info about a particular flow integrator
LatticeGPU.YM.wfl_euler
— Function wfl_euler(::Type{T}, eps::T, tol::T)
Euler scheme integrator for the Wilson Flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
LatticeGPU.YM.zfl_euler
— Function zfl_euler(::Type{T}, eps::T, tol::T)
Euler scheme integrator for the Zeuthen flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
LatticeGPU.YM.wfl_rk2
— Function wfl_rk2(::Type{T}, eps::T, tol::T)
Second order Runge-Kutta integrator for the Wilson flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
LatticeGPU.YM.zfl_rk2
— Function zfl_rk2(::Type{T}, eps::T, tol::T)
Second order Runge-Kutta integrator for the Zeuthen flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
LatticeGPU.YM.wfl_rk3
— Function wfl_rk3(::Type{T}, eps::T, tol::T)
Third order Runge-Kutta integrator for the Wilson flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
LatticeGPU.YM.zfl_rk3
— Function Zfl_rk3(::Type{T}, eps::T, tol::T)
Third order Runge-Kutta integrator for the Zeuthen flow. The fixed step size is given by eps
and the tolerance for the adaptive integrators by tol
.
Integrating the flow equations
LatticeGPU.YM.flw
— Function function flw(U, int::FlowIntr{NI,T}, ns::Int64, gp::GaugeParm, lp::SpaceParm, ymws::YMworkspace)
Integrates the flow equations with the integration scheme defined by int
performing ns
steps with fixed step size. The configuration U
is overwritten.
LatticeGPU.YM.flw_adapt
— Function function flw_adapt(U, int::FlowIntr{NI,T}, tend::T, gp::GaugeParm, lp::SpaceParm, ymws::YMworkspace)
Integrates the flow equations with the integration scheme defined by int
using the adaptive step size integrator up to tend
with the tolerance defined in int
. The configuration U
is overwritten.
Observables
LatticeGPU.YM.Eoft_plaq
— Functionfunction Eoft_plaq([Eslc,] U, gp::GaugeParm, lp::SpaceParm, ymws::YMworkspace)
Measure the action density E(t)
using the plaquette discretization. If the argument Eslc
is given the contribution for each Euclidean time slice and plane are returned.
LatticeGPU.YM.Eoft_clover
— Functionfunction Eoft_clover([Eslc,] U, gp::GaugeParm, lp::SpaceParm, ymws::YMworkspace)
Measure the action density E(t)
using the clover discretization. If the argument Eslc
the contribution for each Euclidean time slice and plane are returned.
LatticeGPU.YM.Qtop
— FunctionQtop([Qslc,] U, gp::GaugeParm, lp::SpaceParm, ymws::YMworkspace)
Measure the topological charge Q
of the configuration U
using the clover definition of the field strength tensor. If the argument Qslc
is present the contribution for each Euclidean time slice are returned. Only wors in 4D.