hoomd.update

Overview

BoxResize

Resizes the box between an initial and final box.

CustomUpdater

User-defined updater.

FilterUpdater

Update sets of particles associated with a filter.

RemoveDrift

Remove the average drift from a system restrained on a lattice.

Details

Updaters.

class hoomd.update.BoxResize(trigger, box1, box2, variant, filter=hoomd._hoomd.ParticleFilter)

Resizes the box between an initial and final box.

When part of a hoomd.Simulation updater list, this object will resize the box between the initial and final boxes passed. The behavior is a linear interpolation between the initial and final boxes where the minimum of the variant is tagged to box1 and the maximum is tagged to box2. All values between the minimum and maximum result in a box that is the interpolation of the three lengths and tilt factors of the initial and final boxes.

Note

The passed Variant must be bounded (i.e. it must have a true minimum and maximum) or the behavior of the updater is undefined.

Note

Currently for MPI simulations the rescaling of particles does not work properly in HPMC.

Parameters
box1

The box associated with the minimum of the passed variant.

Type

hoomd.Box

box2

The box associated with the maximum of the passed variant.

Type

hoomd.Box

variant

A variant used to interpolate between the two boxes.

Type

hoomd.variant.Variant

trigger

The trigger to activate this updater.

Type

hoomd.trigger.Trigger

filter

The subset of particles to update.

Type

hoomd.filter.ParticleFilter

get_box(timestep)

Get the box for a given timestep.

Parameters

timestep (int) – The timestep to use for determining the resized box.

Returns

The box used at the given timestep. None before the first call to Simulation.run.

Return type

Box

static update(state, box, filter=hoomd._hoomd.ParticleFilter)

Immediately scale the particle in the system state to the given box.

Parameters
class hoomd.update.CustomUpdater(trigger, action)

User-defined updater.

Parameters

CustomUpdater is a hoomd.operation.Updater that wraps a user-defined hoomd.custom.Action object so the action can be added to a hoomd.Operations instance for use with hoomd.Simulation objects.

Updaters modify the system state.

class hoomd.update.FilterUpdater(trigger, filters)

Update sets of particles associated with a filter.

HOOMD caches the particles selected by hoomd.filter.ParticleFilter instances to avoid the cost of re-running the filter on every time step. This means that unless the particles selected by a filter are recomputed the set of particles an operation works on is static. This class provides a mechanism to update the cached list of particles periodically. For example, use it to update the particles operated on by an MD integration method so that the integration method applies to particles in a given region of space.

Note

If needed to improve performance, use a hoomd.trigger.Trigger subclass, to update only when there is a known change to the particles that a filter would select.

Note

Some actions automatically recompute all filter particles such as adding or removing particles.

Parameters
trigger

The trigger associated with the updater.

Type

hoomd.trigger.Trigger

__eq__(other)

Return whether two objects are equivalent.

property filters

filters to update select particles.

Type

list[hoomd.filter.ParticleFilter]

class hoomd.update.RemoveDrift(reference_positions, trigger=1)

Remove the average drift from a system restrained on a lattice.

Parameters
  • reference_positions ((N_particles, 3) – float): the reference positions of the lattice \([\mathrm{length}]\).

  • trigger (hoomd.trigger.Trigger) – Select the timesteps to remove drift.

During the time steps specified by trigger, the mean drift \(\Delta\vec{r}\) from the reference_positions (\(\vec{r}_{ref, i}\)) is substracted from the current particle positions (\(\vec{r}_i\)). The drift is then given is given by:

\[\Delta\vec{r} = \frac{1}{\mathrm{N_{particles}}} \sum_{i=1}^\mathrm{N_{particles}} \mathrm{min\_image}(\vec{r}_i - \vec{r}_{ref,i})\]