hpmc.integrate

Overview

hpmc.integrate.convex_polygon HPMC integration for convex polygons (2D).
hpmc.integrate.convex_polyhedron HPMC integration for convex polyhedra (3D).
hpmc.integrate.convex_spheropolygon HPMC integration for convex spheropolygons (2D).
hpmc.integrate.convex_spheropolyhedron HPMC integration for spheropolyhedra (3D).
hpmc.integrate.ellipsoid HPMC integration for ellipsoids (2D/3D).
hpmc.integrate.faceted_sphere HPMC integration for faceted spheres (3D).
hpmc.integrate.interaction_matrix Define pairwise interaction matrix
hpmc.integrate.mode_hpmc Base class HPMC integrator.
hpmc.integrate.polyhedron HPMC integration for general polyhedra (3D).
hpmc.integrate.simple_polygon HPMC integration for simple polygons (2D).
hpmc.integrate.sphere HPMC integration for spheres (2D/3D).
hpmc.integrate.sphere_union HPMC integration for unions of spheres (3D).
hpmc.integrate.sphinx HPMC integration for sphinx particles (3D).

Details

class hoomd.hpmc.integrate.convex_polygon(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4)

HPMC integration for convex polygons (2D).

Parameters:
  • seed (int) – Random number seed
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.

Note

For concave polygons, use simple_polygon.

Convex polygon parameters:

  • vertices (required) - vertices of the polygon as is a list of (x,y) tuples of numbers (distance units)

    • Vertices MUST be specified in a counter-clockwise order.
    • The origin MUST be contained within the vertices.
    • Points inside the polygon MUST NOT be included.
    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance (e.g. don’t put the origin right next to an edge).
  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

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

Examples:

