sisl.geom.NeighborFinder

class sisl.geom.NeighborFinder(geometry: Geometry, R: float | ndarray | None = None, overlap: bool = False, bin_size: float | Tuple[float, float, float] = 2)[source]

Bases: object

Efficient linear scaling finding of neighbors.

Once this class is instantiated, a table is build. Then, the neighbor finder can be queried as many times as wished as long as the geometry doesn’t change.

Note that the radius (R) is used to build the table. Therefore, if one wants to look for neighbors using a different R, one needs to create another finder or call setup.

Parameters:
  • geometry (Geometry) – the geometry to find neighbors in

  • R (float or array-like of shape (geometry.na), optional) –

    The radius to consider two atoms neighbors. If it is a single float, the same radius is used for all atoms. If it is an array, it should contain the radius for each atom.

    If not provided, an array is constructed, where the radius for each atom is its maxR.

  • overlap (bool, optional) –

    If True, two atoms are considered neighbors if their spheres overlap. If False, two atoms are considered neighbors if the second atom is within the sphere of the first atom. Note that this implies that atom i might be atom j’s neighbor while the opposite is not true.

    If not provided, it will be True if R is an array and False if it is a single float.

  • bin_size (float or tuple of float, optional) – the factor for the radius to determine how large the bins are, optionally along each lattice vector. It can minimally be 2, meaning that the maximum radius to consider is twice the radius considered. For larger values, more atoms will be in each bin (and thus fewer bins). Hence, this value can be used to fine-tune the memory requirement by decreasing number of bins, at the cost of a bit more run-time searching bins.

Methods

assert_consistency()

Asserts that the data structure (self._list, self._heads, self._counts) is consistent.

find_close(xyz[, as_pairs])

Find all atoms that are close to some coordinates in space.

find_neighbors([atoms, as_pairs, ...])

Find neighbors as specified in the finder.

find_unique_pairs([self_interaction])

Find unique neighbor pairs within the geometry.

setup([geometry, R, overlap, bin_size])

Prepares everything for neighbor finding.

memory

Memory control of the finder

geometry

nbins

total_nbins

R

R: ndarray
__init__(geometry: Geometry, R: float | ndarray | None = None, overlap: bool = False, bin_size: float | Tuple[float, float, float] = 2)[source]
Parameters:
assert_consistency()[source]

Asserts that the data structure (self._list, self._heads, self._counts) is consistent.

It also stores that the shape is consistent with the stored geometry and the store total_nbins.

find_close(xyz: Sequence, as_pairs: bool = False)[source]

Find all atoms that are close to some coordinates in space.

This routine only executes the action of finding neighbors, the parameters of the search (geometry, R, overlap…) are defined when the finder is initialized or by calling setup.

Parameters:
  • xyz (array-like of shape ([npoints], 3)) – the coordinates for which neighbors are desired.

  • as_pairs (bool, optional) – If True pairs are returned, where the first item is the index of the point and the second one is the atom index. Otherwise a list containing the neighbors for each point is returned.

Returns:

If as_pairs=True:

A np.ndarray of shape (n_pairs, 5) is returned. Each pair ij means that j is a neighbor of i. The three extra columns are the supercell indices of atom j.

If as_pairs=False:

A list containing a numpy array of shape (n_neighs, 4) for each atom.

Return type:

np.ndarray or list

find_neighbors(atoms: sisl.typing.AtomsIndex = None, as_pairs: bool = False, self_interaction: bool = False)[source]

Find neighbors as specified in the finder.

This routine only executes the action of finding neighbors, the parameters of the search (geometry, R, overlap…) are defined when the finder is initialized or by calling setup.

Parameters:
  • atoms (optional) –

    the atoms for which neighbors are desired. Anything that can be sanitized by sisl.Geometry is a valid value.

    If not provided, neighbors for all atoms are searched.

  • as_pairs (bool, optional) – If True pairs of atoms are returned. Otherwise a list containing the neighbors for each atom is returned.

  • self_interaction (bool, optional) – whether to consider an atom a neighbor of itself.

Returns:

If as_pairs=True:

A np.ndarray of shape (n_pairs, 5) is returned. Each pair ij means that j is a neighbor of i. The three extra columns are the supercell indices of atom j.

If as_pairs=False:

A list containing a numpy array of shape (n_neighs, 4) for each atom.

Return type:

np.ndarray or list

find_unique_pairs(self_interaction: bool = False)[source]

Find unique neighbor pairs within the geometry.

A pair of atoms is considered unique if atoms have the same index and correspond to the same unit cell. As an example, the connection atom 0 (unit cell) to atom 5 (1, 0, 0) is not the same as the connection atom 5 (unit cell) to atom 0 (-1, 0, 0).

Note that this routine can not be called if overlap is set to False and the radius is not a single float. In that case, there is no way of defining “uniqueness” since pair ij can have a different threshold radius than pair ji.

Parameters:
  • atoms (optional) –

    the atoms for which neighbors are desired. Anything that can be sanitized by sisl.Geometry is a valid value.

    If not provided, neighbors for all atoms are searched.

  • as_pairs (bool, optional) – If True pairs of atoms are returned. Otherwise a list containing the neighbors for each atom is returned.

  • self_interaction (bool, optional) – whether to consider an atom a neighbor of itself.

Returns:

Each pair ij means that j is a neighbor of i. The three extra columns are the supercell indices of atom j.

Return type:

np.ndarray of shape (n_pairs, 5)

geometry: Geometry
memory: Tuple[str, float] = ('200MB', 1.5)

Memory control of the finder

nbins: Tuple[int, int, int]
setup(geometry: Geometry | None = None, R: float | ndarray | None = None, overlap: bool = None, bin_size: float | Tuple[float, float, float] = 2)[source]

Prepares everything for neighbor finding.

Parameters:
  • geometry (Geometry, optional) –

    the geometry to find neighbors in.

    If not provided, the current geometry is kept.

  • R (float or array-like of shape (geometry.na), optional) –

    The radius to consider two atoms neighbors. If it is a single float, the same radius is used for all atoms. If it is an array, it should contain the radius for each atom.

    If not provided, an array is constructed, where the radius for each atom is its maxR.

    Note that negative R values are allowed

  • overlap (bool) – If True, two atoms are considered neighbors if their spheres overlap. If False, two atoms are considered neighbors if the second atom is within the sphere of the first atom. Note that this implies that atom i might be atom j’s neighbor while the opposite is not true.

  • bin_size (float | Tuple[float, float, float]) – the factor for the radius to determine how large the bins are, optionally along each lattice vector. It can minimally be 2, meaning that the maximum radius to consider is twice the radius considered. For larger values, more atoms will be in each bin (and thus fewer bins). Hence, this value can be used to fine-tune the memory requirement by decreasing number of bins, at the cost of a bit more run-time searching bins.

total_nbins: int