hoomd#

Overview

Box

Define box dimensions.

Operations

A mutable collection of operations which act on a Simulation.

Simulation

Define a simulation.

Snapshot

Self-contained copy of the simulation State.

State

The state of a Simulation object.

Details

HOOMD-blue python package.

hoomd is the top level HOOMD-blue Python package. It consists of the common code shared among all types of HOOMD-blue simulations. The core data structures used to construct a simulation include:

See the table of contents or the modules section for a full list of classes, methods, and variables in the API.

hoomd also contains subpackages that implement specific types of simulations:

See also

Tutorial: Introducing HOOMD-blue

Signal handling

On import, hoomd installs a SIGTERM signal handler that calls sys.exit so that open gsd files have a chance to flush write buffers (hoomd.write.GSD.flush) when a user’s process is terminated. Use signal.signal to adjust this behavior as needed.

class hoomd.Box(Lx, Ly, Lz=0, xy=0, xz=0, yz=0)#

Define box dimensions.

Parameters:
  • Lx (float) – box extent in the x direction \([\mathrm{length}]\).

  • Ly (float) – box extent in the y direction \([\mathrm{length}]\).

  • Lz (float) – box extent in the z direction \([\mathrm{length}]\).

  • xy (float) – tilt factor xy \([\mathrm{dimensionless}]\).

  • xz (float) – tilt factor xz \([\mathrm{dimensionless}]\).

  • yz (float) – tilt factor yz \([\mathrm{dimensionless}]\).

Example simulation box labelled with lengths and vectors.

Particles in a simulation exist in a triclinic box with periodic boundary conditions. A triclinic box is defined by six values: the extents \(L_x\), \(L_y\) and \(L_z\) of the box in the three directions, and three tilt factors \(xy\), \(xz\) and \(yz\).

The parameter matrix is defined in terms of the lattice vectors \(\vec a_1\), \(\vec a_2\) and \(\vec a_3\):

\[\left( \vec a_1, \vec a_2, \vec a_3 \right)\]

The first lattice vector \(\vec a_1\) is parallel to the unit vector \(\vec e_x = (1,0,0)\). The tilt factor \(xy\) indicates how the second lattice vector \(\vec a_2\) is tilted with respect to the first one. Similarly, \(xz\) and \(yz\) indicate the tilt of the third lattice vector \(\vec a_3\) with respect to the first and second lattice vector.

The full cell parameter matrix is:

\[\begin{split}\left(\begin{array}{ccc} L_x & xy L_y & xz L_z \\ 0 & L_y & yz L_z \\ 0 & 0 & L_z \\ \end{array}\right)\end{split}\]

The tilt factors \(xy\), \(xz\) and \(yz\) are dimensionless. The relationships between the tilt factors and the box angles \(\alpha\), \(\beta\) and \(\gamma\) are as follows:

\[\begin{split}\cos\gamma &= \cos(\angle\vec a_1, \vec a_2) &=& \frac{xy}{\sqrt{1+xy^2}}\\ \cos\beta &= \cos(\angle\vec a_1, \vec a_3) &=& \frac{xz}{\sqrt{1+xz^2+yz^2}}\\ \cos\alpha &= \cos(\angle\vec a_2, \vec a_3) &=& \frac{xy \cdot xz + yz}{\sqrt{1+xy^2} \sqrt{1+xz^2+yz^2}}\end{split}\]

Given an arbitrarily oriented lattice with box vectors \(\vec v_1, \vec v_2, \vec v_3\), the parameters for the rotated box can be found as follows:

\[\begin{split}L_x &= v_1\\ a_{2x} &= \frac{\vec v_1 \cdot \vec v_2}{v_1}\\ L_y &= \sqrt{v_2^2 - a_{2x}^2}\\ xy &= \frac{a_{2x}}{L_y}\\ L_z &= \vec v_3 \cdot \frac{\vec v_1 \times \vec v_2}{\left| \vec v_1 \times \vec v_2 \right|}\\ a_{3x} &= \frac{\vec v_1 \cdot \vec v_3}{v_1}\\ xz &= \frac{a_{3x}}{L_z}\\ yz &= \frac{\vec v_2 \cdot \vec v_3 - a_{2x}a_{3x}}{L_y L_z}\end{split}\]

Box images

HOOMD-blue always stores particle positions \(\vec{r}\) inside the primary box image which includes the origin at the center. The primary box image include the left, bottom, and back face while excluding the right, top, and front face. In cubic boxes, this implies that the particle coordinates in the primary box image are in the interval \(\left[ -\frac{L}{2},\frac{L}{2} \right)\).

Unless otherwise noted in the documentation, operations apply the minimum image convention when computing pairwise interactions between particles:

\[\vec{r}_{ij} = \mathrm{minimum\_image}(\vec{r}_j - \vec{r}_i)\]

When running simulations with a fixed box size, use the particle images \(\vec{n}\) to compute the unwrapped coordinates:

\[\vec{r}_\mathrm{unwrapped} = \vec{r} + n_x \vec{a}_1 + n_y \vec{a}_2 + n_z \vec{a}_3\]

Two dimensional systems

Set Lz == 0 to make the box 2D. 2D boxes ignore xz and yz. Changing the box dimensionality from 2D to 3D (or from 3D to 2D) during a simulation will result in undefined behavior.

In 2D boxes, volume is in units of \([\mathrm{length}]^2\).

Factory Methods

Box has factory methods to enable easier creation of boxes: cube, square, from_matrix, and from_box. See each method’s documentation for more details.

Example:

box = hoomd.Box(Lx=10, Ly=20, Lz=30, xy=0.5, xz=0.2, yz=0.1)
property L#

The box lengths, [Lx, Ly, Lz] \([\mathrm{length}]\).

