Strain Data in PESummary

The pesummary.gw.file.strain.StrainData class

pesummary can read in frame files through the universal read function. pesummary stores strain data in the pesummary.gw.file.strain.StrainData class. This object is inherited from the GWPy TimeSeries object and therefore all GWPy TimeSeries methods can be used with this class. Strain data can be read in directly with the StrainData class through the read class method. For example,

>>> from pesummary.gw.file.strain import StrainData
>>> f = StrainData.read("frame_file.gwf", channel="channel")

The StrainData object extends gwpy’s compatible frame file formats. For instance the StrainData object can read in Bilby pickle files which contain the gravitational wave strain. For example,

>>> f = StrainData.read("bilby_strain_data.pickle")

The StrainData class also offers the fetch_open_frame method which allows the user to fetch frame files from GWOSC for a given event,

>>> from pesummary.gw.file.strain import StrainData
>>> f = StrainData.fetch_open_frame(
...    "GW190412", IFO="L1", sampling_rate=4096., duration=32,
...    channel="L1:GWOSC-4KHZ_R1_STRAIN"
... )
class pesummary.gw.file.strain.StrainData(*args, IFO='H1', **kwargs)[source]

Class to extend the gwpy.timeseries.TimeSeries plotting functions to include the pesummary plots

Parameters:

IFO (str, optional) – IFO for which the strain data corresponds too. This is used to determine the color on plots. Default ‘H1’

gwpy

original gwpy TimeSeries object

Type:

gwpy.timeseries.TimeSeries

IFO

IFO for which the strain data corresponds too

Type:

str

strain_dict

dictionary of strain data

Type:

dict

plot:

Generate a plot based on the stored data

classmethod fetch_open_data(*args, **kwargs)[source]

Fetch open-access data from the LIGO Open Science Center

Parameters:
  • ifo (str) – the two-character prefix of the IFO in which you are interested, e.g. ‘L1’

  • start (~gwpy.time.LIGOTimeGPS, float, str, optional) – GPS start time of required data, defaults to start of data found; any input parseable by ~gwpy.time.to_gps is fine

  • end (~gwpy.time.LIGOTimeGPS, float, str, optional) – GPS end time of required data, defaults to end of data found; any input parseable by ~gwpy.time.to_gps is fine

  • sample_rate (float, optional,) – the sample rate of desired data; most data are stored by GWOSC at 4096 Hz, however there may be event-related data releases with a 16384 Hz rate, default: 4096

  • version (int, optional) – version of files to download, defaults to highest discovered version

  • format (str, optional) –

    the data format to download and parse, default: 'h5py'

  • host (str, optional) – HTTP host name of GWOSC server to access

  • verbose (bool, optional, default: False) – print verbose output while fetching data

  • cache (bool, optional) – save/read a local copy of the remote URL, default: False; useful if the same remote data are to be accessed multiple times. Set GWPY_CACHE=1 in the environment to auto-cache.

  • **kwargs – any other keyword arguments are passed to the TimeSeries.read method that parses the file that was downloaded

Examples

>>> from gwpy.timeseries import (TimeSeries, StateVector)
>>> print(TimeSeries.fetch_open_data('H1', 1126259446, 1126259478))
TimeSeries([  2.17704028e-19,  2.08763900e-19,  2.39681183e-19,
            ...,   3.55365541e-20,  6.33533516e-20,
              7.58121195e-20]
           unit: Unit(dimensionless),
           t0: 1126259446.0 s,
           dt: 0.000244140625 s,
           name: Strain,
           channel: None)
>>> print(StateVector.fetch_open_data('H1', 1126259446, 1126259478))
StateVector([127,127,127,127,127,127,127,127,127,127,127,127,
             127,127,127,127,127,127,127,127,127,127,127,127,
             127,127,127,127,127,127,127,127]
            unit: Unit(dimensionless),
            t0: 1126259446.0 s,
            dt: 1.0 s,
            name: Data quality,
            channel: None,
            bits: Bits(0: data present
                       1: passes cbc CAT1 test
                       2: passes cbc CAT2 test
                       3: passes cbc CAT3 test
                       4: passes burst CAT1 test
                       5: passes burst CAT2 test
                       6: passes burst CAT3 test,
                       channel=None,
                       epoch=1126259446.0))

For the StateVector, the naming of the bits will be format-dependent, because they are recorded differently by GWOSC in different formats.

Notes

StateVector data are not available in txt.gz format.

classmethod fetch_open_frame(event, **kwargs)[source]

Fetch open frame files for a given event

Parameters:
  • sampling_rate (int, optional) – sampling rate of strain data you wish to download. Default 16384

  • format (str, optional) – format of strain data you wish to download. Default “gwf”

  • duration (int, optional) – duration of strain data you wish to download. Default 32

  • IFO (str, optional) – detector strain data you wish to download. Default ‘L1’

  • **kwargs (dict, optional) – all additional kwargs passed to StrainData.read

plot(*args, type='td', **kwargs)[source]

Generate a plot displaying the gravitational wave strain data

