hoomd.hpmc.pair#

Overview

AngularStep

Angular-step pair potential (HPMC).

LennardJones

Lennard-Jones pair potential (HPMC).

Pair

Pair potential base class (HPMC).

Step

Step function pair potential (HPMC).

Union

Treat particles as extended bodies.

Details

Pair Potentials for Monte Carlo.

Define \(U_{\mathrm{pair},ij}\) for use with HPMCIntegrator, which will sum all the energy from all Pair potential instances in the pair_potentials list.

Example:

simulation.operations.integrator.pair_potentials = [pair]
class hoomd.hpmc.pair.AngularStep(isotropic_potential)#

Bases: Pair

Angular-step pair potential (HPMC).

Parameters:

isotropic_potential (hoomd.hpmc.pair.Pair) – the isotropic part of the interaction between the particles.

AngularStep computes the the given isotropic potential multiplied by a step function that is dependent on the relative orientation between any two patches:

\[U(\vec{r}_{ij}, \mathbf{q}_i, \mathbf{q}_j) = U_\mathrm{isotropic}(\vec{r}_{ij}) \cdot \max \left(1, \sum_{m=1}^{N_{\mathrm{patches},i}} \sum_{m=1}^{N_{\mathrm{patches},j}} f(\mathbf{q}_i \vec{d}_{n,i} \mathbf{q}_i^*, \mathbf{q}_j \vec{d}_{m,j} \mathbf{q}_j^*, \delta_{n,i}, \delta_{m,j}) \right)\]

where \(U_\mathrm{isotropic}\) is the isotropic potential. For a given particle \(i\), \(N_{\mathrm{patches},i}\) is the number of patches, \(\vec{d}_{n,i}\) is the n-th director, and \(\delta_{n,i}\) is the n-th delta. \(f_{ij}(\vec{a}, \vec{b}, \delta_a, \delta_b)\) is an orientational masking function given by:

\[\begin{split}f(\vec{a}, \vec{b}, \delta_a, \delta_b) = \begin{cases} 1 & \hat{a} \cdot \hat{r}_{ij} \ge \cos \delta_{a} \land \hat{b} \cdot \hat{r}_{ji} \ge \cos \delta_{b} \\ 0 & \text{otherwise} \\ \end{cases}\end{split}\]

One example of this form of potential is the Kern-Frenkel model that is composed of a square well potential and an orientational masking function.

Example

angular_step = hoomd.hpmc.pair.AngularStep(
               isotropic_potential=square_well)
angular_step.mask['A'] = dict(directors=[(1.0, 0, 0)], deltas=[0.1])
simulation.operations.integrator.pair_potentials = [angular_step]

Set the patch directors \(\vec{d}_m\) and delta \(\delta_m\) values for each particle type. Patch directors are the directional unit vectors that represent the patch locations on a particle, and deltas are the half opening angles of the patch in radian.

mask#

The mask definition.

The mask describes the distribution of patches on the particle’s surface and the masking function determines the interaction scale factor as a function of two interacting particle’s masks.

The dictionary has the following keys:

  • directors (list [tuple [float, float, float]]): List of directional vectors of the patches on a particle.

  • deltas (list [float]): List of delta values (the half opening angle of the patch in radian) of the patches.

property isotropic_potential#

Get the isotropic part of the interactions between patchy particles.

This property returns the isotropic component of pairwise interaction potentials for patchy particle systems.

Example

angular_step.isotropic_potential
class hoomd.hpmc.pair.LennardJones(default_r_cut=None, default_r_on=0.0, mode='none')#

Bases: Pair

Lennard-Jones pair potential (HPMC).

Parameters:
  • default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).

  • default_r_on (float) – Default XPLOR on radius \([\mathrm{length}]\).

  • mode (str) – Energy shifting/smoothing mode.

LennardJones computes the Lennard-Jones pair potential between every pair of particles in the simulation state. The functional form of the potential, including its behavior under shifting modes, is identical to that in the MD pair potential hoomd.md.pair.LJ.

Example

lennard_jones =  hoomd.hpmc.pair.LennardJones()
lennard_jones.params[('A', 'A')] = dict(epsilon=1, sigma=1, r_cut=2.5)
simulation.operations.integrator.pair_potentials = [lennard_jones]
params#

The potential parameters. The dictionary has the following keys:

  • epsilon (float, required) - Energy well depth \(\varepsilon\) \([\mathrm{energy}]\).

  • sigma (float, required) - Characteristic length scale \(\sigma\) \([\mathrm{length}]\).

  • r_cut (float): Cutoff radius \([\mathrm{length}]\). Defaults to the value given in default_r_cut on construction.

  • r_on (float): XPLOR on radius \([\mathrm{length}]\). Defaults to the value given in default_r_on on construction.

Type: TypeParameter [tuple [particle_type, particle_type], dict]

mode#

The energy shifting/smoothing mode: Possible values are: "none", "shift", and "xplor".

Example

lennard_jones.mode = 'shift'

Type: str

class hoomd.hpmc.pair.Pair#

Bases:

Pair potential base class (HPMC).

Pair potentials define energetic interactions between pairs of particles in hoomd.hpmc.integrate.HPMCIntegrator. Particles within a cutoff distance interact with an energy that is a function the type and orientation of the particles and the vector pointing from the i particle to the j particle center.

Note

The base class Pair implements common attributes (energy, for example) and may be used in for isinstance or issubclass checks. Pair should not be instantiated directly by users.

property energy#

Potential energy contributed by this potential \([\mathrm{energy}]\).

Typically:

\[U = \sum_{i=0}^\mathrm{N_particles-1} \sum_{j=i+1}^\mathrm{N_particles-1} U_{\mathrm{pair},ij}\]

See hoomd.hpmc.integrate for the full expression which includes the evaluation over multiple images when the simulation box is small.

Example

logger.add(obj=pair, quantities=['energy'])

(Loggable: category=”scalar”)

Type:

float

class hoomd.hpmc.pair.Step#

Bases: Pair

Step function pair potential (HPMC).

Step computes a user-defined step function pair potential between every pair of particles in the simulation state. The functional form of the potential is:

\[\begin{split}U(r) = \begin{cases} \varepsilon_0 & r < r_0 \\ \varepsilon_k & r_{k-1} \le r < r_{k}, k \in [1,n-1] \\ 0 & r \ge r_{n-1} \\ \end{cases}\end{split}\]

Where \(\varepsilon_k\) is the element \(k\) of the epsilon list and \(r_k\) is the element \(k\) of the r list.

Example

step =  hoomd.hpmc.pair.Step()
step.params[('A', 'A')] = dict(epsilon=[1, -1], r=[0.5, 1.5])
simulation.operations.integrator.pair_potentials = [step]
params#

The potential parameters. The dictionary has the following keys:

  • epsilon (list [float], required) - Function values \(\varepsilon\) \([\mathrm{energy}]\).

  • r (list [float], required) - Points at which function values are defined \([\mathrm{length}]\). The values of r must be listed in monotonically increasing order.

Type: TypeParameter [tuple [particle_type, particle_type], dict]

class hoomd.hpmc.pair.Union(constituent_potential, leaf_capacity=0)#

Bases: Pair

Treat particles as extended bodies.

Parameters:
  • constituent_potential (hoomd.hpmc.pair.Pair) – The pair potential to apply between constituent points.

  • leaf_capacity (int) – Maximum number of leaf nodes in the tree data structure used by this class. The default leaf_capacity=0 uses an all N*M code path.

Union computes the potential energy between sets of constituent points that rigidly transform about each particle. The union potential between a pair of particles is:

\[U(\vec{r}_{ij}, \mathbf{q}_i, \mathbf{q}_j)) = \sum_{a=1}^{N_{\mathrm{constituents},i}} \sum_{b=1}^{N_{\mathrm{constituents},j}} U_\mathrm{constituent}(\mathbf{q}_j \vec{P}_{j,b} \mathbf{q}_j^* - \vec{P}_{i,a}, \mathbf{q}_i \mathbf{Q}_{i,a}, \mathbf{q}_j \mathbf{Q}_{j,b})\]

where \(N_{\mathrm{constituents},i}\) is the number of constituents on the \(i\) particle and \(U_\mathrm{constituent}\) is the potential evaluated by the given constituent_potential. \(\vec{P}_{i,a}\) and \(\mathbf{Q}_{i,a}\) are the constituent postitions and orientations with index \(a\) on the \(i\) particle. \(U_\mathrm{constituent}\) also depends on the constituent particle types and charges (not shown in the equation).

See also

hoomd.md.constrain.Rigid implements a similar evaluation for MD simulations.

Important

Unlike hoomd.md.constrain.Rigid, Union does not automatically include the central particle in the evaluation. You must add a constituent particle with position 0,0,0 when desired.

Note

constituent_potential may be isotropic or have orientation dependence, but it may not be another Union potential.

Tip

The default leaf_capacity=0 performs best for unions with small numbers of constituents. leaf_capacity=4 (or any other non-zero value) activates a tree algorithm that may perform better when you have many constituents.

Example

union = hoomd.hpmc.pair.Union(constituent_potential=lennard_jones)
union.body['R'] = dict(types=['A', 'A', 'A'],
                       positions=[(-1,0,0), (0,0,0), (1,0,0)])
union.body['A'] = None

simulation.operations.integrator.pair_potentials = [union]

The particle types used as constituents must be particle types present in the system state, even when there are no actual particles of that type. As shown above, set the body for constituent types to None (which is equivalent to dict(types=[], positions=[])).

body#

The body definition.

Define the position and orientation of each constituent point relative to the position and orientation of the particle (i.e. in the particle reference frame). Set a particle type name for each constituent point which will be used to determine the constituent potential parameters.

The dictionary has the following keys:

  • types (list [str]): List of types of constituent points.

  • positions (list [tuple [float, float, float]]): List of relative positions of constituent points.

  • orientations (list [tuple [float, float, float, float]]): List of orientations (as quaternions) of constituent points (optional, defaults to [(1,0,0,0)] * len(positions)).

  • charges (list [float]): List of charges of constituent points (optional, defaults to [0] * len(positions)).

Type: TypeParameter [particle_type, dict] or None

leaf_capacity#

Maximum number of leaf nodes in the tree data structure used by this class. Set leaf_capacity=0 to use an all N*M code path.

Example

union.leaf_capacity = 4
Type:

int

property constituent_potential#

Interactions between constituent points.

Example

union.constituent_potential
Type:

hpmc.pair.Pair

Modules