TimeSeriesBase#

class gwpy.timeseries.TimeSeriesBase(
data: ArrayLike1D,
unit: UnitLike = None,
t0: SupportsToGps | None = None,
dt: float | Quantity | None = None,
sample_rate: float | Quantity | None = None,
times: ArrayLike1D | None = None,
channel: Channel | str | None = None,
name: str | None = None,
**kwargs,
)[source]#

Bases: Series

An Array with time-domain metadata.

Parameters:
valuearray_like

Input data array.

unitUnit, optional

Physical unit of these data.

t0LIGOTimeGPS, float, str, optional

GPS epoch associated with these data, any input parsable by to_gps is fine.

dtfloat, Quantity, optional

Time between successive samples (seconds), can also be given inversely via sample_rate.

sample_ratefloat, Quantity, optional

The rate of samples per second (Hertz), can also be given inversely via dt.

timesarray-like

The complete array of GPS times accompanying the data for this series. This argument takes precedence over t0 and dt so should be given in place of these if relevant, not alongside.

namestr, optional

Descriptive title for this array.

channelChannel, str, optional

Source data stream for these data.

dtypedtype, optional

Input data type.

copybool, optional, default: False

Choose to copy the input data to new memory.

subokbool, optional, default: True

Allow passing of sub-classes by the array generator.

Attributes Summary

dt

Time (seconds) between successive samples.

duration

Duration of this series in seconds.

epoch

GPS epoch for these data.

get

Get data from any source.

read

Read data into a TimeSeriesBase.

sample_rate

Data rate for this TimeSeries in samples per second (Hertz).

span

Time (seconds) spanned by this series.

t0

GPS start time of this series.

times

Array of GPS times for each sample.

write

Write data from a TimeSeriesBase.

Methods Summary

fetch(channel, start, end, *[, host, port, ...])

Fetch data from NDS.

fetch_open_data(ifo, start, end[, ...])

Fetch open-access data from GWOSC.

find(channel, start, end, *[, observatory, ...])

Find and return data for multiple channels using GWDataFind.

from_arrakis(series, *[, copy])

Construct a new series from an arrakis.Series object.

from_lal(lalts, *[, copy])

Generate a new TimeSeries from a LAL TimeSeries of any type.

from_nds2_buffer(buffer, *[, scaled, copy])

Construct a new series from an nds2.buffer object.

from_pycbc(pycbcseries, *[, copy])

Convert a pycbc.types.timeseries.TimeSeries into a TimeSeries.

plot([method, figsize, xscale])

Plot the data for this timeseries.

to_lal()

Convert this TimeSeries into a LAL TimeSeries.

to_pycbc(*[, copy])

Convert this TimeSeries into a PyCBC TimeSeries.

Attributes Documentation

dt[source]#

Time (seconds) between successive samples.

duration[source]#

Duration of this series in seconds.

Type:

Quantity scalar

epoch[source]#

GPS epoch for these data.

This attribute is stored internally by the t0 attribute.

get#

Get data from any source.

read#

Read data into a TimeSeriesBase.

sample_rate[source]#

Data rate for this TimeSeries in samples per second (Hertz).

This attribute is stored internally by the dx attribute

span[source]#

Time (seconds) spanned by this series.

t0[source]#

GPS start time of this series.

times[source]#

Array of GPS times for each sample.

write#

Write data from a TimeSeriesBase.

Methods Documentation

classmethod fetch(
channel: str | Channel,
start: SupportsToGps,
end: SupportsToGps,
*,
host: str | None = None,
port: int | None = None,
verbose: bool | str = False,
connection: nds2.connection | None = None,
verify: bool = False,
pad: float | None = None,
allow_tape: bool | None = None,
scaled: bool | None = None,
type: int | str | None = None,
dtype: int | str | None = None,
) Self[source]#

Fetch data from NDS.

Parameters:
channelstr, Channel

The name (or representation) of the data channel to fetch.

startLIGOTimeGPS, float, str

GPS start time of required data, any input parseable by to_gps is fine

endLIGOTimeGPS, float, str, optional

GPS end time of required data, defaults to end of data found; any input parseable by to_gps is fine

hoststr, optional

URL of NDS server to use, if blank will try any server (in a relatively sensible order) to get the data

One of connection or host must be given.

portint, optional

Port number for NDS server query, must be given with host.

verifybool, optional

Check channels exist in database before asking for data. Default is True.

verbosebool, optional

This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.

connectionnds2.connection, optional

Open NDS connection to use. Default is to open a new connection using host and port arguments.

One of connection or host must be given.

padfloat, optional

Float value to insert between gaps. Default behaviour is to raise an exception when any gaps are found.

scaledbool, optional

Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.

allow_tapebool, optional

Allow data access from slow tapes. If host or connection is given, the default is to do whatever the server default is, otherwise servers will be searched with allow_tape=False first, then allow_tape=True if that fails.

typeint, str, optional

NDS2 channel type integer or string name to match. Default is to search for any channel type.

dtypenumpy.dtype, str, type, or dict, optional

NDS2 data type to match. Default is to search for any data type.

classmethod fetch_open_data(
ifo: str,
start: SupportsToGps,
end: SupportsToGps,
sample_rate: float = 4096,
version: int | None = None,
format: Literal['gwf', 'hdf5'] = 'hdf5',
host: str = 'https://gwosc.org',
*,
verbose: bool | None = None,
cache: bool | None = None,
**kwargs,
) Self[source]#

Fetch open-access data from GWOSC.

This is just a shim around TimeSeries.get(..., source='gwosc').

Parameters:
ifostr

