BoxResize¶
- class hoomd.update.BoxResize(trigger, box, filter=hoomd._hoomd.ParticleFilter)¶
Bases:
Updater
Vary the simulation box size as a function of time.
- Parameters:
trigger (hoomd.trigger.trigger_like) – The trigger to activate this updater.
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.
BoxResize
resizes the simulation box as a function of time. For each particle \(i\) matched byfilter
,BoxResize
scales 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
box
evaluated 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,
BoxResize
wraps all particles \(j\) back into the new box:\[\vec{r_j} \leftarrow \mathrm{minimum\_image}_{\vec{a}_k}' (\vec{r}_j)\]Warning
Rescaling particles in HPMC simulations with hard particles may introduce overlaps.
Note
When using rigid bodies, ensure that the
BoxResize
updater is last in the operations updater list. Immediately after theBoxResize
updater triggers, rigid bodies (hoomd.md.constrain.Rigid
) will be temporarily deformed.hoomd.md.Integrator
will run after the last updater and resets the constituent particle positions before computing forces.Example:
box_resize = hoomd.update.BoxResize( trigger=hoomd.trigger.Periodic(10), box=inverse_volume_ramp, ) simulation.operations.updaters.append(box_resize)
Members inherited from
AutotunedObject
:- property kernel_parameters¶
Kernel parameters.
Read more...
- property is_tuning_complete¶
Check if kernel parameter tuning is complete.
Read more...
- tune_kernel_parameters()¶
Start tuning kernel parameters.
Read more...
Members inherited from
Integrator
:- trigger¶
The trigger to activate this operation.
Read more...
Members defined in
BoxResize
:- 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:
- 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 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)