hpmc.field

Overview

hpmc.field.callback Use a python-defined energy function in MC integration
hpmc.field.external_field_composite Manage multiple external fields.
hpmc.field.frenkel_ladd_energy Compute the Frenkel-Ladd Energy of a crystal.
hpmc.field.lattice_field Restrain particles on a lattice
hpmc.field.wall Manage walls (an external field type).

Details

Apply external fields to HPMC simulations.

class hoomd.hpmc.field.callback(mc, energy_function, composite=False)

Use a python-defined energy function in MC integration

Parameters:
  • mc (hoomd.hpmc.integrate) – MC integrator.
  • callback (callable) – A python function to evaluate the energy of a configuration
  • composite (bool) – True if this evaluator is part of a composite external field

Example:

def energy(snapshot):
    # evaluate the energy in a linear potential gradient along the x-axis
    gradient = (5,0,0)
    e = 0
    for p in snap.particles.position:
        e -= numpy.dot(gradient,p)
    return e

mc = hpmc.integrate.sphere(seed=415236);
mc.shape_param.set('A',diameter=1.0)
hpmc.field.callback(mc=mc, energy_function=energy);
run(100)
disable()

Disables the compute.

Examples:

c.disable()

Executing the disable command will remove the compute from the system. Any hoomd.run() command executed after disabling a compute will not be able to log computed values with hoomd.analyze.log.

A disabled compute can be re-enabled with enable().

enable()

Enables the compute.

Examples:

c.enable()

See disable().

restore_state()

Restore the state information from the file used to initialize the simulations

class hoomd.hpmc.field.external_field_composite(mc, fields=None)

Manage multiple external fields.

Parameters:
  • mc (hoomd.hpmc.integrate) – MC integrator (don’t specify a new integrator later, external_field_composite will continue to use the old one)
  • fields (list) – List of external fields to combine together.

external_field_composite allows the user to create and compute multiple external fields. Once created use add_field() to add a new field.

Once initialized, the compute provides a log quantities that other external fields create. See those external fields to find the quantities.

Examples:

mc = hpmc.integrate.shape(...);
walls = hpmc.field.walls(...)
lattice = hpmc.field.lattice(...)
composite_field = hpmc.field.external_field_composite(mc, fields=[walls, lattice])
add_field(fields)

Add an external field to the ensemble.

Parameters:fields (list) – list of fields to add

Example:

mc = hpmc.integrate.shape(...);
composite_field = hpmc.compute.external_field_composite(mc)
walls = hpmc.compute.walls(..., setup=False)
lattice = hpmc.compute.lattice(..., setup=False)
composite_field.add_field(fields=[walls, lattice])
disable()

Disables the compute.

Examples:

c.disable()

Executing the disable command will remove the compute from the system. Any hoomd.run() command executed after disabling a compute will not be able to log computed values with hoomd.analyze.log.

A disabled compute can be re-enabled with enable().

enable()

Enables the compute.

Examples:

c.enable()

See disable().

restore_state()

Restore the state information from the file used to initialize the simulations

class hoomd.hpmc.field.frenkel_ladd_energy(mc, ln_gamma, q_factor, r0, q0, drift_period, symmetry=[])

Compute the Frenkel-Ladd Energy of a crystal.

Parameters:
  • ln_gamma (float) – log of the translational spring constant
  • q_factor (float) – scale factor between the translational spring constant and rotational spring constant
  • r0 (list) – reference lattice positions
  • q0 (list) – reference lattice orientations
  • drift_period (int) – period call the remove drift updater

frenkel_ladd_energy interacts with lattice_field and hoomd.hpmc.update.remove_drift.

Once initialized, the compute provides the log quantities from the lattice_field.

Warning

The lattice energies and standard deviations logged by lattice_field are multiplied by the spring constant. As a result, when computing the free energies from frenkel_ladd_energy class, instead of integrating the free energy over the spring constants, you should integrate over the natural log of the spring constants.

Example:

mc = hpmc.integrate.convex_polyhedron(seed=seed);
mc.shape_param.set("A", vertices=verts)
mc.set_params(d=0.005, a=0.005)
#set the FL parameters
fl = hpmc.compute.frenkel_ladd_energy(mc=mc, ln_gamma=0.0, q_factor=10.0, r0=rs, q0=qs, drift_period=1000)
disable()

Disables the compute.

Examples:

c.disable()

Executing the disable command will remove the compute from the system. Any hoomd.run() command executed after disabling a compute will not be able to log computed values with hoomd.analyze.log.