The two-character prefix of the IFO in which you are interested, e.g. 'L1'.

startLIGOTimeGPS, float, str, optional

GPS start time of required data, defaults to start of data found; any input parseable by to_gps is fine.

endLIGOTimeGPS, float, str, optional

GPS end time of required data, defaults to end of data found; any input parseable by to_gps is fine.

sample_ratefloat, 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.

versionint, optional

Version of files to download, defaults to highest discovered version.

formatstr, optional

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

hoststr, optional

HTTP host name of GWOSC server to access.

verbosebool, optional

This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.

cachebool, 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.

timeoutfloat, optional

The time to wait for a response from the GWOSC server.

kwargs

Any other keyword arguments are passed to the TimeSeries.read method that parses the file that was downloaded.

Notes

StateVector data are not available in txt.gz format.

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: quality/simple,
            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.

classmethod find(
channel: str | Channel,
start: SupportsToGps,
end: SupportsToGps,
*,
observatory: str | None = None,
frametype: str | None = None,
frametype_match: str | re.Pattern | None = None,
host: str | None = None,
urltype: str | None = 'file',
ext: str = 'gwf',
pad: float | None = None,
scaled: bool | None = None,
allow_tape: bool | None = None,
parallel: int = 1,
verbose: bool | str = False,
**readargs,
) Self[source]#

Find and return data for multiple channels using GWDataFind.

This method uses gwdatafind to discover the URLs that provide the requested data, then reads those files using TimeSeriesDict.read().

This is just a shim around TimeSeries.get(..., source='gwdatafind').

Parameters:
channelstr

Name of data channel to find.

startLIGOTimeGPS, float, str

GPS start time of required data, any input parseable by to_gps is fine.

endLIGOTimeGPS, float, str

GPS end time of required data, defaults to end of data found; any input parseable by to_gps is fine.

observatorystr, optional

The observatory to use when searching for data. Default is to use the observatory from the channel name prefix, but this should be specified when searching for data in a multi-observatory dataset (e.g. observatory='HLV').

frametypestr, optional

Name of frametype (dataset) in which this channel is stored. Default is to search all available datasets for a match, which can be very slow.

frametype_matchstr, optional

Regular expression to use for frametype matching.

hoststr, optional

Name of the GWDataFind server to use. Default is set by gwdatafind.utils.get_default_host.

urltypestr, optional

The URL type to use. Default is “file” to use paths available on the file system.

extstr, optional

The file extension for which to search. “gwf” is the only file extension supported, but this may be extended in the future.

padfloat, optional

Value with which to fill gaps in the source data, by default gaps will result in a ValueError.

scaledbool, optional

Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.

parallelint, optional

Number of parallel processes to use.

allow_tapebool, optional

Allow reading from frame files on (slow) magnetic tape.

verbosebool, optional

This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.

readargs

Any other keyword arguments to be passed to read().

classmethod from_arrakis(
series: arrakis.block.Series,
*,
copy: bool = True,
**metadata,
) Self[source]#

Construct a new series from an arrakis.Series object.

Parameters:
seriesarrakis.Series

The input Arrakis data series to read.

copybool, optional

If True, copy the contained data array to new to a new array.

metadata

Any other metadata keyword arguments to pass to the TimeSeries constructor.

Returns:
timeseriesTimeSeries

A new TimeSeries containing the data from the arrakis.Series and the appropriate metadata.

classmethod from_lal(lalts: LALTimeSeriesType, *, copy: bool = True) Self[source]#

Generate a new TimeSeries from a LAL TimeSeries of any type.

classmethod from_nds2_buffer(
buffer: nds2.buffer,
*,
scaled: bool | None = None,
copy: bool = True,
**metadata,
) Self[source]#

Construct a new series from an nds2.buffer object.

Requires: NDS2

Parameters:
buffernds2.buffer

The input NDS2-client buffer to read.

scaledbool, optional

Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.

copybool, optional

Tf True, copy the contained data array to new to a new array.

metadata

Any other metadata keyword arguments to pass to the TimeSeries constructor.

Returns:
timeseriesTimeSeries

A new TimeSeries containing the data from the nds2.buffer, and the appropriate metadata.

classmethod from_pycbc(
pycbcseries: pycbc.types.TimeSeries,
*,
copy: bool = True,
) Self[source]#

Convert a pycbc.types.timeseries.TimeSeries into a TimeSeries.

Parameters:
pycbcseriespycbc.types.timeseries.TimeSeries

The input PyCBC TimeSeries array.

copybool, optional

If True, copy these data to a new array.

Returns:
timeseriesTimeSeries

A GWpy version of the input timeseries.

plot(
method: str = 'plot',
figsize: tuple[int, int] = (12, 4),
xscale: str = 'auto-gps',
**kwargs,
) Plot[source]#

Plot the data for this timeseries.

Returns:
figureFigure

The newly created figure, with populated Axes.

See also

matplotlib.pyplot.figure

For documentation of keyword arguments used to create the figure.

matplotlib.figure.Figure.add_subplot

For documentation of keyword arguments used to create the axes.

matplotlib.axes.Axes.plot

For documentation of keyword arguments used in rendering the data.

to_lal() LALTimeSeriesType[source]#

Convert this TimeSeries into a LAL TimeSeries.

Note

This operation always copies data to new memory.

to_pycbc(*, copy: bool = True) pycbc.types.TimeSeries[source]#

Convert this TimeSeries into a PyCBC TimeSeries.

Parameters:
copybool, optional, default: True

If True, copy these data to a new array.

Returns:
timeseriesTimeSeries

A PyCBC representation of this TimeSeries.