ALJ

class hoomd.md.pair.aniso.ALJ(nlist, default_r_cut=None)

Bases: AnisotropicPair

Anistropic LJ force.

Parameters:

ALJ computes the Lennard-Jones force between anisotropic particles as described in Ramasubramani, V. et al. 2020, using the formula:

\[U(r, r_c) = U_0(r) + U_c(r_c)\]

The first term is the central interaction \(U_0\), the standard center-center interaction between two Lennard-Jones particles with center-center distance \(r\). The second term is the contact interaction \(U_c\), computed from the smallest distance between the surfaces of the two shapes \(r_c\). The central and contact interactions are defined as follows:

\[ \begin{align}\begin{aligned}&U_0(r) = 4 \varepsilon \left[ \left( \frac{\sigma}{r} \right)^{12} - \left( \frac{\sigma}{r} \right)^{6} \right]\\&U_c(r_c) = 4 \varepsilon_c(\varepsilon) \left[ \left( \frac{\sigma_c}{r_c} \right)^{12} - \left( \frac{\sigma_c}{r_c} \right)^{6} \right]\end{aligned}\end{align} \]

where \(\varepsilon\) (epsilon) affects strength of both the central and contact interactions, \(\varepsilon_c\) is an energy coefficient set proportional to \(\varepsilon\) to preserve the shape, \(\sigma\) is the interaction distance of the central term computed as the average of \(\sigma_i\) (sigma_i) and \(\sigma_j\) (sigma_j). Lastly, ALJ uses the contact ratios \(\beta_i\) (contact_ratio_i) and \(\beta_j\) (contact_ratio_j) to compute the contact sigma \(\sigma_c\) as follows:

\[ \begin{align}\begin{aligned}\sigma_c &= \frac{1}{2} \left[\sigma_{ci} + \sigma_{cj} \right]\\\sigma_{ci} &= \beta_i \cdot \sigma_i\\\sigma_{cj} &= \beta_j \cdot \sigma_j\end{aligned}\end{align} \]

The total potential energy is therefore the sum of two interactions, a central Lennard-Jones potential and a radially-shifted Lennard-Jones potential where the shift is anisotropic and depends on the extent of the shape in each direction.

Each term has an independent cutoff at which the energy is set to zero. The behavior of these cutoffs is dependent on whether a user requires LJ or Weeks-Chandler-Anderson (WCA)-like (repulsive-only) behavior. This behavior is controlled using the alpha parameter, which can take on the following values:

Set alpha based on range of the center-center and contact-contact interactions.

center-center repulsive only

center-center full-range

contact-contact repulsive only

alpha = 0

alpha = 1

contact-contact full-range

alpha = 2

alpha = 3

For polytopes, computing interactions using a single contact point leads to significant instabilities in the torques because the contact point can jump from one end of a face to another in an arbitrarily small time interval. To ameliorate this, the ALJ potential performs a local averaging over all the features associated with the closest simplices on two polytopes. This averaging can be turned off by setting the average_simplices key for the type pair to False.

Specifying only rounding_radii creates an ellipsoid, while specifying only vertices creates a convex polytope (set vertices and faces to empty lists to create the ellipsoid).

Important

The repulsive part of the contact interaction \(U_c(r_c)\) prevents two ALJ particles from approaching closely, effectively rounding the shape by a radius \(\sigma_c\). For this reason, the shape written by type_shapes includes the rounding due to rounding_radii and that due to \(\sigma_c\).

Choosing r_cut:

Set r_cut for each pair of particle types so that ALJ can compute interactions for all possible relative placements and orientations of the particles. The farthest apart two particles can be while still interacting depends on the value of alpha.

In the following list, the first argument to the \(\max\) function is for the center-center interaction. The second argument is for the contact-contact interaction, where \(R_i\) is the radius of the shape’s minimal origin-centered bounding sphere of the particle with type \(i\).

Let \(\lambda_{min} = 2^{1/6}\) be the position of the potential energy minimum of the Lennard-Jones potential and \(\lambda_{cut}^{attractive}\) be a larger value, such as 2.5 (typically used in isotropic LJ systems).

  • For alpha=0:

    \[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{min}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{min}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]
  • For alpha=1:

    \[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{cut}^{attractive}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j}+ \frac{\lambda_{min}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]
  • For alpha=2:

    \[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{min}}{2} (\sigma_i + \sigma_j)), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{cut}^{attractive}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]
  • For alpha=3:

    \[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{cut}^{attractive}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{cut}^{attractive}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]

Warning

Changing dimension in a simulation will invalidate this force and will lead to error or unrealistic behavior.

Note

While the evaluation of the potential is symmetric with respect to the potential parameter labels i and j, the parameters which physically represent a specific particle type must appear in all sets of pair parameters which include that particle type.

Example:

nl = hoomd.md.nlist.Cell(buffer=0.4)
alj = hoomd.md.pair.aniso.ALJ(nl, r_cut=2.5)