A disabled compute can be re-enabled with enable().

enable()

Enables the compute.

Examples:

c.enable()

See disable().

reset_statistics()

Reset the statistics counters.

Example:

mc = hpmc.integrate.sphere(seed=415236);
fl = hpmc.compute.frenkel_ladd_energy(mc=mc, ln_gamma=0.0, q_factor=10.0, r0=rs, q0=qs, drift_period=1000)
ks = np.linspace(1000, 0.01, 100);
for k in ks:
  fl.set_params(ln_gamma=math.log(k), q_factor=10.0);
  fl.reset_statistics();
  run(1000)
restore_state()

Restore the state information from the file used to initialize the simulations

set_params(ln_gamma=None, q_factor=None)

Set the Frenkel-Ladd parameters.

Parameters:
  • ln_gamma (float) – log of the translational spring constant
  • q_factor (float) – scale factor between the translational spring constant and rotational spring constant

Example:

mc = hpmc.integrate.sphere(seed=415236);
fl = hpmc.compute.frenkel_ladd_energy(mc=mc, ln_gamma=0.0, q_factor=10.0, r0=rs, q0=qs, drift_period=1000)
ks = np.linspace(1000, 0.01, 100);
for k in ks:
  fl.set_params(ln_gamma=math.log(k), q_factor=10.0);
  fl.reset_statistics();
  run(1000)
class hoomd.hpmc.field.lattice_field(mc, position=[], orientation=[], k=0.0, q=0.0, symmetry=[], composite=False)

Restrain particles on a lattice

Parameters:
  • mc (hoomd.hpmc.integrate) – MC integrator.
  • position (list) – list of positions to restrain each particle (distance units).
  • orientation (list) – list of orientations to restrain each particle (quaternions).
  • k (float) – translational spring constant.
  • q (float) – rotational spring constant.
  • symmetry (list) – list of equivalent quaternions for the shape.
  • composite (bool) – Set this to True when this field is part of a external_field_composite.

lattice_field specifies that a harmonic spring is added to every particle:

\[\begin{split}V_{i}(r) = k_r*(r_i-r_{oi})^2 \\ V_{i}(q) = k_q*(q_i-q_{oi})^2\end{split}\]

Note

1/2 is not included in the formulas, specify your spring constants accordingly.

  • \(k_r\) - translational spring constant.
  • \(r_{o}\) - lattice positions (in distance units).
  • \(k_q\) - rotational spring constant.
  • \(q_{o}\) - lattice orientations (quaternion)

Once initialized, the compute provides the following log quantities that can be logged via analyze.log:

  • lattice_energy – total lattice energy
  • lattice_energy_pp_avg – average lattice energy per particle multiplied by the spring constant
  • lattice_energy_pp_sigma – standard deviation of the lattice energy per particle multiplied by the spring constant
  • lattice_translational_spring_constant – translational spring constant
  • lattice_rotational_spring_constant – rotational spring constant
  • lattice_num_samples – number of samples used to compute the average and standard deviation

Warning

The lattice energies and standard deviations logged by this class are multiplied by the spring constant.

Example:

mc = hpmc.integrate.sphere(seed=415236);
hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=1000.0);
log = analyze.log(quantities=['lattice_energy'], period=100, filename='log.dat', overwrite=True);
disable()

Disables the compute.

Examples:

c.disable()

Executing the disable command will remove the compute from the system. Any hoomd.run() command executed after disabling a compute will not be able to log computed values with hoomd.analyze.log.

A disabled compute can be re-enabled with enable().

enable()

Enables the compute.

Examples:

c.enable()

See disable().

get_average_energy()
Get the average energy per particle of the lattice field.
This is a collective call and must be called on all ranks.
Example::
mc = hpmc.integrate.sphere(seed=415236); lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=exp(15)); run(20000) avg_eng = lattice.get_average_energy() // should be about 1.5kT
get_energy()
Get the current energy of the lattice field.
This is a collective call and must be called on all ranks.
Example::
mc = hpmc.integrate.sphere(seed=415236); lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=1000.0); run(20000) eng = lattice.get_energy()
get_sigma_energy()
Gives the standard deviation of the average energy per particle of the lattice field.
This is a collective call and must be called on all ranks.
Example::
mc = hpmc.integrate.sphere(seed=415236); lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=exp(15)); run(20000) sig_eng = lattice.get_sigma_energy()
reset(timestep=None)

Reset the statistics counters.

Parameters:timestep (int) – the timestep to pass into the reset function.

