Spectrogram#

class gwpy.spectrogram.Spectrogram(
data: ArrayLike,
unit: UnitLike = None,
t0: SupportsToGps | None = None,
dt: Quantity | float | None = None,
f0: Quantity | float | None = None,
df: Quantity | float | None = None,
times: ArrayLike1D | None = None,
frequencies: ArrayLike1D | None = None,
name: str | None = None,
channel: Channel | str | None = None,
**kwargs,
)[source]#

Bases: Array2D

A 2D array holding a spectrogram of time-frequency data.

Parameters:
valuearray_like

Input data array.

unitUnit, optional

Physical unit of these data.

epochLIGOTimeGPS, float, str, optional

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

sample_ratefloat, Quantity, optional

The rate of samples per second (Hertz).

timesarray-like

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

f0float, Quantity, optional

Starting frequency for these data.

dffloat, Quantity, optional

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 Spectrogram.

write

Write this Spectrogram to a file.

plot([method, figsize, xscale])

Plot the data for this Spectrogram.

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

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

Attributes Summary

band

Frequency band described by these data.

df

Frequency spacing for these data.

dt

Time (seconds) between successive bins.

epoch

GPS epoch for these data.

f0

Starting frequency for these data.

frequencies

Series of frequencies for these data.

read

Read data into a Spectrogram.

span

GPS [start, stop) span for these data.

t0

GPS time of first time bin.

times

Series of GPS times for each sample

write

Write this Spectrogram to a file.

Methods Summary

crop_frequencies([low, high, copy])

Crop this Spectrogram to the specified frequencies.

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

Apply the given filter to this Spectrogram.

from_spectra(*spectra, **kwargs)

Build a new Spectrogram from a list of spectra.

percentile(percentile)

Calculate a given spectral percentile for this Spectrogram.

plot([method, figsize, xscale])

Plot the data for this Spectrogram.

ratio(operand)

Calculate the ratio of this Spectrogram against a reference.

variance([bins, low, high, nbins, log, ...])

Calculate the SpectralVariance of this Spectrogram.

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

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

Attributes Documentation

band[source]#

Frequency band described by these data.

df[source]#

Frequency spacing for these data.

dt[source]#

Time (seconds) between successive bins.

epoch[source]#

GPS epoch for these data.

f0[source]#

Starting frequency for these data.

frequencies[source]#

Series of frequencies for these data.

read#

Read data into a Spectrogram.

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, list

Source of data, any of the following:

  • str path of single data file,

  • str path of LAL-format cache file,

  • list of paths.

*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

span[source]#

GPS [start, stop) span for these data.

t0[source]#

GPS time of first time bin.

times[source]#

Series of GPS times for each sample

write#

Write this Spectrogram 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

crop_frequencies(
low: float | Quantity | None = None,
high: float | Quantity | None = None,
*,
copy: bool = False,
) Spectrogram[source]#

Crop this Spectrogram to the specified frequencies.

Parameters:
lowfloat, optional

Lower frequency bound for cropped Spectrogram.

highfloat, optional

Upper frequency bound for cropped Spectrogram.

copybool, optional

If False return a view of the original data, otherwise create a fresh memory copy.

Returns:
specSpectrogram

A new Spectrogram with a subset of data from the frequency axis

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

Apply the given filter to this Spectrogram.

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

If True, filter definition will be converted from Hertz to Z-domain digital representation, default: False.

sample_ratefloat, Quantity, optional

Sample rate of data (in Hertz), used to apply a digital filter. Defaults to the last frequency value of this Spectrogram (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, default: False.

kwargs

Additional keyword arguments passed to the filter function.

Returns:
resultSpectrogram

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

Raises:
ValueError

If filt arguments cannot be interpreted properly.

classmethod from_spectra(
*spectra: FrequencySeries,
**kwargs,
) Spectrogram[source]#

Build a new Spectrogram from a list of spectra.

Parameters:
*spectraFrequencySeries

One or more frequency series to stack.

dtfloat, Quantity, optional

Stride between given spectra.

kwargs

Other keyword arguments to pass to the constructor.

Returns:
Spectrogram

A new Spectrogram from a vertical stacking of the spectra The new object takes the metadata from the first given FrequencySeries if not given explicitly.

Notes

Each FrequencySeries passed to this constructor must be the same length.

percentile(
percentile: float,
) FrequencySeries[source]#

Calculate a given spectral percentile for this Spectrogram.

Parameters:
percentilefloat

percentile (0 - 100) of the bins to compute

Returns:
spectrumFrequencySeries

the given percentile FrequencySeries calculated from this SpectralVaraicence

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

Plot the data for this Spectrogram.

Parameters:
methodstr, optional

The plotting method to use to render this spectrogram, either 'pcolormesh' (default) or 'imshow'.

figsizetuple of float, optional

(width, height) (inches) of the output figure.

xscalestr, optional

The X-axis scale.

kwargs

All keyword arguments are passed along to underlying functions, see below for references.

Returns:
plotPlot

The Plot containing the data.

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.

gwpy.plot.Axes.imshow
gwpy.plot.Axes.pcolormesh

For documentation of keyword arguments used in rendering the Spectrogram data.

ratio(
operand: FrequencySeries | Quantity | Literal['mean', 'median'],
) Spectrogram[source]#

Calculate the ratio of this Spectrogram against a reference.

Parameters:
operandstr, FrequencySeries, Quantity

A FrequencySeries or Quantity to weight against, or one of

'mean'

Weight against the mean of each spectrum in this Spectrogram.

'median'

Weight against the median of each spectrum in this Spectrogram.

Returns:
spectrogramSpectrogram

A new Spectrogram.

Raises:
ValueError

If operand is given as a str that isn’t supported.

variance(
bins: ArrayLike1D | None = None,
low: float | None = None,
high: float | None = None,
nbins: int = 500,
*,
log: bool = False,
norm: bool = False,
density: bool = False,
) SpectralVariance[source]#

Calculate the SpectralVariance of this Spectrogram.

Parameters:
binsndarray, optional, default None

array of histogram bin edges, including the rightmost edge

lowfloat, optional, default: None

left edge of lowest amplitude bin, only read if bins is not given

highfloat, optional, default: None

right edge of highest amplitude bin, only read if bins is not given

nbinsint, optional, default: 500

number of bins to generate, only read if bins is not given

logbool, optional, default: False

calculate amplitude bins over a logarithmic scale, only read if bins is not given

normbool, optional, default: False

normalise bin counts to a unit sum

densitybool, optional, default: False

normalise bin counts to a unit integral

Returns:
specvarSpectralVariance

2D-array of spectral frequency-amplitude counts

See also

numpy.histogram

for details on specifying bins and weights

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

Filter this Spectrogram 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.

sample_ratefloat, Quantity, optional

Sample rate of data (in Hertz), used to apply a digital filter. Defaults to the last frequency value of this Spectrogram (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:
specgramSpectrogram

The frequency-domain filtered version of the input data.

See also

Spectrogram.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)