NeighborList¶
- class hoomd.md.nlist.NeighborList(buffer, exclusions, rebuild_check_delay, check_dist, mesh, default_r_cut)¶
Bases:
Compute
Base class neighbor list.
NeighborList
is the base class for all neighbor lists.Warning
Users should not instantiate this class directly. The class can be used for
isinstance
orissubclass
checks.
Members inherited from
AutotunedObject
:- property kernel_parameters¶
Kernel parameters.
Read more...
- property is_tuning_complete¶
Check if kernel parameter tuning is complete.
Read more...
- tune_kernel_parameters()¶
Start tuning kernel parameters.
Read more...
Members defined in
NeighborList
:- exclusions¶
Defines which particles to exclude from the neighbor list, as described in
hoomd.md.nlist
.
- r_cut¶
Base cutoff radius for neighbor list queries. \([\mathrm{length}]\). Optional: defaults to the value
default_r_cut
specified on construction.Type:
TypeParameter
[tuple
[particle_type
,particle_type
],float
])
- property cpu_local_nlist_arrays¶
Local nlist arrays on the CPU.
Provides direct access to the neighbor list arrays on the cpu. All data is MPI rank-local.
The
hoomd.md.data.NeighborListLocalAccess
object exposes the internal data representation to efficiently iterate over the neighbor list.Note
The local arrays are read only.
Examples:
with self.cpu_local_nlist_arrays as arrays: nlist_iter = zip(arrays.head_list, arrays.n_neigh) for i, (head, nn) in enumerate(nlist_iter): for j_idx in range(head, head + nn): j = arrays.nlist[j_idx] # i and j are "neighbor" indices
- property gpu_local_nlist_arrays¶
Expose nlist arrays on the GPU.
Provides direct access to the neighbor list arrays on the gpu. All data is MPI rank-local.
The
hoomd.md.data.NeighborListLocalAccessGPU
object exposes the internal data representation to efficiently iterate over the neighbor list.Note
The local arrays are read only.
See also
Examples:
get_local_pairs = cupy.RawKernel(r''' extern "C" __global__ void get_local_pairs( const unsigned int N, const unsigned long* heads, const unsigned int* nns, const unsigned int* nlist, const unsigned int* tags, const unsigned long* offsets, unsigned long* pairs) { unsigned int i = (unsigned int) (blockDim.x * blockIdx.x + threadIdx.x); if (i >= N) return; uint2* pair = (uint2*)pairs; unsigned long head = heads[i]; unsigned int nn = nns[i]; unsigned long offset = offsets[i]; unsigned int tag_i = tags[i]; for (unsigned int idx = 0; idx < nn; idx++) { unsigned int j = nlist[head + idx]; pair[offset + idx] = make_uint2(tag_i, tags[j]); } } ''', 'get_local_pairs') with nlist.gpu_local_nlist_arrays as data: with sim.state.gpu_local_snapshot as snap_data: tags = snap_data.particles.tag_with_ghost tags = tags._coerce_to_ndarray() head_list = data.head_list._coerce_to_ndarray() n_neigh = data.n_neigh._coerce_to_ndarray() raw_nlist = data.nlist._coerce_to_ndarray() N = int(head_list.size) n_pairs = int(cupy.sum(n_neigh)) offsets = cupy.cumsum(n_neigh.astype(cupy.uint64) offsets -= n_neigh[0] device_local_pairs = cupy.zeros( (n_pairs, 2), dtype=cupy.uint32) block = 256 n_grid = (N + 255) // 256 get_local_pairs( (n_grid,), (block,), ( N, head_list, n_neigh, raw_nlist, tags, offsets, device_local_pairs ))
Note
GPU local nlist data is not available if the chosen device for the simulation is
hoomd.device.CPU
.
- property local_pair_list¶
Local pair list.
Note
The local pair list returns rank-local indices, not particle tags.
- Type:
(N_pairs, 2)
numpy.ndarray
ofnumpy.uint32
- property num_builds¶
The number of neighbor list builds.
num_builds
is the number of neighbor list rebuilds performed since the last call toSimulation.run
.(
Loggable
: category=”scalar”, default=False)- Type:
- property pair_list¶
Global pair list.
Note
The pair list returns particle tags, not rank-local indices.
Attention
In MPI parallel execution, the array is available on rank 0 only.
pair_list
isNone
on ranks >= 1.- Type:
(N_pairs, 2)
numpy.ndarray
ofnumpy.uint32
- property shortest_rebuild¶
The shortest period between neighbor list rebuilds.
shortest_rebuild
is the smallest number of time steps between neighbor list rebuilds since the last call toSimulation.run
.(
Loggable
: category=”scalar”)- Type: