TimeSeriesBaseDict#
- class gwpy.timeseries.TimeSeriesBaseDict[source]#
Bases:
dict[str|Channel,_V],Generic[_V]Key-value mapping of named
TimeSeriesBaseobjects.This object is designed to hold data for many different sources (channels) for a single time span. Dictionary keys are ordered by insertion order.
The main entry points for this object are the
read()andfetch()data access methods.Attributes Summary
Retrieve data for multiple names from any data source.
Read data into a
TimeSeriesBaseDict.The GPS
[start, stop)extent of data in thisdict.Write data from a
TimeSeriesBaseDict.Methods Summary
append(other, *[, copy])Append the dict
otherto this one.copy()Return a copy of this dict with each value copied to new memory.
crop([start, end, copy])Crop each entry of this
dict.fetch(channels, start, end, *[, host, port, ...])Fetch data from NDS for a number of channels.
fetch_open_data(detectors, start, end, *[, ...])Fetch open-access data from the LIGO Open Science Center.
find(channels, start, end, *[, observatory, ...])Find and read data from frames for a number of channels.
from_arrakis(block, *[, copy])Construct a new dict from an
arrakis.SeriesBlock.from_nds2_buffers(buffers, *[, scaled, copy])Construct a new dict from a list of
nds2.bufferobjects.plot([label, method, figsize, xscale])Plot the data for this
TimeSeriesBaseDict.prepend(other, **kwargs)Prepend the dict
otherto this one.resample(rate, **kwargs)Resample items in this dict.
step([label, where, figsize, xscale])Create a step plot of this dict.
Attributes Documentation
- get#
Retrieve data for multiple names from any data source.
- read#
Read data into a
TimeSeriesBaseDict.
- write#
Write data from a
TimeSeriesBaseDict.
Methods Documentation
- append( ) Self[source]#
Append the dict
otherto this one.- Parameters:
- other
dictofTimeSeries The container to append to this one.
- copy
bool, optional If
Truecopy data fromotherbefore storing, only affects those keys inotherthat aren’t inself.- **kwargs
Other keyword arguments to send to
TimeSeries.append.
- other
See also
TimeSeries.appendFor details of the underlying series append operation.
- crop( ) Self[source]#
Crop each entry of this
dict.This method calls the
crop()method of all entries and modifies this dict in place.- Parameters:
- start
LIGOTimeGPS,float,str GPS start time of required data, any input parseable by
to_gpsis fine- end
LIGOTimeGPS,float,str, optional GPS end time of required data, defaults to end of data found; any input parseable by
to_gpsis fine- copy
bool, optional, default:False If
Truecopy the data for each entry to fresh memory, otherwise return a view.
- start
See also
TimeSeries.cropfor more details
- classmethod fetch(
- channels: list[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,
Fetch data from NDS for a number of channels.
This is just a shim around
TimeSeriesDict.get(..., source='nds2').- Parameters:
- channels
str,Channel List of names of data channels to find.
- start
LIGOTimeGPS,float,str GPS start time of required data, any input parseable by
to_gpsis fine- end
LIGOTimeGPS,float,str, optional GPS end time of required data, defaults to end of data found; any input parseable by
to_gpsis fine- host
str, optional URL of NDS server to use, if blank will try any server (in a relatively sensible order) to get the data
One of
connectionorhostmust be given.- port
int, optional Port number for NDS server query, must be given with
host.- verify
bool, optional Check channels exist in database before asking for data. Default is
True.- verbose
bool, optional This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.
- connection
nds2.connection, optional Open NDS connection to use. Default is to open a new connection using
hostandportarguments.One of
connectionorhostmust be given.- pad
float, optional Float value to insert between gaps. Default behaviour is to raise an exception when any gaps are found.
- scaled
bool, optional Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.
- allow_tape
bool, optional Allow data access from slow tapes. If
hostorconnectionis given, the default is to do whatever the server default is, otherwise servers will be searched withallow_tape=Falsefirst, thenallow_tape=Trueif that fails.- type
int,str, optional NDS2 channel type integer or string name to match. Default is to search for any channel type.
- dtype
numpy.dtype,str,type, ordict, optional NDS2 data type to match. Default is to search for any data type.
- channels
- Returns:
- data
TimeSeriesBaseDict A new
TimeSeriesBaseDictof (str,TimeSeries) pairs fetched from NDS.
- data
- classmethod fetch_open_data(
- detectors: str,
- start: SupportsToGps,
- end: SupportsToGps,
- *,
- sample_rate: float = 4096,
- version: int | None = None,
- format: str = 'hdf5',
- host: str = 'https://gwosc.org',
- verbose: bool | None = None,
- cache: bool | None = None,
- parallel: int = 1,
- **kwargs,
Fetch open-access data from the LIGO Open Science Center.
This is just a shim around
TimeSeriesDict.get(..., source='gwosc').- Parameters:
- detectors
listofstr List of two-character prefices of the IFOs in which you are interested, e.g.
['H1', 'L1'].- start
LIGOTimeGPS,float,str GPS start time of required data, any input parseable by
to_gpsis fine.- end
LIGOTimeGPS,float,str GPS end time of required data, any input parseable by
to_gpsis fine.- sample_rate
float,Quantity, The sample rate (Hertz) 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.
- version
int Version of files to download, defaults to highest discovered version.
- format
str The data format to download and parse. One of
- “hdf5”
HDF5 data files, read using
h5py.- “gwf”
Gravitational-Wave Frame files, requires
LDAStools.frameCPP.
- host
str Host name of GWOSC server to access.
- verbose
bool This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.
- cache
bool Save/read a local copy of the remote URL, default:
False; useful if the same remote data are to be accessed multiple times. SetGWPY_CACHE=1in the environment to auto-cache.- parallel
int Number of parallel threads to use when downloading data for multiple detectors. Default is
1.- kwargs
Any other keyword arguments are passed to the
TimeSeries.readmethod that parses the file that was downloaded.
- detectors
See also
TimeSeries.fetch_open_dataFor more examples.
TimeSeries.readFor details of how files are read.
Examples
>>> from gwpy.timeseries import TimeSeriesDict >>> print(TimeSeriesDict.fetch_open_data(['H1', 'L1'], 1126259446, 1126259478)) TimeSeriesDict({'H1': <TimeSeries([2.17704028e-19, 2.08763900e-19, 2.39681183e-19, ..., 3.55365541e-20, 6.33533516e-20, 7.58121195e-20] unit=Unit(dimensionless), t0=<Quantity 1.12625945e+09 s>, dt=<Quantity 0.00024414 s>, name='Strain', channel=None)>, 'L1': <TimeSeries([-1.04289994e-18, -1.03586274e-18, -9.89322445e-19, ..., -1.01767748e-18, -9.82876816e-19, -9.59276974e-19] unit=Unit(dimensionless), t0=<Quantity 1.12625945e+09 s>, dt=<Quantity 0.00024414 s>, name='Strain', channel=None)>})
- classmethod find(
- channels: list[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,
Find and read data from frames for a number of channels.
This method uses
gwdatafindto discover the (file://) URLs that provide the requested data, then reads those files usingTimeSeriesDict.read().This is just a shim around
TimeSeriesDict.get(..., source="gwdatafind").- Parameters:
- channels
list List of names of data channels to find.
- start
LIGOTimeGPS,float,str GPS start time of required data, any input parseable by
to_gpsis fine- end
LIGOTimeGPS,float,str GPS end time of required data, defaults to end of data found; any input parseable by
to_gpsis fine- observatory
str, 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').- frametype
str, 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_match
str, optional Regular expression to use for frametype matching.
- host
str, optional Name of the GWDataFind server to use. Default is set by
gwdatafind.utils.get_default_host.- urltype
str, optional The URL type to use. Default is “file” to use paths available on the file system.
- ext
str, optional The file extension for which to search. “gwf” is the only file extension supported, but this may be extended in the future.
- pad
float, optional Value with which to fill gaps in the source data, by default gaps will result in a
ValueError.- scaled
bool, optional Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.
- parallel
int, optional Number of parallel threads to use when reading data.
- allow_tape
bool, optional Allow reading from frame files on (slow) magnetic tape.
- verbose
bool, 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().
- channels
- Raises:
requests.exceptions.HTTPErrorIf the GWDataFind query fails for any reason.
RuntimeErrorIf no files are found to read, or if the read operation fails.
- classmethod from_arrakis(
- block: arrakis.SeriesBlock,
- *,
- copy: bool = True,
- **metadata,
Construct a new dict from an
arrakis.SeriesBlock.- Parameters:
- block
arrakis.SeriesBlock The input Arrakis data to read.
- copy
bool, optional If
True, copy the contained data array to new to a new array.- metadata
Any other metadata keyword arguments to pass to the
TimeSeriesconstructor.
- block
- Returns:
- dict
TimeSeriesDict A new
TimeSeriesDictcontaining the data from the Arrakis block.
- dict
- classmethod from_nds2_buffers( ) Self[source]#
Construct a new dict from a list of
nds2.bufferobjects.Requires:
NDS2- Parameters:
- buffers
listofnds2.buffer The input NDS2-client buffers to read.
- scaled
bool, optional Apply slope and bias calibration to ADC data, for non-ADC data this option has no effect.
- copy
bool, optional If
True, copy the contained data array to new to a new array.- metadata
Any other metadata keyword arguments to pass to the
TimeSeriesconstructor.
- buffers
- Returns:
- dict
TimeSeriesDict A new
TimeSeriesDictcontaining the data from the given buffers.
- dict
- plot(
- label: str = 'key',
- method: str = 'plot',
- figsize: tuple[float, float] = (12, 4),
- xscale: str = 'auto-gps',
- **kwargs,
Plot the data for this
TimeSeriesBaseDict.- Parameters:
- label
str, optional Labelling system to use, or fixed label for all elements Special values include
'key': use the key of theTimeSeriesBaseDict,'name': use thenameof each element
If anything else, that fixed label will be used for all lines.
- method
str, optional The plotting method to use. This can be any method supported by the underlying plotting library (e.g., Matplotlib).
- figsize
tuple[float, float], optional The size of the figure to create, in inches.
- xscale
str, optional The scale of the x-axis. This can be one of
'linear': linear scale'log': logarithmic scale'auto-gps': automatically determine scale based on GPS time
- kwargs
All other keyword arguments are passed to the plotter as appropriate.
- label
- prepend(other: Mapping[str | Channel, _V], **kwargs) Self[source]#
Prepend the dict
otherto this one.- Parameters:
- other
dictofTimeSeries The container to prepend to this one.
- copy
bool, optional If
Truecopy data fromotherbefore storing, only affects those keys inotherthat aren’t inself.- kwargs
Other keyword arguments to send to
TimeSeries.prepend.
- other
See also
TimeSeries.prependFor details of the underlying series prepend operation.
- step(
- label: str = 'key',
- where: Literal['pre', 'post', 'mid'] = 'post',
- figsize: tuple[float, float] = (12, 4),
- xscale: str = 'auto-gps',
- **kwargs,
Create a step plot of this dict.
- Parameters:
- label
str, optional Labelling system to use, or fixed label for all elements. Special values include
'key'Use the key of the
TimeSeriesBaseDict'name'Use the
nameof each element
If anything else, that fixed label will be used for all lines.
- where
str, optional The location of the step change. This can be one of
'pre': the step change occurs before the x value'post': the step change occurs after the x value'mid': the step change occurs at the midpoint of the x value
- figsize
tuple[float, float], optional The size of the figure to create, in inches.
- xscale
str, optional The scale of the x-axis. This can be one of
'linear': linear scale'log': logarithmic scale'auto-gps': automatically determine scale based on GPS time
- kwargs
All other keyword arguments are passed to the plotter as appropriate.
- label