md.constrain

Overview

md.constrain.distance Constrain pairwise particle distances.
md.constrain.rigid Constrain particles in rigid bodies.
md.constrain.sphere Constrain particles to the surface of a sphere.
md.constrain.oneD Constrain particles to move along a specific direction only

Details

Constraints.

Constraint forces constrain a given set of particle to a given surface, to have some relative orientation, or impose some other type of constraint. For example, a group of particles can be constrained to the surface of a sphere with sphere.

As with other force commands in hoomd, multiple constrain commands can be issued to specify multiple constraints, which are additively applied.

Warning

Constraints will be invalidated if two separate constraint commands apply to the same particle.

The degrees of freedom removed from the system by constraints are correctly taken into account when computing the temperature for thermostatting and logging.

class hoomd.md.constrain.distance

Constrain pairwise particle distances.

distance specifies that forces will be applied to all particles pairs for which constraints have been defined.

The constraint algorithm implemented is described in:

  • [1] M. Yoneya, H. J. C. Berendsen, and K. Hirasawa, “A Non-Iterative Matrix Method for Constraint Molecular Dynamics Simulations,” Mol. Simul., vol. 13, no. 6, pp. 395–405, 1994.
  • [2] M. Yoneya, “A Generalized Non-iterative Matrix Method for Constraint Molecular Dynamics Simulations,” J. Comput. Phys., vol. 172, no. 1, pp. 188–197, Sep. 2001.

In brief, the second derivative of the Lagrange multipliers with respect to time is set to zero, such that both the distance constraints and their time derivatives are conserved within the accuracy of the Velocity Verlet scheme, i.e. within \(\Delta t^2\). The corresponding linear system of equations is solved. Because constraints are satisfied at \(t + 2 \Delta t\), the scheme is self-correcting and drifts are avoided.

Warning

In MPI simulations, all particles connected through constraints will be communicated between processors as ghost particles. Therefore, it is an error when molecules defined by constraints extend over more than half the local domain size.

Caution

constrain.distance() does not currently interoperate with integrate.brownian() or integrate.langevin()

Example:

constrain.distance()
disable()

Disable the force.

Example:

force.disable()

Executing the disable command removes the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable()

enable()

Enable the force.

Example:

force.enable()

See disable().

set_params(rel_tol=None)

Set parameters for constraint computation.

Parameters:rel_tol (float) – The relative tolerance with which constraint violations are detected (optional).

Example:

dist = constrain.distance()
dist.set_params(rel_tol=0.0001)
class hoomd.md.constrain.oneD(group, constraint_vector=[0, 0, 1])

Constrain particles to move along a specific direction only

Parameters:
  • group (hoomd.group) – Group on which to apply the constraint.
  • constraint_vector (list) – [x,y,z] list indicating the direction that the particles are restricted to

oneD specifies that forces will be applied to all particles in the given group to constrain them to only move along a given vector.

Example:

constrain.oneD(group=groupA, constraint_vector=[1,0,0])

New in version 2.1.

disable()

Disable the force.

Example:

force.disable()

Executing the disable command removes the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable()

enable()

Enable the force.

Example:

force.enable()

See disable().

class hoomd.md.constrain.rigid

Constrain particles in rigid bodies.

Overview

Rigid bodies are defined by a single central particle and a number of constituent particles. All of these are particles in the HOOMD system configuration and can interact with other particles via force computes. The mass and moment of inertia of the central particle set the full mass and moment of inertia of the rigid body (constituent particle mass is ignored).

The central particle is at the center of mass of the rigid body and the orientation quaternion defines the rotation from the body space into the simulation box. In body space, the center of mass of the body is at 0,0,0 and the moment of inertia is diagonal. You specify the constituent particles to rigid for each type of body in body coordinates. Then, rigid takes control of those particles, and sets their position and orientation in the simulation box relative to the position and orientation of the central particle. rigid also transfers forces and torques from constituent particles to the central particle. Then, MD integrators can use these forces and torques to integrate the equations of motion of the central particles (representing the whole rigid body) forward in time.

Defining bodies

rigid accepts one local body environment per body type. The type of a body is the particle type of the central particle in that body. In this way, each particle of type R in the system configuration defines a body of type R.

As a convenience, you do not need to create placeholder entries for all of the constituent particles in your initial configuration. You only need to specify the positions and orientations of all the central particles. When you call create_bodies(), it will create all constituent particles that do not exist. (those that already exist e.g. in a restart file are left unchanged).

Example that creates rigid rods:

