md.force

Overview

Force

Defines a force for molecular dynamics simulations.

Custom

Custom forces implemented in python.

Active

Active force.

ActiveOnManifold

Active force on a manifold.

Details

Apply forces to particles.

class hoomd.md.force.Force

Defines a force for molecular dynamics simulations.

A Force class computes the force and torque on each particle in the simulation state \(\vec{F}_i\) and \(\vec{\tau}_i\). With a few exceptions (noted in the documentation of the specific force classes), Force subclasses also compute the contribution to the system’s potential energy \(U\) and the the virial tensor \(W\). Force breaks the computation of the total system \(U\) and \(W\) into per-particle and additional terms as detailed in the documentation for each specific Force subclass.

\[\begin{split}U & = U_\mathrm{additional} + \sum_{i=0}^{N_\mathrm{particles}-1} U_i \\ W & = W_\mathrm{additional} + \sum_{i=0}^{N_\mathrm{particles}-1} W_i\end{split}\]

Force represents virial tensors as six element arrays listing the components of the tensor in this order:

\[(W^{xx}, W^{xy}, W^{xz}, W^{yy}, W^{yz}, W^{zz}).\]

The components of the virial tensor for a force on a single particle are:

\[W^{kl}_i = F^k \cdot r_i^l\]

where the superscripts select the x,y, and z components of the vectors. To properly account for periodic boundary conditions, pairwise interactions evaluate the virial:

\[W^{kl}_i = \sum_j F^k_{ij} \cdot \mathrm{minimum\_image}(\vec{r}_j - \vec{r}_i)^l\]

Note

Force is the base class for all molecular dynamics forces and provides common methods. Users should not instantiate this class directly.

property additional_energy

Additional energy term \(U_\mathrm{additional}\) \([\mathrm{energy}]\).

(Loggable: category=”scalar”)

Type

float

property additional_virial

Additional virial tensor term \(W_\mathrm{additional}\) \([\mathrm{energy}]\).

(Loggable: category=”sequence”)

Type

(1, 6) numpy.ndarray of float

property cpu_local_force_arrays

Expose force arrays on the CPU.

Provides direct access to the force, potential energy, torque, and virial data of the particles in the system on the cpu through a context manager. All data is MPI rank-local.

The hoomd.md.data.ForceLocalAccess object returned by this property has four arrays through which one can modify the force data:

Note

The local arrays are read only for built-in forces. Use Custom to implement custom forces.

Examples:

with self.cpu_local_force_arrays as arrays:
    arrays.force[:] = ...
    arrays.potential_energy[:] = ...
    arrays.torque[:] = ...
    arrays.virial[:] = ...
Type

hoomd.md.data.ForceLocalAccess

property energies

Energy contribution \(U_i\) from each particle \([\mathrm{energy}]\).

Attention

In MPI parallel execution, the array is available on rank 0 only. energies is None on ranks >= 1.

(Loggable: category=”particle”)

Type

(N_particles, ) numpy.ndarray of float

property energy

The potential energy \(U\) of the system from this force \([\mathrm{energy}]\).

(Loggable: category=”scalar”)

Type

float

property forces

The force \(\vec{F}_i\) applied to each particle \([\mathrm{force}]\).

Attention

In MPI parallel execution, the array is available on rank 0 only. forces is None on ranks >= 1.

(Loggable: category=”particle”)

Type

(N_particles, 3) numpy.ndarray of float

property gpu_local_force_arrays

Expose force arrays on the GPU.

Provides direct access to the force, potential energy, torque, and virial data of the particles in the system on the gpu through a context manager. All data is MPI rank-local.

The hoomd.md.data.ForceLocalAccessGPU object returned by this property has four arrays through which one can modify the force data:

Note

The local arrays are read only for built-in forces. Use Custom to implement custom forces.

Examples:

with self.gpu_local_force_arrays as arrays:
    arrays.force[:] = ...
    arrays.potential_energy[:] = ...
    arrays.torque[:] = ...
    arrays.virial[:] = ...

Note

GPU local force data is not available if the chosen device for the simulation is hoomd.device.CPU.

Type

hoomd.md.data.ForceLocalAccessGPU

property torques

The torque \(\vec{\tau}_i\) applied to each particle \([\mathrm{force} \cdot \mathrm{length}]\).

Attention

In MPI parallel execution, the array is available on rank 0 only. torques is None on ranks >= 1.

(Loggable: category=”particle”)

Type

(N_particles, 3) numpy.ndarray of float

property virials

Virial tensor contribution \(W_i\) from each particle \([\mathrm{energy}]\).

Attention

To improve performance Force objects only compute virials when needed. When not computed, virials is None. Virials are computed on every step when using a md.methods.NPT or md.methods.NPH integrator, on steps where a writer is triggered (such as write.GSD which may log pressure or virials), or when Simulation.always_compute_pressure is True.

Attention

In MPI parallel execution, the array is available on rank 0 only. virials is None on ranks >= 1.

(Loggable: category=”particle”)

Type

(N_particles, 6) numpy.ndarray of float

class hoomd.md.force.Custom(aniso=False)

Custom forces implemented in python.

Derive a custom force class from Custom, and override the set_forces method to compute forces on particles. Users have direct, zero-copy access to the C++ managed buffers via either the cpu_local_force_arrays or gpu_local_force_arrays property. Choose the property that corresponds to the device you wish to alter the data on. In addition to zero-copy access to force buffers, custom forces have access to the local snapshot API via the _state.cpu_local_snapshot or the _state.gpu_local_snapshot property.