Can be set with a float which sets all lengths, or a length 3 vector.

Example:

box.L = (15, 30, 60)
Type:

(3, ) numpy.ndarray of float

property Lx#

The length of the box in the x dimension \([\mathrm{length}]\).

Example:

box.Lx = 15
Type:

float

property Ly#

The length of the box in the y dimension \([\mathrm{length}]\).

Example:

box.Ly = 30
Type:

float

property Lz#

The length of the box in the z dimension \([\mathrm{length}]\).

Example:

box.Lz = 60
Type:

float

__eq__(other)#

Test if boxes are equal.

__neq__(other)#

Test if boxes are not equal.

__reduce__()#

Reduce values to picklable format.

__repr__()#

Executable representation of the object.

classmethod cube(L)#

Create a cube with side lengths L.

Parameters:

L (float) – The box side length \([\mathrm{length}]\).

Returns:

The created 3D box.

Return type:

hoomd.Box

Example:

box = hoomd.Box.cube(L=13)
property dimensions#

The dimensionality of the box.

If Lz == 0, the box is treated as 2D, otherwise it is 3D. This property is not settable.

Example:

if box.dimensions == 2:
    pass
Type:

int

classmethod from_box(box)#

Initialize a Box instance from a box-like object.

Parameters:

box (box_like) – A box-like object.

Note

If all values are provided, a triclinic box will be constructed. If only Lx, Ly, Lz are provided, an orthorhombic box will be constructed. If only Lx, Ly are provided, a rectangular (2D) box will be constructed.

Returns:

The resulting box object.

Return type:

hoomd.Box

Example:

box = hoomd.Box.from_box(box=[10, 20, 30, 0.5, 0.2, 0.1])
classmethod from_matrix(box_matrix)#

Create a box from an upper triangular matrix.

Parameters:

box_matrix ((3, 3) numpy.ndarray of float) –

An upper triangular matrix representing a box. The values for Lx, Ly, Lz, xy, xz, and yz are related to the matrix:

\[\begin{split}\begin{bmatrix} L_x & L_y \cdot xy & L_z \cdot xz \\ 0 & L_y & L_z \cdot yz \\ 0 & 0 & L_z \end{bmatrix}\end{split}\]

Returns:

The created box.

Return type:

hoomd.Box

Example:

box = hoomd.Box.from_matrix(
    box_matrix = [[10, 12, 14],
                  [0, 8, 16],
                  [0, 0, 18]])
property is2D#

Flag whether the box is 2D.

If Lz == 0, the box is treated as 2D, otherwise it is 3D. This property is not settable.

Example:

if box.is2D:
    pass
Type:

bool

property periodic#

The periodicity of each dimension.

periodic is always (True, True, True) for the box associated with the simulation State. Some components of periodic may be False in the hoomd.data.LocalSnapshot box attribute in MPI domain decomposition simulations. This indicates which box directions are communicated with neighboring ranks (False) and which are not (True).

Type:

(3, ) numpy.ndarray of bool

scale(s)#

Scale box dimensions.

Scales the box in place by the given scale factors. Tilt factors are not modified.

Parameters:

s (float or list[float]) – scale factors in each dimension. If a single float is given then scale all dimensions by s; otherwise, s must be a sequence of 3 values used to scale each dimension.

Returns:

self

Examples:

box.scale(2)
box.scale((1, 2, 4))
classmethod square(L)#

Create a square with side lengths L.

Parameters:

L (float) – The box side length \([\mathrm{length}]\).

Returns:

The created 2D box.

Return type:

hoomd.Box

box = hoomd.Box.square(L=128)
property tilts#

The box tilts, [xy, xz, yz].

Can be set using one tilt for all axes or three tilts. If the box is 2D xz and yz will automatically be set to zero.

Example:

box.tilts = (1.1, 0.8, 0.2)
Type:

(3, ) numpy.ndarray of float

to_matrix()#

(3, 3) numpy.ndarray float: The upper triangular matrix that defines the box.

[[Lx, Ly * xy, Lz * xz],
 [0,  Ly,      Lz * yz],
 [0,  0,       Lz]]

Example:

matrix = box.to_matrix()
property volume#

Volume of the box.

\([\mathrm{length}]^{2}\) in 2D and \([\mathrm{length}]^{3}\) in 3D.

When setting volume the aspect ratio of the box is maintained while the lengths are changed.

Example:

box.volume = 2000
Type:

float

property xy#

The tilt for the xy plane.

Example:

box.xy = 1.1
Type:

float

property xz#

The tilt for the xz plane.

Example:

box.xz = 0.8
Type:

float

property yz#

The tilt for the yz plane.

Example:

box.yz = 0.2
Type:

float

class hoomd.Operations#

A mutable collection of operations which act on a Simulation.

An Operations class instance contains all the operations acting on a simulation. These operations are classes that perform various actions on a hoomd.Simulation. Operations can be added and removed at any point from a hoomd.Operations instance. The class provides the interface defined by collections.abc.Collection. Other methods for manipulating instances mimic Python objects where possible, but the class is not simply a mutable list or set. Operations objects manage multiple independent sequences described below.

The types of operations which can be added to an Operations object are tuners, updaters, integrators, writers, and computes. An Operations instance can have zero or one integrator and any number of tuners, updaters, writers, or computes. To see examples of these types of operations see hoomd.tune (tuners), hoomd.update (updaters), hoomd.hpmc.integrate or hoomd.md.Integrator (integrators), hoomd.write (writers), and hoomd.md.compute.ThermodynamicQuantities (computes).

A given instance of an operation class can only be added to a single Operations container. Likewise, a single instance cannot be added to the same Operations container more than once.

All Operations instances start with a hoomd.tune.ParticleSorter instance in their tuners attribute. This increases simulation performance. However, users can choose to modify or remove this tuner if desired.