Example:

mc = hpmc.integrate.sphere(seed=415236);
lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=1000.0);
ks = np.linspace(1000, 0.01, 100);
for k in ks:
  lattice.set_params(k=k, q=0.0);
  lattice.reset();
  run(1000)
restore_state()

Restore the state information from the file used to initialize the simulations

set_params(k, q)

Set the translational and rotational spring constants.

Parameters:
  • k (float) – translational spring constant.
  • q (float) – rotational spring constant.

Example:

mc = hpmc.integrate.sphere(seed=415236);
lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=1000.0);
ks = np.linspace(1000, 0.01, 100);
for k in ks:
  lattice.set_params(k=k, q=0.0);
  run(1000)
set_references(position=[], orientation=[])

Reset the reference positions or reference orientations.

Parameters:
  • position (list) – list of positions to restrain each particle.
  • orientation (list) – list of orientations to restrain each particle.

Example:

mc = hpmc.integrate.sphere(seed=415236);
lattice = hpmc.field.lattice_field(mc=mc, position=fcc_lattice, k=1000.0);
lattice.set_references(position=bcc_lattice)
class hoomd.hpmc.field.wall(mc, composite=False)

Manage walls (an external field type).

Parameters:

wall allows the user to implement one or more walls. If multiple walls are added, then particles are confined by the INTERSECTION of all of these walls. In other words, particles are confined by all walls if they independently satisfy the confinement condition associated with each separate wall. Once you’ve created an instance of this class, use add_sphere_wall() to add a new spherical wall, add_cylinder_wall() to add a new cylindrical wall, or add_plane_wall() to add a new plane wall.

Specialized overlap checks have been written for supported combinations of wall types and particle shapes. These combinations are: * Sphere particles: sphere walls, cylinder walls, plane walls * Convex polyhedron particles: sphere walls, cylinder walls, plane walls * Convex spheropolyhedron particles: sphere walls

Once initialized, the compute provides the following log quantities that can be logged via hoomd.analyze.log:

  • hpmc_wall_volume : the volume associated with the intersection of implemented walls. This number is only meaningful if the user has initially provided it through set_volume(). It will subsequently change when the box is resized and walls are scaled appropriately.
  • hpmc_wall_sph_rsq-i : the squared radius of the spherical wall indexed by i, beginning at 0 in the order the sphere walls were added to the system.
  • hpmc_wall_cyl_rsq-i : the squared radius of the cylindrical wall indexed by i, beginning at 0 in the order the cylinder walls were added to the system.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.set_volume(4./3.*np.pi);
log = analyze.log(quantities=['hpmc_wall_volume','hpmc_wall_sph_rsq-0'], period=100, filename='log.dat', overwrite=True);
add_cylinder_wall(radius, origin, orientation, inside=True)

Add a cylindrical wall to the simulation.

Parameters:
  • radius (float) – radius of cylindrical wall
  • origin (tuple) – origin (center) of cylindrical wall
  • orientation (tuple) – vector that defines the direction of the long axis of the cylinder. will be normalized automatically by hpmc.
  • inside (bool) – When True, then particles are CONFINED by the wall if they exist entirely inside the cylinder (in the portion of connected space that contains the origin). When False, then particles are CONFINED by the wall if they exist entirely outside the cylinder (in the portion of connected space that does not contain the origin). DEFAULTS to True.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_cylinder_wall(radius = 1.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
add_plane_wall(normal, origin)

Add a plane wall to the simulation.

Parameters:
  • normal (tuple) – vector normal to the plane. this, in combination with a point on the plane, defines the plane entirely. It will be normalized automatically by hpmc. The direction of the normal vector defines the confinement condition associated with the plane wall. If every part of a particle exists in the halfspace into which the normal points, then that particle is CONFINED by the plane wall.
  • origin (tuple) – a point on the plane wall. this, in combination with the normal vector, defines the plane entirely.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_plane_wall(normal = [0, 0, 1], origin = [0, 0, 0]);
add_sphere_wall(radius, origin, inside=True)

Add a spherical wall to the simulation.

Parameters:
  • radius (float) – radius of spherical wall
  • origin (tuple) – origin (center) of spherical wall.
  • inside (bool) – When True, particles are CONFINED by the wall if they exist entirely inside the sphere (in the portion of connected space that contains the origin). When False, then particles are CONFINED by the wall if they exist entirely outside the sphere (in the portion of connected space that does not contain the origin).

Quick Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
count_overlaps(exit_early=False)

Count the overlaps associated with the walls.

Parameters:exit_early (bool) – When True, stop counting overlaps after the first one is found.
Returns:The number of overlaps associated with the walls

A particle “overlaps” with a wall if it fails to meet the confinement condition associated with the wall.

Example

mc = hpmc.integrate.sphere(seed = 415236); ext_wall = hpmc.compute.wall(mc); ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True); run(100) num_overlaps = ext_wall.count_overlaps();

