Locality Module

Overview

freud.locality.NeighborQuery

Class representing a set of points along with the ability to query for neighbors of these points.

freud.locality.NeighborQueryResult

Class encapsulating the output of queries of NeighborQuery objects.

freud.locality.NeighborList

Class representing a certain number of “bonds” between particles.

freud.locality.IteratorLinkCell

Iterates over the particles in a cell.

freud.locality.LinkCell

Supports efficiently finding all points in a set within a certain distance from a given point.

freud.locality.NearestNeighbors

Supports efficiently finding the \(N\) nearest neighbors of each point in a set for some fixed integer \(N\).

freud.locality.AABBQuery

Use an AABB tree to find neighbors.

Details

The freud.locality module contains data structures to efficiently locate points based on their proximity to other points.

Neighbor List

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>

New in version 0.6.4.

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 as freud.locality.LinkCell, freud.locality.NearestNeighbors, or freud.voronoi.Voronoi.

Variables
  • index_i (np.ndarray) – The reference point indices from the last set of points this object was evaluated with. This array is read-only to prevent breakage of find_first_index().

  • index_j (np.ndarray) – The reference point indices from the last set of points this object was evaluated with. This array is read-only to prevent breakage of find_first_index().

  • weights ((\(N_{bonds}\)) np.ndarray) – The per-bond weights from the last set of points this object was evaluated with.

  • segments ((\(N_{ref\_points}\)) np.ndarray) – 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 this object was evaluated with.

  • neighbor_counts ((\(N_{ref\_points}\)) np.ndarray) – 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 this object was evaluated with.

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)

The NeighborList can be indexed to access bond particle indices. Example:

for i, j in lc.nlist[:]:
    print(i, j)
copy

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).

filter

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.index_i] != types[nlist.index_j])
filter_r

Removes bonds that are outside of a given radius range.

Parameters
  • box (freud.box.Box) – Simulation box.

  • ref_points ((\(N_{particles}\), 3) numpy.ndarray) – Reference points to use for filtering.

  • points ((\(N_{particles}\), 3) numpy.ndarray) – Target points to use for filtering.

  • rmax (float) – Maximum bond distance in the resulting neighbor list.

  • rmin (float, optional) – Minimum bond distance in the resulting neighbor list (Default value = 0).

find_first_index

Returns the lowest bond index corresponding to a reference particle with an index \(\geq i\).

Parameters

i (unsigned int) – The particle index.

classmethod from_arrays(type cls, Nref, Ntarget, index_i, index_j, weights=None)

Create a NeighborList from a set of bond information arrays.

Parameters
  • Nref (int) – Number of reference points (corresponding to index_i).

  • Ntarget (int) – Number of target points (corresponding to index_j).

  • index_i (np.ndarray) – Array of integers corresponding to indices in the set of reference points.

  • index_j (np.ndarray) – Array of integers corresponding to indices in the set of target points.

  • weights (np.ndarray, optional) – Array of per-bond weights (if None is given, use a value of 1 for each weight) (Default value = None).

Neighbor Querying

class freud.locality.NeighborQuery

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 and LinkCell 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 the NeighborQuery is through the query() and queryBall() 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.

Module author: Vyas Ramasubramani <vramasub@umich.edu>

New in version 1.1.0.

Parameters
Variables
  • box (freud.box.Box) – The box object used by this data structure.

  • points (np.ndarray) – The array of points in this data structure.

query

Query for nearest neighbors of the provided point.

Parameters
  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • k (int) – The number of nearest neighbors to find.

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices to those in self.points should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

queryBall

Query for all points within a distance r of the provided point(s).

Parameters
  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • r (float) – The distance within which to find neighbors

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices to those in self.points should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

class freud.locality.NeighborQueryResult

Class encapsulating the output of queries of NeighborQuery objects.

Warning

This class should not be instantiated directly, it is the return value of all query* functions of NeighborQuery. The class provides a convenient interface for iterating over query results, and can be transparently converted into a list or a NeighborList 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.

Module author: Vyas Ramasubramani <vramasub@umich.edu>

New in version 1.1.0.

toNList

Convert query result to a freud NeighborList.

Returns

A freud NeighborList containing all neighbor pairs found by the query generating this result object.

Return type

NeighborList

class freud.locality.AABBQuery

Use an AABB tree to find neighbors.

Module author: Bradley Dice <bdice@bradleydice.com>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

New in version 1.1.0.

Variables
  • box (freud.locality.Box) – The simulation box.

  • points (np.ndarray) – The points associated with this class.

query

Query for nearest neighbors of the provided point.

This method has a slightly different signature from the parent method to support querying based on a specified guessed rcut and scaling.

Parameters
  • box (freud.box.Box) – Simulation box.

  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • k (int) – The number of nearest neighbors to find.

  • r (float) – The initial guess of a distance to search to find N neighbors.

  • scale (float) – Multiplier by which to increase r if not enough neighbors are found.

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

queryBall

Query for all points within a distance r of the provided point(s).