cube_verts = [(-0.5, -0.5, -0.5),
              (-0.5, -0.5, 0.5),
              (-0.5, 0.5, -0.5),
              (-0.5, 0.5, 0.5),
              (0.5, -0.5, -0.5),
              (0.5, -0.5, 0.5),
              (0.5, 0.5, -0.5),
              (0.5, 0.5, 0.5)];

cube_faces = [[0, 2, 6],
              [6, 4, 0],
              [5, 0, 4],
              [5,1,0],
              [5,4,6],
              [5,6,7],
              [3,2,0],
              [3,0,1],
              [3,6,2],
              [3,7,6],
              [3,1,5],
              [3,5,7]]

alj.params[("A", "A")] = dict(epsilon=2.0,
                              sigma_i=1.0,
                              sigma_j=1.0,
                              alpha=1,
                              )
alj.shape["A"] = dict(vertices=cube_verts,
                      faces=cube_faces)

The following example shows how to easily get the faces, with vertex indices properly ordered, for a shape with known vertices by using the coxeter package:

Example:

import coxeter

nl = hoomd.md.nlist.Cell(buffer=0.4)
alj = hoomd.md.pair.aniso.ALJ(nl, r_cut=2.5)

cube_verts = [[-0.5, -0.5, -0.5],
              [-0.5, -0.5, 0.5],
              [-0.5, 0.5, -0.5],
              [-0.5, 0.5, 0.5],
              [0.5, -0.5, -0.5],
              [0.5, -0.5, 0.5],
              [0.5, 0.5, -0.5],
              [0.5, 0.5, 0.5]]

cube = coxeter.shapes.ConvexPolyhedron(cube_verts)

alj.params[("A", "A")] = dict(epsilon=2.0,
                              sigma_i=1.0,
                              sigma_j=1.0,
                              alpha=1,
                              )
alj.shape["A"] = dict(vertices=cube.vertices,
                      faces=cube.faces)

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 Force:

additional_energy

Additional energy term. Read more...

additional_virial

Additional virial tensor term \(W_\mathrm{additional}\). Read more...

cpu_local_force_arrays

Local force arrays on the CPU. Read more...

energies

Energy contribution \(U_i\) from each particle. Read more...

energy

The potential energy \(U\) of the system from this force. Read more...

forces

The force \(\vec{F}_i\) applied to each particle. Read more...

gpu_local_force_arrays

Local force arrays on the GPU. Read more...

torques

The torque \(\vec{\tau}_i\) applied to each particle. Read more...

virials

Virial tensor contribution \(W_i\) from each particle. Read more...


Members inherited from Pair:

nlist

Neighbor list used to compute the pair force. Read more...

mode

Energy smoothing/cutoff mode. Read more...

r_cut

Cuttoff radius beyond which the energy and force are 0. Read more...

r_on

Radius at which the XPLOR smoothing function starts. Read more...

compute_energy()

Compute the energy between two sets of particles. Read more...


Members defined in ALJ:

params

The ALJ potential parameters. The dictionary has the following keys:

  • epsilon (float, required) - base energy scale \(\varepsilon\) \([energy]\).

  • sigma_i (float, required) - the insphere diameter of the first particle type, \(\sigma_i\) \([length]\).

  • sigma_j (float, required) - the insphere diameter of the second particle type, \(\sigma_j\) \([length]\).

  • alpha (int, required) - Integer 0-3 indicating whether or not to include the attractive component of the interaction (see above for details).

  • contact_ratio_i (float, optional) - \(\beta_i\), the ratio of the contact sphere diameter of the first type with sigma_i. Defaults to 0.15.

  • contact_ratio_j (float, optional) - \(\beta_j\), the ratio of the contact sphere diameter of the second type with sigma_j. Defaults to 0.15.

  • average_simplices (bool, optional) - Whether to average over simplices. Defaults to True. See class documentation for more information.

Type: hoomd.data.typeparam.TypeParameter [tuple [particle_types, particle_types], dict]

shape

The shape of a given type. The dictionary has the following keys per type:

  • vertices (list [tuple [float, float, float]], required) - The vertices of a convex polytope in 2 or 3 dimensions. The third dimension in 2D is ignored.

  • rounding_radii (tuple [float, float, float] or float) - The semimajor axes of a rounding ellipsoid \(R_{\mathrm{rounding},i}\). If a single value is specified, the rounding ellipsoid is a sphere. Defaults to (0.0, 0.0, 0.0).

  • faces (list [list [int]], required) - The faces of the polyhedron specified as a list of list of integers. The indices corresponding to the vertices must be ordered counterclockwise with respect to the face normal vector pointing outward from the origin.

Type: hoomd.data.typeparam.TypeParameter [particle_types, dict]

property type_shapes

The shape specification for use with GSD files for visualization.

This is not meant to be used for access to shape information in Python. See the attribute shape for programmatic assess. Use this property to log shape for visualization and storage through the GSD file type.

(Loggable: category=”object”)

Type:

list [dict [str, any]]