Locality Module
Overview
Use an Axis-Aligned Bounding Box (AABB) tree [HAN+16] to find neighbors. |
|
Filter an Existing |
|
Filter a |
|
Supports efficiently finding all points in a set within a certain distance from a given point. |
|
Class representing bonds between two sets of points. |
|
Class representing a set of points along with the ability to query for neighbors of these points. |
|
Class encapsulating the output of queries of NeighborQuery objects. |
|
Replicate periodic images of points inside a box. |
|
Computes Voronoi diagrams using voro++. |
Details
The freud.locality
module contains data structures to efficiently
locate points based on their proximity to other points.
- class freud.locality.AABBQuery
Bases:
NeighborQuery
Use an Axis-Aligned Bounding Box (AABB) tree [HAN+16] to find neighbors.
Also available as
freud.AABBQuery
.- Parameters
box (
freud.box.Box
) – Simulation box.points ((\(N\), 3)
numpy.ndarray
) – The points to use to build the tree.
- box
The box object used by this data structure.
- Type
- classmethod from_system(type cls, system, dimensions=None)
Create a
NeighborQuery
from any system-like object.The standard concept of a system in freud is any object that provides a way to access a box-like object (anything that can be coerced to a box by
freud.box.Box.from_box()
) and an array-like object (according to NumPy’s definition) of particle positions that turns into a \(N\times 3\) array.Supported types for
system
include:A sequence of
(box, points)
wherebox
is aBox
andpoints
is anumpy.ndarray
.Objects with attributes
box
andpoints
.hoomd.data
snapshot
- Parameters
system (system-like object) – Any object that can be converted to a
NeighborQuery
.dimensions (int) – Whether the object is 2 or 3 dimensional. It may be inferred if not provided, but in some cases inference is not possible, in which case it will default to 3 (Default value = None).
- Returns
The same
NeighborQuery
object if one is given, or an instance ofNeighborQuery
built from an inferredbox
andpoints
.- Return type
- plot(self, ax=None, title=None, *args, **kwargs)
Plot system box and points.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. The axis projection (2D or 3D) must match the dimensionality of the system (Default value =None
).title (str) – Title of the plot (Default value =
None
).*args – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.**kwargs – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.
- Returns
Axis and point data for the plot.
- Return type
tuple (
matplotlib.axes.Axes
,matplotlib.collections.PathCollection
)
- points
The array of points in this data structure.
- Type
np.ndarray
- query(self, query_points, query_args)
Query for nearest neighbors of the provided point.
- Parameters
query_points ((\(N\), 3)
numpy.ndarray
) – Points to query for.query_args (dict) – Query arguments determining how to find neighbors. For information on valid query argument, see the Query API.
- Returns
Results object containing the output of this query.
- Return type
- class freud.locality.Filter
Bases:
_PairCompute
Filter an Existing
NeighborList
.This class serves as the base class for all NeighborList filtering methods in freud. Filtering a
NeighborList
requires first computing the unfilteredNeighborList
from a system and a set of query arguments. Then, based on the arrangement of particles, their shapes, and other criteria determined by the derived class, some of the neighbors are removed from the unfilteredNeighborList
.The compute method of each
Filter
class will take a system object along with a neighbors dictionary specifying query arguments. Theneighbors
dictionary along with the system object will be used to build the unfiltered neighborlist, which will then be filtered according to the filter class. After the calculation, the filtered neighborlist will be available as the propertyfiltered_nlist
.Warning
This class is abstract and should not be instantiated directly.
- compute(self, system, neighbors=None, query_points=None)
Filter a
Neighborlist
.- Parameters
system – Any object that is a valid argument to
freud.locality.NeighborQuery.from_system
.neighbors (
freud.locality.NeighborList
or dict, optional) – Either aNeighborList
of neighbor pairs to use for the unfiltered neighbor list, or a dictionary of query arguments. IfNone
, an unfiltered neighborlist will be created such that all pairs of particles are neighbors viaNeighborList.all_pairs()
(Default value =None
).query_points ((\(N_{query\_points}\), 3)
np.ndarray
, optional) – Query points used to calculate the unfiltered neighborlist. Uses the system’s points ifNone
(Default value =None
).
- default_query_args
No default query arguments.
- property filtered_nlist
The filtered neighbor list.
- Type
- property unfiltered_nlist
The unfiltered neighbor list.
- Type
- class freud.locality.FilterSANN
Bases:
Filter
Filter a
NeighborList
via the SANN method.The Solid Angle Nearest Neighbor (SANN) method [vMFVF12] is a parameter-free algorithm for the identification of nearest neighbors. The SANN method attributes to each possible neighbor a solid angle and determines the cutoff radius by the requirement that the sum of the solid angles is 4π.
For performance considerations, SANN is implemented as a way of filtering a pre-existing set of neighbors due to the high performance cost of sorting all \(N^2\) particle pairs by distance. For a more in-depth explanation of the neighborlist filter concept in freud, see
Filter
.Warning
Due to the above design decision, it is possible that the unfiltered neighborlist will not contain enough neighbors to completely fill the neighbor shell of some particles in the system. The
allow_incomplete_shell
argument toFilterSANN
’s constructor controls whether a warning or exception is raised in these cases.Note
The
filtered_nlist
computed by this class will be sorted by distance.- Parameters
allow_incomplete_shell (bool) – Whether particles with incomplete neighbor shells are allowed in the filtered neighborlist. If True, a warning will be raised if there are particles with incomplete neighbors shells in the filtered neighborlist. If False, an exception will be raised in the same case (Default value =
False
).
- compute(self, system, neighbors=None, query_points=None)
Filter a
Neighborlist
.- Parameters
system – Any object that is a valid argument to
freud.locality.NeighborQuery.from_system
.neighbors (
freud.locality.NeighborList
or dict, optional) – Either aNeighborList
of neighbor pairs to use for the unfiltered neighbor list, or a dictionary of query arguments. IfNone
, an unfiltered neighborlist will be created such that all pairs of particles are neighbors viaNeighborList.all_pairs()
(Default value =None
).query_points ((\(N_{query\_points}\), 3)
np.ndarray
, optional) – Query points used to calculate the unfiltered neighborlist. Uses the system’s points ifNone
(Default value =None
).
- default_query_args
No default query arguments.
- property filtered_nlist
The filtered neighbor list.
- Type
- property unfiltered_nlist
The unfiltered neighbor list.
- Type
- class freud.locality.LinkCell
Bases:
NeighborQuery
Supports efficiently finding all points in a set within a certain distance from a given point.
Also available as
freud.LinkCell
.- Parameters
box (
freud.box.Box
) – Simulation box.points ((\(N\), 3)
numpy.ndarray
) – The points to bin into the cell list.cell_width (float, optional) – Width of cells. If not provided,
LinkCell
will estimate a cell width based on the number of points and the box size, assuming a constant density of points in the box.
- box
The box object used by this data structure.
- Type
- classmethod from_system(type cls, system, dimensions=None)
Create a
NeighborQuery
from any system-like object.The standard concept of a system in freud is any object that provides a way to access a box-like object (anything that can be coerced to a box by
freud.box.Box.from_box()
) and an array-like object (according to NumPy’s definition) of particle positions that turns into a \(N\times 3\) array.Supported types for
system
include:A sequence of
(box, points)
wherebox
is aBox
andpoints
is anumpy.ndarray
.Objects with attributes
box
andpoints
.hoomd.data
snapshot
- Parameters
system (system-like object) – Any object that can be converted to a
NeighborQuery
.dimensions (int) – Whether the object is 2 or 3 dimensional. It may be inferred if not provided, but in some cases inference is not possible, in which case it will default to 3 (Default value = None).
- Returns
The same
NeighborQuery
object if one is given, or an instance ofNeighborQuery
built from an inferredbox
andpoints
.- Return type
- plot(self, ax=None, title=None, *args, **kwargs)
Plot system box and points.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. The axis projection (2D or 3D) must match the dimensionality of the system (Default value =None
).title (str) – Title of the plot (Default value =
None
).*args – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.**kwargs – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.
- Returns
Axis and point data for the plot.
- Return type
tuple (
matplotlib.axes.Axes
,matplotlib.collections.PathCollection
)
- points
The array of points in this data structure.
- Type
np.ndarray
- query(self, query_points, query_args)
Query for nearest neighbors of the provided point.
- Parameters
query_points ((\(N\), 3)
numpy.ndarray
) – Points to query for.query_args (dict) –
Query arguments determining how to find neighbors. For information on valid query argument, see the Query API.
- Returns
Results object containing the output of this query.
- Return type
- class freud.locality.NeighborList
Bases:
object
Class representing bonds between two sets of points.
Compute classes contain a set of bonds between two sets of position arrays (“query points” and “points”) and hold a list of index pairs \(\left(i, j\right)\) where \(i < N_{query\_points}, j < N_{points}\) corresponding to neighbor pairs between the two sets.
For efficiency, all bonds must be sorted by the query point index, from least to greatest. Bonds have an query point index \(i\) and a point index \(j\). The first bond index corresponding to a given query point can be found in \(\log(N_{bonds})\) time using
find_first_index()
, because bonds are ordered by the query point index.Note
Typically, there is no need to instantiate this class directly. In most cases, users should manipulate
freud.locality.NeighborList
objects received from a neighbor search algorithm, such asfreud.locality.LinkCell
,freud.locality.AABBQuery
, orfreud.locality.Voronoi
.Also available as
freud.NeighborList
.Example:
# Assume we have position as Nx3 array aq = freud.locality.AABBQuery(box, positions) nlist = aq.query(positions, {'r_max': 3}).toNeighborList() # Get all vectors from central particles to their neighbors rijs = (positions[nlist.point_indices] - positions[nlist.query_point_indices]) rijs = box.wrap(rijs)
The NeighborList can be indexed to access bond particle indices. Example:
for i, j in nlist[:]: print(i, j)
- classmethod all_pairs(type cls, system, query_points=None, exclude_ii=True)
Create a NeighborList where all pairs of points are neighbors.
More explicitly, this method returns a NeighborList in which all pairs of points \(i\), \(j\) are neighbors. Pairs such that \(i = j\) can also be excluded using the
exclude_ii
option. The weight of all neighbors pairs in the returned list will be 1.- Parameters
system – Any object that is valid argument to
freud.locality.NeighborQuery.from_system
.query_points ((\(N_{query\_points}\), 3)
np.ndarray
, optional) – Query points used to create neighbor pairs. Uses the system’s points ifNone
(Default value =None
).exclude_ii (bool) – Whether to exclude pairs of particles with the same point index in the output neighborlist (Default value =
True
).
- copy(self, other=None)
Create a copy. If other is given, copy its contents into this object. Otherwise, return a copy of this object.
- Parameters
other (
freud.locality.NeighborList
, optional) – A NeighborList to copy into this object (Default value =None
).
- distances
The distances for each bond.
- Type
(\(N_{bonds}\))
np.ndarray
- filter(self, filt)
Removes bonds that satisfy a boolean criterion.
- Parameters
filt (
np.ndarray
) – Boolean-like array of bonds to keep (True means the bond will not be removed).
Note
This method modifies this object in-place.
Example:
# Keep only the bonds between particles of type A and type B nlist.filter(types[nlist.query_point_indices] != types[nlist.point_indices])
- filter_r(self, float r_max, float r_min=0)
Removes bonds that are outside of a given radius range.
- find_first_index(self, unsigned int i)
Returns the lowest bond index corresponding to a query particle with an index \(\geq i\).
- Parameters
i (unsigned int) – The particle index.
- classmethod from_arrays(type cls, num_query_points, num_points, query_point_indices, point_indices, distances, weights=None)
Create a NeighborList from a set of bond information arrays.
Example:
import freud import numpy as np box = freud.box.Box(2, 3, 4, 0, 0, 0) query_points = np.array([[0, 0, 0], [0, 0, 1]]) points = np.array([[0, 0, -1], [0.5, -1, 0]]) query_point_indices = np.array([0, 0, 1]) point_indices = np.array([0, 1, 1]) distances = box.compute_distances( query_points[query_point_indices], points[point_indices]) num_query_points = len(query_points) num_points = len(points) nlist = freud.locality.NeighborList.from_arrays( num_query_points, num_points, query_point_indices, point_indices, distances)
- Parameters
num_query_points (int) – Number of query points (corresponding to
query_point_indices
).num_points (int) – Number of points (corresponding to
point_indices
).query_point_indices (
np.ndarray
) – Array of integers corresponding to indices in the set of query points.point_indices (
np.ndarray
) – Array of integers corresponding to indices in the set of points.distances (
np.ndarray
) – Array of distances between corresponding query points and points.weights (
np.ndarray
, optional) – Array of per-bond weights (ifNone
is given, use a value of 1 for each weight) (Default value =None
).
- neighbor_counts
A neighbor count array indicating the number of neighbors for each query point.
- Type
(\(N_{query\_points}\))
np.ndarray
- num_points
The number of points.
All point indices are less than this value.
- Type
unsigned int
- num_query_points
The number of query points.
All query point indices are less than this value.
- Type
unsigned int
- point_indices
The point indices for each bond. This array is read-only to prevent breakage of
find_first_index()
. Equivalent to indexing with[:, 1]
.- Type
(\(N_{bonds}\))
np.ndarray
- query_point_indices
The query point indices for each bond. This array is read-only to prevent breakage of
find_first_index()
. Equivalent to indexing with[:, 0]
.- Type
(\(N_{bonds}\))
np.ndarray
- segments
A segment array indicating the first bond index for each query point.
- Type
(\(N_{query\_points}\))
np.ndarray
- sort(self, bool by_distance=False)
Sort the entries in the neighborlist.
- Parameters
by_distance (bool) – If
True
, this method sorts the neighborlist entries byquery_point_index
, thendistance
, thenpoint_index
. IfFalse
, this method sorts the NeighborList entries byquery_point_index
, thenpoint_index
, thendistance
(Default value =False
).
- weights
The weights for each bond. By default, bonds have a weight of 1.
- Type
(\(N_{bonds}\))
np.ndarray
- class freud.locality.NeighborQuery
Bases:
object
Class representing a set of points along with the ability to query for neighbors of these points.
Warning
This class should not be instantiated directly. The subclasses
AABBQuery
andLinkCell
provide the intended interfaces.The
NeighborQuery
class represents the abstract interface for neighbor finding. The class contains a set of points and a simulation box, the latter of which is used to define the system and the periodic boundary conditions required for finding neighbors of these points. The primary mode of interacting with theNeighborQuery
is through thequery()
andqueryBall()
functions, which enable finding either the nearest neighbors of a point or all points within a distance cutoff, respectively. Subclasses of NeighborQuery implement these methods based on the nature of the underlying data structure.- Parameters
box (
freud.box.Box
) – Simulation box.points ((\(N\), 3)
numpy.ndarray
) – Point coordinates to build the structure.
- box
The box object used by this data structure.
- Type
- classmethod from_system(type cls, system, dimensions=None)
Create a
NeighborQuery
from any system-like object.The standard concept of a system in freud is any object that provides a way to access a box-like object (anything that can be coerced to a box by
freud.box.Box.from_box()
) and an array-like object (according to NumPy’s definition) of particle positions that turns into a \(N\times 3\) array.Supported types for
system
include:A sequence of
(box, points)
wherebox
is aBox
andpoints
is anumpy.ndarray
.Objects with attributes
box
andpoints
.hoomd.data
snapshot
- Parameters
system (system-like object) – Any object that can be converted to a
NeighborQuery
.dimensions (int) – Whether the object is 2 or 3 dimensional. It may be inferred if not provided, but in some cases inference is not possible, in which case it will default to 3 (Default value = None).
- Returns
The same
NeighborQuery
object if one is given, or an instance ofNeighborQuery
built from an inferredbox
andpoints
.- Return type
- plot(self, ax=None, title=None, *args, **kwargs)
Plot system box and points.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. The axis projection (2D or 3D) must match the dimensionality of the system (Default value =None
).title (str) – Title of the plot (Default value =
None
).*args – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.**kwargs – Passed on to
mpl_toolkits.mplot3d.Axes3D.plot()
ormatplotlib.axes.Axes.plot()
.
- Returns
Axis and point data for the plot.
- Return type
tuple (
matplotlib.axes.Axes
,matplotlib.collections.PathCollection
)
- points
The array of points in this data structure.
- Type
np.ndarray
- query(self, query_points, query_args)
Query for nearest neighbors of the provided point.
- Parameters
query_points ((\(N\), 3)
numpy.ndarray
) – Points to query for.query_args (dict) –
Query arguments determining how to find neighbors. For information on valid query argument, see the Query API.
- Returns
Results object containing the output of this query.
- Return type
- class freud.locality.NeighborQueryResult
Bases:
object
Class encapsulating the output of queries of NeighborQuery objects.
Warning
This class should not be instantiated directly, it is the return value of the
query()
method ofNeighborQuery
. The class provides a convenient interface for iterating over query results, and can be transparently converted into a list or aNeighborList
object.The
NeighborQueryResult
makes it easy to work with the results of queries and convert them to various natural objects. Additionally, the result is a generator, making it easy for users to lazily iterate over the object.- toNeighborList(self, sort_by_distance=False)
Convert query result to a freud
NeighborList
.- Parameters
sort_by_distance (bool) – If
True
, sort neighboring bonds by distance. IfFalse
, sort neighboring bonds by point index (Default value =False
).- Returns
A
NeighborList
containing all neighbor pairs found by the query generating this result object.- Return type
- class freud.locality.PeriodicBuffer
Bases:
_Compute
Replicate periodic images of points inside a box.
- property buffer_box
The buffer box, expanded to hold the replicated points.
- Type
- property buffer_ids
The buffer point ids.
- Type
\(\left(N_{buffer}\right)\)
numpy.ndarray
- property buffer_points
The buffer point positions.
- Type
\(\left(N_{buffer}, 3\right)\)
numpy.ndarray
- compute(self, system, buffer, bool images=False, include_input_points=False)
Compute the periodic buffer.
- Parameters
system – Any object that is a valid argument to
freud.locality.NeighborQuery.from_system
.buffer (float or list of 3 floats) – Buffer distance for replication outside the box.
images (bool, optional) – If
False
,buffer
is a distance. IfTrue
,buffer
is a number of images to replicate in each dimension. Note that one image adds half of a box length to each side, meaning that one image doubles the box side lengths, two images triples the box side lengths, and so on. (Default value =False
).include_input_points (bool, optional) – Whether the original points provided by
system
are included in the buffer, (Default value =False
).
- class freud.locality.Voronoi
Bases:
_Compute
Computes Voronoi diagrams using voro++.
Voronoi diagrams (Wikipedia) are composed of convex polytopes (polyhedra in 3D, polygons in 2D) called cells, corresponding to each input point. The cells bound a region of Euclidean space for which all contained points are closer to a corresponding input point than any other input point. A ridge is defined as a boundary between cells, which contains points equally close to two or more input points.
The voro++ library [Ryc09] is used for fast computations of the Voronoi diagram.
- compute(self, system)
Compute Voronoi diagram.
- Parameters
system – Any object that is a valid argument to
freud.locality.NeighborQuery.from_system
.
- property nlist
Returns the computed
NeighborList
.The
NeighborList
computed by this class is weighted. In 2D systems, the bond weight is the length of the ridge (boundary line) between the neighboring points’ Voronoi cells. In 3D systems, the bond weight is the area of the ridge (boundary polygon) between the neighboring points’ Voronoi cells. The weights are not normalized, and the weights for each query point sum to the surface area (perimeter in 2D) of the polytope.It is possible for pairs of points to appear multiple times in the neighbor list. For example, in a small unit cell, points may neighbor one another on multiple sides because of periodic boundary conditions.
- Returns
Neighbor list.
- Return type
- plot(self, ax=None, color_by_sides=True, cmap=None)
Plot Voronoi diagram.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)color_by_sides (bool) – If
True
, color cells by the number of sides. IfFalse
, random colors are used for each cell. (Default value =True
)cmap (str) – Colormap name to use (Default value =
None
).
- Returns
Axis with the plot.
- Return type
- property polytopes
A list of
numpy.ndarray
defining Voronoi polytope vertices for each cell.- Type
list[
numpy.ndarray
]
- property volumes
Returns an array of Voronoi cell volumes (areas in 2D).
- Type
\(\left(N_{points} \right)\)
numpy.ndarray