Storing Particle Shape#

Overview#

Questions#

  • How can I store particle shape for use with visualization tools?

Objectives#

  • Demonstrate logging type_shapes to a GSD file.

  • Explain that OVITO can read this information.

Boilerplate code#

[1]:
import gsd.hoomd
import hoomd

Particle Shape#

HPMC integrators and some anisotropic MD pair potentials model particles that have a well defined shape. You can save this shape definition to a GSD file for use in analysis and visualization workflows. In particular, OVITO will read this shape information and render particles appropriately.

Define the Simulation#

This section executes the hard particle simulation from a previous tutorial. See Introducing HOOMD-blue for a complete description of this code.

[3]:
cpu = hoomd.device.CPU()
sim = hoomd.Simulation(device=cpu, seed=2)
mc = hoomd.hpmc.integrate.ConvexPolyhedron()
mc.shape['octahedron'] = dict(vertices=[
    (-0.5, 0, 0),
    (0.5, 0, 0),
    (0, -0.5, 0),
    (0, 0.5, 0),
    (0, 0, -0.5),
    (0, 0, 0.5),
])
sim.operations.integrator = mc
sim.create_state_from_gsd(
    filename='../00-Introducing-HOOMD-blue/compressed.gsd')
sim.run(0)

Logging particle shape to a GSD file#

The type_shapes loggable quantity is a representation of the particle shape for each type following the type_shapes specification for the GSD file format. In HPMC simulations, the integrator provides type_shapes:

[4]:
mc.loggables
[4]:
{'state': 'state',
 'map_overlaps': 'sequence',
 'overlaps': 'scalar',
 'translate_moves': 'sequence',
 'rotate_moves': 'sequence',
 'mps': 'scalar',
 'type_shapes': 'object'}
[5]:
mc.type_shapes
[5]:
[{'type': 'ConvexPolyhedron',
  'rounding_radius': 0,
  'vertices': [[-0.5, 0, 0],
   [0.5, 0, 0],
   [0, -0.5, 0],
   [0, 0.5, 0],
   [0, 0, -0.5],
   [0, 0, 0.5]]}]

Add the type_shapes quantity to a Logger.

[6]:
logger = hoomd.logging.Logger()
logger.add(mc, quantities=['type_shapes'])

Write the simulation trajectory to a GSD file along with the logged quantities:

[7]:
gsd_writer = hoomd.write.GSD(filename='trajectory.gsd',
                             trigger=hoomd.trigger.Periodic(10000),
                             mode='xb',
                             filter=hoomd.filter.All(),
                             log=logger)
sim.operations.writers.append(gsd_writer)

Run the simulation:

[8]:
sim.run(20000)

As discussed in a previous section, delete the simulation so it is safe to open the GSD file for reading in the same process.

[9]:
del sim, gsd_writer, logger, mc, cpu

Reading logged shapes from a GSD file#

You can access the shape from scripts using the gsd package:

[10]:
traj = gsd.hoomd.open('trajectory.gsd', 'rb')

type_shapes is a special quantity available via particles.type_shapes rather than the log dictionary:

[11]:
traj[0].particles.type_shapes
[11]:
[{'type': 'ConvexPolyhedron',
  'rounding_radius': 0,
  'vertices': [[-0.5, 0, 0],
   [0.5, 0, 0],
   [0, -0.5, 0],
   [0, 0.5, 0],
   [0, 0, -0.5],
   [0, 0, 0.5]]}]

Open the file in OVITO and it will read the shape definition and render particles appropriately.

In this section, you have logged particle shape to a GSD file during a simulation so that visualization and analysis tools can access it. The next section shows how to write formatted output.