Note

An Operations object is created by default when a new simulation is created.

__contains__(operation)#

Whether an operation is contained in this container.

Parameters:

operation – Returns whether this exact operation is contained in the collection.

Example:

operation in simulation.operations
__getstate__()#

Get the current state of the operations container for pickling.

__iadd__(operation)#

Works the same as Operations.add.

Parameters:

operation (hoomd.operation.Operation) – A HOOMD-blue tuner, updater, integrator, writer, or compute to add to the object.

Example:

simulation.operations += operation
__isub__(operation)#

Works the same as Operations.remove.

Parameters:

operation (hoomd.operation.Operation) – A HOOMD-blue integrator, tuner, updater, integrator, analyzer, or compute to remove from the collection.

Example:

simulation.operations -= operation
__iter__()#

Iterates through all contained operations.

Example:

for operation in simulation.operations:
    pass
__len__()#

Return the number of operations contained in this collection.

Example:

len(simulation.operations)
add(operation)#

Add an operation to this container.

Adds the provided operation to the appropriate attribute of the Operations instance.

Parameters:

operation (hoomd.operation.Operation) – A HOOMD-blue tuner, updater, integrator, writer, or compute to add to the collection.

Raises:

TypeError – If operation is not of a valid type.

Note

Since only one integrator can be associated with an Operations object at a time, this removes the current integrator when called with an integrator operation. Also, the integrator property cannot be set to None using this function. Use operations.integrator = None explicitly for this.

Example:

simulation.operations.add(operation)
property computes#

A list of compute operations.

Holds the list of computes associated with this collection. The list can be modified as a standard Python list.

Type:

list[hoomd.operation.Compute]

property integrator#

An MD or HPMC integrator object.

Operations objects have an initial integrator property of None. Can be set to MD or HPMC integrators. The property can also be set to None.

Examples:

simulation.operations.integrator = hoomd.md.Integrator(dt=0.001)
simulation.operations.integrator = None
Type:

hoomd.operation.Integrator

property is_tuning_complete#

Check whether all children have completed tuning.

True when is_tuning_complete is True for all children.

Note

In MPI parallel execution, is_tuning_complete is True only when all children on all ranks have completed tuning.

Example:

while (not simulation.operations.is_tuning_complete):
    simulation.run(1000)
Type:

bool

remove(operation)#

Remove an operation from the Operations object.

Remove the item from the collection whose Python object id is the same as operation.

Parameters:

operation (hoomd.operation.Operation) – A HOOMD-blue integrator, tuner, updater, integrator, or compute to remove from the container.

Raises:
  • ValueError – If operation is not found in this container.

  • TypeError – If operation is not of a valid type.

Example:

simulation.operations.remove(operation)
tune_kernel_parameters()#

Start tuning kernel parameters in all children.

Example:

simulation.operations.tune_kernel_parameters()
property tuners#

A list of tuner operations.

Holds the list of tuners associated with this collection. The list can be modified as a standard Python list.

Type:

list[hoomd.operation.Tuner]

property updaters#

A list of updater operations.

Holds the list of updaters associated with this collection. The list can be modified as a standard Python list.

Type:

list[hoomd.operation.Updater]

property writers#

A list of writer operations.

Holds the list of writers associated with this collection. The list can be modified as a standard Python list.

Type:

list[hoomd.operation.Writer]

class hoomd.Simulation(device, seed=None)#

Define a simulation.

Parameters:

Simulation is the central class that defines a simulation, including the state of the system, the operations that apply to the state during a simulation run, and the device to use when executing the simulation.

seed sets the seed for the random number generator used by all operations added to this Simulation.

Newly initialized Simulation objects have no state. Call create_state_from_gsd or create_state_from_snapshot to initialize the simulation’s state.

Example:

simulation = hoomd.Simulation(device=hoomd.device.CPU(), seed=1)
__del__()#

Clean up dangling references to simulation.

property always_compute_pressure#

Always compute the virial and pressure (defaults to False).

By default, HOOMD only computes the virial and pressure on timesteps where it is needed (when hoomd.write.GSD writes log data to a file or when using an NPT integrator). Set always_compute_pressure to True to make the per particle virial, net virial, and system pressure available to query any time by property or through the hoomd.logging.Logger interface.

Note

Enabling this flag will result in a moderate performance penalty when using MD pair potentials.

Example:

simulation.always_compute_pressure = True
Type:

bool

create_state_from_gsd(filename, frame=-1, domain_decomposition=(None, None, None))#

Create the simulation state from a GSD file.

Parameters:
  • filename (str) – GSD file to read

  • frame (int) – Index of the frame to read from the file. Negative values index back from the last frame in the file.

  • domain_decomposition (tuple) – Choose how to distribute the state across MPI ranks with domain decomposition. Provide a tuple of 3 integers indicating the number of evenly spaced domains in the x, y, and z directions (e.g. (8,4,2)). Provide a tuple of 3 lists of floats to set the fraction of the simulation box to include in each domain. The sum of each list of floats must be 1.0 (e.g. ([0.25, 0.75], [0.2, 0.8], [1.0])).

When timestep is None before calling, create_state_from_gsd sets timestep to the value in the selected GSD frame in the file.

Note

Set any or all of the domain_decomposition tuple elements to None and create_state_from_gsd will select a value that minimizes the surface area between the domains (e.g. (2,None,None)). The domains are spaced evenly along each automatically selected direction. The default value of (None, None, None) will automatically select the number of domains in all directions.

Example:

simulation.create_state_from_gsd(filename=gsd_filename)
create_state_from_snapshot(snapshot, domain_decomposition=(None, None, None))#

Create the simulation state from a Snapshot.

