hoomd.hpmc.integrate#

Overview

ConvexPolygon

Convex polygon hard particle Monte Carlo integrator.

ConvexPolyhedron

Convex polyhedron hard particle Monte Carlo integrator.

ConvexSpheropolygon

Convex spheropolygon hard particle Monte Carlo integrator.

ConvexSpheropolyhedron

Convex spheropolyhedron hard particle Monte Carlo integrator.

ConvexSpheropolyhedronUnion

Convex spheropolyhedron union hard particle Monte Carlo integrator.

Ellipsoid

Ellipsoid hard particle Monte Carlo integrator.

FacetedEllipsoid

Faceted ellipsoid hard particle Monte Carlo integrator.

FacetedEllipsoidUnion

Faceted ellispod union hard particle Monte Carlo integrator.

HPMCIntegrator

Base class hard particle Monte Carlo integrator.

Polyhedron

Polyhedron hard particle Monte Carlo integrator.

SimplePolygon

Simple polygon hard particle Monte Carlo integrator.

Sphere

Sphere hard particle Monte Carlo integrator.

SphereUnion

Sphere union hard particle Monte Carlo integrator.

Sphinx

Sphinx hard particle Monte Carlo integrator.

Details

Hard particle Monte Carlo integrators.

Metropolis Monte Carlo

The hard particle Monte Carlo (HPMC) integrator HPMCIntegrator samples equilibrium system states using the Metropolis Monte Carlo method. In this method, HPMCIntegrator takes the existing system state in the configuration \(C = (\vec{r}_0, \vec{r}_1, \ldots \vec{r}_{N_\mathrm{particles}-1}, \mathbf{q}_0, \mathbf{q}_2, \ldots \mathbf{q}_{N_\mathrm{particles}-1})\) with potential energy \(U\) and perturbs it to a trial configuration \(C^t\) with potential energy \(U^t\) leading to an energy difference \(\Delta U = U^t - U\). The trial move is accepted with the probability:

\[\begin{split}p_\mathrm{accept} = \begin{cases} \exp(-\beta \Delta U) & \Delta U > 0 \\ 1 & \Delta U \le 0 \\ \end{cases}\end{split}\]

When the trial move is accepted, the system state is set to the the trial configuration. When it is not accepted, the move is rejected and the state is not modified.

Temperature

HPMC assumes that \(\beta = \frac{1}{kT} = 1\). This is not relevant to systems of purely hard particles where \(\Delta U\) is either 0 or \(\infty\). To adjust the effective temperature in systems with finite interactions (see Energy evaluation below), scale the magnitude of the energetic interactions accordingly.

Local trial moves

HPMCIntegrator generates local trial moves for a single particle \(i\) at a time. The move is either a translation move or a rotation move, selected randomly with the probability of a translation move set by HPMCIntegrator.translation_move_probability (\(p_\mathrm{translation}\)).

The form of the trial move depends on the dimensionality of the system. Let \(u\) be a random value in the interval \([0,1]\), \(\vec{v}\) be a random vector uniformly distributed within the ball of radius 1, and \(\mathbf{w}\) be a random unit quaternion from the set of uniformly distributed rotations. Then the 3D trial move for particle \(i\) is:

\[\begin{split}\begin{cases} \left( \begin{array}{l} \vec{r}^t_i = \vec{r}_i + d_i \vec{v}, \\ \mathbf{q}^t_i = \mathbf{q}_i \end{array} \right) & u \le p_\mathrm{translation} \\ \left( \begin{array}{l} \vec{r}^t_i = \vec{r}_i, \\ \mathbf{q}^t_i = \frac{\mathbf{q}_i + a_i \mathbf{w}} {\vert \mathbf{q}_i + a_i \mathbf{w} \vert} \end{array} \right) & u > p_\mathrm{translation} \\ \end{cases}\end{split}\]

where \(d_i\) is the translation move size for particle \(i\) (set by particle type with HPMCIntegrator.d) and \(a_i\) is the rotation move size (set by particle type with HPMCIntegrator.a).

In 2D boxes, let \(\vec{v}\) be a random vector uniformly distributed within the disk of radius 1 in the x,y plane and \(\alpha\) be a random angle in radians in the interval \([-a_i,a_i]\). Form a quaternion that rotates about the z axis by \(\alpha\): \(\mathbf{w} = (\cos(\alpha/2), 0, 0, \sin(\alpha/2))\). The 2D trial move for particle \(i\) is:

\[\begin{split}\begin{cases} \left( \begin{array}{l} \vec{r}^t_i = \vec{r}_i + d_i \vec{v}, \\ \mathbf{q}^t_i = \mathbf{q}_i \end{array} \right) & u \le p_\mathrm{translation} \\ \left( \begin{array}{l} \vec{r}^t_i = \vec{r}_i, \\ \mathbf{q}^t_i = \frac{\mathbf{q}_i \cdot \mathbf{w}} {\vert \mathbf{q}_i \cdot \mathbf{w} \vert} \end{array} \right) & u > p_\mathrm{translation} \\ \end{cases}\end{split}\]

Note

For non-orientable spheres, \(p_\mathrm{translation} = 1\).

Timesteps

In the serial CPU implementation, HPMCIntegrator performs nselect trial moves per particle in each timestep (which defaults to 4). To achieve detailed balance at the level of a timestep, HPMCIntegrator randomly chooses with equal probability to loop through particles in forward index or reverse index order (random selection severely degrades performance due to cache incoherency). In the GPU and MPI implementations, trial moves are performed in parallel for particles in active domains while leaving particles on the border fixed (see Anderson 2016 for a full description). As a consequence, a single timestep may perform more or less than nselect trial moves per particle when using the parallel code paths. Monitor the number of trial moves performed with HPMCIntegrator.translate_moves and HPMCIntegrator.rotate_moves.