Parameters
  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • r (float) – The distance within which to find neighbors

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices to those in self.points should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

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>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • cell_width (float) – Maximum distance to find particles within.

  • points (np.ndarray, optional) – The points associated with this class, if used as a NeighborQuery object, i.e. built on one set of points that can then be queried against. (Defaults to None).

Variables
  • box (freud.box.Box) – Simulation box.

  • num_cells (unsigned int) – The number of cells in the box.

  • nlist (freud.locality.NeighborList) – The neighbor list stored by this object, generated by compute().

Note

2D: freud.locality.LinkCell properly handles 2D boxes. The points must be passed in as [x, y, 0]. Failing to set z=0 will lead to undefined behavior.

Example:

# Assume positions are an Nx3 array
lc = LinkCell(box, 1.5)
lc.compute(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)
compute

Update the data structure for the given set of points and compute a NeighborList.

Parameters
  • box (freud.box.Box) – Simulation box.

  • ref_points ((\(N_{particles}\), 3) numpy.ndarray) – Reference point coordinates.

  • points ((\(N_{particles}\), 3) numpy.ndarray, optional) – Point coordinates (Default value = None).

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

getCell

Returns the index of the cell containing the given point.

Parameters

point (\(\left(3\right)\) numpy.ndarray) – Point coordinates \(\left(x,y,z\right)\).

Returns

Cell index.

Return type

unsigned int

getCellNeighbors

Returns the neighboring cell indices of the given cell.

Parameters

cell (unsigned int) – Cell index.

Returns

Array of cell neighbors.

Return type

\(\left(N_{neighbors}\right)\) numpy.ndarray

itercell

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

query

Query for nearest neighbors of the provided point.

Parameters
  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • k (int) – The number of nearest neighbors to find.

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices to those in self.points should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

queryBall

Query for all points within a distance r of the provided point(s).

Parameters
  • points ((\(N\), 3) numpy.ndarray) – Points to query for.

  • r (float) – The distance within which to find neighbors

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices to those in self.points should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

Returns

Results object containing the output of this query.

Return type

NeighborQueryResult

class freud.locality.IteratorLinkCell

Iterates over the particles in a cell.

Module author: Joshua Anderson <joaander@umich.edu>

Example:

# Grab particles in cell 0
for j in linkcell.itercell(0):
    print(positions[j])
next

Implements iterator interface

Nearest Neighbors

class freud.locality.NearestNeighbors(rmax, n_neigh, scale=1.1, strict_cut=False)

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 (default): rmax will be expanded to find the requested number of neighbors. If rmax increases to the point that a cell list cannot be constructed, a warning will be raised and the neighbors already 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 if the 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, default is False.

Variables
  • UINTMAX (unsigned int) – Value of C++ UINTMAX used to pad the arrays.

  • box (freud.box.Box) – Simulation box.

  • num_neighbors (unsigned int) – The number of neighbors this object will find.

  • n_ref (unsigned int) – The number of particles this object found neighbors of.

  • r_max (float) – Current nearest neighbors search radius guess.

  • wrapped_vectors (\(\left(N_{particles}\right)\) numpy.ndarray) – The wrapped vectors padded with -1 for empty neighbors.

  • r_sq_list (\(\left(N_{particles}, N_{neighbors}\right)\) numpy.ndarray) – The Rsq values list.

  • nlist (freud.locality.NeighborList) – The neighbor list stored by this object, generated by compute().

Example:

nn = NearestNeighbors(2, 6)
nn.compute(box, positions, positions)
hexatic = order.HexOrderParameter(2)
hexatic.compute(box, positions, nlist=nn.nlist)
compute

Update the data structure for the given set of points.

Parameters
  • box (freud.box.Box) – Simulation box.

  • ref_points ((\(N_{particles}\), 3) numpy.ndarray) – Reference point coordinates.

  • points ((\(N_{particles}\), 3) numpy.ndarray, optional) – Point coordinates. Defaults to ref_points if not provided or None.

  • exclude_ii (bool, optional) – Set this to True if pairs of points with identical indices should be excluded. If this is None, it will be treated as True if points is None or the same object as ref_points (Defaults to None).

getNeighborList

Return the entire neighbor list.

Returns

Indices of up to \(N_{neighbors}\) points that are neighbors of the \(N_{particles}\) reference points, padded with UINTMAX if fewer neighbors than requested were found.

Return type

\(\left(N_{particles}, N_{neighbors}\right)\) numpy.ndarray

getNeighbors

Return the \(N\) nearest neighbors of the reference point with index \(i\).

Parameters

i (unsigned int) – Index of the reference point whose neighbors will be returned.

Returns

Indices of points that are neighbors of reference point \(i\), padded with UINTMAX if fewer neighbors than requested were found.

Return type

\(\left(N_{neighbors}\right)\) numpy.ndarray

getRsq

Return the squared distances to 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 to the \(N\) nearest neighbors.

Return type

\(\left(N_{particles}\right)\) numpy.ndarray