13.2.4. Neighbor search library — MDAnalysis.lib.nsgrid

13.2.4.1. About the code

This Neighbor search library is a serialized Cython version greatly inspired by the NS grid search implemented in GROMACS .

GROMACS 4.x code (more precisely nsgrid.c and ns.c ) was used as reference to write this file.

GROMACS 4.x code is released under the GNU Public Licence v2.

13.2.4.2. About the algorithm

The neighbor search implemented here is based on cell lists which allow computation of pairs [1] with a cost of \(O(N)\), instead of \(O(N^2)\). The basic algorithm is described in Appendix F, Page 552 of Understanding Molecular Dynamics: From Algorithm to Applications by Frenkel and Smit.

In brief, the algorithm divides the domain into smaller subdomains called cells and distributes every particle to these cells based on their positions. Subsequently, any distance based query first identifies the corresponding cell position in the domain followed by distance evaluations within the identified cell and neighboring cells only. Care must be taken to ensure that cellsize is greater than the desired search distance, otherwise all of the neighbours might not reflect in the results.

[1]a pair correspond to two particles that are considered as neighbors .

New in version 0.19.0.

13.2.4.3. Classes

class MDAnalysis.lib.nsgrid.FastNS(cutoff, coords, box, max_gridsize=5000, pbc=True)

Grid based search between two group of atoms

Instantiates a class object which uses _PBCBox and _NSGrid to construct a cuboidal grid in an orthogonal brick shaped box.

Minimum image convention is used for distance evaluations if pbc is set to True.

Initialize the grid and sort the coordinates in respective cells by shifting the coordinates in a brick shaped box. The brick shaped box is defined by _PBCBox and cuboidal grid is initialize by _NSGrid. If box is supplied, periodic shifts along box vectors are used to contain all the coordinates inside the brick shaped box. If box is not supplied, the range of coordinates i.e. [xmax, ymax, zmax] - [xmin, ymin, zmin] should be used to construct a pseudo box. Subsequently, the origin should also be shifted to [xmin, ymin, zmin]. These arguments must be provided to the function.

Parameters:
  • cutoff (float) – Desired cutoff distance
  • coords (numpy.ndarray) – atom coordinates of shape (N, 3) for N atoms. dtype=numpy.float32. For Non-PBC calculations, all the coords must be within the bounding box specified by box
  • box (numpy.ndarray) – Box dimension of shape (6, ). The dimensions must be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions: [lx, ly, lz, alpha, beta, gamma]. For non-PBC evaluations, provide an orthogonal bounding box (dtype = numpy.float32)
  • max_gridsize (int) – maximum number of cells in the grid. This parameter can be tuned for superior performance.
  • pbc (boolean) – Handle to switch periodic boundary conditions on/off [True]

Note

  • pbc=False Only works for orthogonal boxes.
  • Care must be taken such that all particles are inside the bounding box as defined by the box argument for non-PBC calculations.
  • In case of Non-PBC calculations, a bounding box must be provided to encompass all the coordinates as well as the search coordinates. The dimension should be similar to box argument but for an orthogonal box. For instance, one valid set of argument for box for the case of no PBC could be [10, 10, 10, 90, 90, 90]
  • Following operations are advisable for non-PBC calculations
lmax = all_coords.max(axis=0)
lmin = all_coords.min(axis=0)
pseudobox[:3] = 1.1*(lmax - lmin)
pseudobox[3:] = 90.
shift = all_coords.copy()
shift -= lmin
gridsearch = FastNS(max_cutoff, shift, box=pseudobox, pbc=False)
search(self, search_coords)

Search a group of atoms against initialized coordinates

Creates a new grid with the query atoms and searches against the initialized coordinates. The search is exclusive i.e. only the pairs (i, j) such that atom[i] from query atoms and atom[j] from the initialized set of coordinates is stored as neighbors.

PBC-aware/non PBC-aware calculations are automatically enabled during the instantiation of :class:FastNS.

Parameters:search_coords (numpy.ndarray) – Query coordinates of shape (N, 3) where N is the number of queries
Returns:results – An NSResults object holding neighbor search results, which can be accessed by its methods get_indices(), get_distances(), get_pairs(), and get_pair_distances().
Return type:NSResults

Note

For non-PBC aware calculations, the current implementation doesn’t work if any of the query coordinates lies outside the box supplied to FastNS.

Searches all the pairs within the initialized coordinates

All the pairs among the initialized coordinates are registered in hald the time. Although the algorithm is still the same, but the distance checks can be reduced to half in this particular case as every pair need not be evaluated twice.

Returns:results – An NSResults object holding neighbor search results, which can be accessed by its methods get_indices(), get_distances(), get_pairs(), and get_pair_distances().
Return type:NSResults
class MDAnalysis.lib.nsgrid.NSResults(dreal cutoff, real[:, ::1] coords, real[:, ::1] searchcoords)

Class to store the results

All outputs from FastNS are stored in an instance of this class. All methods of FastNS return an instance of this class, which can be used to generate the desired results on demand.

Parameters:
  • cutoff (float) – Specified cutoff distance
  • coords (numpy.ndarray) – Array with coordinates of atoms of shape (N, 3) for N particles. dtype must be numpy.float32
  • searchcoords (numpy.ndarray) – Array with query coordinates. Shape must be (M, 3) for M queries. dtype must be numpy.float32
get_distances(self)

Distance corresponding to individual neighbors of query atom

For every queried atom i, a list of all the distances from its neighboring atoms can be obtained from get_distances()[i]. Every distance[i][j] will correspond to the distance between atoms atom[i] from the query atoms and atom[indices[j]] from the initialized set of coordinates, where indices can be obtained by get_indices()

Returns:distances – Every element i.e. distances[i] will be an array of shape m where m is the number of neighbours of query atom[i].
results = NSResults()
distances = results.get_distances()
Return type:list

See also

get_indices()

get_indices(self)

Individual neighbours of query atom

For every queried atom i, an array of all its neighbors indices can be obtained from get_indices()[i]

Returns:indices – Indices of neighboring atoms. Every element i.e. indices[i] will be a list of size m where m is the number of neighbours of query atom[i].
results = NSResults()
indices = results.get_indices()

indices[i] will be a list of neighboring atoms of atom[i] from query atoms atom. indices[i][j] will give the atom-id of initial coordinates such that initial_atom[indices[i][j]] is a neighbor of atom[i].

Return type:list
get_pair_distances(self)

Returns all the distances corresponding to each pair of neighbors

Returns an array of shape N where N is the number of pairs among the query atoms and initial atoms within a specified distance. Every element [i] corresponds to the distance between pairs[i, 0] and pairs[i, 1], where pairs is the array obtained from get_pairs()

Returns:distances – distances between pairs of query and initial atom coordinates of shape N
Return type:numpy.ndarray

See also

get_pairs()

get_pairs(self)

Returns all the pairs within the desired cutoff distance

Returns an array of shape (N, 2), where N is the number of pairs between reference and configuration within the specified distance. For every pair (i, j), reference[i] and configuration[j] are atom positions such that reference is the position of query atoms while configuration coontains the position of group of atoms used to search against the query atoms.

Returns:pairs – pairs of atom indices of neighbors from query and initial atom coordinates of shape (N, 2)
Return type:numpy.ndarray