hoomd.md#

Details

Molecular dynamics.

In molecular dynamics simulations, HOOMD-blue numerically integrates the degrees of freedom in the system as a function of time under the influence of forces. To perform MD simulations, assign a MD Integrator to the hoomd.Simulation operations. Provide the Integrator with lists of integration methods, forces, and constraints to apply during the integration. Use hoomd.md.minimize.FIRE to perform energy minimization.

MD updaters (hoomd.md.update) perform additional operations during the simulation, including rotational diffusion and establishing shear flow. Use MD computes (hoomd.md.compute) to compute the thermodynamic properties of the system state.

See also

Tutorial: Introducing Molecular Dynamics

class hoomd.md.HalfStepHook(*args: Any, **kwargs: Any)#

HalfStepHook base class.

HalfStepHook provides an interface to perform computations during the half-step of a hoomd.md.Integrator.

update(timestep)#

Called during the half-step of a hoomd.md.Integrator.

This method should provide the implementation of any computation that the user wants to execute at each timestep in the middle of the integration routine.

class hoomd.md.Integrator(dt, integrate_rotational_dof=False, forces=None, constraints=None, methods=None, rigid=None, half_step_hook=None)#

Molecular dynamics time integration.

Parameters:
  • dt (float) – Integrator time step size \([\mathrm{time}]\).

  • methods (Sequence[hoomd.md.methods.Method]) – Sequence of integration methods. The default value of None initializes an empty list.

  • forces (Sequence[hoomd.md.force.Force]) – Sequence of forces applied to the particles in the system. The default value of None initializes an empty list.

  • integrate_rotational_dof (bool) – When True, integrate rotational degrees of freedom.

  • constraints (Sequence[hoomd.md.constrain.Constraint]) – Sequence of constraint forces applied to the particles in the system. The default value of None initializes an empty list. Rigid body objects (i.e. hoomd.md.constrain.Rigid) are not allowed in the list.

  • rigid (hoomd.md.constrain.Rigid) – An object defining the rigid bodies in the simulation.

  • half_step_hook (hoomd.md.HalfStepHook) – Enables the user to perform arbitrary computations during the half-step of the integration.

Integrator is the top level class that orchestrates the time integration step in molecular dynamics simulations. The integration methods define the equations of motion to integrate under the influence of the given forces and constraints.

Each method applies the given equations of motion to a subset of particles in the simulation state. See the documentation for each method for details on what equations of motion it solves. The intersection of the subsets must be null.

Integrator computes the net force, torque, energy, and virial on each particle as a sum of those applied by hoomd.md.force.Force objects in the forces and constraints lists:

\[\begin{split}\vec{F}_{\mathrm{net},i} &= \sum_{f \in \mathrm{forces}} \vec{F}_i^f \\ \vec{\tau}_{\mathrm{net},i} &= \sum_{f \in \mathrm{forces}} \vec{\tau}_i^f \\ U_{\mathrm{net},i} &= \sum_{f \in \mathrm{forces}} U_i^f \\ W_{\mathrm{net},i} &= \sum_{f \in \mathrm{forces}} W_i^f \\\end{split}\]

Integrator also computes the net additional energy and virial

\[\begin{split}U_{\mathrm{net},\mathrm{additional}} &= \sum_{f \in \mathrm{forces}} U_\mathrm{additional}^f \\ W_{\mathrm{net},\mathrm{additional}} &= \sum_{f \in \mathrm{forces}} W_\mathrm{additional}^f \\\end{split}\]

See md.force.Force for definitions of these terms. Constraints are a special type of force used to enforce specific constraints on the system state, such as distances between particles with hoomd.md.constrain.Distance. Integrator handles rigid bodies as a special case, as it only integrates the degrees of freedom of each body’s center of mass. See hoomd.md.constrain.Rigid for details.

Degrees of freedom

Integrator always integrates the translational degrees of freedom. It optionally integrates one or more rotational degrees of freedom for a given particle i when all the following conditions are met:

  • The intergration method supports rotational degrees of freedom.

  • integrate_rotational_dof is True.

  • The moment of inertia is non-zero \(I^d_i > 0\).

Each particle may have zero, one, two, or three rotational degrees of freedom.

Note

By default, integrate_rotational_dof is False. gsd and hoomd.Snapshot also set particle moments of inertia to 0 by default.

Classes

Classes of the following modules can be used as elements in methods:

The classes of following modules can be used as elements in forces:

The classes of the following module can be used as elements in constraints:

Examples:

nlist = hoomd.md.nlist.Cell()
lj = hoomd.md.pair.LJ(nlist=nlist)
lj.params.default = dict(epsilon=1.0, sigma=1.0)
lj.r_cut[('A', 'A')] = 2**(1/6)
nve = hoomd.md.methods.NVE(filter=hoomd.filter.All())
integrator = hoomd.md.Integrator(dt=0.001, methods=[nve], forces=[lj])
sim.operations.integrator = integrator
dt#

Integrator time step size \([\mathrm{time}]\).

Type:

float

methods#

List of integration methods.

Type:

list[hoomd.md.methods.Method]

forces#

List of forces applied to the particles in the system.

Type:

list[hoomd.md.force.Force]

integrate_rotational_dof#

When True, integrate rotational degrees of freedom.

Type:

bool

constraints#

List of constraint forces applied to the particles in the system.

Type:

list[hoomd.md.constrain.Constraint]

rigid#

The rigid body definition for the simulation associated with the integrator.

Type:

hoomd.md.constrain.Rigid

half_step_hook#

User defined implementation to perform computations during the half-step of the integration.

Type:

hoomd.md.HalfStepHook

__setattr__(attr, value)#

Hande group DOF update when setting integrate_rotational_dof.

property linear_momentum#

The linear momentum vector of the system \([\mathrm{mass} \cdot \mathrm{velocity}]\).

\[\vec{p} = \sum_{i=0}^\mathrm{N_particles-1} m_i \vec{v}_i\]

(Loggable: category=”sequence”)

Type:

tuple(float,float,float)

Modules