HOOMD examples#

gsd.hoomd provides high-level access to HOOMD schema GSD files.

View the page source to find unformatted example code.

Import the module#

In [1]: import gsd.hoomd

Define a frame#

In [2]: frame = gsd.hoomd.Frame()

In [3]: frame.particles.N = 4

In [4]: frame.particles.types = ['A', 'B']

In [5]: frame.particles.typeid = [0,0,1,1]

In [6]: frame.particles.position = [[0,0,0],[1,1,1], [-1,-1,-1], [1,-1,-1]]

In [7]: frame.configuration.box = [3, 3, 3, 0, 0, 0]

gsd.hoomd.Frame stores the state of a single system configuration, or frame, in the file. Instantiate this class to create a system configuration. All fields default to None. Each field is written to the file when not None and when the data does not match the data in the first frame or defaults specified in the schema.

Create a hoomd gsd file#

In [8]: f = gsd.hoomd.open(name='file.gsd', mode='w')

Use gsd.hoomd.open to open a GSD file as a gsd.hoomd.HOOMDTrajectory instance.

Write frames to a gsd file#

In [9]: def create_frame(i):
   ...:     frame = gsd.hoomd.Frame()
   ...:     frame.configuration.step = i
   ...:     frame.particles.N = 4+i
   ...:     frame.particles.position = numpy.random.random(size=(4+i,3))
   ...:     return frame
   ...: 

In [10]: f = gsd.hoomd.open(name='example.gsd', mode='w')

In [11]: f.extend( (create_frame(i) for i in range(10)) )

In [12]: f.append( create_frame(10) )

In [13]: len(f)
Out[13]: 11

gsd.hoomd.HOOMDTrajectory is similar to a sequence of gsd.hoomd.Frame objects. The append and extend methods add frames to the trajectory.

Tip

When using extend, pass in a generator or generator expression to avoid storing the entire trajectory in memory before writing it out.

Randomly index frames#

In [14]: f = gsd.hoomd.open(name='example.gsd', mode='r')

In [15]: frame = f[5]

In [16]: frame.configuration.step
Out[16]: 5

In [17]: frame.particles.N
Out[17]: 9

In [18]: frame.particles.position
Out[18]: 
array([[0.37598178, 0.5648982 , 0.25179958],
       [0.8603388 , 0.8267959 , 0.05918318],
       [0.27709654, 0.7600235 , 0.3987156 ],
       [0.09090893, 0.74983144, 0.9653515 ],
       [0.17442188, 0.02355867, 0.20670001],
       [0.35886744, 0.8444683 , 0.71650344],
       [0.882093  , 0.33218572, 0.7594378 ],
       [0.829573  , 0.68029374, 0.648241  ],
       [0.16908108, 0.3694354 , 0.7125407 ]], dtype=float32)

gsd.hoomd.HOOMDTrajectory supports random indexing of frames in the file. Indexing into a trajectory returns a gsd.hoomd.Frame.

Slicing and selection#

Use the slicing operator to select individual frames or a subset of a trajectory.

In [19]: f = gsd.hoomd.open(name='example.gsd', mode='r')

In [20]: for frame in f[5:-2]:
   ....:     print(frame.configuration.step, end=' ')
   ....: 
5 6 7 8 
In [21]: every_2nd_frame = f[::2]  # create a view of a trajectory subset

In [22]: for frame in every_2nd_frame[:4]:
   ....:     print(frame.configuration.step, end=' ')
   ....: 
0 2 4 6 

Slicing a trajectory creates a trajectory view, which can then be queried for length or sliced again.

Pure python reader#

In [23]: f = gsd.pygsd.GSDFile(open('example.gsd', 'rb'))

In [24]: trajectory = gsd.hoomd.HOOMDTrajectory(f);

In [25]: trajectory[3].particles.position
Out[25]: 
array([[0.2050102 , 0.32875705, 0.549524  ],
       [0.6593699 , 0.5504239 , 0.8610174 ],
       [0.06630748, 0.70967233, 0.01630579],
       [0.5533299 , 0.866326  , 0.84913945],
       [0.74628955, 0.47101706, 0.7816548 ],
       [0.58378327, 0.01991365, 0.60347694],
       [0.55813897, 0.062085  , 0.26027662]], dtype=float32)