Parameters:
  • snapshot (Snapshot or gsd.hoomd.Frame) – Snapshot to initialize the state from. A gsd.hoomd.Frame will first be converted to a hoomd.Snapshot.

  • domain_decomposition (tuple) – Choose how to distribute the state across MPI ranks with domain decomposition. Provide a tuple of 3 integers indicating the number of evenly spaced domains in the x, y, and z directions (e.g. (8,4,2)). Provide a tuple of 3 lists of floats to set the fraction of the simulation box to include in each domain. The sum of each list of floats must be 1.0 (e.g. ([0.25, 0.75], [0.2, 0.8], [1.0])).

When timestep is None before calling, create_state_from_snapshot sets timestep to 0.

Note

Set any or all of the domain_decomposition tuple elements to None and create_state_from_snapshot will select a value that minimizes the surface area between the domains (e.g. (2,None,None)). The domains are spaced evenly along each automatically selected direction. The default value of (None, None, None) will automatically select the number of domains in all directions.

Example:

simulation.create_state_from_snapshot(snapshot=snapshot)
property device#

Device used to execute the simulation.

Type:

hoomd.device.Device

property final_timestep#

run will end at this timestep.

final_timestep is the timestep on which the currently executing run will complete.

Example:

logger.add(obj=simulation, quantities=['final_timestep'])

(Loggable: category=”scalar”)

Type:

float

property initial_timestep#

run started at this timestep.

initial_timestep is the timestep on which the currently executing run started.

Example:

logger.add(obj=simulation, quantities=['initial_timestep'])

(Loggable: category=”scalar”)

Type:

float

property loggables#

Name, category mapping of loggable quantities.

Type:

dict[str, str]

property operations#

The operations that apply to the state.

The operations apply to the state during the simulation run when scheduled.

See also

run

Type:

hoomd.Operations

run(steps, write_at_start=False)#

Advance the simulation a number of steps.

Parameters:
  • steps (int) – Number of steps to advance the simulation.

  • write_at_start (bool) – When True, writers with triggers that evaluate True for the initial step will be executed before the time step loop.

Note

Initialize the simulation’s state before calling run.

Simulation applies its operations to the state during each time step in the order: tuners, updaters, integrator, then writers following the logic in this pseudocode:

if write_at_start:
    for writer in operations.writers:
        if writer.trigger(timestep):
            writer.write(timestep)

end_step = timestep + steps
while timestep < end_step:
    for tuner in operations.tuners:
        if tuner.trigger(timestep):
            tuner.tune(timestep)

    for updater in operations.updaters:
        if updater.trigger(timestep):
            updater.update(timestep)

    if operations.integrator is not None:
        operations.integrator(timestep)

    timestep += 1

    for writer in operations.writers:
        if writer.trigger(timestep):
            writer.write(timestep)

This order of operations ensures that writers (such as hoomd.write.GSD) capture the final output of the last step of the run loop. For example, a writer with a trigger hoomd.trigger.Periodic(period=100, phase=0) active during a run(500) would write on steps 100, 200, 300, 400, and 500. Set write_at_start=True on the first call to run to also obtain output at step 0.

Warning

Using write_at_start=True in subsequent calls to run will result in duplicate output frames.

Example:

simulation.run(1_000)
property seed#

Random number seed.

Seeds are in the range [0, 65535]. When set, seed will take only the lowest 16 bits of the given value.

HOOMD-blue uses a deterministic counter based pseudorandom number generator. Any time a random value is needed, HOOMD-blue computes it as a function of the user provided seed seed (16 bits), the current timestep (lower 40 bits), particle identifiers, MPI ranks, and other unique identifying values as needed to sample uncorrelated values: random_value = f(seed, timestep, ...)

Example:

simulation.seed = 2

(Loggable: category=”scalar”)

Type:

int

property state#

The current simulation state.

Type:

hoomd.State

property timestep#

The current simulation time step.

timestep is read only after creating the simulation state.

Note

Methods like create_state_from_gsd will set the initial timestep from the input. Set timestep before creating the simulation state to override values from create_ methods:

sim.timestep = 5000
sim.create_state_from_gsd('gsd_at_step_10000000.gsd')
assert sim.timestep == 5000

(Loggable: category=”scalar”)

Type:

int

property tps#

The average number of time steps per second.

tps is the number of steps executed divided by the elapsed walltime in seconds. It is updated during the run loop and remains fixed after run completes.

Warning

The elapsed walltime and timestep are reset at the beginning of each call to run. Thus, tps will provide noisy estimates of performance at the start and stable long term averages after many timesteps.

Tip

Use the total elapsed wall time and timestep to average the timesteps executed per second at desired intervals.

Example:

logger.add(obj=simulation, quantities=['tps'])

(Loggable: category=”scalar”)

Type:

float

property walltime#

The walltime spent during the last call to run.

walltime is the number seconds that the last call to run took to complete. It is updated during the run loop and remains fixed after run completes.

Note

walltime resets to 0 at the beginning of each call to run.

Example:

logger.add(obj=simulation, quantities=['walltime'])

(Loggable: category=”scalar”)

Type:

float

class hoomd.Snapshot(communicator=None)#

Self-contained copy of the simulation State.

Parameters:

communicator (Communicator) – MPI communicator to be used when accessing the snapshot.

See State and gsd.hoomd.Frame for detailed documentation on the components of Snapshot.

Note

Snapshot is duck-type compatible with gsd.hoomd.Frame except that arrays in Snapshot are not assignable. You can edit their contents: e.g. snapshot.particles.typeid[:] == 0.

Warning

Data is only present on the root rank:

if snapshot.communicator.rank == 0:
    pos = snapshot.particles.position[0]

Example:

snapshot = hoomd.Snapshot()
communicator#

MPI communicator.

Type:

Communicator

property angles#

Angles.

angles.N#

