Posterior Samples I/O (ligo.skymap.io.hdf5)

Read HDF5 posterior sample chain HDF5 files.

ligo.skymap.io.hdf5.read_samples(filename, path=None, tablename='posterior_samples')[source] [edit on github]

Read an HDF5 sample chain file.

Parameters:
filenamestr

The path of the HDF5 file on the filesystem.

pathstr, optional

The path of the dataset within the HDF5 file.

tablenamestr, optional

The name of table to search for recursively within the HDF5 file. By default, search for ‘posterior_samples’.

Returns:
chainastropy.table.Table

The sample chain as an Astropy table.

Examples

Test reading a file written using the Python API:

>>> import os.path
>>> import tempfile
>>> table = Table([
...     Column(np.ones(10), name='foo', meta={'vary': FIXED}),
...     Column(np.arange(10), name='bar', meta={'vary': LINEAR}),
...     Column(np.arange(10) * np.pi, name='bat', meta={'vary': CIRCULAR}),
...     Column(np.arange(10), name='baz', meta={'vary': OUTPUT})
... ])
>>> with tempfile.TemporaryDirectory() as dir:
...     filename = os.path.join(dir, 'test.hdf5')
...     write_samples(table, filename, path='foo/bar/posterior_samples')
...     len(read_samples(filename))
10

Test reading a file that was written using the LAL HDF5 C API:

>>> from importlib.resources import files
>>> with files('ligo.skymap.io.tests.data').joinpath(
...         'test.hdf5').open('rb') as f:
...     table = read_samples(f)
>>> table.colnames
['uvw', 'opq', 'lmn', 'ijk', 'def', 'abc', 'ghi', 'rst']
ligo.skymap.io.hdf5.write_samples(table, filename, metadata=None, **kwargs)[source] [edit on github]

Write an HDF5 sample chain file.

Parameters:
tableastropy.table.Table

The sample chain as an Astropy table.

filenamestr

The path of the HDF5 file on the filesystem.

metadata: dict (optional)

Dictionary of (path, value) pairs of metadata attributes to add to the output file

kwargs: dict

Any keyword arguments for astropy.table.Table.write.

Examples

Check that we catch columns that are supposed to be FIXED but are not:

>>> table = Table([
...     Column(np.arange(10), name='foo', meta={'vary': FIXED})
... ])
>>> write_samples(table, 'bar.hdf5', 'bat/baz')
Traceback (most recent call last):
    ...
AssertionError: 
Arrays are not equal
Column foo is a fixed column, but its values are not identical
...

And now try writing an arbitrary example to a temporary file.

>>> import os.path
>>> import tempfile
>>> table = Table([
...     Column(np.ones(10), name='foo', meta={'vary': FIXED}),
...     Column(np.arange(10), name='bar', meta={'vary': LINEAR}),
...     Column(np.arange(10) * np.pi, name='bat', meta={'vary': CIRCULAR}),
...     Column(np.arange(10), name='baz', meta={'vary': OUTPUT}),
...     Column(np.ones(10), name='plugh'),
...     Column(np.arange(10), name='xyzzy')
... ])
>>> with tempfile.TemporaryDirectory() as dir:
...     write_samples(
...         table, os.path.join(dir, 'test.hdf5'), path='bat/baz',
...         metadata={'bat/baz': {'widget': 'shoephone'}})