FrequencySeries#

class gwpy.frequencyseries.FrequencySeries(
data: ArrayLike1D,
unit: UnitLike = None,
f0: Quantity | float | None = None,
df: Quantity | float | None = None,
frequencies: ArrayLike1D | None = None,
name: str | None = None,
epoch: SupportsToGps | None = None,
channel: Channel | str | None = None,
**kwargs,
)[source]#

Bases: Series

A data array holding some metadata to represent a frequency series.

Parameters:
valuearray_like

Input data array.

unitUnit, optional

Physical unit of these data.

f0float, Quantity, optional, default: 0

Starting frequency for these data.

dffloat, Quantity, optional, default: 1

Frequency resolution for these data.

frequenciesarray-like

The complete array of frequencies indexing the data. This argument takes precedence over f0 and df so should be given in place of these if relevant, not alongside.

epochLIGOTimeGPS, float, str, optional

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

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.

Notes

Key methods:

read

Read data into a FrequencySeries.

write

Write this FrequencySeries to a file.

plot([method, xscale])

Plot the data for this FrequencySeries.

zpk(zeros, poles, gain, *[, analog, ...])

Filter this FrequencySeries by applying a zero-pole-gain filter.

Attributes Summary

df

Frequency spacing of this FrequencySeries

f0

Starting frequency for this FrequencySeries

frequencies

Series of frequencies for each sample

read

Read data into a FrequencySeries.

write

Write this FrequencySeries to a file.

Methods Summary

filter(filt, *[, analog, sample_rate, unit, ...])

Apply a filter to this FrequencySeries.

from_lal(lalfs, *[, copy])

Generate a new FrequencySeries from a LAL FrequencySeries.

from_pycbc(fs, *[, copy])

Convert a pycbc.types.frequencyseries.FrequencySeries.

ifft()

Compute the one-dimensional discrete inverse Fourier transform.

interpolate(df)

Interpolate this FrequencySeries to a new resolution.

plot([method, xscale])

Plot the data for this FrequencySeries.

to_lal()

Convert this FrequencySeries into a LAL FrequencySeries.

to_pycbc(*[, copy])

Convert this FrequencySeries into a PyCBC FrequencySeries.

zpk(zeros, poles, gain, *[, analog, ...])

Filter this FrequencySeries by applying a zero-pole-gain filter.

Attributes Documentation

df[source]#

Frequency spacing of this FrequencySeries

f0[source]#

Starting frequency for this FrequencySeries

frequencies[source]#

Series of frequencies for each sample

read#

Read data into a FrequencySeries.

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:
sourcestr, os.PathLike, file, or list of these

Source of data, any of the following:

  • Path of a single data file

  • List of data file paths

  • Path of LAL-format cache file

args

Other arguments are (in general) specific to the given format.

formatstr, optional

Source format identifier. If not given, the format will be detected if possible. See below for list of acceptable formats.

kwargs

Other keywords are (in general) specific to the given format.

Raises:
IndexError

If source is an empty list.

write#

Write this FrequencySeries to a file.

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:
targetstr

Output filename.

formatstr, optional

Output format identifier. If not given, the format will be detected if possible. See below for list of acceptable formats.

Methods Documentation

filter(
filt: FilterCompatible,
*,
analog: bool = False,
sample_rate: QuantityLike | None = None,
unit: str = 'rad/s',
normalize_gain: bool = False,
inplace: bool = False,
) Self[source]#

Apply a filter to this FrequencySeries.

The input filter argument is designed to accept any filter created by the scipy.signal filter design functions, and operates on the conventions of that module.

Parameters:
filtnumpy.ndarray or tuple

The filter to be applied. This can be specified in any of the following forms, with the appropriate number of elements in the tuple:

  • numpy.ndarray - 1D array of FIR filter coefficients.

  • tuple[numpy.ndarray, numpy.ndarray] - numerator/demoinator polynomials of the transfer function.

  • numpy.ndarray - 2D array of SOS coefficients.

  • tuple[numpy.ndarray, numpy.ndarray, float] - zero-pole-gain representation.

analogbool, optional

Type of ZPK being applied, if analog=True all parameters will be converted in the Z-domain for digital filtering via the bilinear transform.

sample_ratefloat, Quantity, optional

Sample rate of data (in Hertz), used to apply a digital filter. Defaults to the last frequency value of this FrequencySeries (i.e. the Nyquist frequency).

unitstr, optional

For analogue ZPK filters, the units in which the zeros and poles are specified. Either 'Hz' or 'rad/s' (default).

normalize_gainbool, optional