Random numbers

HPMCIntegrator uses a pseudorandom number stream to generate the trial moves. Set the seed using hoomd.Simulation.seed. Given the same seed, the same initial configuration, and the same execution configuration (device and MPI configuration), HPMCIntegrator, will produce exactly the same trajectory.

Note

Due to limited floating point precision, full trajectory reproducibility is only possible with the same binary installation running on the same hardware device. Compiler optimizations, changes to the HOOMD source code, and machine specific code paths may lead to different trajectories.

Energy evaluation

HPMCIntegrator evaluates the energy of a configuration from a number of terms:

\[U = U_{\mathrm{pair}} + U_{\mathrm{shape}} + U_{\mathrm{external}}\]

To enable simulations of small systems, the pair and shape energies evaluate interactions between pairs of particles in multiple box images:

\[\begin{split}U_{\mathrm{pair}} = & \sum_{i=0}^{N_\mathrm{particles}-1} \sum_{j=i+1}^{N_\mathrm{particles}-1} U_{\mathrm{pair},ij}(\vec{r}_j - \vec{r}_i, \mathbf{q}_i, \mathbf{q}_j) \\ + & \sum_{i=0}^{N_\mathrm{particles}-1} \sum_{j=i}^{N_\mathrm{particles}-1} \sum_{\vec{A} \in B_\mathrm{images}, \vec{A} \ne \vec{0}} U_{\mathrm{pair},ij}(\vec{r}_j - (\vec{r}_i + \vec{A}), \mathbf{q}_i, \mathbf{q}_j)\end{split}\]

where \(\vec{A} = h\vec{a}_1 + k\vec{a}_2 + l\vec{a}_3\) is a vector that translates by periodic box images and the set of box images includes all image vectors necessary to find interactions between particles in the primary image with particles in periodic images The first sum evaluates interactions between particle \(i\) with other particles (not itself) in the primary box image. The second sum evaluates interactions between particle \(i\) and all potentially interacting periodic images of all particles (including itself). HPMCIntegrator computes \(U_{\mathrm{shape}}\) similarly (see below).

External potentials apply to each particle individually:

\[U_\mathrm{external} = \sum_{i=0}^\mathrm{N_particles-1} U_{\mathrm{external},i}(\vec{r}_i, \mathbf{q}_i)\]

Potential classes in hoomd.hpmc.pair evaluate \(U_{\mathrm{pair},ij}\). HPMC sums all the Pair potentials in pair_potentials and the CPPPotentialBase potential in pair_potential during integration.

Similarly, potential classes in hoomd.hpmc.external evaluate \(U_{\mathrm{external},i}\). Assign a class instance to HPMCIntegrator.external_potential to apply it during integration.

Shape overlap tests

HPMCIntegrator performs shape overlap tests to evaluate \(U_{\mathrm{shape}}\). Let \(S\) be the set of all points inside the shape in the local coordinate system of the shape:

\[S = \{ \vec{a} \in \mathbb{R}^3 : \vec{a} \enspace \mathrm{inside\ the\ shape} \}\]

See the subclasses of HPMCIntegrator for formal definitions of the shapes, whose parameters are set by particle type. Let \(S_i\) refer specifically to the shape for particle \(i\).

The quaternion \(\mathbf{q}\) represents a rotation of the shape from its local coordinate system to the given orientation:

\[S(\mathbf{q}) = \{ \mathbf{q}\vec{a}\mathbf{q}^* : \vec{a} \in S \}\]

The full transformation from the local shape coordinate system to the simulation box coordinate system includes a rotation and translation:

\[S(\mathbf{q}, \vec{r}) = \{ \mathbf{q}\vec{a}\mathbf{q}^* + \vec{r} : \vec{a} \in S \}\]

HPMCIntegrator defines the shape overlap test for two shapes:

\[\mathrm{overlap}(S_1, S_2) = S_1 \bigcap S_2\]

To check for overlaps between two particles in the box, rotating both shapes from their local frame to the box frame, and translate \(S_2\) relative to particle 1:

\[\mathrm{overlap}(S_1(\mathbf{q}_1), S_2(\mathbf{q}_2, \vec{r}_2 - \vec{r}_1))\]

The complete hard shape interaction energy for a given configuration is:

\[\begin{split}U_\mathrm{shape} = \quad & \infty \cdot \sum_{i=0}^{N_\mathrm{particles}-1} \sum_{j=i+1}^{N_\mathrm{particles}-1} \left[ \mathrm{overlap}\left( S_i(\mathbf{q}_i), S_j(\mathbf{q}_j, \vec{r}_j - \vec{r}_i) \right) \ne \emptyset \right] \\ + & \infty \cdot \sum_{i=0}^{N_\mathrm{particles}-1} \sum_{j=i}^{N_\mathrm{particles}-1} \sum_{\vec{A} \in B_\mathrm{images}, \vec{A} \ne \vec{0}} \left[ \mathrm{overlap}\left( S_i(\mathbf{q}_i), S_j(\mathbf{q}_j, \vec{r}_j - (\vec{r}_i + \vec{A})) \right) \ne \emptyset \right]\end{split}\]

where the square brackets denote the Iverson bracket.