Number of angles.

Type:

int

angles.types#

Names of the angle types

Type:

list[str]

angles.typeid#

Angle type id.

Type:

(N,) numpy.ndarray of uint32

angles.group#

Tags of the particles in the angle.

Type:

(N, 3) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.angles.N = 1
    snapshot.angles.group[:] = [[0, 1, 2]]
    snapshot.angles.types = ['A-B-B']
    snapshot.angles.typeid[:] = [0]
property bonds#

Bonds.

bonds.N#

Number of bonds.

Type:

int

bonds.types#

Names of the bond types

Type:

list[str]

bonds.typeid#

Bond type id.

Type:

(N,) numpy.ndarray of uint32

bonds.group#

Tags of the particles in the bond.

Type:

(N, 2) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.bonds.N = 2
    snapshot.bonds.group[:] = [[0, 1], [2, 3]]
    snapshot.bonds.types = ['A-B']
    snapshot.bonds.typeid[:] = [0, 0]
property configuration#

Snapshot box configuration.

dimensions#

Number of dimensions

Type:

int

box#

Simulation box parameters [Lx, Ly, Lz, xy, xz, yz].

Type:

tuple[float, float, float, float, float, float]

Note

box accepts any values that Box.from_box allows when setting.

See also

Box

Example:

snapshot.configuration.box = [10, 20, 30, 0.1, 0.2, 0.3]
property constraints#

Constraints.

constraints.N#

Number of constraints.

Type:

int

constraints.value#

Constraint length.

Type:

(N, ) numpy.ndarray of float

constraints.group#

Tags of the particles in the constraint.

Type:

(N, 2) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.constraints.N = 2
    snapshot.constraints.group[:] = [[0, 1], [2, 3]]
    snapshot.constraints.value[:] = [1, 1]
property dihedrals#

Dihedrals.

dihedrals.N#

Number of dihedrals.

Type:

int

dihedrals.types#

Names of the dihedral types

Type:

list[str]

dihedrals.typeid#

Dihedral type id.

Type:

(N,) numpy.ndarray of uint32

dihedrals.group#

Tags of the particles in the dihedral.

Type:

(N, 4) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.dihedrals.N = 1
    snapshot.dihedrals.group[:] = [[0, 1, 2, 3]]
    snapshot.dihedrals.types = ['A-B-B-A']
    snapshot.dihedrals.typeid[:] = [0]
classmethod from_gsd_frame(gsd_snap, communicator)#

Constructs a hoomd.Snapshot from a gsd.hoomd.Frame object.

Parameters:

Tip

Use Simulation.create_state_from_gsd to efficiently initialize the system state from a GSD file.

Note

from_gsd_frame only accesses the gsd_snap argument on rank 0. In MPI simulations, avoid duplicating memory and file reads by reading GSD files only on rank 0 and passing gsd_snap=None on other ranks.

classmethod from_gsd_snapshot(gsd_snap, communicator)#

Constructs a hoomd.Snapshot from a gsd.hoomd.Snapshot object.

Deprecated since version 4.0.0: Use from_gsd_frame.

property impropers#

Impropers.

impropers.N#

Number of impropers.

Type:

int

impropers.types#

Names of the improper types

Type:

list[str]

impropers.typeid#

Improper type id.

Type:

(N,) numpy.ndarray of uint32

impropers.group#

Tags of the particles in the improper.

Type:

(N, 4) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

if snapshot.communicator.rank == 0:
    snapshot.impropers.N = 1
    snapshot.impropers.group[:] = [[0, 1, 2, 3]]
    snapshot.impropers.types = ['A-B-B-A']
    snapshot.impropers.typeid[:] = [0]
property mpcd#

MPCD data.

mpcd.N#

Number of MPCD particles.

Type:

int

mpcd.position#

Particle position \([\mathrm{length}]\).

Type:

(N, 3) numpy.ndarray of float

mpcd.velocity#

Particle velocity \([\mathrm{velocity}]\).

Type:

(N, 3) numpy.ndarray of float

mpcd.types#

Names of the particle types.

Type:

list[str]

mpcd.typeid#

Particle type id.

Type:

(N, ) numpy.ndarray of uint32

mpcd.mass#

Particle mass.

Type:

float

Note

Set N to change the size of the arrays.

Note

This attribute is only available when HOOMD-blue is built with the MPCD component.

property pairs#

Special pairs.

pairs.N#

Number of special pairs.

Type:

int

pairs.types#

Names of the special pair types

Type:

list[str]

pairs.typeid#

Special pair type id.

Type:

(N,) numpy.ndarray of uint32

pairs.group#

Tags of the particles in the special pair.

Type:

(N, 2) numpy.ndarray of uint32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.pairs.N = 2
    snapshot.pairs.group[:] = [[0, 1], [2, 3]]
    snapshot.pairs.types = ['A-B']
    snapshot.pairs.typeid[:] = [0, 0]
property particles#

Particles.

particles.N#

Number of particles in the snapshot.

Type:

int

particles.types#

Names of the particle types.

Type:

list[str]

particles.position#

Particle position \([\mathrm{length}]\).

Type:

(N, 3) numpy.ndarray of float

particles.orientation#

Particle orientation.

Type:

(N, 4) numpy.ndarray of float

particles.typeid#

Particle type id.

Type:

(N, ) numpy.ndarray of uint32

particles.mass#

Particle mass \([\mathrm{mass}]\).

Type:

(N, ) numpy.ndarray of float

particles.charge#

Particle charge \([\mathrm{charge}]\).

Type:

(N, ) numpy.ndarray of float

particles.diameter#

Particle diameter \([\mathrm{length}]\).

Type:

(N, ) numpy.ndarray of float

particles.body#

Particle body.

Type:

(N, ) numpy.ndarray of int32