disable()

Disables the compute.

Examples:

c.disable()

Executing the disable command will remove the compute from the system. Any hoomd.run() command executed after disabling a compute will not be able to log computed values with hoomd.analyze.log.

A disabled compute can be re-enabled with enable().

enable()

Enables the compute.

Examples:

c.enable()

See disable().

get_curr_box()

Get the simulation box that the wall class is currently storing.

Returns:The boxdim object that the wall class is currently storing.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.set_volume(4./3.*np.pi);
run(100)
curr_box = ext_wall.get_curr_box();
get_cylinder_wall_param(index, param)

Access a parameter associated with a particular cylinder wall.

Parameters:
  • index (int) – index of the cylinder wall to be accessed. indices begin at 0 in the order the cylinder walls were added to the system.
  • param (str) – name of parameter to be accessed. options are “rsq” (squared radius of cylinder wall), “origin” (origin of cylinder wall), “orientation” (orientation of cylinder wall), and “inside” (confinement condition associated with cylinder wall).
Returns:

Value of queried parameter.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_cylinder_wall(radius = 1.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
rsq = ext_wall.get_cylinder_wall_param(index = 0, param = "rsq");
get_num_cylinder_walls()

Get the current number of cylinder walls in the simulation.

Returns:The current number of cylinder walls in the simulation.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_cylinder_wall(radius = 1.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
num_cyl_walls = ext_wall.get_num_cylinder_walls();
get_num_plane_walls()

Get the current number of plane walls in the simulation.

Returns:The current number of plane walls in the simulation.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_plane_wall(normal = [0, 0, 1], origin = [0, 0, 0]);
num_plane_walls = ext_wall.get_num_plane_walls();
get_num_sphere_walls()

Get the current number of sphere walls in the simulation.

Returns: the current number of sphere walls in the simulation

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
num_sph_walls = ext_wall.get_num_sphere_walls();
get_plane_wall_param(index, param)

Access a parameter associated with a particular plane wall.

Parameters:
  • index (int) – index of the plane wall to be accessed. indices begin at 0 in the order the plane walls were added to the system.
  • param (str) – name of parameter to be accessed. options are “normal” (vector normal to the plane wall), and “origin” (point on the plane wall)
Returns:

Value of queried parameter.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_plane_wall(normal = [0, 0, 1], origin = [0, 0, 0]);
n = ext_wall.get_plane_wall_param(index = 0, param = "normal");
get_sphere_wall_param(index, param)

Access a parameter associated with a particular sphere wall.

Parameters:
  • index (int) – index of the sphere wall to be accessed. indices begin at 0 in the order the sphere walls were added to the system.
  • param (str) – name of parameter to be accessed. options are “rsq” (squared radius of sphere wall), “origin” (origin of sphere wall), and “inside” (confinement condition associated with sphere wall)
Returns:

Value of queried parameter.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
rsq = ext_wall.get_sphere_wall_param(index = 0, param = "rsq");
get_volume()

Get the current volume associated with the intersection of all walls in the system.

If this quantity has not previously been set by the user, this returns a meaningless value.

Returns:The current volume associated with the intersection of all walls in the system.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.set_volume(4./3.*np.pi);
run(100)
curr_vol = ext_wall.get_volume();
remove_cylinder_wall(index)

Remove a particular cylinder wall from the simulation.

Parameters:index (int) – index of the cylinder wall to be removed. indices begin at 0 in the order the cylinder walls were added to the system.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_cylinder_wall(radius = 1.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
ext_wall.remove_cylinder_wall(index = 0);
remove_plane_wall(index)

Remove a particular plane wall from the simulation.

Parameters:index (int) – index of the plane wall to be removed. indices begin at 0 in the order the plane walls were added to the system.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_plane_wall(normal = [0, 0, 1], origin = [0, 0, 0]);
ext_wall.remove_plane_wall(index = 0);
remove_sphere_wall(index)

Remove a particular sphere wall from the simulation.

Parameters:index (int) – index of the sphere wall to be removed. indices begin at 0 in the order the sphere walls were added to the system.

Quick Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.remove_sphere_wall(index = 0);
restore_state()

Restore the state information from the file used to initialize the simulations

set_curr_box(Lx=None, Ly=None, Lz=None, xy=None, xz=None, yz=None)

Set the simulation box that the wall class is currently storing.

You may want to set this independently so that you can cleverly control whether or not the walls actually scale in case you manually resize your simulation box. The walls scale automatically when they get the signal that the global box, associated with the system definition, has scaled. They do so, however, with a scale factor associated with the ratio of the volume of the global box to the volume of the box that the walls class is currently storing. (After the scaling the box that the walls class is currently storing is updated appropriately.) If you want to change the simulation box WITHOUT scaling the walls, then, you must first update the simulation box that the walls class is storing, THEN update the global box associated with the system definition.

Example:

init_box = hoomd.data.boxdim(L=10, dimensions=3);
snap = hoomd.data.make_snapshot(N=1, box=init_box, particle_types=['A']);
system = hoomd.init.read_snapshot(snap);
system.particles[0].position = [0,0,0];
system.particles[0].type = 'A';
mc = hpmc.integrate.sphere(seed = 415236);
mc.shape_param.set('A', diameter = 2.0);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 3.0, origin = [0, 0, 0], inside = True);
ext_wall.set_curr_box(Lx=2.0*init_box.Lx, Ly=2.0*init_box.Ly, Lz=2.0*init_box.Lz, xy=init_box.xy, xz=init_box.xz, yz=init_box.yz);
system.sysdef.getParticleData().setGlobalBox(ext_wall.get_curr_box()._getBoxDim())
set_cylinder_wall(index, radius, origin, orientation, inside=True)

Change the parameters associated with a particular cylinder wall.

Parameters:
  • index (int) – index of the cylinder wall to be modified. indices begin at 0 in the order the cylinder walls were added to the system.
  • radius (float) – New radius of cylindrical wall
  • origin (tuple) – New origin (center) of cylindrical wall
  • orientation (tuple) – New vector that defines the direction of the long axis of the cylinder. will be normalized automatically by hpmc.
  • inside (bool) – New confinement condition. When True, then particles are CONFINED by the wall if they exist entirely inside the cylinder (in the portion of connected space that contains the origin). When False, then particles are CONFINED by the wall if they exist entirely outside the cylinder (in the portion of connected space that does not contain the origin). DEFAULTS to True.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_cylinder_wall(radius = 1.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
ext_wall.set_cylinder_wall(index = 0, radius = 3.0, origin = [0, 0, 0], orientation = [0, 0, 1], inside = True);
set_plane_wall(index, normal, origin)

Change the parameters associated with a particular plane wall.

Parameters:
  • index (int) – index of the plane wall to be modified. indices begin at 0 in the order the plane walls were added to the system.
  • normal (tuple) – new vector normal to the plane. this, in combination with a point on the plane, defines the plane entirely. It will be normalized automatically by hpmc. The direction of the normal vector defines the confinement condition associated with the plane wall. If every part of a particle exists in the halfspace into which the normal points, then that particle is CONFINED by the plane wall.
  • origin (tuple) – new point on the plane wall. this, in combination with the normal vector, defines the plane entirely.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_plane_wall(normal = [0, 0, 1], origin = [0, 0, 0]);
ext_wall.set_plane_wall(index = 0, normal = [0, 0, 1], origin = [0, 0, 1]);
set_sphere_wall(index, radius, origin, inside=True)

Change the parameters associated with a particular sphere wall.

Parameters:
  • index (int) – index of the sphere wall to be modified. indices begin at 0 in the order the sphere walls were added to the system.
  • radius (float) – New radius of spherical wall
  • origin (tuple) – New origin (center) of spherical wall.
  • inside (bool) – New confinement condition. When True, particles are CONFINED by the wall if they exist entirely inside the sphere (in the portion of connected space that contains the origin). When False, then particles are CONFINED by the wall if they exist entirely outside the sphere (in the portion of connected space that does not contain the origin).

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.set_sphere_wall(index = 0, radius = 3.0, origin = [0, 0, 0], inside = True);
set_volume(volume)

Set the volume associated with the intersection of all walls in the system.

This number will subsequently change when the box is resized and walls are scaled appropriately.

Example:

mc = hpmc.integrate.sphere(seed = 415236);
ext_wall = hpmc.compute.wall(mc);
ext_wall.add_sphere_wall(radius = 1.0, origin = [0, 0, 0], inside = True);
ext_wall.set_volume(4./3.*np.pi);