Ewald summation for Biot–Savart
The main originality of the VortexPasta solver is that it adapts the Ewald summation method to accelerate the computation of the Biot–Savart law along vortex filaments. See for example Arnold and Holm (2005) for a nice introduction to Ewald methods applied to the electrostatic interaction between point charges.
The adaptation of these methods to the vortex filament model is described in Polanco (2025). That paper also explains the role of the different parameters entering the method on accuracy and performance.
Splitting the Biot–Savart integral
The basic idea of the method is to split the Biot–Savart integral into short- and long-range parts:
where the short-range and long-range scalar functions
for all ; decays exponentially with , so that long-range interactions can be neglected beyond some cut-off distance when computing the short-range velocity ; is non-singular and smooth at , which implies that must quickly tend to 0 as . In periodic domains, this enables the use of the fast Fourier transform (FFT) to efficiently estimate long-range interactions.
The traditional choice of splitting functions satisfying all these properties is:
Here
The two splitting functions are plotted below. In the horizontal axis, the scale
Code for this figure
using SpecialFunctions: erf, erfc
using CairoMakie
rs = 2.0.^(-4:0.1:3)
gs(αr) = erfc(αr) + 2αr / sqrt(π) * exp(-αr^2) # short-range
gl(αr) = erf(αr) - 2αr / sqrt(π) * exp(-αr^2) # long-range
xticks = LogTicks(-4:1:3)
yticks = LogTicks(-16:4:0)
fig = Figure(size = (600, 400), fontsize = 18)
ax = Axis(fig[1, 1]; xticks, yticks, xscale = log2, yscale = log10, xlabel = L"αr", ylabel = "Splitting function")
ylims!(ax, 1e-17, 4)
ls = lines!(ax, rs, gs.(rs); label = L"g^<(r)")
ll = lines!(ax, rs, gl.(rs); label = L"g^>(r)")
lines!(ax, rs, erfc.(rs); label = L"\mathrm{erfc}(αr)", linestyle = :dot, color = ls.color)
lines!(ax, rs, erf.(rs); label = L"\mathrm{erf}(αr)", linestyle = :dot, color = ll.color)
let rs = 2.0.^(range(-4, -1; length = 3)), color = :grey20 # plot ~r^3 slope
ys = @. 0.2 * rs^3
lines!(ax, rs, ys; linestyle = :dash, color)
text!(ax, rs[2], ys[2]; text = L"r^3", align = (:left, :top), color)
end
let rs = 2.0.^(range(-4, -2; length = 3)), color = :grey20 # plot ~r^1 slope
ys = @. 0.5 * (2 / sqrt(π)) * rs
lines!(ax, rs, ys; linestyle = :dash, color)
text!(ax, rs[2], ys[2]; text = L"r", align = (:left, :top), color)
end
axislegend(ax; position = (0, 0), labelsize = 20)
save("splitting_functions.svg", fig)For small
Short-range velocity
As seen in the figure above, the short-range splitting function is dominant for small
Therefore, the short-range velocity is obtained as
i.e. the integral is performed over vortex points which are sufficiently close to the point of interest
Note that one can use a cell lists algorithm to further speed-up the search for nearby points.
Long-range velocity
The long-range velocity
More precisely, the long-range vorticity is given by the convolution
1. Estimating the vorticity in Fourier space
The idea is to first expand the (actual) vorticity field in Fourier series:
Here
To estimate this integral, we discretise the curves defining the support
The weights are related to the length of the discrete segments and to the quadrature rule used to estimate the integrals over segments (see the Numerical integration section). The important point here is that, since the discrete points are generally not on an equispaced grid, one cannot directly use the fast Fourier transform (FFT) to efficiently evaluate these coefficients. Nevertheless, they can be efficiently and accurately estimated using the non-uniform FFT (NUFFT) algorithm. More precisely, this corresponds to a type-1 NUFFT, which converts from non-uniform sources in physical space to uniform wavenumbers
2. Coarse-grained vorticity and velocity in Fourier space
Once we have obtained the
Similarly, the curl operator can be easily inverted in Fourier space to get the coarse-grained velocity:
Zero mean vorticity condition
The velocity is well defined only if
This condition is automatically satisfied when dealing with closed vortex filaments. This may however not be the case for infinite filaments.
To avoid this issue, we always set the mean vorticity to zero. In the case of infinite unclosed filaments, this effectively means that a non-zero total filament charge will be compensated by a uniform background vorticity BiotSavart.background_vorticity_correction!.
3. Notes on required resolution
Above we have assumed that we can evaluate Fourier coefficients for any wavenumber
Similarly to the cut-off distance in physical space, one can expect that the appropriate value of
Code for this figure
using CairoMakie
ks_α = range(0, 12.2; step = 0.1)
xticks = 0:12
yticks = LogTicks(-16:2:0)
ys = @. exp(-ks_α^2 / 4)
lines(
ks_α, ys;
axis = (yscale = log10, xticks, yticks, xlabel = L"k/α", ylabel = L"\exp \, \left[ -k^2 / 4α^2 \right]",),
figure = (fontsize = 20, size = (600, 400),),
)
save("gaussian_kalpha.svg", current_figure())4. Physical velocity at filament locations
The last step is to evaluate, from the coarse-grained velocity
for a set of locations
This operation can be efficiently computed using a type-2 NUFFT (from uniform wavenumbers
From the above steps, we can directly write the large-scale velocity at
which is analogous to the term for long-range electrostatic forces in e.g. Arnold and Holm (2005), Eq. (17).
Desingularisation
As discussed in the VFM page, the Biot–Savart integral must be desingularised if the velocity is to be evaluated on a point
In this case, the actual velocity that we want to compute is
which basically means that both short-range and long-range components should skip the computation of their respective integrals in the vicinity of
Note
The terms local/non-local and short/long-range are completely orthogonal and should not be confused! The first pair of terms refer to the VFM while the second pair refers to Ewald summation. For example, as discussed right below, it is possible (and it makes sense) to compute the local long-range velocity component.
For the short-range component, it is straightforward to skip the local part of the modified Biot–Savart integral: one just integrates over segments within the cut-off distance, but excluding the two segments which are in direct contact with
In the case of the long-range component special care is needed, as the procedure detailed above computes the full (coarse-grained) Biot–Savart integral, including the "singular" region. As discussed above, the modified long-range integral is actually not singular, so this is not a problem. But it means that we now need to subtract the local long-range velocity:
which can be numerically computed using the same method as for the non-local short-range velocity.
To summarise, in practice the velocity on a vortex position
where the first term is computed using an analytical expression based on Taylor expansions, the second and third terms are estimated by numerical integration (using quadratures), and the fourth term is indirectly computed using (non-uniform) fast Fourier transforms.
Numerical integration
Both the short-range and long-range computations use quadrature rules to approximate line integrals using numerical integration. The general strategy that has been implemented is to perform integrals segment-by-segment, and to use a quadrature rule of few points for each segment.
To be more explicit, say that we have a single closed vortex line (for example the trefoil in the VFM page) and that we want to compute the short-range non-local velocity induced by the vortex on
Here
where
Parameter selection
As detailed above, this method introduces a few parameters which must be tuned for accuracy and performance. In fact, most of these parameters are related and should not be treated independently. For instance, the physical- and Fourier-space cut-offs
In practice, one recommended way of setting the parameters is as follows:
Start by setting the physical domain period
. It is convenient and standard practice to choose , which means that the corresponding wavenumbers in long-range computations will be integers (in general, with ). But in principle one can choose any positive value of . Set the number of Fourier modes
in each direction. Choosing also sets the maximum resolved wavenumber as well as the physical grid spacing associated to the long-range fields. The value of should be tuned to have a good balance between the time spent on short-range and long-range computations. Now set the non-dimensional parameter
to the desired accuracy. A reasonable choice is , which roughly gives relative accuracy. From
and , one obtains the remaining parameters and .
For simplicity here we have assumed that the domain size
Size of FFTs
It may also be a good idea to "round" the chosen value of NonuniformFFTsBackend in VortexPasta uses an oversampling of
Another parameter to choose is the size of the quadrature rules for numerical integration. Using Gauss–Legendre quadratures, integrals seem to converge quite fast using a small number of quadrature nodes per filament segment. Typically, using 3 nodes seems to be enough when using quintic splines to describe filaments.
Ewald method for the streamfunction
As in the case of the velocity, the streamfunction integral decays slowly with the distance
In the case of the streamfunction, the idea of Ewald summation is to split the singular and slowly-decaying Green's function
where
This leads to modified integrals for the short-range and long-range streamfunction. For example, the short-range integral is given by:
Taking the curl of those integrals eventually leads to the