particles.moment_inertia#

Particle moment of inertia \([\mathrm{mass} \cdot \mathrm{length}^2]\).

Type:

(N, 3) numpy.ndarray of float

particles.velocity#

Particle velocity \([\mathrm{velocity}]\).

Type:

(N, 3) numpy.ndarray of float

particles.angmom#

Particle angular momentum \([\mathrm{mass} \cdot \mathrm{velocity} \cdot \mathrm{length}]\).

Type:

(N, 4) numpy.ndarray of float

particles.image#

Particle image.

Type:

(N, 3) numpy.ndarray of int32

Note

Set N to change the size of the arrays.

Example:

if snapshot.communicator.rank == 0:
    snapshot.particles.N = 4
    snapshot.particles.position[:] = [[0, 0, 0],
                                    [1, 0, 0],
                                    [0, 1, 0],
                                    [0, 0, 1]]
    snapshot.particles.types = ['A', 'B']
    snapshot.particles.typeid[:] = [0, 1, 1, 0]
replicate(nx, ny, nz=1)#

Replicate the snapshot along the periodic box directions.

Parameters:
  • nx (int) – Number of times to replicate in the x direction.

  • ny (int) – Number of times to replicate in the y direction.

  • nz (int) – Number of times to replicate in the z direction.

Performs the same operation as State.replicate on a Snapshot.

Returns:

self

Example:

snapshot.replicate(nx=2, ny=2, nz=2)
wrap()#

Wrap particles into the snapshot box.

Returns:

self

Example:

snapshot.wrap()
class hoomd.State(simulation, snapshot, domain_decomposition)#

The state of a Simulation object.

Note

This object cannot be directly instantiated. Use Simulation.create_state_from_gsd and Simulation.create_state_from_snapshot to instantiate a State object as part of a simulation.

Overview

State stores the data that describes the thermodynamic microstate of a Simulation object. This data consists of the box, particles, bonds, angles, dihedrals, impropers, special pairs, and constraints.

Box

The simulation box describes the space that contains the particles as a Box object.

Particles

The state contains N_particles particles. Each particle has a position, orientation, type id, body, mass, moment of inertia, charge, diameter, velocity, angular momentum, image, and tag:

  • \(\vec{r}\): position \([\mathrm{length}]\) - X,Y,Z cartesian coordinates defining the position of the particle in the box.

  • \(\mathbf{q}\): orientation \([\mathrm{dimensionless}]\) - Unit quaternion defining the rotation from the particle’s local reference frame to the box reference frame. The four components are in the order \(s\), \(a_x\), \(a_y\), \(a_z\) for the in complex notation \(s + a_x i + a_y j + a_z k\).

  • particle_typeid: type id \([\mathrm{dimensionless}]\) - An integer in the interval [0,len(particle_types)) that identifies the particle’s type. particle_types maps type ids to names with: name = particle_types[particle_typeid].

  • particle_body: body id \([\mathrm{dimensionless}]\) - An integer that identifies the particle’s rigid body. A value of -1 indicates that this particle does not belong to a body. A positive value indicates that the particle belongs to the body particle_body. This particle is the central particle of a body when the body id is equal to the tag \(\mathrm{particle\_body} = \mathrm{particle\_tag}\). (used by hoomd.md.constrain.Rigid)

  • \(m\): mass \([\mathrm{mass}]\) - The particle’s mass.

  • \(I\): moment of inertia \([\mathrm{mass} \cdot \mathrm{length}^2]\) - \(I_{xx}\), \(I_{yy}\), \(I_{zz}\) elements of the diagonal moment of inertia tensor in the particle’s local reference frame. The off-diagonal elements are 0.

  • \(q\): charge \([\mathrm{charge}]\) - The particle’s charge.

  • \(d\): diameter \([\mathrm{length}]\) - Deprecated in v3.0.0. HOOMD-blue reads and writes particle diameters, but does not use them in any computations.

  • \(\vec{v}\): velocity \([\mathrm{velocity}]\) - X,Y,Z components of the particle’s velocity in the box’s reference frame.

  • \(\mathbf{P_S}\): angular momentum \([\mathrm{mass} \cdot \mathrm{velocity} \cdot \mathrm{length}]\) - Quaternion defining the particle’s angular momentum (see note).

  • \(\vec{n}\) : image \([\mathrm{dimensionless}]\) - Integers x,y,z that record how many times the particle has crossed each of the periodic box boundaries.

  • particle_tag : tag \([\mathrm{dimensionless}]\) - An integer that uniquely identifies a given particle. The particles are stored in tag order when writing and initializing to/from a GSD file or snapshot: \(\mathrm{particle\_tag}_i = i\). When accessing data in local snapshots, particles may be in any order.

Note

HOOMD stores angular momentum as a quaternion because that is the form used when integrating the equations of motion (see Kamberaj 2005). The angular momentum quaternion \(\mathbf{P_S}\) is defined with respect to the orientation quaternion of the particle \(\mathbf{q}\) and the vector angular momentum of the particle, lifted into pure imaginary quaternion form \(\mathbf{S}^{(4)}\) as:

\[\mathbf{P_S} = 2 \mathbf{q} \times \mathbf{S}^{(4)}\]

. Following this, the angular momentum vector \(\vec{S}\) in the particle’s local reference frame is:

\[\vec{S} = \frac{1}{2}im(\mathbf{q}^* \times \mathbf{P_S})\]

Bonded groups

