The universal read function

Builtin functions

pesummay offers a read function which allows for configuration files, injection files, prior files and nearly result files formats to be read in. For a full list of result files that can be read in with the core module see here. For a full list of result files that can be read in with the gw module see here. Below we show a few examples. of how the read function works.

First we import the universal read function,

>>> from pesummary.io import read

A dat file containing a set of posterior samples can then be read in with,

>>> data = read("example.dat")

Of course, if you would like, you may specify the package that you wish to use when reading the file. This is done with the package kwarg,

>>> data = read("example.dat", package="core")

For details about the returned posterior samples object, see the core package tutorial or the gw package tutorial.

As with the above, a config file can be read in with,

>>> config = read("config.ini")

This will then return a dictionary containing the section headers and various settings. If your config file does not have a ini extension, this can still be read in with,

>>> config = read("config.txt", file_format="ini")

A skymap may also be read in with the following,

>>> skymap = read("skymap.fits")

Of course, if your skymap does not have the fits extension, this skymap can still be read in with,

>>> skymap = read("skymap.gz", skymap=True)

A frame file containing gravitational wave strain data can be read in with the following,

>>> strain = read("frame_file.gwf", channel="channel")

All files with a gwf or lcf extension are treated as frame files. If your frame file does not have this extension, this frame file can still be read in with,

>>> strain = read("frame_file.gwf", channel="channel", strain=True)

strain is a pesummary.gw.file.strain.StrainData object. For details about this object see Strain Data in PESummary.

Custom functions

Of course, you might have a file in a format which pesummary is unable to read in with the inbuilt functions. As a result of the modularity of pesummary, we may define a class which is capable of reading in this custom file format and pass it as a kwarg to the universal read function. Below we show a couple of examples,

# Here we show how to read in a result file which cannot be read in using
# the pesummary builtin functions. This result file only contains a single
# analysis
import numpy as np
import h5py
from pesummary.core.file.formats.base_read import SingleAnalysisRead
from pesummary.io import read

# First let us create a dictionary containing some posterior samples
data = {
    "a": np.random.uniform(1, 5, 1000),
    "b": np.random.uniform(1, 5, 1000)
}

# Next, lets make a file which cannot be read in with pesummary
f = h5py.File("example.h5", "w")
samples = f.create_group("my_posterior_samples")
a = samples.create_group("a")
b = samples.create_group("b")
a.create_dataset("combined", data=data["a"])
b.create_dataset("combined", data=data["b"])
f.close()

# We now show that it cannot be read in using the pesummary inbuilt functions
try:
    f = read("example.h5")
    raise
except Exception:
    print("Failed to read in 'example.h5' with the pesummary inbuilt functions")

# Next, lets define a class which inherits from the SingleAnalysisRead class
# that is able to extract the posterior samples from our custom result file


class CustomReadClass(SingleAnalysisRead):
    """Class to read in our custom file

    Parameters
    ----------
    path_to_results_file: str
        path to the result file you wish to read in
    """
    def __init__(self, path_to_results_file, **kwargs):
        super(CustomReadClass, self).__init__(path_to_results_file, **kwargs)
        self.load(self.custom_load_function)

    def custom_load_function(self, path, **kwargs):
        """Function to load data from a custom hdf5 file
        """
        import h5py

        f = h5py.File(path, 'r')
        parameters = list(f["my_posterior_samples"].keys())
        samples = np.array([
            f["my_posterior_samples"][param]["combined"] for param in
            parameters
        ]).T
        f.close()
        # Return a dictionary of data
        return {
            "parameters": parameters, "samples": samples
        }


# Now, lets read in the result file
f = read("example.h5", cls=CustomReadClass)
print("Read in with class: {}".format(f.__class__))

# Now lets confirm that we have extracted the samples correctly
np.testing.assert_almost_equal(data["a"], f.samples_dict["a"])
np.testing.assert_almost_equal(data["b"], f.samples_dict["b"])

# We can now use the builtin pesummary plotting methods
samples = f.samples_dict
fig = samples.plot("a", type="hist")
fig.savefig("test.png")
# Here we show how to read in a result file which cannot be read in using
# the pesummary builtin functions. This result file contains multiple
# analyses
import numpy as np
import h5py
from pesummary.core.file.formats.base_read import MultiAnalysisRead
from pesummary.io import read

# First let us create a dictionary containing some posterior samples
data = {
    "run_1": {
        "a": np.random.uniform(1, 5, 1000),
        "b": np.random.uniform(1, 5, 1000)
    }, "run_2": {
        "a": np.random.uniform(1, 2, 1000),
        "b": np.random.uniform(1, 2, 1000)
    }
}

# Next, lets make a file which cannot be read in with pesummary
f = h5py.File("example.h5", "w")
samples = f.create_group("my_posterior_samples")
a = samples.create_group("a")
b = samples.create_group("b")
a.create_dataset("run_1", data=data["run_1"]["a"])
a.create_dataset("run_2", data=data["run_2"]["a"])
b.create_dataset("run_1", data=data["run_1"]["b"])
b.create_dataset("run_2", data=data["run_2"]["b"])
f.close()

# We now show that it cannot be read in using the pesummary inbuilt functions
try:
    f = read("example.h5")
    raise
except Exception:
    print("Failed to read in 'example.h5' with the pesummary inbuilt functions")

# Next, lets define a class which inherits from the SingleAnalysisRead class
# that is able to extract the posterior samples from our custom result file


class CustomReadClass(MultiAnalysisRead):
    """Class to read in our custom file

    Parameters
    ----------
    path_to_results_file: str
        path to the result file you wish to read in
    """
    def __init__(self, path_to_results_file, **kwargs):
        super(CustomReadClass, self).__init__(path_to_results_file, **kwargs)
        self.load(self.custom_load_function)

    def custom_load_function(self, path, **kwargs):
        """Function to load data from a custom hdf5 file
        """
        import h5py

        f = h5py.File(path, 'r')
        parameters = list(f["my_posterior_samples"].keys())
        labels = list(f["my_posterior_samples"][parameters[0]].keys())
        samples = [
            np.array([
                f["my_posterior_samples"][param][label] for param in
                parameters
            ]).T for label in labels
        ]
        f.close()
        # Return a dictionary of data
        return {
            "parameters": [parameters] * 2, "samples": samples,
            "labels": labels
        }


# Now, lets read in the result file
f = read("example.h5", cls=CustomReadClass)
print("Read in with class: {}".format(f.__class__))

# Now lets confirm that we have extracted the samples correctly
np.testing.assert_almost_equal(data["run_1"]["a"], f.samples_dict["run_1"]["a"])
np.testing.assert_almost_equal(data["run_1"]["b"], f.samples_dict["run_1"]["b"])
np.testing.assert_almost_equal(data["run_2"]["a"], f.samples_dict["run_2"]["a"])
np.testing.assert_almost_equal(data["run_2"]["b"], f.samples_dict["run_2"]["b"])

# We can now use the builtin pesummary plotting methods
samples = f.samples_dict
fig, _, _, _ = samples.plot(["a", "b"], type="reverse_triangle")
fig.savefig("test.png")