Whether to normalize the gain when converting from Hz to rad/s.

  • False (default): Multiply zeros/poles by -2π but leave gain unchanged. This matches the LIGO GDS ‘f’ plane convention (plane='f' in s2z()).

  • True: Normalize gain to preserve frequency response magnitude. Gain is scaled by \(|∏p_i/∏z_i| · (2π)^{(n_p - n_z)}\). Use this when your filter was designed with the transfer function \(H(f) = k·∏(f-z_i)/∏(f-p_i)\) in Hz. This matches the LIGO GDS ‘n’ plane convention (plane='n' in s2z()).

Only used for analogue filters in Hz (analog=True, unit="Hz").

inplacebool, optional

If True, this array will be overwritten with the filtered version.

Returns:
resultFrequencySeries

The filtered version of the input FrequencySeries, if inplace=True was given, this is just a reference to the modified input array.

Raises:
ValueError

If filt arguments cannot be interpreted properly.

See also

FrequencySeries.zpk

For applying a zero-pole-gain filter, including in other units (e.g. poles and zeros specified in Hertz).

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

Generate a new FrequencySeries from a LAL FrequencySeries.

Any type of LAL FrequencySeries is supported.

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

Convert a pycbc.types.frequencyseries.FrequencySeries.

Parameters:
fspycbc.types.frequencyseries.FrequencySeries

The input PyCBC FrequencySeries array.

copybool, optional

If True, copy these data to a new array.

Returns:
spectrumFrequencySeries

A GWpy version of the input frequency series.

ifft() TimeSeries[source]#

Compute the one-dimensional discrete inverse Fourier transform.

Returns:
outTimeSeries

The normalised, real-valued TimeSeries.

See also

numpy.fft.irfft

The inverse (real) FFT function.

Notes

This method applies the necessary normalisation such that the condition holds:

>>> timeseries = TimeSeries([1.0, 0.0, -1.0, 0.0], sample_rate=1.0)
>>> timeseries.fft().ifft() == timeseries
interpolate(df: float) Self[source]#

Interpolate this FrequencySeries to a new resolution.

Parameters:
dffloat

Desired frequency resolution of the interpolated FrequencySeries, in Hz.

Returns:
outFrequencySeries

The interpolated version of the input FrequencySeries.

See also

numpy.interp

For the underlying 1-D linear interpolation scheme.

plot(
method: str = 'plot',
xscale: str = 'log',
**kwargs,
) Plot[source]#

Plot the data for this FrequencySeries.

to_lal() LALFrequencySeriesType[source]#

Convert this FrequencySeries into a LAL FrequencySeries.

Returns:
lalspecFrequencySeries

An XLAL-format FrequencySeries of a given type, e.g. REAL8FrequencySeries.

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

Convert this FrequencySeries into a PyCBC FrequencySeries.

Parameters:
copybool, optional

If True, copy these data to a new array.

Returns:
frequencyseriespycbc.types.frequencyseries.FrequencySeries

A PyCBC representation of this FrequencySeries.

zpk(
zeros: ArrayLike1D,
poles: ArrayLike1D,
gain: float,
*,
analog: bool = True,
sample_rate: QuantityLike | None = None,
unit: str = 'rad/s',
normalize_gain: bool = False,
) Self[source]#

Filter this FrequencySeries by applying a zero-pole-gain filter.

Parameters:
zerosarray-like

List of zero frequencies (in Hertz).

polesarray-like

List of pole frequencies (in Hertz).

gainfloat

DC gain of filter.

analogbool, optional

Type of ZPK being applied, if analog=True all parameters will be converted in the Z-domain for digital filtering via the bilinear transform.

sample_ratefloat, Quantity, optional

Sample rate of data (in Hertz), used to apply a digital filter. Defaults to the last frequency value of this FrequencySeries (i.e. the Nyquist frequency).

unitstr, optional

For analogue ZPK filters, the units in which the zeros and poles are specified. Either 'Hz' or 'rad/s' (default).

normalize_gainbool, optional

Whether to normalize the gain when converting from Hz to rad/s.

  • False (default): Multiply zeros/poles by -2π but leave gain unchanged. This matches the LIGO GDS ‘f’ plane convention (plane='f' in s2z()).

  • True: Normalize gain to preserve frequency response magnitude. Gain is scaled by \(|∏p_i/∏z_i| · (2π)^{(n_p - n_z)}\). Use this when your filter was designed with the transfer function \(H(f) = k·∏(f-z_i)/∏(f-p_i)\) in Hz. This matches the LIGO GDS ‘n’ plane convention (plane='n' in s2z()).

Only used for analogue filters in Hz (analog=True, unit="Hz").

Returns:
spectrumFrequencySeries

The frequency-domain filtered version of the input data.

See also

FrequencySeries.filter

For details on how a digital ZPK-format filter is applied.

Examples

To apply a zpk filter with file poles at 100 Hz, and five zeros at 1 Hz (giving an overall DC gain of 1e-10):

>>> data2 = data.zpk([100]*5, [1]*5, 1e-10)