hoomd.update¶
Overview
Vary the simulation box size as a function of time. |
|
User-defined updater. |
|
Update sets of particles associated with a filter. |
|
Remove the average drift from a restrained system. |
Details
Updaters.
Updater operations modify the simulation state when they act.
- class hoomd.update.BoxResize(trigger, box1=None, box2=None, variant=None, filter=hoomd._hoomd.ParticleFilter, box=None)¶
Bases:
UpdaterVary the simulation box size as a function of time.
- Parameters:
trigger (hoomd.trigger.trigger_like) – The trigger to activate this updater.
box1 (hoomd.box.box_like) – The box associated with the minimum of the passed variant.
box2 (hoomd.box.box_like) – The box associated with the maximum of the passed variant.
variant (hoomd.variant.variant_like) – A variant used to interpolate between the two boxes.
filter (hoomd.filter.filter_like) – The subset of particle positions to update (defaults to
hoomd.filter.All).box (hoomd.variant.box.BoxVariant) – Box as a function of time.
BoxResizeresizes the simulation box as a function of time. For each particle \(i\) matched byfilter,BoxResizescales the particle to fit in the new box:\[\vec{r}_i \leftarrow s_x \vec{a}_1' + s_y \vec{a}_2' + s_z \vec{a}_3' - \frac{\vec{a}_1' + \vec{a}_2' + \vec{a}_3'}{2}\]where \(\vec{a}_k'\) are the new box vectors determined by
boxevaluated at the current time step and the scale factors are determined by the current particle position \(\vec{r}_i\) and the previous box vectors \(\vec{a}_k\):\[\vec{r}_i = s_x \vec{a}_1 + s_y \vec{a}_2 + s_z \vec{a}_3 - \frac{\vec{a}_1 + \vec{a}_2 + \vec{a}_3}{2}\]After scaling particles that match the filter,
BoxResizewraps all particles \(j\) back into the new box:\[\vec{r_j} \leftarrow \mathrm{minimum\_image}_{\vec{a}_k}' (\vec{r}_j)\]Note
For backward compatibility, you may set
box1,box2, andvariantwhich is equivalent to:box = hoomd.variant.box.Interpolate(initial_box=box1, final_box=box2, variant=variant)
Warning
Rescaling particles in HPMC simulations with hard particles may introduce overlaps.
Note
When using rigid bodies, ensure that the
BoxResizeupdater is last in the operations updater list. Immediately after theBoxResizeupdater triggers, rigid bodies (hoomd.md.constrain.Rigid) will be temporarily deformed.hoomd.md.Integratorwill run after the last updater and resets the constituent particle positions before computing forces.Deprecated since version 4.6.0: The arguments
box1,box2, andvariantare deprecated. Usebox.Example:
box_resize = hoomd.update.BoxResize(trigger=hoomd.trigger.Periodic(10), box=inverse_volume_ramp) simulation.operations.updaters.append(box_resize)
- box¶
The box as a function of time.
Example:
box_resize.box = inverse_volume_ramp
- filter¶
The subset of particles to update.
Example:
filter_ = box_resize.filter
- Type:
- property box1¶
The box associated with the minimum of the passed variant.
Deprecated since version 4.6.0: Use
box.- Type:
- property box2¶
The box associated with the maximum of the passed variant.
Deprecated since version 4.6.0: Use
box.- Type:
- 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.
Nonebefore the first call toSimulation.run.- Return type:
Example:
box = box_resize.get_box(1_000_000)
- static update(state, box, filter=hoomd._hoomd.ParticleFilter)¶
Immediately scale the particle in the system state to the given box.
- Parameters:
state (State) – System state to scale.
box (hoomd.box.box_like) – New box.
filter (hoomd.filter.filter_like) – The subset of particles to update (defaults to
hoomd.filter.All).
Example:
hoomd.update.BoxResize.update(state=simulation.state, box=box)
- class hoomd.update.CustomUpdater(trigger, action)¶
Bases:
CustomOperation,UpdaterUser-defined updater.
- Parameters:
action (hoomd.custom.Action) – The action to call.
trigger (hoomd.trigger.trigger_like) – Select the timesteps to call the action.
CustomUpdateris ahoomd.operation.Updaterthat wraps a user-definedhoomd.custom.Actionobject so the action can be added to ahoomd.Operationsinstance for use withhoomd.Simulationobjects.Updaters modify the system state.
Example:
custom_updater = hoomd.update.CustomUpdater( action=custom_action, trigger=hoomd.trigger.Periodic(1000)) simulation.operations.updaters.append(custom_updater)
See also
The base class
hoomd.custom.CustomOperation.
- class hoomd.update.FilterUpdater(trigger, filters)¶
Bases:
UpdaterUpdate sets of particles associated with a filter.
- Parameters:
trigger (hoomd.trigger.trigger_like) – A trigger to use for determining when to update particles associated with a filter.
filters (list[hoomd.filter.filter_like]) – A list of
hoomd.filter.filter_likeobjects to update.
hoomd.Simulationcaches the particles selected byhoomd.filter.filter_likeobjects to avoid the cost of re-running the filter on every time step. The particles selected by a filter will remain static until recomputed. This class provides a mechanism to update the cached list of particles. For example, periodically update a MD integration method’s group so that the integration method applies to particles in a given region of space.Tip
To improve performance, use a
hoomd.trigger.Triggersubclass, 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 to the
hoomd.Simulation.state.Example:
filter_updater = hoomd.update.FilterUpdater( trigger=hoomd.trigger.Periodic(1_000), filters=[filter1, filter2])
- __eq__(other)¶
Return whether two objects are equivalent.
- property filters¶
filters to update select particles.
Example:
filter_updater.filters = [filter1, filter2]
- Type:
- class hoomd.update.RemoveDrift(reference_positions, trigger=1)¶
Bases:
UpdaterRemove the average drift from a restrained system.
- Parameters:
reference_positions ((N_particles, 3)
numpy.ndarrayoffloat) – the reference positions \([\mathrm{length}]\).trigger (hoomd.trigger.trigger_like) – Select the timesteps to remove drift.
RemoveDriftcomputes the mean drift \(\vec{D}\) from the givenreference_positions(\(\vec{r}_{ref, i}\)):\[\vec{D} = \frac{1}{\mathrm{N_{particles}}} \sum_{i=0}^\mathrm{N_{particles-1}} \mathrm{minimum\_image}(\vec{r}_i - \vec{r}_{ref,i})\]RemoveDriftthen shifts all particles in the system by \(-\vec{D}\):\[\vec{r}_i \leftarrow \mathrm{minimum\_image}(\vec{r}_i - \vec{D})\]Tip
Use
RemoveDriftwithhoomd.hpmc.external.field.Harmonicto improve the accuracy of Frenkel-Ladd calculations.Example:
remove_drift = hoomd.update.RemoveDrift( reference_positions=[(0,0,0), (1,0,0)]) simulation.operations.updaters.append(remove_drift)
- reference_positions¶
the reference positions \([\mathrm{length}]\).
Example:
remove_drift.reference_positions = [(0,0,0), (1,0,0)]
- Type:
(N_particles, 3)
numpy.ndarrayoffloat