mc = hpmc.integrate.convex_polygon(seed=415236)
mc = hpmc.integrate.convex_polygon(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', vertices=[(-0.5, -0.5), (0.5, -0.5), (0.5, 0.5), (-0.5, 0.5)]);
print('vertices = ', mc.shape_param['A'].vertices)
get_type_shapes()

Get all the types of shapes in the current simulation :returns: A list of dictionaries, one for each particle type in the system. Currently assumes that all 3D shapes are convex.

class hoomd.hpmc.integrate.convex_polyhedron(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False, max_verts=8)

HPMC integration for convex polyhedra (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – (Override the automatic choice for the number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.
  • max_verts (int) – Set the maximum number of vertices in a polyhedron.

Convex polyhedron parameters:

  • vertices (required) - vertices of the polyhedron as is a list of (x,y,z) tuples of numbers (distance units)

    • The origin MUST be contained within the vertices.
    • The origin centered circle that encloses all verticies should be of minimal size for optimal performance (e.g. don’t put the origin right next to a face).
  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

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

Example:

mc = hpmc.integrate.convex_polyhedron(seed=415236)
mc = hpmc.integrate.convex_polyhedron(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', 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_param['A'].vertices)

Depletants Example:

mc = hpmc.integrate.convex_polyhedron(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=3,depletant_type='B')
mc.shape_param.set('A', 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_param.set('B', vertices=[(0.05, 0.05, 0.05), (0.05, -0.05, -0.05), (-0.05, 0.05, -0.05), (-0.05, -0.05, 0.05)]);
get_type_shapes()

Get all the types of shapes in the current simulation :returns: A list of dictionaries, one for each particle type in the system. Currently assumes that all 3D shapes are convex.

class hoomd.hpmc.integrate.convex_spheropolygon(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4)

HPMC integration for convex spheropolygons (2D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.

Spheropolygon parameters:

  • vertices (required) - vertices of the polygon as is a list of (x,y) tuples of numbers (distance units)

    • The origin MUST be contained within the shape.
    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance (e.g. don’t put the origin right next to an edge).
  • sweep_radius (default: 0.0) - the radius of the sphere swept around the edges of the polygon (distance units) - optional

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Useful cases:

  • A 1-vertex spheropolygon is a disk.
  • A 2-vertex spheropolygon is a spherocylinder.

Warning

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

Examples:

mc = hpmc.integrate.convex_spheropolygon(seed=415236)
mc = hpmc.integrate.convex_spheropolygon(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', vertices=[(-0.5, -0.5), (0.5, -0.5), (0.5, 0.5), (-0.5, 0.5)], sweep_radius=0.1, ignore_statistics=False);
mc.shape_param.set('A', vertices=[(0,0)], sweep_radius=0.5, ignore_statistics=True);
print('vertices = ', mc.shape_param['A'].vertices)
get_type_shapes()

Get all the types of shapes in the current simulation :returns: A list of dictionaries, one for each particle type in the system. Currently assumes that all 3D shapes are convex.

class hoomd.hpmc.integrate.convex_spheropolyhedron(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False, max_verts=8)

HPMC integration for spheropolyhedra (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.
  • max_verts (int) – Set the maximum number of vertices in a polyhedron.

A sperholpolyhedron can also represent spheres (0 or 1 vertices), and spherocylinders (2 vertices).

Spheropolyhedron parameters:

  • vertices (required) - vertices of the polyhedron as is a list of (x,y,z) tuples of numbers (distance units)

    • The origin MUST be contained within the vertices.
    • The origin centered sphere that encloses all verticies should be of minimal size for optimal performance (e.g. don’t put the origin right next to a face).
    • A sphere can be represented by specifying zero vertices (i.e. vertices=[]) and a non-zero radius R
    • Two vertices and a non-zero radius R define a prolate spherocylinder.
  • sweep_radius (default: 0.0) - the radius of the sphere swept around the edges of the polygon (distance units) - optional

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

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

Example:

mc = hpmc.integrate.convex_spheropolyhedron(seed=415236)
mc = hpmc.integrate.convex_spheropolyhedron(seed=415236, d=0.3, a=0.4)
mc.shape_param['tetrahedron'].set(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_param['A'].vertices)
mc.shape_param['SphericalDepletant'].set(vertices=[], sweep_radius=0.1, ignore_statistics=True);

Depletants example:

mc = hpmc.integrate.convex_spheropolyhedron(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nR=3,depletant_type='SphericalDepletant')
mc.shape_param['tetrahedron'].set(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_param['SphericalDepletant'].set(vertices=[], sweep_radius=0.1);
class hoomd.hpmc.integrate.ellipsoid(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False)

HPMC integration for ellipsoids (2D/3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.

Ellipsoid parameters:

  • a (required) - principle axis a of the ellipsoid (radius in the x direction) (distance units)

  • b (required) - principle axis b of the ellipsoid (radius in the y direction) (distance units)

  • c (required) - principle axis c of the ellipsoid (radius in the z direction) (distance units)

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Example:

mc = hpmc.integrate.ellipsoid(seed=415236)
mc = hpmc.integrate.ellipsoid(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', a=0.5, b=0.25, c=0.125);
print('ellipsoids parameters (a,b,c) = ', mc.shape_param['A'].a, mc.shape_param['A'].b, mc.shape_param['A'].c)

Depletants Example:

mc = hpmc.integrate.ellipsoid(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=50,depletant_type='B')
mc.shape_param.set('A', a=0.5, b=0.25, c=0.125);
mc.shape_param.set('B', a=0.05, b=0.05, c=0.05);
class hoomd.hpmc.integrate.faceted_sphere(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False)

HPMC integration for faceted spheres (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.

A faceted sphere is a sphere interesected with halfspaces. The equation defining each halfspace is given by:

\[n_i\cdot r + b_i \le 0\]

where \(n_i\) is the face normal, and \(b_i\) is the offset.

Warning

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

Faceted sphere parameters:

  • normals (required) - list of (x,y,z) tuples defining the facet normals (distance units)

  • offsets (required) - list of offsets (distance unit^2)

  • diameter (required) - diameter of sphere

  • vertices (required) - list of vertices for intersection polyhedron

  • origin (required) - origin vector

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

Planes must not be coplanar.

Example:

mc = hpmc.integrate.faceted_sphere(seed=415236)
mc = hpmc.integrate.faceted_sphere(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', normals=[(-1,0,0),(1,0,0),(0,-1,0),(0,1,0),(0,0,-1),(0,0,1)],diameter=1.0);
print('diameter = ', mc.shape_param['A'].diameter)

Depletants Example:

mc = hpmc.integrate.pathcy_sphere(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=3,depletant_type='B')
mc.shape_param.set('A', normals=[(-1,0,0),(1,0,0),(0,-1,0),(0,1,0),(0,0,-1),(0,0,1)],diameter=1.0);
mc.shape_param.set('B', normals=[],diameter=0.1);
class hoomd.hpmc.integrate.interaction_matrix

Define pairwise interaction matrix

All shapes use interaction_matrix to define the interaction matrix between different pairs of particles indexed by type. The set of pair coefficients is a symmetric matrix defined over all possible pairs of particle types.

By default, all elements of the interaction matrix are 1, that means that overlaps are checked between all pairs of types. To disable overlap checking for a specific type pair, set the coefficient for that pair to 0.

Access the interaction matrix with a saved integrator object like so:

from hoomd import hpmc

mc = hpmc.integrate.some_shape(arguments...)
mv.overlap_checks.set('A', 'A', enable=False)
mc.overlap_checks.set('A', 'B', enable=True)
mc.overlap_checks.set('B', 'B', enable=False)

New in version 2.1.

set(a, b, enable)

Sets parameters for one type pair.

Parameters:
  • a (str) – First particle type in the pair (or a list of type names)
  • b (str) – Second particle type in the pair (or a list of type names)
  • enable – Set to True to enable overlap checks for this pair, False otherwise

By default, all interaction matrix elements are set to ‘True’.

It is not an error, to specify matrix elements for particle types that do not exist in the simulation.

There is no need to specify matrix elements for both pairs ‘A’, ‘B’ and ‘B’, ‘A’. Specifying only one is sufficient.

To set the same elements between many particle types, provide a list of type names instead of a single one. All pairs between the two lists will be set to the same parameters.

Examples:

interaction_matrix.set('A', 'A', False);
interaction_matrix.set('B', 'B', False);
interaction_matrix.set('A', 'B', True);
interaction_matrix.set(['A', 'B', 'C', 'D'], 'F', True);
interaction_matrix.set(['A', 'B', 'C', 'D'], ['A', 'B', 'C', 'D'], False);
class hoomd.hpmc.integrate.mode_hpmc(implicit)

Base class HPMC integrator.

mode_hpmc is the base class for all HPMC integrators. It provides common interface elements. Users should not instantiate this class directly. Methods documented here are available to all hpmc integrators.

count_overlaps()

Count the number of overlaps.

Returns:The number of overlaps in the current system configuration

Example:

mc = hpmc.integrate.shape(..);
mc.shape_param.set(....);
run(100)
num_overlaps = mc.count_overlaps();
get_a(type=None)

Get the maximum trial rotation.

Parameters:type (str) – Type name to query.
Returns:The current value of the ‘a’ parameter of the integrator.
get_configurational_bias_ratio()

Get the average ratio of configurational bias attempts to depletant insertion moves.

Returns:The average configurational bias ratio during the last hoomd.run().

Example:

mc = hpmc.integrate.shape(..,implicit=True);
mc.shape_param.set(....);
run(100)
cb_ratio = mc.get_configurational_bias_ratio();
get_counters()

Get all trial move counters.

Returns:A dictionary containing all trial moves counted during the last hoomd.run().

The dictionary contains the entries:

  • translate_accept_count - count of the number of accepted translate moves
  • translate_reject_count - count of the number of rejected translate moves
  • rotate_accept_count - count of the number of accepted rotate moves
  • rotate_reject_count - count of the number of rejected rotate moves
  • overlap_checks - estimate of the number of overlap checks performed
  • translate_acceptance - Average translate acceptance ratio over the run
  • rotate_acceptance - Average rotate acceptance ratio over the run
  • move_count - Count of the number of trial moves during the run
get_d(type=None)

Get the maximum trial displacement.

Parameters:type (str) – Type name to query.
Returns:The current value of the ‘d’ parameter of the integrator.
get_move_ratio()

Get the current probability of attempting translation moves.

Returns: The current value of the ‘move_ratio’ parameter of the integrator.

get_mps()

Get the number of trial moves per second.

Returns:The number of trial moves per second performed during the last hoomd.run().
get_nselect()

Get nselect parameter.

Returns:The current value of the ‘nselect’ parameter of the integrator.
get_ntrial()

Get ntrial parameter.

Returns:The current value of the ‘ntrial’ parameter of the integrator.
get_rotate_acceptance()

Get the average acceptance ratio for rotate moves.

Returns:The average rotate accept ratio during the last hoomd.run().

Example:

mc = hpmc.integrate.shape(..);
mc.shape_param.set(....);
run(100)
t_accept = mc.get_rotate_acceptance();
get_translate_acceptance()

Get the average acceptance ratio for translate moves.

Returns:The average translate accept ratio during the last hoomd.run().

Example:

mc = hpmc.integrate.shape(..);
mc.shape_param.set(....);
run(100)
t_accept = mc.get_translate_acceptance();
get_type_shapes()

Get all the types of shapes in the current simulation

Since this behaves differently for different types of shapes, the default behavior just raises an exception. Subclasses can override this to properly return

map_overlaps()

Build an overlap map of the system

Returns:List of tuples. True/false value of the i,j entry indicates overlap/non-overlap of the ith and jth particles (by tag)

Example

mc = hpmc.integrate.shape(…) mc.shape_param.set(…) overlap_map = np.asarray(mc.map_overlaps())

set_params(d=None, a=None, move_ratio=None, nselect=None, nR=None, depletant_type=None, ntrial=None)

Changes parameters of an existing integration mode.

Parameters:
  • d (float) – (if set) Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – (if set) Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – (if set) New value for the move ratio.
  • nselect (int) – (if set) New value for the number of particles to select for trial moves in one cell.
  • nR (int) – (if set) Implicit depletants only: Number density of implicit depletants in free volume.
  • depletant_type (str) – (if set) Implicit depletants only: Particle type to use as implicit depletant.
  • ntrial (int) – (if set) Implicit depletants only: Number of re-insertion attempts per overlapping depletant.
setup_pos_writer(pos, colors={})

Set pos_writer definitions for specified shape parameters.

Parameters:

setup_pos_writer() uses the shape_param settings to specify the shape definitions (via set_def) to the provided pos file writer. This overrides any previous values specified to hoomd.deprecated.dump.pos.set_def().

colors allows you to set per-type colors for particles. Specify colors as strings in the injavis format. When colors is not specified for a type, all colors default to 005984FF.

Examples:

mc = hpmc.integrate.shape(...);
mc.shape_param.set(....);
pos = pos_writer.dumpy.pos("dump.pos", period=100);
mc.setup_pos_writer(pos, colors=dict(A='005984FF'));
class hoomd.hpmc.integrate.polyhedron(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False)

HPMC integration for general polyhedra (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.

Polyhedron parameters:

  • vertices (required) - vertices of the polyhedron as is a list of (x,y,z) tuples of numbers (distance units)

    • The origin MUST strictly be contained in the generally nonconvex volume defined by the vertices and faces
    • The origin centered circle that encloses all verticies should be of minimal size for optimal performance (e.g. don’t put the origin right next to a face).
  • faces (required) - a list of vertex indices for every face

  • sweep_radius (default: 0.0) - rounding radius applied to polyhedron

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

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

Example:

mc = hpmc.integrate.polyhedron(seed=415236)
mc = hpmc.integrate.polyhedron(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', 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 = [(7, 3, 1, 5), (7, 5, 4, 6), (7, 6, 2, 3), (3, 2, 0, 1), (0, 2, 6, 4), (1, 0, 4, 5)]);
print('vertices = ', mc.shape_param['A'].vertices)
print('faces = ', mc.shape_param['A'].faces)

Depletants Example:

mc = hpmc.integrate.polyhedron(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=3,depletant_type='B')
faces = [(7, 3, 1, 5), (7, 5, 4, 6), (7, 6, 2, 3), (3, 2, 0, 1), (0, 2, 6, 4), (1, 0, 4, 5)];
mc.shape_param.set('A', 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 = faces);
mc.shape_param.set('B', vertices=[(-0.05, -0.05, -0.05), (-0.05, -0.05, 0.05), (-0.05, 0.05, -0.05), (-0.05, 0.05, 0.05), \
    (0.05, -0.05, -0.05), (0.05, -0.05, 0.05), (0.05, 0.05, -0.05), (0.05, 0.05, 0.05)], faces = faces);
class hoomd.hpmc.integrate.simple_polygon(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4)

HPMC integration for simple polygons (2D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.

Note

For simple polygons that are not concave, use convex_polygon, it will execute much faster than simple_polygon.

Simple polygon parameters:

  • vertices (required) - vertices of the polygon as is a list of (x,y) tuples of numbers (distance units)

    • Vertices MUST be specified in a counter-clockwise order.
    • The polygon may be concave, but edges must not cross.
    • The origin doesn’t necessarily need to be inside the shape.
    • The origin centered circle that encloses all vertices should be of minimal size for optimal performance.
  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Warning

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

Examples:

mc = hpmc.integrate.simple_polygon(seed=415236)
mc = hpmc.integrate.simple_polygon(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', vertices=[(0, 0.5), (-0.5, -0.5), (0, 0), (0.5, -0.5)]);
print('vertices = ', mc.shape_param['A'].vertices)
get_type_shapes()

Get all the types of shapes in the current simulation :returns: A list of dictionaries, one for each particle type in the system. Currently assumes that all 3D shapes are convex.

class hoomd.hpmc.integrate.sphere(seed, d=0.1, nselect=4, implicit=False)

HPMC integration for spheres (2D/3D).

Parameters:
  • seed (int) – Random number seed
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.

Hard particle Monte Carlo integration method for spheres.

Sphere parameters:

  • diameter (required) - diameter of the sphere (distance units)

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Examples:

mc = hpmc.integrate.sphere(seed=415236)
mc = hpmc.integrate.sphere(seed=415236, d=0.3)
mc.shape_param.set('A', diameter=1.0)
mc.shape_param.set('B', diameter=2.0)
print('diameter = ', mc.shape_param['A'].diameter)

Depletants Example:

mc = hpmc.integrate.sphere(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=8,nR=3,depletant_type='B')
mc.shape_param.set('A', diameter=1.0)
mc.shape_param.set('B', diameter=.1)
get_type_shapes()

Get all the types of shapes in the current simulation :returns: A list of dictionaries, one for each particle type in the system. Currently assumes that all 3D shapes are convex.

class hoomd.hpmc.integrate.sphere_union(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False, max_members=8)

HPMC integration for unions of spheres (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.
  • max_members (int) – Set the maximum number of members in the sphere union * .. versionadded:: 2.1

Sphere union parameters:

  • diameters (required) - list of diameters of the spheres (distance units).

  • centers (required) - list of centers of constituent spheres in particle coordinates.

  • overlap (default: 1 for all spheres) - only check overlap between constituent particles for which overlap [i] & overlap[j] is !=0, where ‘&’ is the bitwise AND operator.

    • New in version 2.1.

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking.

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Example:

mc = hpmc.integrate.sphere_union(seed=415236)
mc = hpmc.integrate.sphere_union(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', diameters=[1.0, 1.0], centers=[(-0.25, 0.0, 0.0), (0.25, 0.0, 0.0)]);
print('diameter of the first sphere = ', mc.shape_param['A'].members[0].diameter)
print('center of the first sphere = ', mc.shape_param['A'].centers[0])

Depletants Example:

mc = hpmc.integrate.sphere_union(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=50,depletant_type='B')
mc.shape_param.set('A', diameters=[1.0, 1.0], centers=[(-0.25, 0.0, 0.0), (0.25, 0.0, 0.0)]);
mc.shape_param.set('B', diameters=[0.05], centers=[(0.0, 0.0, 0.0)]);
class hoomd.hpmc.integrate.sphinx(seed, d=0.1, a=0.1, move_ratio=0.5, nselect=4, implicit=False)

HPMC integration for sphinx particles (3D).

Parameters:
  • seed (int) – Random number seed.
  • d (float) – Maximum move displacement, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • a (float) – Maximum rotation move, Scalar to set for all types, or a dict containing {type:size} to set by type.
  • move_ratio (float) – Ratio of translation moves to rotation moves.
  • nselect (int) – The number of trial moves to perform in each cell.
  • implicit (bool) – Flag to enable implicit depletants.

Sphinx particles are dimpled spheres (spheres with ‘positive’ and ‘negative’ volumes).

Sphinx parameters:

  • diameters - diameters of spheres (positive OR negative real numbers)

  • centers - centers of spheres in local coordinate frame

  • ignore_statistics (default: False) - set to True to disable ignore for statistics tracking

  • ignore_overlaps (default: False) - set to True to disable overlap checks between this and other types with ignore_overlaps=True

Quick Example:

mc = hpmc.integrate.sphinx(seed=415236)
mc = hpmc.integrate.sphinx(seed=415236, d=0.3, a=0.4)
mc.shape_param.set('A', centers=[(0,0,0),(1,0,0)], diameters=[1,.25])
print('diameters = ', mc.shape_param['A'].diameters)

Depletants Example:

mc = hpmc.integrate.sphinx(seed=415236, d=0.3, a=0.4, implicit=True)
mc.set_param(nselect=1,nR=3,depletant_type='B')
mc.shape_param.set('A', centers=[(0,0,0),(1,0,0)], diameters=[1,-.25])
mc.shape_param.set('B', centers=[(0,0,0)], diameters=[.15])