*args: tuple

all arguments are passed to the plotting function

type: str

name of the plot you wish to make

**kwargs: dict

all additional kwargs are passed to the plotting function

Subfunctions:

pesummary.gw.plots.detchar.spectrogram

Generate a spectrogram from the timeseries

strain: dict

dictionary of gw.py timeseries objects containing the strain data for each IFO

vmin: float, optional

minimum for the colormap

vmax: float, optional

maximum for the colormap

cmap: str, optional

cmap for the plot. See matplotlib.pyplot.colormaps() for options

ylim: list, optional

list to give the lower and upper bound of the plot

pesummary.gw.plots.detchar.omegascan

Generate an omegascan from the timeseries

strain: dict

dictionary of gw.py timeseries objects containing the strain data for each IFO

gps: float

gps time you wish to center your omegascan around

window: float, optional

window around gps time to generate omagescan for. Default 4s

vmin: float, optional

minimum for the colormap

vmax: float, optional

maximum for the colormap

cmap: str, optional

cmap for the plot. See matplotlib.pyplot.colormaps() for options

ylim: list, optional

list to give the lower and upper bound of the plot

pesummary.gw.plots.detchar.time_domain_strain_data

Plot the strain data in the time domain. Code based on the GW150914

tutorial provided by gwpy: https://gwpy.github.io/docs/latest/examples/signal/gw150914.html

pesummary.gw.plots.detchar.frequency_domain_strain_data

Plot the strain data in the frequency domain

strain: dict

dictionary of gw.py timeseries objects containing the strain data for each IFO

window: Bool, optional

if True, apply a window to the data before applying FFT to the data. Default True

window_kwargs: dict, optional

optional kwargs for the window function

resolution: float, optional

resolution to downsample the frequency domain data. Default 1./512

fmin: float, optional

lowest frequency to start plotting the data

fmax: float, optional

highest frequency to stop plotting the data

asd: dict, optional

dictionary containing the ASDs for each detector to plot ontop of the detector data

classmethod read(*args, IFO='H1', **kwargs)[source]

Read data into a TimeSeries

Arguments and keywords depend on the output format, see the online documentation for full details for each format, the parameters below are common to most formats.

Parameters:
  • source (str, list) –

    Source of data, any of the following:

    • str path of single data file,

    • str path of LAL-format cache file,

    • list of paths.

  • name (str, ~gwpy.detector.Channel) – the name of the channel to read, or a Channel object.

  • start (~gwpy.time.LIGOTimeGPS, float, str, optional) – GPS start time of required data, defaults to start of data found; any input parseable by ~gwpy.time.to_gps is fine

  • end (~gwpy.time.LIGOTimeGPS, float, str, optional) – GPS end time of required data, defaults to end of data found; any input parseable by ~gwpy.time.to_gps is fine

  • format (str, optional) – source format identifier. If not given, the format will be detected if possible. See below for list of acceptable formats.

  • nproc (int, optional) – number of parallel processes to use, serial process by default.

  • pad (float, optional) – value with which to fill gaps in the source data, by default gaps will result in a ValueError.

Raises:

IndexError – if source is an empty list

Notes

The available built-in formats are:

Format

Read

Write

Auto-identify

csv

Yes

Yes

Yes

gwf

Yes

Yes

Yes

gwf.framecpp

Yes

Yes

No

gwf.framel

Yes

Yes

No

gwf.lalframe

Yes

Yes

No

hdf5

Yes

Yes

Yes

hdf5.gwosc

Yes

No

No

txt

Yes

Yes

Yes

wav

Yes

No

No

The pesummary.gw.file.strain.StrainDataDict class

If you wish to load numerous frame files from different detectors, pesummary offers the pesummary.gw.file.strain.StrainDataDict class to read in these files. As with the StrainData class, StrainDataDict offers a read class method to load a dictionary of frame files.

>>> from pesummary.gw.file.strain import StrainDataDict
>>> data = {
...     "H1": "./H-H1_LOSC_4_V2-1126257414-4096.gwf",
...     "L1": "./L-L1_LOSC_4_V2-1126257414-4096.gwf"
... }
>>> channels = {"H1": "H1:LOSC-STRAIN", "L1": "L1:LOSC-STRAIN"}
>>> strain = StrainDataDict.read(data, channels=channels)

The output is a dictionary keyed by the IFO with each value being a StrainData object.

class pesummary.gw.file.strain.StrainDataDict(*args)[source]

Class to store multiple StrainData objects from different IFOs

Parameters:

data (dict) – dict keyed by IFO and values StrainData objects

Examples

>>> from pesummary.gw.file.strain import StrainDataDict
>>> data = {
...     "H1": "./H-H1_LOSC_4_V2-1126257414-4096.gwf",
...     "L1": "./L-L1_LOSC_4_V2-1126257414-4096.gwf"
... }
>>> channels = {"H1": "H1:LOSC-STRAIN", "L1": "L1:LOSC-STRAIN"}
>>> strain = StrainDataDict.read(data, channels=channels)