You can use GSD without needing to compile C code to read GSD files using gsd.pygsd.GSDFile in combination with gsd.hoomd.HOOMDTrajectory. It only supports the rb mode and does not read files as fast as the C implementation. It takes in a python file-like object, so it can be used with in-memory IO classes, and grid file classes that access data over the internet.

Warning

gsd.pygsd is slow. Use gsd.hoomd.open whenever possible.

Access logged data#

In [26]: with gsd.hoomd.open(name='log-example.gsd', mode='w') as f:
   ....:     frame = gsd.hoomd.Frame()
   ....:     frame.particles.N = 4
   ....:     for i in range(10):
   ....:         frame.configuration.step = i*100
   ....:         frame.log['particles/net_force'] = numpy.array([[-1,2,-3+i],
   ....:                                                            [0,2,-4],
   ....:                                                            [-3,2,1],
   ....:                                                            [1,2,3]],
   ....:                                                           dtype=numpy.float32)
   ....:         frame.log['value/potential_energy'] = 1.5+i
   ....:         f.append(frame)
   ....: 

Logged data is stored in the log dictionary as numpy arrays. Place data into this dictionary directly without the 'log/' prefix and gsd will include it in the output. Store per-particle quantities with the prefix particles/. Choose another prefix for other quantities.

In [27]: log = gsd.hoomd.read_log(name='log-example.gsd', scalar_only=True)

In [28]: list(log.keys())
Out[28]: ['configuration/step', 'log/value/potential_energy']

In [29]: log['log/value/potential_energy']
Out[29]: array([ 1.5,  2.5,  3.5,  4.5,  5.5,  6.5,  7.5,  8.5,  9.5, 10.5])

In [30]: log['configuration/step']
Out[30]: array([  0, 100, 200, 300, 400, 500, 600, 700, 800, 900], dtype=uint64)

Read logged data from the log dictionary.

Note

Logged data must be a convertible to a numpy array of a supported type.

In [31]: with gsd.hoomd.open(name='example.gsd', mode='w') as f:
   ....:     frame = gsd.hoomd.Frame()
   ....:     frame.particles.N = 4
   ....:     frame.log['invalid'] = dict(a=1, b=5)
   ....:     f.append(frame)
   ....: 
---------------------------------------------------------------------------
ValueError                                Traceback (most recent call last)
Cell In[31], line 5
      3 frame.particles.N = 4
      4 frame.log['invalid'] = dict(a=1, b=5)
----> 5 f.append(frame)

File ~/checkouts/readthedocs.org/user_builds/gsd/envs/stable/lib/python3.11/site-packages/gsd/hoomd.py:786, in HOOMDTrajectory.append(self, frame)
    784 # write log data
    785 for log, data in frame.log.items():
--> 786     self.file.write_chunk('log/' + log, data)
    788 self.file.end_frame()

File ~/checkouts/readthedocs.org/user_builds/gsd/envs/stable/lib/python3.11/site-packages/gsd/fl.pyx:611, in gsd.fl.GSDFile.write_chunk()

ValueError: invalid type for chunk: log/invalid

Use multiprocessing#

import multiprocessing

def count_particles(args):
   t, frame_idx = args
   return len(t[frame_idx].particles.position)

with gsd.hoomd.open(name='example.gsd', mode='r') as t:
   with multiprocessing.Pool(processes=multiprocessing.cpu_count()) as pool:
      result = pool.map(count_particles, [(t, frame_idx) for frame_idx in range(len(t))])

 result

gsd.hoomd.HOOMDTrajectory can be pickled when in read mode to allow for multiprocessing through Python’s multiprocessing library. Here count_particles finds the number of particles in each frame and appends it to a list.

Using the command line#

The GSD library provides a command line interface for reading files with first-class support for reading HOOMD GSD files. The CLI opens a Python interpreter with a file opened in a specified mode.

$ gsd read -s hoomd 'example.gsd'
...
File: example.gsd
Number of frames: 11

The GSD file handle is available via the "handle" variable.
For supported schema, you may access the trajectory using the "traj" variable.
Type "help(handle)" or "help(traj)" for more information.
The gsd and gsd.fl packages are always loaded.
Schema-specific modules (e.g. gsd.hoomd) are loaded if available.

>>> len(traj)
11
>>> traj[0].particles.position.shape == (4, 3)
True
>>> handle.read_chunk(0, 'particles/N')
array([4], dtype=uint32)