The state contains N_bonds bonds, N_angles angles, N_dihedrals dihedrals, N_impropers impropers, and N_special_pairs special pairs. Each of these data structures is similar, differing in the number of particles in the group and what operations use them. Bonds, angles, dihedrals, and impropers contain 2, 3, 4, and 4 particles per group respectively. Bonds specify the toplogy used when computing energies and forces in md.bond, angles define the same for md.angle, dihedrals for md.dihedral and impropers for md.improper. These collectively implement bonding potentials used in molecular dynamics force fields. Like bonds, special pairs define connections between two particles, but special pairs are intended to adjust the 1-4 pairwise interactions in some molecular dynamics force fields: see md.special_pair. Each bonded group is defined by a type id, the group members, and a tag.

  • bond_typeid: type id \([\mathrm{dimensionless}]\) - An integer in the interval [0,len(bond_types)) that identifies the bond’s type. bond_types maps type ids to names with: name = bond_types[bond_typeid]. Similarly, angle_types lists the angle types, dihedral_types lists the dihedral types, improper_types lists the improper types, and special_pair_types lists the special pair types.

  • bond_group: A list of integers in the interval \([0, \max(\mathrm{particle\_tag})]\) that defines the tags of the particles in the bond (2), angle in angle_group (3), dihedral in dihedral_group (4), improper in improper_group (4), or special pair in pair_group (2).

  • bond_tag : tag \([\mathrm{dimensionless}]\) - An integer that uniquely identifies a given bond. The bonds are in tag order when writing and initializing to/from a GSD file or snapshot \(\mathrm{bond\_tag}_i = i\). When accessing data in local snapshots, bonds may be in any order. The same applies to angles with angle_tag, dihedrals with dihedral_tag, impropers with improper_tag, and special pairs with pair_tag.

Constraints

The state contains N_constraints distance constraints between particles. These constraints are used by hoomd.md.constrain.Distance. Each distance constraint consists of a distance value and the group members.

  • constraint_group: A list of 2 integers in the interval \([0, \max(\mathrm{particle\_tag})]\) that identifies the tags of the particles in the constraint.

  • \(d\): constraint value \([\mathrm{length}]\) - The distance between particles in the constraint.

MPI domain decomposition

When running in serial or on 1 MPI rank, the entire simulation state is stored in that process. When using more than 1 MPI rank, HOOMD-blue employs a domain decomposition approach to split the simulation box an integer number of times in the x, y, and z directions (domain_decomposition). Each MPI rank stores and operates on the particles local to that rank. Local particles are those contained within the region defined by the split planes (domain_decomposition_split_fractions). Each MPI rank communicates with its neighbors to obtain the properties of particles near the boundary between ranks (ghost particles) so that it can compute interactions across the boundary.

Accessing Data

Two complementary APIs provide access to the state data: local snapshots that access data directly available on the local MPI rank (including the local and ghost particles) and global snapshots that collect the entire state on rank 0. See State.cpu_local_snapshot, State.gpu_local_snapshot, get_snapshot, and set_snapshot for information about these data access patterns.

See also

To write the simulation to disk, use write.GSD.

property N_angles#

The number of angles in the simulation state.

Example:

N_angles = simulation.state.N_angles
Type:

int

property N_bonds#

The number of bonds in the simulation state.

Example:

N_bonds = simulation.state.N_bonds
Type:

int

property N_constraints#

The number of constraints in the simulation state.

Example:

N_constraints = simulation.state.N_constraints
Type:

int

property N_dihedrals#

The number of dihedrals in the simulation state.

Example:

N_dihedrals = simulation.state.N_dihedrals
Type:

int

property N_impropers#

The number of impropers in the simulation state.

Example:

N_impropers = simulation.state.N_impropers
Type:

int

property N_particles#

The number of particles in the simulation state.

Example:

N_particles = simulation.state.N_particles
Type:

int

property N_special_pairs#

The number of special pairs in the simulation state.

Type:

int

property angle_types#

List of all angle types in the simulation state.

Example:

angle_types = simulation.state.angle_types
Type:

list[str]

property bond_types#

List of all bond types in the simulation state.

Example:

bond_types = simulation.state.bond_types
Type:

list[str]

property box#

A copy of the current simulation box.

Note

The box property cannot be set. Call set_box to set a new simulation box.

Example:

box = simulation.state.box
Type:

hoomd.Box

property cpu_local_snapshot#

Expose simulation data on the CPU.

Provides access directly to the system state’s particle, bond, angle, dihedral, improper, constaint, and pair data through a context manager. Data in State.cpu_local_snapshot is MPI rank local, and the hoomd.data.LocalSnapshot object is only usable within a context manager (i.e. with sim.state.cpu_local_snapshot as data:). Attempts to assess data outside the context manager will result in errors. The local snapshot interface is similar to that of Snapshot.

The hoomd.data.LocalSnapshot data access is mediated through hoomd.data.array.HOOMDArray objects. This lets us ensure memory safety when directly accessing HOOMD-blue’s data. The interface provides zero-copy access (zero-copy is guaranteed on CPU, access may be zero-copy if running on GPU).

Changing the data in the buffers exposed by the local snapshot will change the data across the HOOMD-blue simulation. For a trivial example, this example would set all particle z-axis positions to 0.

with simulation.state.cpu_local_snapshot as local_snapshot:
    local_snapshot.particles.position[:, 2] = 0

Note

The state’s box and the number of particles, bonds, angles, dihedrals, impropers, constaints, and pairs cannot change within the context manager.

Note

Getting a local snapshot object is order \(O(1)\) and setting a single value is of order \(O(1)\).

Type:

hoomd.data.LocalSnapshot

property dihedral_types#

List of all dihedral types in the simulation state.

Example:

angle_types = simulation.state.angle_types
Type:

list[str]

property domain_decomposition#

Number of domains in the x, y, and z directions.

Type:

tuple(int, int, int)

property domain_decomposition_split_fractions#

Box fractions of the domain split planes in the x, y, and z directions.

Type:

tuple(list[float], list[float], list[float])

get_snapshot()#

Make a copy of the simulation current state.