See also

See the documentation in hoomd.State for more information on the local snapshot API.

Examples:

class MyCustomForce(hoomd.force.Custom):
    def __init__(self):
        super().__init__(aniso=True)

    def set_forces(self, timestep):
        with self.cpu_local_force_arrays as arrays:
            arrays.force[:] = -5
            arrays.torque[:] = 3
            arrays.potential_energy[:] = 27
            arrays.virial[:] = np.arange(6)[None, :]

In addition, since data is MPI rank-local, there may be ghost particle data associated with each rank. To access this read-only ghost data, access the property name with either the prefix ghost_ of the suffix _with_ghost.

Note

Pass aniso=True to the md.force.Custom constructor if your custom force produces non-zero torques on particles.

Examples:

class MyCustomForce(hoomd.force.Custom):
    def __init__(self):
        super().__init__()

    def set_forces(self, timestep):
        with self.cpu_local_force_arrays as arrays:
            # access only the ghost particle forces
            ghost_force_data = arrays.ghost_force

            # access torque data on this rank and ghost torque data
            torque_data = arrays.torque_with_ghost

Note

When accessing the local force arrays, always use a context manager.

Note

The shape of the exposed arrays cannot change while in the context manager.

Note

All force data buffers are MPI rank local, so in simulations with MPI, only the data for a single rank is available.

Note

Access to the force buffers is constant (O(1)) time.

abstract set_forces(timestep)

Set the forces in the simulation loop.

Parameters

timestep (int) – The current timestep in the simulation.

class hoomd.md.force.Active(filter)

Bases: Force

Active force.

Parameters

filter (hoomd.filter) – Subset of particles on which to apply active forces.

Active computes an active force and torque on all particles selected by the filter:

\[\begin{split}\vec{F}_i = \mathbf{q}_i \vec{f}_i \mathbf{q}_i^* \\ \vec{\tau}_i = \mathbf{q}_i \vec{u}_i \mathbf{q}_i^*,\end{split}\]

where \(\vec{f}_i\) is the active force in the local particle coordinate system (set by type active_force) and \(\vec{u}_i\) is the active torque in the local particle coordinate system (set by type in active_torque.

Note

To introduce rotational diffusion to the particle orientations, use create_diffusion_updater.

Examples:

all = hoomd.filter.All()
active = hoomd.md.force.Active(
    filter=hoomd.filter.All()
    )
active.active_force['A','B'] = (1,0,0)
active.active_torque['A','B'] = (0,0,0)
rotational_diffusion_updater = active.create_diffusion_updater(
    trigger=10)
sim.operations += rotational_diffusion_updater

Note

The energy and virial associated with the active force are 0.

filter

Subset of particles on which to apply active forces.

Type

hoomd.filter

active_force

Active force vector in the local reference frame of the particle \([\mathrm{force}]\). It is defined per particle type and stays constant during the simulation.

Type: TypeParameter [particle_type, tuple [float, float, float]]

active_torque

Active torque vector in the local reference frame of the particle \([\mathrm{force} \cdot \mathrm{length}]\). It is defined per particle type and stays constant during the simulation.

Type: TypeParameter [particle_type, tuple [float, float, float]]

create_diffusion_updater(trigger, rotational_diffusion)

Create a rotational diffusion updater for this active force.

Parameters
Returns

The rotational diffusion updater.

Return type

hoomd.md.update.ActiveRotationalDiffusion

class hoomd.md.force.ActiveOnManifold(filter, manifold_constraint)

Bases: Active

Active force on a manifold.

Parameters

ActiveOnManifold computes a constrained active force and torque on all particles selected by the filter similar to Active. ActiveOnManifold restricts the forces to the local tangent plane of the manifold constraint. For more information see Active.

Hint

Use ActiveOnManifold with a md.methods.rattle integration method with the same manifold constraint.

Note

To introduce rotational diffusion to the particle orientations, use create_diffusion_updater. The rotational diffusion occurs in the local tangent plane of the manifold.

Examples:

all = filter.All()
sphere = hoomd.md.manifold.Sphere(r=10)
active = hoomd.md.force.ActiveOnManifold(
    filter=hoomd.filter.All(), rotation_diff=0.01,
    manifold_constraint = sphere
    )
active.active_force['A','B'] = (1,0,0)
active.active_torque['A','B'] = (0,0,0)
filter

Subset of particles on which to apply active forces.

Type

hoomd.filter.ParticleFilter

manifold_constraint

Manifold constraint.

Type

hoomd.md.manifold.Manifold

active_force

Active force vector in the local reference frame of the particle \([\mathrm{force}]\). It is defined per particle type and stays constant during the simulation.

Type: TypeParameter [particle_type, tuple [float, float, float]]

active_torque

Active torque vector in local reference frame of the particle \([\mathrm{force} \cdot \mathrm{length}]\). It is defined per particle type and stays constant during the simulation.

Type: TypeParameter [particle_type, tuple [float, float, float]]

create_diffusion_updater(trigger, rotational_diffusion)

Create a rotational diffusion updater for this active force.

Parameters
Returns

The rotational diffusion updater.

Return type

hoomd.md.update.ActiveRotationalDiffusion