md.force

Overview

md.force.active Active force.
md.force.constant Constant force.
md.force.dipole Treat particles as dipoles in an electric field.

Details

Apply forces to particles.

class hoomd.md.force.active(seed, f_lst, group, orientation_link=True, orientation_reverse_link=False, rotation_diff=0, constraint=None)

Active force.

Parameters:
  • seed (int) – required user-specified seed number for random number generator.
  • f_list (list) – An array of (x,y,z) tuples for the active force vector for each individual particle.
  • group (hoomd.group) – Group for which the force will be set
  • orientation_link (bool) – When True, particle orientation is coupled to the active force vector. Only relevant for non-point-like anisotropic particles.
  • orientation_reverse_link (bool) – When True, the active force vector is coupled to particle orientation. Useful for for using a particle’s orientation to log the active force vector. Quaternion rotation assumes base vector of (0,0,1).
  • rotation_diff (float) – rotational diffusion constant, \(D_r\), for all particles in the group.
  • constraint (hoomd.md.update.constraint_ellipsoid) – such as update.constraint_ellipsoid.

active specifies that an active force should be added to all particles. Obeys \(\delta {\bf r}_i = \delta t v_0 \hat{p}_i\), where \(v_0\) is the active velocity. In 2D \(\hat{p}_i = (\cos \theta_i, \sin \theta_i)\) is the active force vector for particle \(i\); and the diffusion of the active force vector follows \(\delta \theta / \delta t = \sqrt{2 D_r / \delta t} \Gamma\), where \(D_r\) is the rotational diffusion constant, and the gamma function is a unit-variance random variable, whose components are uncorrelated in time, space, and between particles. In 3D, \(\hat{p}_i\) is a unit vector in 3D space, and diffusion follows \(\delta \hat{p}_i / \delta t = \sqrt{2 D_r / \delta t} \Gamma (\hat{p}_i (\cos \theta - 1) + \hat{p}_r \sin \theta)\), where \(\hat{p}_r\) is an uncorrelated random unit vector. The persistence length of an active particle’s path is \(v_0 / D_r\).

Attention

active() does not support MPI execution.

Examples:

force.active( seed=13, f_list=[tuple(3,0,0) for i in range(N)])

ellipsoid = update.constraint_ellipsoid(group=groupA, P=(0,0,0), rx=3, ry=4, rz=5)
force.active( seed=7, f_list=[tuple(1,2,3) for i in range(N)], orientation_link=False, rotation_diff=100, constraint=ellipsoid)
disable(log=False)

Disable the force.

Parameters:log (bool) – Set to True if you plan to continue logging the potential energy associated with this force.

Examples:

force.disable()
force.disable(log=True)

Executing the disable command will remove the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable().

By setting log to True, the values of the force can be logged even though the forces are not applied in the simulation. For forces that use cutoff radii, setting log=True will cause the correct r_cut values to be used throughout the simulation, and therefore possibly drive the neighbor list size larger than it otherwise would be. If log is left False, the potential energy associated with this force will not be available for logging.

enable()

Enable the force.

Examples:

force.enable()

See disable().

get_energy(group)

Get the energy of a particle group.

Parameters:group (hoomd.group) – The particle group to query the energy for.
Returns:The last computed energy for the members in the group.

Examples:

g = group.all()
energy = force.get_energy(g)
class hoomd.md.force.constant(fx, fy, fz, group=None)

Constant force.

Parameters:
  • fx (float) – x-component of the force (in force units).
  • fy (float) – y-component of the force (in force units).
  • fz (float) – z-component of the force (in force units).
  • group (hoomd.group) – Group for which the force will be set.

constant specifies that a constant force should be added to every particle in the simulation or optionally to all particles in a group.

Examples:

force.constant(fx=1.0, fy=0.5, fz=0.25)
const = force.constant(fx=0.4, fy=1.0, fz=0.5)
const = force.constant(fx=0.4, fy=1.0, fz=0.5,group=fluid)
disable(log=False)

Disable the force.

Parameters:log (bool) – Set to True if you plan to continue logging the potential energy associated with this force.

Examples:

force.disable()
force.disable(log=True)

Executing the disable command will remove the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable().

By setting log to True, the values of the force can be logged even though the forces are not applied in the simulation. For forces that use cutoff radii, setting log=True will cause the correct r_cut values to be used throughout the simulation, and therefore possibly drive the neighbor list size larger than it otherwise would be. If log is left False, the potential energy associated with this force will not be available for logging.

enable()

Enable the force.

Examples:

force.enable()

See disable().

get_energy(group)

Get the energy of a particle group.

Parameters:group (hoomd.group) – The particle group to query the energy for.
Returns:The last computed energy for the members in the group.

Examples:

g = group.all()
energy = force.get_energy(g)
class hoomd.md.force.dipole(field_x, field_y, field_z, p)

Treat particles as dipoles in an electric field.

Parameters:
  • field_x (float) – x-component of the field (units?)
  • field_y (float) – y-component of the field (units?)
  • field_z (float) – z-component of the field (units?)
  • p (float) – magnitude of the particles’ dipole moment in the local z direction

Examples:

force.external_field_dipole(field_x=0.0, field_y=1.0 ,field_z=0.5, p=1.0)
const_ext_f_dipole = force.external_field_dipole(field_x=0.0, field_y=1.0 ,field_z=0.5, p=1.0)
disable(log=False)

Disable the force.

Parameters:log (bool) – Set to True if you plan to continue logging the potential energy associated with this force.

Examples:

force.disable()
force.disable(log=True)

Executing the disable command will remove the force from the simulation. Any hoomd.run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable().

By setting log to True, the values of the force can be logged even though the forces are not applied in the simulation. For forces that use cutoff radii, setting log=True will cause the correct r_cut values to be used throughout the simulation, and therefore possibly drive the neighbor list size larger than it otherwise would be. If log is left False, the potential energy associated with this force will not be available for logging.

enable()

Enable the force.

Examples:

force.enable()

See disable().

get_energy(group)

Get the energy of a particle group.

Parameters:group (hoomd.group) – The particle group to query the energy for.
Returns:The last computed energy for the members in the group.

Examples:

g = group.all()
energy = force.get_energy(g)
set_params(field_x, field_y, field_z, p)

Change the constant field and dipole moment.

Parameters:
  • field_x (float) – x-component of the field (units?)
  • field_y (float) – y-component of the field (units?)
  • field_z (float) – z-component of the field (units?)
  • p (float) – magnitude of the particles’ dipole moment in the local z direction

Examples:

const_ext_f_dipole = force.external_field_dipole(field_x=0.0, field_y=1.0 ,field_z=0.5, p=1.0)
const_ext_f_dipole.setParams(field_x=0.1, field_y=0.1, field_z=0.0, p=1.0))