md.data#
Overview
Access HOOMD-Blue force data buffers on the CPU. |
|
Access HOOMD-Blue force data buffers on the GPU. |
|
Access HOOMD-Blue neighbor list data buffers on the CPU. |
|
Access HOOMD-Blue neighbor list data buffers on the GPU. |
Details
Force data local access.
ForceLocalAccess, ForceLocalAccessGPU, and related classes provide direct
access to the data buffers managed by hoomd.md.force.Force. This means that
MPI rank locality must be considered in accessing the arrays in a multi-rank
simulation.
- class hoomd.md.data.ForceLocalAccess(force_obj, state)#
Access HOOMD-Blue force data buffers on the CPU.
- force#
Local force data. \([\mathrm{force}]\)
- Type:
(N_particles, 3)
hoomd.data.arrayoffloat
- potential_energy#
Local potential energy data. \([\mathrm{energy}]\)
- Type:
(N_particles,)
hoomd.data.arrayoffloat
- torque#
Local torque data. \([\mathrm{force} \cdot \mathrm{length}]\)
- Type:
(N_particles, 3)
hoomd.data.arrayoffloat
- virial#
Local virial data. \([\mathrm{energy}]\)
- Type:
(N_particles, 6)
hoomd.data.arrayoffloat
- class hoomd.md.data.ForceLocalAccessGPU(*args, **kwargs)#
Access HOOMD-Blue force data buffers on the GPU.
- force#
Local force data. \([\mathrm{force}]\)
- Type:
(N_particles, 3)
hoomd.data.arrayoffloat
- potential_energy#
Local potential energy data. \([\mathrm{energy}]\)
- Type:
(N_particles,)
hoomd.data.arrayoffloat
- torque#
Local torque data. \([\mathrm{force} \cdot \mathrm{length}]\)
- Type:
(N_particles, 3)
hoomd.data.arrayoffloat
- virial#
Local virial data. \([\mathrm{energy}]\)
- Type:
(N_particles, 6)
hoomd.data.arrayoffloat
- class hoomd.md.data.NeighborListLocalAccess(nlist_obj, state)#
Access HOOMD-Blue neighbor list data buffers on the CPU.
The internal
NeighborListimplementation of HOOMD is comprised of essentially three array buffers. The buffers are:nlist: Ragged array of neighbor data.head_list: Indexes for particles to read from the neighbor list.n_neigh: Number of neighbors for each particle.
The neighbor indices of particle \(i\) are stored in the slice
nlist[head_list[i]:head_list[i]+n_neigh[i]]. The result of access outside of these bounds is undefined. Thehalf_nlistproperty is used to query whether the neighbor list stores a single copy for each pair (True), or two copies for each pair (False). Under MPI, pairs that cross domains are stored twice, once in each domain rank.- head_list#
Local head list.
- Type:
(N_particles,)
hoomd.data.arrayofunsigned long
- n_neigh#
Number of neighbors.
- Type:
(N_particles,)
hoomd.data.arrayofunsigned int
- nlist#
Raw neighbor list data.
- Type:
(…)
hoomd.data.arrayofunsigned int
- half_nlist#
Convenience property to check if the storage mode is ‘half’.
- Type:
bool
- class hoomd.md.data.NeighborListLocalAccessGPU(*args, **kwargs)#
Access HOOMD-Blue neighbor list data buffers on the GPU.
The internal
NeighborListimplementation of HOOMD is comprised of essentially three array buffers. The buffers are:nlist: Ragged array of neighbor data.head_list: Indexes for particles to read from the neighbor list.n_neigh: Number of neighbors for each particle.
The neighbor indices of particle \(i\) are stored in the slice
nlist[head_list[i]:head_list[i]+n_neigh[i]]. The result of access outside of these bounds is undefined. Thehalf_nlistproperty is used to query whether the neighbor list stores a single copy for each pair (True), or two copies for each pair (False). Under MPI, pairs that cross domains are stored twice, once in each domain rank.- head_list#
Local head list.
- Type:
(N_particles,)
hoomd.data.arrayofunsigned long
- n_neigh#
Number of neighbors.
- Type:
(N_particles,)
hoomd.data.arrayofunsigned int
- nlist#
Raw neighbor list data.
- Type:
(…)
hoomd.data.arrayofunsigned int
- half_nlist#
Convenience property to check if the storage mode is ‘half’.
- Type:
bool
Modules