State.get_snapshot makes a copy of the simulation state and makes it available in a single object. set_snapshot resets the internal state to that in the given snapshot. Use these methods to implement techniques like hybrid MD/MC or umbrella sampling where entire system configurations need to be reset to a previous one after a rejected move.

Note

Data across all MPI ranks and from GPUs is gathered on the root MPI rank’s memory. When accessing data in MPI simulations, use a if snapshot.communicator.rank == 0: conditional to access data arrays only on the root rank.

Note

State.get_snapshot is an order \(O(N_{particles} + N_{bonds} + \ldots)\) operation.

See also

set_snapshot

Returns:

The current simulation state

Return type:

Snapshot

Example:

snapshot = simulation.state.get_snapshot()
property gpu_local_snapshot#

Expose simulation data on the GPU.

Provides access directly to the system state’s particle, bond, angle, dihedral, improper, constaint, and pair data through a context manager. Data in State.gpu_local_snapshot is GPU local, and the hoomd.data.LocalSnapshotGPU object is only usable within a context manager (i.e. with sim.state.gpu_local_snapshot as data:). Attempts to assess data outside the context manager will result in errors. The local snapshot interface is similar to that of Snapshot.

The hoomd.data.LocalSnapshotGPU data access is mediated through hoomd.data.array.HOOMDGPUArray objects. This helps us maintain memory safety when directly accessing HOOMD-blue’s data. The interface provides zero-copy access on the GPU (assuming data was last accessed on the GPU).

Changing the data in the buffers exposed by the local snapshot will change the data across the HOOMD-blue simulation. For a trivial example, this example would set all particle z-axis positions to 0.

with simulation.state.gpu_local_snapshot as local_snapshot:
    local_snapshot.particles.position[:, 2] = 0

Warning

This property is only available when running on a GPU (or multiple GPUs).

Note

The state’s box and the number of particles, bonds, angles, dihedrals, impropers, constaints, and pairs cannot change within the context manager.

Note

Getting a local snapshot object is order \(O(1)\) and setting a single value is of order \(O(1)\).

Type:

hoomd.data.LocalSnapshotGPU

property improper_types#

List of all improper types in the simulation state.

Example:

improper_types = simulation.state.improper_types
Type:

list[str]

property particle_types#

List of all particle types in the simulation state.

Example:

particle_types = simulation.state.particle_types
Type:

list[str]

replicate(nx, ny, nz=1)#

Replicate the state of the system along the periodic box directions.

Parameters:
  • nx (int) – Number of times to replicate in the x direction.

  • ny (int) – Number of times to replicate in the y direction.

  • nz (int) – Number of times to replicate in the z direction.

replicate makes the system state nx * ny * nz times larger. In each of the new periodic box images, it places a copy of the initial state with the particle positions offset to locate them in the image and the bond, angle, dihedral, improper, and pair group tags offset to apply to the copied particles. All other particle properties (mass, typeid, velocity, charge, …) are copied to the new particles without change.

After placing the particles, replicate expands the simulation box by a factor of nx, ny, and nz in the direction of the first, second, and third box lattice vectors respectively and adjusts the particle positions to center them in the new box.

Example:

simulation.state.replicate(nx=2, ny=2, nz=2)
set_box(box)#

Set a new simulation box.

Parameters:

box (hoomd.box.box_like) – New simulation box.

Note

All particles must be inside the new box. set_box does not change any particle properties.

Example:

simulation.state.set_box(box)
set_snapshot(snapshot)#

Restore the state of the simulation from a snapshot.

Also calls update_group_dof to count the number of degrees in the system with the new state.

Parameters:

snapshot (Snapshot) – Snapshot of the system from get_snapshot

Warning

set_snapshot can only make limited changes to the simulation state. While it can change the number of particles/bonds/etc… or their properties, it cannot change the number or names of the particle/bond/etc.. types.

Note

set_snapshot is an order \(O(N_{particles} + N_{bonds} + \ldots)\) operation and is very expensive when the simulation device is a GPU.

Example:

simulation.state.set_snapshot(snapshot)
property special_pair_types#

List of all special pair types in the simulation state.

Example:

special_pair_types = simulation.state.special_pair_types
Type:

list[str]

thermalize_particle_momenta(filter, kT)#

Assign random values to particle momenta.

Parameters:

thermalize_particle_momenta assigns the selected particle’s velocities and angular momentum to random values drawn from a Gaussian distribution consistent with the given thermal energy kT.

Velocity

thermalize_particle_momenta assigns random velocities to the x and y components of each particle’s velocity. When the simulation box is 3D, it also assigns a random velocity to the z component. When the simulation box is 2D, it sets the z component to 0. Finally, sets the center of mass velocity of the selected particles to 0.

Angular momentum

thermalize_particle_momenta assigns random angular momenta to each rotational degree of freedom that has a non-zero moment of intertia. Each particle can have 0, 1, 2, or 3 rotational degrees of freedom as determine by its moment of inertia.

Example:

simulation.state.thermalize_particle_momenta(
    filter=hoomd.filter.All(),
    kT=1.5)
property types#

dictionary of all types in the state.

Combines the data from particle_types, bond_types, angle_types, dihedral_types, improper_types, and special_pair_types into a dictionary with keys matching the property names.

Type:

dict[str, list[str]]

update_group_dof()#

Schedule an update to the number of degrees of freedom in each group.

update_group_dof requests that Simulation update the degrees of freedom provided to each group by the Integrator. Simulation will perform this update at the start of Simulation.run or at the start of the next timestep during an ongoing call to Simulation.run.

This method is called automatically when:

Call update_group_dof manually to force an update, such as when you modify particle moments of inertia with cpu_local_snapshot.

Example:

box = simulation.state.update_group_dof()

Modules