Note

While this notation is written in as sums over all particles HPMCIntegrator uses spatial data structures to evaluate these calculations efficiently. Similarly, while the overlap test is notated as a set intersection, HPMCIntegrator employs efficient computational geometry algorithms to determine whether there is or is not an overlap.

Implicit depletants

Set HPMCIntegrator.depletant_fugacity to activate the implicit depletant code path. This inerts depletant particles during every trial move and modifies the acceptance criterion accordingly. See Glaser 2015 for details.

Deprecated since version 4.4.0: depletant_fugacity > 0 is deprecated.

class hoomd.hpmc.integrate.ConvexPolygon(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Convex polygon hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of convex polygons. The shape \(S\) of a convex polygon includes the points inside and on the surface of the convex hull of the vertices (see shape). For example:

Example of a convex polygon with vertex labels.

Important

ConvexPolygon simulations must be performed in 2D systems.

See also

Use SimplePolygon for concave polygons.

Wall support.

ConvexPolygon supports no hoomd.wall geometries.

Examples:

mc = hoomd.hpmc.integrate.ConvexPolygon(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(vertices=[(-0.5, -0.5),
                               (0.5, -0.5),
                               (0.5, 0.5),
                               (-0.5, 0.5)]);
print('vertices = ', mc.shape["A"]["vertices"])
shape#

The shape parameters for each particle type. The dictionary has the following keys.

  • vertices (list [tuple [float, float]], required) - vertices of the polygon \([\mathrm{length}]\).

    • Vertices MUST be specified in a counter-clockwise order.

    • The origin MUST be contained within the polygon.

    • Points inside the polygon MUST NOT be included.

    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (float, default: 0.0) - Ignored, but present because ConvexPolygon shares data structures with ConvexSpheropolygon \([\mathrm{length}]\).

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'Polygon', 'sweep_radius': 0,
  'vertices': [[-0.5, -0.5], [0.5, -0.5], [0.5, 0.5], [-0.5, 0.5]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.ConvexPolyhedron(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Convex polyhedron hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of convex polyhedra. The shape \(S\) of a convex polyhedron includes the points inside and on the surface of the convex hull of the vertices (see shape). For example:

Example of a convex polyhedron with vertex labels.

See also

Use Polyhedron for concave polyhedra.

Wall support.

ConvexPolyhedron supports all hoomd.wall geometries.

Example:

mc = hpmc.integrate.ConvexPolyhedron(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(vertices=[(0.5, 0.5, 0.5),
                               (0.5, -0.5, -0.5),
                               (-0.5, 0.5, -0.5),
                               (-0.5, -0.5, 0.5)]);
print('vertices = ', mc.shape["A"]["vertices"])
shape#

The shape parameters for each particle type. The dictionary has the following keys.

  • vertices (list [tuple [float, float, float]], required) - vertices of the polyhedron \([\mathrm{length}]\).

    • The origin MUST be contained within the polyhedron.

    • The origin centered sphere that encloses all vertices should be of minimal size for optimal performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (float, default: 0.0) - Ignored, but present because ConvexPolyhedron shares data structures with ConvexSpheropolyhedron \([\mathrm{length}]\).

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'ConvexPolyhedron', 'sweep_radius': 0,
  'vertices': [[0.5, 0.5, 0.5], [0.5, -0.5, -0.5],
               [-0.5, 0.5, -0.5], [-0.5, -0.5, 0.5]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.ConvexSpheropolygon(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Convex spheropolygon hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of convex spheropolygons. The shape \(S\) of a convex spheropolygon includes the points inside and on the surface of the convex hull of the vertices plus a disk (with radius sweep_radius)swept along the perimeter (see shape). For example:

Example of a convex spheropolygon with vertex and sweep labels.

Important

ConvexSpheropolygon simulations must be performed in 2D systems.

Tip

To model mixtures of convex polygons and convex spheropolygons, use ConvexSpheropolygon and set the sweep radius to 0 for some shape types.

Tip

A 1-vertex spheropolygon is a disk and a 2-vertex spheropolygon is a rectangle with half disk caps.

Wall support.

ConvexSpheropolygon supports no hoomd.wall geometries.

Examples:

mc = hoomd.hpmc.integrate.ConvexSpheropolygon(default_d=0.3,
                                             default_a=0.4)
mc.shape["A"] = dict(vertices=[(-0.5, -0.5),
                               (0.5, -0.5),
                               (0.5, 0.5),
                               (-0.5, 0.5)],
                     sweep_radius=0.1);

mc.shape["A"] = dict(vertices=[(0,0)],
                     sweep_radius=0.5,
                     ignore_statistics=True);

print('vertices = ', mc.shape["A"]["vertices"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • vertices (list [tuple [float, float]], required) - vertices of the polygon \([\mathrm{length}]\).

    • The origin MUST be contained within the spheropolygon.

    • Points inside the polygon should not be included.

    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (default: 0.0) - radius of the disk swept around the edges of the polygon \([\mathrm{length}]\). Set a non-zero sweep_radius to create a spheropolygon.

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'Polygon', 'sweep_radius': 0.1,
  'vertices': [[-0.5, -0.5], [0.5, -0.5], [0.5, 0.5], [-0.5, 0.5]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.ConvexSpheropolyhedron(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Convex spheropolyhedron hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of convex spheropolyhedra. The shape \(S\) of a convex spheropolyhedron includes the points inside and on the surface of the convex hull of the vertices plus a sphere (with radius sweep_radius) swept along the perimeter (see shape). See ConvexSpheropolygon for a visual example in 2D.

Tip

A 1-vertex spheropolygon is a sphere and a 2-vertex spheropolygon is a spherocylinder.

Wall support.

ConvexSpheropolyhedron supports the hoomd.wall.Sphere and hoomd.wall.Plane geometries.

Example:

mc = hpmc.integrate.ConvexSpheropolyhedron(default_d=0.3, default_a=0.4)
mc.shape['tetrahedron'] = dict(vertices=[(0.5, 0.5, 0.5),
                                         (0.5, -0.5, -0.5),
                                         (-0.5, 0.5, -0.5),
                                         (-0.5, -0.5, 0.5)]);
print('vertices = ', mc.shape['tetrahedron']["vertices"])

mc.shape['SphericalDepletant'] = dict(vertices=[],
                                      sweep_radius=0.1,
                                      ignore_statistics=True);

Depletants example:

mc = hpmc.integrate.ConvexSpheropolyhedron(default_d=0.3, default_a=0.4)
mc.shape["tetrahedron"] = dict(vertices=[(0.5, 0.5, 0.5),
                                         (0.5, -0.5, -0.5),
                                         (-0.5, 0.5, -0.5),
                                         (-0.5, -0.5, 0.5)]);
mc.shape["SphericalDepletant"] = dict(vertices=[], sweep_radius=0.1);
mc.depletant_fugacity["SphericalDepletant"] = 3.0
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • vertices (list [tuple [float, float, float]], required) - vertices of the polyhedron \([\mathrm{length}]\).

    • The origin MUST be contained within the polyhedron.

    • The origin centered sphere that encloses all vertices should be of minimal size for optimal performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (float, default: 0.0) - radius of the sphere swept around the surface of the polyhedron \([\mathrm{length}]\). Set a non-zero sweep_radius to create a spheropolyhedron.

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'ConvexPolyhedron', 'sweep_radius': 0.1,
  'vertices': [[0.5, 0.5, 0.5], [0.5, -0.5, -0.5],
               [-0.5, 0.5, -0.5], [-0.5, -0.5, 0.5]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.ConvexSpheropolyhedronUnion(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Convex spheropolyhedron union hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of unions of convex sphereopolyhedra. The union shape \(S\) is the set union of the given convex spheropolyhedra:

\[S = \bigcup_k S_k(\mathbf{q}_k, \vec{r}_k)\]

Each constituent shape in the union has its own shape parameters \(S_k\), position \(\vec{r}_k\), and orientation \(\mathbf{q}_k\) (see shape).

Note

This shape uses an internal OBB tree for fast collision queries. Depending on the number of constituent spheropolyhedra in the tree, different values of the number of spheropolyhedra per leaf node may yield different performance. The capacity of leaf nodes is configurable.

Wall support.

ConvexSpheropolyhedronUnion supports no hoomd.wall geometries.

Example:

mc = hoomd.hpmc.integrate.ConvexSpheropolyhedronUnion(default_d=0.3,
                                                      default_a=0.4)
cube_verts = [[-1,-1,-1],
              [-1,-1,1],
              [-1,1,1],
              [-1,1,-1],
              [1,-1,-1],
              [1,-1,1],
              [1,1,1],
              [1,1,-1]]
mc.shape["A"] = dict(shapes=[cube_verts, cube_verts],
                     positions=[(0, 0, 0), (0, 0, 1)],
                     orientations=[(1, 0, 0, 0), (1, 0, 0, 0)],
                     overlap=[1, 1]);
print('vertices of the first cube = ',
      mc.shape["A"]["shapes"][0]["vertices"])
print('center of the first cube = ', mc.shape["A"]["positions"][0])
print('orientation of the first cube = ',
      mc.shape_param["A"]["orientations"][0])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • shapes (list [dict], required) - Shape parameters for each spheropolyhedron in the union. See ConvexSpheropolyhedron.shape for the accepted parameters.

  • positions (list [tuple [float, float, float]], required) - Position of each spheropolyhedron in the union. \([\mathrm{length}]\)

  • orientations (list [ tuple [float, float, float, float]], default: None) - Orientation of each spheropolyhedron in the union. When not None, orientations must have a length equal to that of positions. When None (the default), orientations is initialized with all [1,0,0,0]’s.

  • overlap (list [int], default: None) - Check for overlaps between constituent particles when overlap [i] & overlap[j] is nonzero (& is the bitwise AND operator). When not None, overlap must have a length equal to that of positions. When None (the default), overlap is initialized with all 1’s.

  • capacity (int, default: 4) - set the maximum number of particles per leaf node to adjust performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Type:

TypeParameter [particle type, dict]

class hoomd.hpmc.integrate.Ellipsoid(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Ellipsoid hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of ellipsoids. The shape \(S\) includes all points inside and on the surface of an ellipsoid:

\[S = \left \{ \vec{r} : \frac{r_x^2}{a^2} + \frac{r_y^2}{b^2} + \frac{r_z^2}{c^2} \le 1 \right\}\]

where \(r_x\), \(r_y\), \(r_z\) are the components of \(\vec{r}\), and the parameters \(a\), \(b\), and \(c\) are the half axes of the ellipsoid set in shape.

Wall support.

Ellipsoid supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.Ellipsoid(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(a=0.5, b=0.25, c=0.125);
print('ellipsoids parameters (a,b,c) = ',
      mc.shape["A"]["a"],
      mc.shape["A"]["b"],
      mc.shape["A"]["c"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • a (float, required) - half axis of ellipsoid in the x direction \([\mathrm{length}]\)

  • b (float, required) - half axis of ellipsoid in the y direction \([\mathrm{length}]\)

  • c (float, required) - half axis of ellipsoid in the z direction \([\mathrm{length}]\)

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'Ellipsoid', 'a': 1.0, 'b': 1.5, 'c': 1}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.FacetedEllipsoid(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Faceted ellipsoid hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of faceted ellipsoids. The shape \(S\) of a faceted ellipsoid is the intersection of an ellipsoid with a convex polyhedron defined through halfspaces (see shape). The equation defining each halfspace is given by:

\[\vec{n}_i\cdot \vec{r} + b_i \le 0\]

where \(\vec{n}_i\) is the face normal, and \(b_i\) is the offset.

Wall support.

FacetedEllipsoid supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.FacetedEllipsoid(default_d=0.3, default_a=0.4)

# half-space intersection
slab_normals = [(-1,0,0),(1,0,0),(0,-1,0),(0,1,0),(0,0,-1),(0,0,1)]
slab_offsets = [-0.1,-1,-.5,-.5,-.5,-.5]

# polyedron vertices
slab_verts = [[-.1,-.5,-.5],
              [-.1,-.5,.5],
              [-.1,.5,.5],
              [-.1,.5,-.5],
              [1,-.5,-.5],
              [1,-.5,.5],
              [1,.5,.5],
              [1,.5,-.5]]

mc.shape["A"] = dict(normals=slab_normals,
                     offsets=slab_offsets,
                     vertices=slab_verts,
                     a=1.0,
                     b=0.5,
                     c=0.5);
print('a = {}, b = {}, c = {}',
      mc.shape["A"]["a"], mc.shape["A"]["b"], mc.shape["A"]["c"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • normals (list [tuple [float, float, float]], required) - facet normals \(\\vec{n}_i\).

  • offsets (list [float], required) - list of offsets \(b_i\) \([\mathrm{length}^2]\)

  • a (float, required) - half axis of ellipsoid in the x direction \([\mathrm{length}]\)

  • b (float, required) - half axis of ellipsoid in the y direction \([\mathrm{length}]\)

  • c (float, required) - half axis of ellipsoid in the z direction \([\mathrm{length}]\)

  • vertices (list [tuple [float, float, float]], default: []) - list of vertices for intersection polyhedron (see note below). \([\mathrm{length}]\)

  • origin (tuple [float, float, float], default: (0,0,0)) - A point inside the shape. \([\mathrm{length}]\)

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Important

The origin must be chosen so as to lie inside the shape, or the overlap check will not work. This condition is not checked.

Warning

Planes must not be coplanar.

Note

For simple intersections with planes that do not intersect within the sphere, the vertices list can be left empty. When specified, the half-space intersection of the normals must match the convex polyhedron defined by the vertices (if non-empty), the half-space intersection is not calculated automatically.

Type:

TypeParameter[particle type, dict]

class hoomd.hpmc.integrate.FacetedEllipsoidUnion(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Faceted ellispod union hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of unions of faceted ellipsoids. The union shape \(S\) is the set union of the given faceted ellipsoids:

\[S = \bigcup_k S_k(\mathbf{q}_k, \vec{r}_k)\]

Each constituent shape in the union has its own shape parameters \(S_k\), position \(\vec{r}_k\), and orientation \(\mathbf{q}_k\) (see shape).

Note

This shape uses an internal OBB tree for fast collision queries. Depending on the number of constituent faceted ellipsoids in the tree, different values of the number of faceted ellipsoids per leaf node may yield different performance. The capacity of leaf nodes is configurable.

Wall support.

FacetedEllipsoidUnion supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.FacetedEllipsoidUnion(default_d=0.3, default_a=0.4)

# make a prolate Janus ellipsoid
# cut away -x halfspace
normals = [(-1,0,0)]
offsets = [0]
slab_normals = [(-1,0,0),(1,0,0),(0,-1,0),(0,1,0),(0,0,-1),(0,0,1)]
slab_offsets = [-0.1,-1,-.5,-.5,-.5,-.5)

# polyedron vertices
slab_verts = [[-.1,-.5,-.5],
              [-.1,-.5,.5],
              [-.1,.5,.5],
              [-.1,.5,-.5],
              [1,-.5,-.5],
              [1,-.5,.5],
              [1,.5,.5],
              [1,.5,-.5]]

faceted_ellipsoid1 = dict(normals=slab_normals,
                          offsets=slab_offsets,
                          vertices=slab_verts,
                          a=1.0,
                          b=0.5,
                          c=0.5);
faceted_ellipsoid2 = dict(normals=slab_normals,
                          offsets=slab_offsets,
                          vertices=slab_verts,
                          a=0.5,
                          b=1,
                          c=1);

mc.shape["A"] = dict(shapes=[faceted_ellipsoid1, faceted_ellipsoid2],
                     positions=[(0, 0, 0), (0, 0, 1)],
                     orientations=[(1, 0, 0, 0), (1, 0, 0, 0)],
                     overlap=[1, 1]);

print('offsets of the first faceted ellipsoid = ',
      mc.shape["A"]["shapes"][0]["offsets"])
print('normals of the first faceted ellipsoid = ',
      mc.shape["A"]["shapes"][0]["normals"])
print('vertices of the first faceted ellipsoid = ',
      mc.shape["A"]["shapes"][0]["vertices"]
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • shapes (list [ dict], required) - Shape parameters for each faceted ellipsoid in the union. See FacetedEllipsoid.shape for the accepted parameters.

  • positions (list [tuple [float, float, float]], required) - Position of each faceted ellipsoid in the union. \([\mathrm{length}]\)

  • orientations (list [tuple [float, float, float, float]], default: None) - Orientation of each faceted ellipsoid in the union. When not None, orientations must have a length equal to that of positions. When None (the default), orientations is initialized with all [1,0,0,0]’s.

  • overlap (list [int], default: None) - Check for overlaps between constituent particles when overlap [i] & overlap[j] is nonzero (& is the bitwise AND operator). When not None, overlap must have a length equal to that of positions. When None (the default), overlap is initialized with all 1’s.

  • capacity (int, default: 4) - set the maximum number of particles per leaf node to adjust performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Type:

TypeParameter [particle type, dict]

class hoomd.hpmc.integrate.HPMCIntegrator(default_d, default_a, translation_move_probability, nselect)#

Bases: Integrator

Base class hard particle Monte Carlo integrator.

HPMCIntegrator is the base class for all HPMC integrators. The attributes documented here are available to all HPMC integrators.

See also

The module level documentation hoomd.hpmc.integrate describes the hard particle Monte Carlo algorithm.

Warning

This class should not be instantiated by users. The class can be used for isinstance or issubclass checks.

Ignoring overlap checks

Set elements of interaction_matrix to False to disable overlap checks between specific pairs of particle types.

Writing type_shapes to GSD files.

Use a Logger in combination with a HPMC integrator and a GSD writer to write type_shapes to the GSD file for use with OVITO. For example:

mc = hoomd.hpmc.integrate.Sphere()
log = hoomd.logging.Logger()
log.add(mc, quantities=['type_shapes'])
gsd = hoomd.write.GSD(
    'trajectory.gsd', hoomd.trigger.Periodic(1000), log=log)

Threading

HPMC integrators use threaded execution on multiple CPU cores only when placing implicit depletants (depletant_fugacity != 0).

Deprecated since version 4.4.0: num_cpu_threads >= 1 is deprecated. Set num_cpu_threads = 1.

Mixed precision

All HPMC integrators use reduced precision floating point arithmetic when checking for particle overlaps in the local particle reference frame.

Parameters

a#

Maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

Type:

TypeParameter [particle type, float]

d#

Maximum size of displacement trial moves \([\mathrm{length}]\).

Type:

TypeParameter [particle type, float]

depletant_fugacity#

Depletant fugacity \([\mathrm{volume}^{-1}]\) (default: 0)

Allows setting the fugacity per particle type, e.g. 'A' refers to a depletant of type A.

Type:

TypeParameter [ particle type, float]

depletant_ntrial#

Multiplicative factor for the number of times a depletant is inserted. This factor is accounted for in the acceptance criterion so that detailed balance is unchanged. Higher values of ntrial (than one) can be used to reduce the variance of the free energy estimate and improve the acceptance rate of the Markov chain.

Type:

TypeParameter [particle type, int]

interaction_matrix#

Set to False for a pair of particle types to disable overlap checks between particles of those types (default: True).

Type:

TypeParameter [ tuple [particle type, particle type], bool]

translation_move_probability#

Fraction of moves to be selected as translation moves.

Type:

float

nselect#

Number of trial moves to perform per particle per timestep.

Type:

int

Attributes

__dir__()#

Expose all attributes for dynamic querying in notebooks and IDEs.

property counters#

Trial move counters.

The counter object has the following attributes:

  • translate: tuple [int, int] - Number of accepted and rejected translate trial moves.

  • rotate: tuple [int, int] - Number of accepted and rejected rotate trial moves.

  • overlap_checks: int - Number of overlap checks performed.

  • overlap_errors: int - Number of overlap checks that were too close to resolve.

Note

The counts are reset to 0 at the start of each hoomd.Simulation.run.

Type:

dict

property external_potential#

The user-defined potential energy field integrator.

Defines the external energy \(U_{\mathrm{external},i}\). Defaults to None. May be set to an object from hoomd.hpmc.external.

property is_tuning_complete#

Check if kernel parameter tuning is complete.

True when tuning is complete and kernel_parameters has locked optimal parameters for all active kernels used by this object.

Example:

while (not operation.is_tuning_complete):
    simulation.run(1000)
Type:

bool

property kernel_parameters#

Kernel parameters.

The dictionary maps GPU kernel names to tuples of integers that control how the kernel executes on the GPU. These values will change during the tuning process and remain static after tuning completes. Set the kernel parameters for one or more kernels to force specific values and stop tuning.

Warning

The keys and valid values in this dictionary depend on the hardware device, the HOOMD-blue version, and the value of class attributes. Keys and/or valid values may change even with new patch releases of HOOMD-blue.

Provided that you use the same HOOMD-blue binary on the same hardware and execute a script with the same parameters, you may save the tuned values from one run and load them in the next.

Examples:

kernel_parameters = operation.kernel_parameters
operation.kernel_parameters = kernel_parameters
Type:

dict[str, tuple[float]]

property loggables#

Name, category mapping of loggable quantities.

Type:

dict[str, str]

property map_overlaps#

List of overlapping particles.

The list contains one entry for each overlapping pair of particles. When a tuple (i,j) is present in the list, there is an overlap between the particles with tags i and j.

Attention

map_overlaps does not support MPI parallel simulations. It returns None when there is more than one MPI rank.

(Loggable: category=”sequence”)

Type:

list[tuple[int, int]]

property mps#

Number of trial moves performed per second.

Note

The count is reset at the start of each hoomd.Simulation.run.

(Loggable: category=”scalar”)

Type:

float

property overlaps#

Number of overlapping particle pairs.

(Loggable: category=”scalar”)

Type:

int

property pair_energy#

Total potential energy contributed by all pair potentials \([\mathrm{energy}]\).

(Loggable: category=”scalar”)

Type:

float

property pair_potential#

The user-defined pair potential.

Defines the pairwise particle interaction energy \(U_{\mathrm{pair},ij}\). Defaults to None. May be set to an object from hoomd.hpmc.pair.

property pair_potentials#

Pair potentials to apply.

Example

Type:

list[hoomd.hpmc.pair.Pair]

property rotate_moves#

Count of the accepted and rejected rotate moves.

Note

The counts are reset to 0 at the start of each hoomd.Simulation.run.

(Loggable: category=”sequence”)

Type:

tuple[int, int]

property translate_moves#

Count of the accepted and rejected translate moves.

Note

The counts are reset to 0 at the start of each hoomd.Simulation.run.

(Loggable: category=”sequence”)

Type:

tuple[int, int]

tune_kernel_parameters()#

Start tuning kernel parameters.

After calling tune_kernel_parameters, AutotunedObject will vary the kernel parameters in subsequent time steps, check the run time of each, and lock to the fastest performing parameters after the scan is complete.

Example:

operation.tune_kernel_parameters()
class hoomd.hpmc.integrate.Polyhedron(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Polyhedron hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of general polyhedra. The shape \(S\) contains the points inside the polyhedron defined by vertices and faces (see shape). Polyhedron supports triangle meshes and spheres only. The mesh must be free of self-intersections.

See also

Use ConvexPolyhedron for faster performance with convex polyhedra.

Note

This shape uses an internal OBB tree for fast collision queries. Depending on the number of constituent faces in the tree, different values of the number of faces per leaf node may yield different optimal performance. The capacity of leaf nodes is configurable.

Wall support.

Polyhedron supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.Polyhedron(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(vertices=[(-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)],
                    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]])
print('vertices = ', mc.shape["A"]["vertices"])
print('faces = ', mc.shape["A"]["faces"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • vertices (list [tuple [float, float, float]], required) - vertices of the polyhedron \([\mathrm{length}]\).

    • The origin MUST strictly be contained in the generally nonconvex volume defined by the vertices and faces.

    • The origin centered sphere that encloses all vertices should be of minimal size for optimal performance.

  • faces (list [tuple [int, int, int], required) - Vertex indices for every triangle in the mesh.

    • For visualization purposes, the faces MUST be defined with a counterclockwise winding order to produce an outward normal.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (float, default: 0.0) - radius of the sphere swept around the surface of the polyhedron \([\mathrm{length}]\). Set a non-zero sweep_radius to create a spheropolyhedron.

  • overlap (list [int], default: None) - Check for overlaps between faces when overlap [i] & overlap[j] is nonzero (& is the bitwise AND operator). When not None, overlap must have a length equal to that of faces. When None (the default), overlap is initialized with all 1’s.

  • capacity (int, default: 4) - set the maximum number of particles per leaf node to adjust performance.

  • origin (tuple [float, float, float], default: (0,0,0)) - a point strictly inside the shape, needed for correctness of overlap checks.

  • hull_only (bool, default: False) - When True, only check for intersections between the convex hulls.

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'Mesh', 'vertices': [[0.5, 0.5, 0.5], [0.5, -0.5, -0.5],
  [-0.5, 0.5, -0.5], [-0.5, -0.5, 0.5]],
  'faces': [[0, 1, 2], [0, 3, 1], [0, 2, 3], [1, 3, 2]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.SimplePolygon(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Simple polygon hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of simple polygons. The shape \(S\) of a simple polygon includes the points inside and on the surface of the simple polygon defined by the vertices (see shape). For example:

Example of a simple polygon with vertex labels.

Important

SimplePolygon simulations must be performed in 2D systems.

See also

Use ConvexPolygon for faster performance with convex polygons.

Wall support.

SimplePolygon supports no hoomd.wall geometries.

Examples:

mc = hpmc.integrate.SimplePolygon(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(vertices=[(0, 0.5),
                               (-0.5, -0.5),
                               (0, 0),
                               (0.5, -0.5)]);
print('vertices = ', mc.shape["A"]["vertices"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • vertices (list [tuple [float, float]], required) - vertices of the polygon \([\mathrm{length}]\).

    • Vertices MUST be specified in a counter-clockwise order.

    • The polygon may be concave, but edges must not cross.

    • The origin may be inside or outside the shape.

    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • sweep_radius (float, default: 0.0) - Ignored, but present because SimplePolygon shares data structures with ConvexSpheropolygon \([\mathrm{length}]\).

Warning

HPMC does not check that all vertex requirements are met. Undefined behavior will result when they are violated.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Example

>>> mc.type_shapes()
[{'type': 'Polygon', 'sweep_radius': 0,
  'vertices': [[-0.5, -0.5], [0.5, -0.5], [0.5, 0.5], [-0.5, 0.5]]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.Sphere(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Sphere hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of spheres. The shape \(S\) includes all points inside and on the surface of a sphere:

\[S = \left \{ \vec{r} : \frac{\vec{r}\cdot\vec{r}}{(d/2)^2} \le 1 \right\}\]

where \(d\), is the diameter set in shape. When the shape parameter orientable is False (the default), Sphere only applies translation trial moves and ignores translation_move_probability.

Tip

Use spheres with diameter=0 in conjunction with pair potentials for Monte Carlo simulations of particles with no hard core.

Tip

Use Sphere in a 2D simulation to perform Monte Carlo on hard disks.

Wall support.

Sphere supports all hoomd.wall geometries.

Examples:

mc = hoomd.hpmc.integrate.Sphere(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(diameter=1.0)
mc.shape["B"] = dict(diameter=2.0)
mc.shape["C"] = dict(diameter=1.0, orientable=True)
print('diameter = ', mc.shape["A"]["diameter"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • diameter (float, required) - Sphere diameter \([\mathrm{length}]\).

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

  • orientable (bool, default: False) - set to True to allow rotation moves on this particle type.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Examples

The types will be ‘Sphere’ regardless of dimensionality.

>>> mc.type_shapes
[{'type': 'Sphere', 'diameter': 1},
 {'type': 'Sphere', 'diameter': 2}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.SphereUnion(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Sphere union hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of unions of spheres. The union shape \(S\) is the set union of the given spheres:

\[S = \bigcup_k S_k(\mathbf{q}_k, \vec{r}_k)\]

Each constituent shape in the union has its own shape parameters \(S_k\), position \(\vec{r}_k\), and orientation \(\mathbf{q}_k\) (see shape).

Note

This shape uses an internal OBB tree for fast collision queries. Depending on the number of constituent spheres in the tree, different values of the number of spheres per leaf node may yield different performance. The capacity of leaf nodes is configurable.

Wall support.

SphereUnion supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.SphereUnion(default_d=0.3, default_a=0.4)
sphere1 = dict(diameter=1)
sphere2 = dict(diameter=2)
mc.shape["A"] = dict(shapes=[sphere1, sphere2],
                     positions=[(0, 0, 0), (0, 0, 1)],
                     orientations=[(1, 0, 0, 0), (1, 0, 0, 0)],
                     overlap=[1, 1])
print('diameter of the first sphere = ',
      mc.shape["A"]["shapes"][0]["diameter"])
print('center of the first sphere = ', mc.shape["A"]["positions"][0])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • shapes (list [dict], required) - Shape parameters for each sphere in the union. See Sphere.shape for the accepted parameters.

  • positions (list [tuple [float, float, float]], required) - Position of each sphere in the union. \([\mathrm{length}]\)

  • orientations (list [tuple [float, float, float, float]], default: None) - Orientation of each sphere in the union. When not None, orientations must have a length equal to that of positions. When None (the default), orientations is initialized with all [1,0,0,0]’s.

  • overlap (list [int], default: None) - Check for overlaps between constituent particles when overlap [i] & overlap[j] is nonzero (& is the bitwise AND operator). When not None, overlap must have a length equal to that of positions. When None (the default), overlap is initialized with all 1’s.

  • capacity (int, default: 4) - set the maximum number of particles per leaf node to adjust performance.

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Type:

TypeParameter [particle type, dict]

property type_shapes#

Description of shapes in type_shapes format.

Examples

The type will be ‘SphereUnion’ regardless of dimensionality.

>>> mc.type_shapes
[{'type': 'SphereUnion',
  'centers': [[0, 0, 0], [0, 0, 1]],
  'diameters': [1, 0.5]},
 {'type': 'SphereUnion',
  'centers': [[1, 2, 3], [4, 5, 6]],
  'diameters': [0.5, 1]}]

(Loggable: category=”object”)

Type:

list[dict]

class hoomd.hpmc.integrate.Sphinx(default_d=0.1, default_a=0.1, translation_move_probability=0.5, nselect=4)#

Bases: HPMCIntegrator

Sphinx hard particle Monte Carlo integrator.

Parameters:
  • default_d (float) – Default maximum size of displacement trial moves \([\mathrm{length}]\).

  • default_a (float) – Default maximum size of rotation trial moves \([\mathrm{dimensionless}]\).

  • translation_move_probability (float) – Fraction of moves that are translation moves.

  • nselect (int) – Number of trial moves to perform per particle per timestep.

Perform hard particle Monte Carlo of sphere unions and differences, depending on the sign of the diameter. The shape \(S\) is:

\[S = \left(\bigcup_{k,d_k\ge 0} S_k((1, 0, 0, 0), \vec{r}_k) \right) \setminus \left(\bigcup_{k,d_k < 0} S_k((1, 0, 0, 0), \vec{r}_k) \right)\]

Where \(d_k\) is the diameter given in shape, \(\vec{r}_k\) is the center given in shape and \(S_k\) is the set of points in a sphere or diameter \(|d_k|\).

Wall support.

Sphinx supports no hoomd.wall geometries.

Example:

mc = hpmc.integrate.Sphinx(default_d=0.3, default_a=0.4)
mc.shape["A"] = dict(centers=[(0,0,0),(1,0,0)], diameters=[1,.25])
print('diameters = ', mc.shape["A"]["diameters"])
shape#

The shape parameters for each particle type. The dictionary has the following keys:

  • diameters (list [float], required) - diameters of spheres (positive OR negative real numbers) \([\mathrm{length}]\).

  • centers (list [tuple [float, float, float], required) - centers of spheres in local coordinate frame \([\mathrm{length}]\).

  • ignore_statistics (bool, default: False) - set to True to ignore tracked statistics.

Type:

TypeParameter [particle type, dict]