Simulation
Details
Base simulation class for flowerMD.
- class flowermd.base.simulation.Simulation(*args: Any, **kwargs: Any)
The simulation context management class.
This class takes the output of the Initialization class and sets up a hoomd-blue simulation.
- Parameters:
initial_state (hoomd.snapshot.Snapshot or str, required) – A snapshot to initialize a simulation from, or a path to a GSD file to initialize a simulation from.
forcefield (List of hoomd.md.force.Force, required) – List of HOOMD force objects to add to the integrator.
reference_values (dict, default {}) – A dictionary of reference values for mass, length, and energy.
dt (float, default 0.0001) – Initial value for dt, the size of simulation timestep.
device (hoomd.device, default hoomd.device.auto_select()) – The CPU or GPU device to use for the simulation.
seed (int, default 42) – Seed passed to integrator when randomizing velocities.
gsd_write (int, default 1e4) – Period to write simulation snapshots to gsd file.
gsd_file_name (str, default "trajectory.gsd") – The file name to use for the GSD file
gsd_max_buffer_size (int, default 64 * 1024 * 1024) – Size (in bytes) to buffer in memory before writing GSD to file.
log_write (int, default 1e3) – Period to write simulation data to the log file.
log_file_name (str, default "sim_data.txt") – The file name to use for the .txt log file
thermostat (flowermd.utils.HOOMDThermostats, default) – HOOMDThermostats.MTTK The thermostat to use for the simulation.
- add_force(hoomd_force)
Add a force to the simulation.
- Parameters:
hoomd_force (hoomd.md.force.Force, required) – The force to add to the simulation.
- add_walls(wall_axis, sigma, epsilon, r_cut, r_extrap=0)
Add hoomd.md.external.wall.LJ forces to the simulation.
- Parameters:
wall_axis (np.ndarray, shape=(3,), dtype=float, required) – The axis of the wall in (x, y, z) order.
sigma (float, required) – The sigma parameter of the LJ wall.
epsilon (float, required) – The epsilon parameter of the LJ wall.
r_cut (float, required) – The cutoff radius of the LJ wall.
r_extrap (float, default 0) – The extrapolation radius of the LJ wall.
- adjust_epsilon(scale_by=None, shift_by=None, type_filter=None)
Adjust the epsilon parameter of the Lennard-Jones pair force.
- Parameters:
scale_by (float, default None) – The factor to scale epsilon by.
shift_by (float, default None) – The amount to shift epsilon by.
type_filter (list of str, default None) – A list of particle pair types to apply the adjustment to.
- adjust_sigma(scale_by=None, shift_by=None, type_filter=None)
Adjust the sigma parameter of the Lennard-Jones pair force.
- Parameters:
scale_by (float, default None) – The factor to scale sigma by.
shift_by (float, default None) – The amount to shift sigma by.
type_filter (list of str, default None) – A list of particle pair types to apply the adjustment to.
- property box_lengths
The simulation box lengths.
If reference length is set, the box lengths are scaled back to the non-reduced values.
- property box_lengths_reduced
The simulation box lengths in reduced units.
- property density
The density of the system.
- property density_reduced
The density of the system in reduced units.
- property dt
The simulation timestep.
- flush_writers()
Flush all write buffers to file.
- property forces
The list of forces in the simulation.
- classmethod from_snapshot_forces(initial_state, forcefield, **kwargs)
Initialize a simulation from an initial state and HOOMD forces.
- Parameters:
initial_state (gsd.hoomd.Snapshot or str) – A snapshot to initialize a simulation from, or a path to a GSD file.
forcefield (List of HOOMD force objects, required) – List of HOOMD force objects to add to the integrator.
- classmethod from_system(system, **kwargs)
Initialize a simulation from a flowermd.base.System object.
- Parameters:
system (flowermd.base.System, required) – A flowermd.base.System object.
- property integrate_group
The group of particles to apply the integrator to.
Default is all particles.
- property mass
The total mass of the system.
If reference mass is set, the mass is scaled back to the non-reduced value.
- property mass_reduced
The total mass of the system in reduced units.
- property method
The integrator method used by the simulation.
- property nlist
The neighbor list used by the Lennard-Jones pair force.
- pickle_forcefield(file_path='forcefield.pickle')
Pickle the list of HOOMD forces.
This method useful for saving the forcefield of a simulation to a file and reusing it for restarting a simulation or running a different simulation.
- Parameters:
file_path (str, default "forcefield.pickle") – The path to save the pickle file to.
Examples
In this example, a simulation is initialized and run for 1000 steps. The forcefield is then pickled and saved to a file. The forcefield is then loaded from the pickle file and used to run a tensile simulation.
from flowermd.base import Pack, Simulation from flowermd.library import PPS, OPLS_AA_PPS, Tensile import pickle pps_mols = PPS(num_mols=10, lengths=5) pps_system = Pack(molecules=[pps_mols], force_field=OPLS_AA_PPS(), r_cut=2.5, density=0.5, auto_scale=True, scale_charges=True) sim = Simulation(initial_state=pps_system.hoomd_snapshot, forcefield=pps_system.hoomd_forcefield) sim.run_NVT(n_steps=1e3, kT=1.0, tau_kt=1.0) sim.pickle_forcefield("pps_forcefield.pickle") with open("pps_forcefield.pickle", "rb") as f: pps_forcefield = pickle.load(f) tensile_sim = Tensile(initial_state=pps_system.hoomd_snapshot, forcefield=pps_forcefield, tensile_axis=(1, 0, 0)) tensile_sim.run_tensile(strain=0.05, kT=2.0, n_steps=1e3, period=10)
- property real_timestep
The simulation timestep in real units.
- property reference_energy
The reference energy for the simulation.
- property reference_length
The reference length for the simulation.
- property reference_mass
The reference mass for the simulation.
- property reference_values
The reference values for the simulation in form of a dictionary.
- remove_force(hoomd_force)
Remove a force from the simulation.
- Parameters:
hoomd_force (hoomd.md.force.Force, required) – The force to remove from the simulation.
- remove_walls(wall_axis)
Remove LJ walls from the simulation.
- Parameters:
wall_axis (np.ndarray, shape=(3,), dtype=float, required) – The axis of the wall in (x, y, z) order.
- run_NPT(n_steps, kT, pressure, tau_kt, tau_pressure, couple='xyz', box_dof=[True, True, True, False, False, False], rescale_all=False, gamma=0.0, thermalize_particles=True, write_at_start=True)
Run the simulation in the NPT ensemble.
- Parameters:
n_steps (int, required) – Number of steps to run the simulation.
kT (int or hoomd.variant.Ramp, required) – The temperature to use during the simulation.
pressure (int or hoomd.variant.Ramp, required) – The pressure to use during the simulation.
tau_kt (float, required) – Thermostat coupling period (in simulation time units).
tau_pressure (float, required) – Barostat coupling period.
couple (str, default "xyz") – Couplings of diagonal elements of the stress tensor/
box_dof (list of bool;) – optional default [True, True, True, False, False, False] Degrees of freedom of the box.
rescale_all (bool, default False) – Rescale all particles, not just those in the group.
gamma (float, default 0.0) – Friction constant for the box degrees of freedom,
thermalize_particles (bool, default True) – When set to True, assigns random velocities to all particles.
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
- run_NVE(n_steps, write_at_start=True)
Run the simulation in the NVE ensemble.
- Parameters:
n_steps (int, required) – Number of steps to run the simulation.
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
- run_NVT(n_steps, kT, tau_kt, thermalize_particles=True, write_at_start=True)
Run the simulation in the NVT ensemble.
- Parameters:
n_steps (int, required) – Number of steps to run the simulation.
kT (int or hoomd.variant.Ramp, required) – The temperature to use during the simulation.
tau_kt (float, required) – Thermostat coupling period (in simulation time units).
thermalize_particles (bool, default True) – When set to True, assigns random velocities to all particles.
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
- run_displacement_cap(n_steps, maximum_displacement=0.001, write_at_start=True)
NVE integrator with a cap on the maximum displacement per time step.
DisplacementCapped method is mostly useful for initially relaxing a system with overlapping particles. Putting a cap on the max particle displacement prevents Hoomd Particle Out of Box execption. Once the system is relaxed, other run methods (NVE, NVT, etc) can be used.
- Parameters:
n_steps (int, required) – Number of steps to run the simulation.
maximum_displacement (float, default 1e-3) – Maximum displacement per step (length)
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
- run_langevin(n_steps, kT, tally_reservoir_energy=False, default_gamma=1.0, default_gamma_r=(1.0, 1.0, 1.0), thermalize_particles=True, write_at_start=True)
Run the simulation using the Langevin dynamics integrator.
- Parameters:
n_steps (int, required) – Number of steps to run the simulation.
kT (int or hoomd.variant.Ramp, required) – The temperature to use during the simulation.
tally_reservoir_energy (bool, default False) –
- When set to True, energy exchange between the thermal reservoir
and the particles is tracked.
default_gamma (float, default 1.0) – The default drag coefficient to use for all particles.
default_gamma_r (tuple of floats, default (1.0, 1.0, 1.0)) – The default rotational drag coefficient to use for all particles.
thermalize_particles (bool, default True) – When set to True, assigns random velocities to all particles.
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
- run_update_volume(final_box_lengths, n_steps, period, kT, tau_kt, thermalize_particles=True, write_at_start=True)
Run an NVT simulation while shrinking or expanding simulation box.
The simulation box is updated using hoomd.update.BoxResize and the final box lengths are set to final_box_lengths.
See flowermd.utils.get_target_volume_mass_density and flowermd.utils.get_target_volume_number_density which are helper functions that can be used to get final_box_lengths.
- Parameters:
final_box_lengths (np.ndarray or unyt.array.unyt_array, shape=(3,), required # noqa: E501) – The final box edge lengths in (x, y, z) order.
n_steps (int, required) – Number of steps to run during volume update.
period (int, required) – The number of steps ran between each box update iteration.
kT (float or hoomd.variant.Ramp, required) – The temperature to use during volume update.
tau_kt (float, required) – Thermostat coupling period (in simulation time units).
write_at_start (bool, default True) – When set to True, triggers writers that evaluate to True for the initial step to execute before the next simulation time step.
Examples
In this example, a low density system is initialized with Pack and a box matching a density of 1.1 g/cm^3 is passed into final_box_lengths.
import unyt from flowermd.base import Pack, Simulation from flowermd.library import PPS, OPLS_AA_PPS pps_mols = PPS(num_mols=20, lengths=15) pps_system = Pack( molecules=[pps_mols], force_field=OPLS_AA_PPS(), r_cut=2.5, density=0.5, auto_scale=True, scale_charges=True ) sim = Simulation( initial_state=pps_system.hoomd_snapshot, forcefield=pps_system.hoomd_forcefield ) target_box = flowermd.utils.get_target_box_mass_density( density=1.1 * unyt.g/unyt.cm**3, mass=sim.mass.to("g") ) sim.run_update_volume( n_steps=1e4, kT=1.0, tau_kt=1.0, final_box_lengths=target_box )
- save_restart_gsd(file_path='restart.gsd')
Save a GSD file of the current simulation state.
This method is useful for saving the state of a simulation to a file and reusing it for restarting a simulation or running a different simulation.
- Parameters:
file_path (str, default "restart.gsd") – The path to save the GSD file to.
Examples
This example is similar to the example in pickle_forcefield. The only difference is that the simulation state is also saved to a GSD file.
from flowermd.base import Pack, Simulation from flowermd.library import PPS, OPLS_AA_PPS, Tensile import pickle pps_mols = PPS(num_mols=10, lengths=5) pps_system = Pack(molecules=[pps_mols], force_field=OPLS_AA_PPS(), r_cut=2.5, density=0.5, auto_scale=True, scale_charges=True) sim = Simulation(initial_state=pps_system.hoomd_snapshot, forcefield=pps_system.hoomd_forcefield) sim.run_NVT(n_steps=1e3, kT=1.0, tau_kt=1.0) sim.pickle_forcefield("pps_forcefield.pickle") sim.save_restart_gsd("pps_restart.gsd") with open("pps_forcefield.pickle", "rb") as f: pps_forcefield = pickle.load(f) tensile_sim = Tensile(initial_state="pps_restart.gsd", forcefield=pps_forcefield, tensile_axis=(1, 0, 0)) tensile_sim.run_tensile(strain=0.05, kT=2.0, n_steps=1e3, period=10)
- set_integrator_method(integrator_method, method_kwargs)
Create an initial (or updates the existing) integrator method.
This doesn’t need to be called directly; instead the various run functions use this method to update the integrator method as needed.
- Parameters:
integrrator_method (hoomd.md.method, required) – Instance of one of the hoomd.md.method options.
method_kwargs (dict, required) – A diction of parameter:value for the integrator method used.
- temperature_ramp(n_steps, kT_start, kT_final)
Create a temperature ramp.
- Parameters:
n_steps (int, required) – The number of steps to ramp the temperature over.
kT_start (float, required) – The starting temperature.
kT_final (float, required) – The final temperature.
- property thermostat
The thermostat used for the simulation.
- property volume
The simulation volume.
- property volume_reduced
The simulation volume in reduced units.