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,
Bases:
SeriesAn
Arraywith time-domain metadata.- Parameters:
- valuearray_like
Input data array.
- unit
Unit, optional Physical unit of these data.
- t0
LIGOTimeGPS,float,str, optional GPS epoch associated with these data, any input parsable by
to_gpsis fine.- dt
float,Quantity, optional Time between successive samples (seconds), can also be given inversely via
sample_rate.- sample_rate
float,Quantity, optional The rate of samples per second (Hertz), can also be given inversely via
dt.- times
array-like The complete array of GPS times accompanying the data for this series. This argument takes precedence over
t0anddtso should be given in place of these if relevant, not alongside.- name
str, optional Descriptive title for this array.
- channel
Channel,str, optional Source data stream for these data.
- dtype
dtype, optional Input data type.
- copy
bool, optional, default:False Choose to copy the input data to new memory.
- subok
bool, optional, default:True Allow passing of sub-classes by the array generator.
Attributes Summary
Time (seconds) between successive samples.
Duration of this series in seconds.
GPS epoch for these data.
Get data from any source.
Read data into a
TimeSeriesBase.Data rate for this
TimeSeriesin samples per second (Hertz).Time (seconds) spanned by this series.
GPS start time of this series.
Array of GPS times for each sample.
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.Seriesobject.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.bufferobject.from_pycbc(pycbcseries, *[, copy])Convert a
pycbc.types.timeseries.TimeSeriesinto aTimeSeries.plot([method, figsize, xscale])Plot the data for this timeseries.
to_lal()Convert this
TimeSeriesinto a LAL TimeSeries.to_pycbc(*[, copy])Convert this
TimeSeriesinto a PyCBCTimeSeries.Attributes Documentation
- get#
Get data from any source.
- read#
Read data into a
TimeSeriesBase.
- sample_rate[source]#
Data rate for this
TimeSeriesin samples per second (Hertz).This attribute is stored internally by the
dxattribute
- 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,
Fetch data from NDS.
- Parameters:
- channel
str,Channel The name (or representation) of the data channel to fetch.
- 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.
- channel
- 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,
Fetch open-access data from GWOSC.
This is just a shim around
TimeSeries.get(..., source='gwosc').- Parameters:
- ifo
str The two-character prefix of the IFO in which you are interested, e.g.
'L1'.- start
LIGOTimeGPS,float,str, optional GPS start time of required data, defaults to start of data found; 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.- 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''hdf5''gwf'- requireslalframe
- host
str, optional HTTP host name of GWOSC server to access.
- verbose
bool, optional This argument is deprecated and will be removed in a future release. Use DEBUG-level logging instead, see Logging with GWpy.
- 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. SetGWPY_CACHE=1in the environment to auto-cache.- timeout
float, optional The time to wait for a response from the GWOSC server.
- kwargs
Any other keyword arguments are passed to the
TimeSeries.readmethod that parses the file that was downloaded.
- ifo
Notes
StateVectordata are not available intxt.gzformat.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 beformat-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,
Find and return data for multiple channels using GWDataFind.
This method uses
gwdatafindto discover the URLs that provide the requested data, then reads those files usingTimeSeriesDict.read().This is just a shim around
TimeSeries.get(..., source='gwdatafind').- Parameters:
- channel
str Name of data channel 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 processes to use.
- 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().
- channel
- classmethod from_arrakis(
- series: arrakis.block.Series,
- *,
- copy: bool = True,
- **metadata,
Construct a new series from an
arrakis.Seriesobject.- Parameters:
- series
arrakis.Series The input Arrakis data series 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.
- series
- Returns:
- timeseries
TimeSeries A new
TimeSeriescontaining the data from thearrakis.Seriesand the appropriate metadata.
- timeseries
- 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( ) Self[source]#
Construct a new series from an
nds2.bufferobject.Requires:
NDS2- Parameters:
- buffer
nds2.buffer The input NDS2-client buffer 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 Tf
True, copy the contained data array to new to a new array.- metadata
Any other metadata keyword arguments to pass to the
TimeSeriesconstructor.
- buffer
- Returns:
- timeseries
TimeSeries A new
TimeSeriescontaining the data from thends2.buffer, and the appropriate metadata.
- timeseries
- classmethod from_pycbc(
- pycbcseries: pycbc.types.TimeSeries,
- *,
- copy: bool = True,
Convert a
pycbc.types.timeseries.TimeSeriesinto aTimeSeries.- Parameters:
- pycbcseries
pycbc.types.timeseries.TimeSeries The input PyCBC
TimeSeriesarray.- copy
bool, optional If
True, copy these data to a new array.
- pycbcseries
- Returns:
- timeseries
TimeSeries A GWpy version of the input timeseries.
- timeseries
- plot( ) Plot[source]#
Plot the data for this timeseries.
- Returns:
- figure
Figure The newly created figure, with populated Axes.
- figure
See also
matplotlib.pyplot.figureFor documentation of keyword arguments used to create the figure.
matplotlib.figure.Figure.add_subplotFor documentation of keyword arguments used to create the axes.
matplotlib.axes.Axes.plotFor documentation of keyword arguments used in rendering the data.
- to_lal() LALTimeSeriesType[source]#
Convert this
TimeSeriesinto a LAL TimeSeries.Note
This operation always copies data to new memory.
- to_pycbc(*, copy: bool = True) pycbc.types.TimeSeries[source]#
Convert this
TimeSeriesinto a PyCBCTimeSeries.- Parameters:
- Returns:
- timeseries
TimeSeries A PyCBC representation of this
TimeSeries.
- timeseries