# Place the type R central particles
uc = hoomd.lattice.unitcell(N = 1,
                            a1 = [10.8, 0,   0],
                            a2 = [0,    1.2, 0],
                            a3 = [0,    0,   1.2],
                            dimensions = 3,
                            position = [[0,0,0]],
                            type_name = ['R'],
                            mass = [1.0],
                            moment_inertia = [[0,
                                               1/12*1.0*8**2,
                                               1/12*1.0*8**2]],
                            orientation = [[1, 0, 0, 0]]);
system = hoomd.init.create_lattice(unitcell=uc, n=[2,18,18]);

# Add constituent particles of type A and create the rods
system.particles.types.add('A');
rigid = hoomd.md.constrain.rigid();
rigid.set_param('R',
                types=['A']*8,
                positions=[(-4,0,0),(-3,0,0),(-2,0,0),(-1,0,0),
                           (1,0,0),(2,0,0),(3,0,0),(4,0,0)]);

rigid.create_bodies()

Danger

Automatic creation of constituent particles can change particle tags. If bonds have been defined between particles in the initial configuration, or bonds connect to constituent particles, rigid bodies should be created manually.

When you create the constituent particles manually (i.e. in an input file or with snapshots), the central particle of a rigid body must have a lower tag than all of its constituent particles. Constituent particles follow in monotonically increasing tag order, corresponding to the order they were defined in the argument to set_param(). The order of central and contiguous particles need not to be contiguous. Additionally, you must set the body field for each of the particles in the rigid body to the tag of the central particle (for both the central and constituent particles). Set body to -1 for particles that do not belong to a rigid body. After setting an initial configuration that contains properly defined bodies and all their constituent particles, call validate_bodies() to verify that the bodies are defined and prepare the constraint.

You must call either create_bodies() or validate_bodies() prior to starting a simulation hoomd.run().

Integrating bodies

Most integrators in HOOMD support the integration of rotational degrees of freedom. When there are rigid bodies present in the system, do not apply integrators to the constituent particles, only the central and non-rigid particles.

Example:

rigid = hoomd.group.rigid_center();
hoomd.md.integrate.langevin(group=rigid, kT=1.0, seed=42);

Thermodynamic quantities of bodies

HOOMD computes thermodynamic quantities (temperature, kinetic energy, etc…) appropriately when there are rigid bodies present in the system. When it does so, it ignores all constituent particles and computes the translational and rotational energies of the central particles, which represent the whole body. hoomd.analyze.log can log the translational and rotational energy terms separately.

Restarting simulations with rigid bodies.

To restart, use hoomd.dump.gsd to write restart files. GSD stores all of the particle data fields needed to reconstruct the state of the system, including the body tag, rotational momentum, and orientation of the body. Restarting from a gsd file is equivalent to manual constituent particle creation. You still need to specify the same local body space environment to rigid as you did in the earlier simulation.

create_bodies(create=True)

Create copies of rigid bodies.

Parameters:create (bool) – When True, create rigid bodies, otherwise validate existing ones.
disable()

Disable the force.

Example:

force.disable()

Executing the disable command removes the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable()

enable()

Enable the force.

Example:

force.enable()

See disable().

set_param(type_name, types, positions, orientations=None, charges=None, diameters=None)

Set constituent particle types and coordinates for a rigid body.

Parameters:
  • type_name (str) – The type of the central particle
  • types (list) – List of types of constituent particles
  • positions (list) – List of relative positions of constituent particles
  • orientations (list) – List of orientations of constituent particles (optional)
  • charge (list) – List of charges of constituent particles (optional)
  • diameters (list) – List of diameters of constituent particles (optional)

Caution

The constituent particle type must be exist. If it does not exist, it can be created on the fly using system.particles.types.add('A_const') (see hoomd.data).

Example:

rigid = constrain.rigid()
rigid.set_param('A', types = ['A_const', 'A_const'], positions = [(0,0,1),(0,0,-1)])
rigid.set_param('B', types = ['B_const', 'B_const'], positions = [(0,0,.5),(0,0,-.5)])
validate_bodies()

Validate that bodies are well defined and prepare for the simulation run.

class hoomd.md.constrain.sphere(group, P, r)

Constrain particles to the surface of a sphere.

Parameters:
  • group (hoomd.group) – Group on which to apply the constraint.
  • P (tuple) – (x,y,z) tuple indicating the position of the center of the sphere (in distance units).
  • r (float) – Radius of the sphere (in distance units).

sphere specifies that forces will be applied to all particles in the given group to constrain them to a sphere. Currently does not work with Brownian or Langevin dynamics (hoomd.md.integrate.brownian and hoomd.md.integrate.langevin).

Example:

constrain.sphere(group=groupA, P=(0,10,2), r=10)
disable()

Disable the force.

Example:

force.disable()

Executing the disable command removes the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable()

enable()

Enable the force.

Example:

force.enable()

See disable().