CellLists

Types

VortexPasta.CellLists.PeriodicCellListType
PeriodicCellList{N, T}
PeriodicCellList(
    ::Type{T}, rs_cut::NTuple{N, Real}, periods::NTuple{N, Real},
    [nsubdiv = static(1)];
    [to_coordinate::Function = identity],
)

Construct a cell list for dealing with pair interactions.

Above, N is the number of spatial dimensions, and T is the type of each element. In the simplest cases, T can simply describe a coordinate in $N$-dimensional space (e.g. T = SVector{N, Float64}). One can also deal with more complicated elements which include more information. As an example, see further below for how to deal with filament segments in 3D space.

The cutoff radii rs_cut (which can be different in each direction) don't need to exactly divide the domain period L into equal pieces, but it's recommended that it does so for performance reasons.

Optionally, one can choose to subdivide each cell (of size ≈ rcut) onto nsubdiv subcells. This can significantly improve performance, since it allows to discard some spurious pair interactions (i.e. beyond the chosen cutoff radius) as described here. In practice, a value of 2 or 3 can significantly improve performance compared to no subdivision (1). For convenience, it can be passed as nsubdiv = Val(M) or as nsubdiv = static(M).

Infinite non-periodic domains (in the sense of period = Infinity()) are not supported.

Dealing with filament segments

One of the possible uses of PeriodicCellList is to classify filament segments (which are typically shorter than the cutoff radius) according to their spatial location. In that case, T is not a simple coordinate, but may contain more information including things like (1) the filament the segment belongs to, and (2) the location of the segment within the filament. As there is no unique way of associating a coordinate to a segment, one should pass the to_coordinate argument which "converts" the segment to a coordinate in space. For instance, the passed to_coordinate function may return the midpoint of the segment, which will be used to determine the cell associated to the segment.

The applied cutoff radius $r_{\text{cut}}$ should be much larger than the maximum segment length $ℓ$, or should at least account for $ℓ$. Basically, if one wants an actual cut-off radius $r₀$, then the applied cutoff radius passed to the constructor should be $r_{\text{cut}} = r₀ + ℓ$. Otherwise, a small amount of interactions within $[r₀ - ℓ, r₀]$ may be missed.

source

Functions

VortexPasta.CellLists.add_element!Function
add_element!(cl::PeriodicCellList{N, T}, el::T, [x⃗])

Add element to the cell list.

Determines the cell associated to the element and then appends the element to that cell.

Optionally, one may pass the coordinate location $x⃗$ associated to the element. Otherwise, it will be obtained from the element according to

x⃗ = to_coordinate(el)

where to_coordinate corresponds to the keyword argument of PeriodicCellList.

source
VortexPasta.CellLists.finalise_cells!Function
CellLists.finalise_cells!(cl::PeriodicCellList)

"Finalise" cells before iterating over its elements.

This function performs operations needed to iterate over elements of the cell lists, such as filling ghost cells. It must be called once after all elements have been added to the cell list using add_element! and before one starts iterating using nearby_elements.

source
Base.empty!Function
Base.empty!(tsf::TimeSeriesFile)

Reset TimeSeriesFile, removing all entries.

source
Base.empty!(cl::PeriodicCellList) -> cl

Remove all elements from the cell list.

source
VortexPasta.CellLists.nearby_elementsFunction
nearby_elements(cl::PeriodicCellList{N}, x⃗)

Return an iterator over the elements that are sufficiently close to the point x⃗.

The iterator returns the elements which are likely to be within the cutoff radius of the point x⃗. More precisely, it returns elements in the same cell as x⃗ as well as in neighbouring cells.

Here x⃗ should be a coordinate, usually represented by an SVector{N} or an NTuple{N}.

source