Locality Module¶
The locality module contains data structures to efficiently locate points based on their proximity to other points.
NeighborList¶
-
class
freud.locality.
NeighborList
¶ Class representing a certain number of “bonds” between particles. Computation methods will iterate over these bonds when searching for neighboring particles.
NeighborList objects are constructed for two sets of position arrays A (alternatively reference points; of length \(n_A\)) and B (alternatively target points; of length \(n_B\)) and hold a set of \(\left(i, j\right): i < n_A, j < n_B\) index pairs corresponding to near-neighbor points in A and B, respectively.
For efficiency, all bonds for a particular reference particle i are contiguous and bonds are stored in order based on reference particle index i. The first bond index corresponding to a given particle can be found in \(log(n_{bonds})\) time using
find_first_index()
.Module author: Matthew Spellings <mspells@umich.edu>
Note
Typically, in python you will only manipulate a
freud.locality.NeighborList
object that you receive from a neighbor search algorithm, such asfreud.locality.LinkCell
andfreud.locality.NearestNeighbors
.Example:
# assume we have position as Nx3 array lc = LinkCell(box, 1.5).compute(box, positions) nlist = lc.nlist # get all vectors from central particles to their neighbors rijs = positions[nlist.index_j] - positions[nlist.index_i] box.wrap(rijs)
-
copy
(self, other=None)¶ Create a copy. If other is given, copy its contents into ourself; otherwise, return a copy of ourself.
-
filter
(self, filt)¶ Removes bonds that satisfy a boolean criterion.
Parameters: filt – Boolean-like array of bonds to keep (True => bond stays) 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.index_i] != types[nlist.index_j])
-
filter_r
(self, box, ref_points, points, float rmax, float rmin=0)¶ Removes bonds that are outside of a given radius range.
Parameters: - ref_points – reference points to use for filtering
- points – target points to use for filtering
- rmax – maximum bond distance in the resulting neighbor list
- rmin – minimum bond distance in the resulting neighbor list
Note
This method modifies this object in-place
-
find_first_index
(self, unsigned int i)¶ Returns the lowest bond index corresponding to a reference particle with index >=i
-
from_arrays
(type cls, Nref, Ntarget, index_i, index_j, weights=None)¶ Create a NeighborList from a set of bond information arrays.
Parameters: - Nref (unsigned int) – Number of reference points (corresponding to index_i)
- Ntarget (unsigned int) – Number of target points (corresponding to index_j)
- index_i (Array-like of unsigned ints, length num_bonds) – Array of integers corresponding to indices in the set of reference points
- index_j (Array-like of unsigned ints, length num_bonds) – Array of integers corresponding to indices in the set of target points
- weights (Array-like of floats, length num_bonds) – Array of per-bond weights (if None is given, use a value of 1 for each weight)
-
index_i
¶ Returns the reference point indices from the last set of points we were evaluated with. This array is read-only to prevent breakage of find_first_index.
-
index_j
¶ Returns the target point indices from the last set of points we were evaluated with. This array is read-only to prevent breakage of find_first_index.
-
neighbor_counts
¶ Returns a neighbor count array, which is an array of length N_ref indicating the number of neighbors for each reference particle from the last set of points we were evaluated with.
-
segments
¶ Returns a segment array, which is an array of length N_ref indicating the first bond index for each reference particle from the last set of points we were evaluated with.
-
weights
¶ Returns the per-bond weights from the last set of points we were evaluated with
-
LinkCell¶
-
class
freud.locality.
LinkCell
(box, cell_width)¶ Supports efficiently finding all points in a set within a certain distance from a given point.
Module author: Joshua Anderson <joaander@umich.edu>
Parameters: - box (
freud.box.Box
) – simulation box - cell_width (float) – Maximum distance to find particles within
Note
freud.locality.LinkCell
supports 2D boxes; in this case, make sure to set the z coordinate of all points to 0.Example:
# assume we have position as Nx3 array lc = LinkCell(box, 1.5) lc.computeCellList(box, positions) for i in range(positions.shape[0]): # cell containing particle i cell = lc.getCell(positions[0]) # list of cell's neighboring cells cellNeighbors = lc.getCellNeighbors(cell) # iterate over neighboring cells (including our own) for neighborCell in cellNeighbors: # iterate over particles in each neighboring cell for neighbor in lc.itercell(neighborCell): pass # do something with neighbor index # using NeighborList API dens = density.LocalDensity(1.5, 1, 1) dens.compute(box, positions, nlist=lc.nlist)
-
box
¶ return – freud Box :rtype:
freud.box.Box
-
compute
(self, box, ref_points, points=None, exclude_ii=None)¶ Update the data structure for the given set of points and compute a NeighborList
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape= \(\left(N_{refpoints}, 3\right)\), dtype=numpy.float32
) – reference point coordinates - points (
numpy.ndarray
, shape= \(\left(N_{points}, 3\right)\), dtype=numpy.float32
) – point coordinates - exlude_ii – True if pairs of points with identical indices should be excluded; if None, is set to True if points is None or the same object as ref_points
- box (
-
computeCellList
(self, box, ref_points, points=None, exclude_ii=None)¶ Update the data structure for the given set of points and compute a NeighborList
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape= \(\left(N_{refpoints}, 3\right)\), dtype=numpy.float32
) – reference point coordinates - points (
numpy.ndarray
, shape= \(\left(N_{points}, 3\right)\), dtype=numpy.float32
) – point coordinates - exlude_ii – True if pairs of points with identical indices should be excluded; if None, is set to True if points is None or the same object as ref_points
- box (
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box
-
getCell
(self, point)¶ Returns the index of the cell containing the given point
Parameters: point ( numpy.ndarray
, shape= \(\left(3\right)\), dtype=numpy.float32
) – point coordinates \(\left(x,y,z\right)\)Returns: cell index Return type: unsigned int
-
getCellNeighbors
(self, cell)¶ Returns the neighboring cell indices of the given cell
Parameters: cell (unsigned int) – Cell index Returns: array of cell neighbors Return type: numpy.ndarray
, shape= \(\left(N_{neighbors}\right)\), dtype=numpy.uint32
-
getNumCells
(self)¶ Returns: the number of cells in this box Return type: unsigned int
-
itercell
(self, unsigned int cell)¶ Return an iterator over all particles in the given cell
Parameters: cell (unsigned int) – Cell index Returns: iterator to particle indices in specified cell Return type: iter
-
nlist
¶ Returns the neighbor list stored by this object, generated by compute.
-
num_cells
¶ return – the number of cells in this box :rtype: unsigned int
- box (
NearestNeighbors¶
-
class
freud.locality.
NearestNeighbors
(rmax, n_neigh)¶ Supports efficiently finding the N nearest neighbors of each point in a set for some fixed integer N.
- strict_cut = True: rmax will be strictly obeyed, and any particle which has fewer than N neighbors will have values of UINT_MAX assigned
- strict_cut = False: rmax will be expanded to find requested number of neighbors. If rmax increases to the point that a cell list cannot be constructed, a warning will be raised and neighbors found will be returned
Module author: Eric Harper <harperic@umich.edu>
Parameters: - rmax (float) – Initial guess of a distance to search within to find N neighbors
- n_neigh (unsigned int) – Number of neighbors to find for each point
- scale (float) – multiplier by which to automatically increase rmax value by if requested number of neighbors is not found. Only utilized if strict_cut is False. Scale must be greater than 1
- strict_cut (bool) – whether to use a strict rmax or allow for automatic expansion
Example:
nn = NearestNeighbors(2, 6) nn.compute(box, positions, positions) hexatic = order.HexOrderParameter(2) hexatic.compute(box, positions, nlist=nn.nlist)
-
UINTMAX
¶ return – value of C++ UINTMAX used to pad the arrays :rtype: unsigned int
-
box
¶ return – freud Box :rtype:
freud.box.Box
-
compute
(self, box, ref_points, points, exclude_ii=None)¶ Update the data structure for the given set of points
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – coordinated of reference points - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – coordinates of points - exlude_ii – True if pairs of points with identical indices should be excluded; if None, is set to True if points is None or the same object as ref_points
- box (
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box
-
getNRef
(self)¶ Returns: the number of particles this object found neighbors of Return type: unsigned int
-
getNeighborList
(self)¶ Return the entire neighbors list
Returns: Neighbor List Return type: numpy.ndarray
, shape= \(\left(N_{particles}, N_{neighbors}\right)\), dtype=numpy.uint32
-
getNeighbors
(self, unsigned int i)¶ Return the N nearest neighbors of the reference point with index i
Parameters: i (unsigned int) – index of the reference point to fetch the neighboring points of
-
getNumNeighbors
(self)¶ Returns: the number of neighbors this object will find Return type: unsigned int
-
getRMax
(self)¶ Return the current neighbor search distance guess :return: nearest neighbors search radius :rtype: float
-
getRsq
(self, unsigned int i)¶ Return the Rsq values for the N nearest neighbors of the reference point with index i
Parameters: i (unsigned int) – index of the reference point of which to fetch the neighboring point distances Returns: squared distances of the N nearest neighbors Returns: Neighbor List Return type: numpy.ndarray
, shape= \(\left(N_{particles}\right)\), dtype=numpy.float32
-
getRsqList
(self)¶ Return the entire Rsq values list
Returns: Rsq list Return type: numpy.ndarray
, shape= \(\left(N_{particles}, N_{neighbors}\right)\), dtype=numpy.float32
-
getUINTMAX
(self)¶ Returns: value of C++ UINTMAX used to pad the arrays Return type: unsigned int
-
getWrappedVectors
(self)¶ Return the wrapped vectors for computed neighbors. Array padded with -1 for empty neighbors
Returns: wrapped vectors Return type: numpy.ndarray
, shape= \(\left(N_{particles}\right)\), dtype=numpy.float32
-
n_ref
¶ return – the number of particles this object found neighbors of :rtype: unsigned int
-
nlist
¶ Returns the neighbor list stored by this object, generated by compute.
-
num_neighbors
¶ return – the number of neighbors this object will find :rtype: unsigned int
-
r_max
¶ Return the current neighbor search distance guess :return: nearest neighbors search radius :rtype: float
-
r_sq_list
¶ Return the entire Rsq values list
Returns: Rsq list Return type: numpy.ndarray
, shape= \(\left(N_{particles}, N_{neighbors}\right)\), dtype=numpy.float32
-
setCutMode
(self, strict_cut)¶ Set mode to handle rmax by Nearest Neighbors.
- strict_cut = True: rmax will be strictly obeyed, and any particle which has fewer than N neighbors will have values of UINT_MAX assigned
- strict_cut = False: rmax will be expanded to find requested number of neighbors. If rmax increases to the point that a cell list cannot be constructed, a warning will be raised and neighbors found will be returned
Parameters: strict_cut (bool) – whether to use a strict rmax or allow for automatic expansion
-
setRMax
(self, float rmax)¶ Update the neighbor search distance guess :param rmax: nearest neighbors search radius :type rmax: float
-
wrapped_vectors
¶ Return the wrapped vectors for computed neighbors. Array padded with -1 for empty neighbors
Returns: wrapped vectors Return type: numpy.ndarray
, shape= \(\left(N_{particles}\right)\), dtype=numpy.float32