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,
Bases:
Array2DA 2D array holding a spectrogram of time-frequency data.
- Parameters:
- valuearray_like
Input data array.
- unit
Unit, optional Physical unit of these data.
- epoch
LIGOTimeGPS,float,str, optional GPS epoch associated with these data, any input parsable by
to_gpsis fine.- sample_rate
float,Quantity, optional The rate of samples per second (Hertz).
- times
array-like The complete array of GPS times accompanying the data for this series. This argument takes precedence over
epochandsample_rateso should be given in place of these if relevant, not alongside.- f0
float,Quantity, optional Starting frequency for these data.
- df
float,Quantity, optional Frequency resolution for these data.
- frequencies
array-like The complete array of frequencies indexing the data. This argument takes precedence over
f0anddfso should be given in place of these if relevant, not alongside.- epoch
LIGOTimeGPS,float,str, optional GPS epoch associated with these data, any input parsable by
to_gpsis fine.- 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.
Notes
Key methods:
Read data into a
Spectrogram.Write this
Spectrogramto a file.plot([method, figsize, xscale])Plot the data for this
Spectrogram.zpk(zeros, poles, gain, *[, analog, ...])Filter this
Spectrogramby applying a zero-pole-gain filter.Attributes Summary
Frequency band described by these data.
Frequency spacing for these data.
Time (seconds) between successive bins.
GPS epoch for these data.
Starting frequency for these data.
Series of frequencies for these data.
Read data into a
Spectrogram.GPS [start, stop) span for these data.
GPS time of first time bin.
Series of GPS times for each sample
Write this
Spectrogramto a file.Methods Summary
crop_frequencies([low, high, copy])Crop this
Spectrogramto the specified frequencies.filter(filt, *[, analog, sample_rate, unit, ...])Apply the given filter to this
Spectrogram.from_spectra(*spectra, **kwargs)Build a new
Spectrogramfrom 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
Spectrogramagainst a reference.variance([bins, low, high, nbins, log, ...])Calculate the
SpectralVarianceof thisSpectrogram.zpk(zeros, poles, gain, *[, analog, ...])Filter this
Spectrogramby applying a zero-pole-gain filter.Attributes Documentation
- 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:
- source
str,list Source of data, any of the following:
- *args
Other arguments are (in general) specific to the given
format.- format
str, 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.
- source
- Raises:
IndexErrorif
sourceis an empty list
- write#
Write this
Spectrogramto 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.
Methods Documentation
- crop_frequencies( ) Spectrogram[source]#
Crop this
Spectrogramto the specified frequencies.- Parameters:
- low
float, optional Lower frequency bound for cropped
Spectrogram.- high
float, optional Upper frequency bound for cropped
Spectrogram.- copy
bool, optional If
Falsereturn a view of the original data, otherwise create a fresh memory copy.
- low
- Returns:
- spec
Spectrogram A new
Spectrogramwith a subset of data from the frequency axis
- spec
- filter(
- filt: FilterCompatible,
- *,
- analog: bool = False,
- sample_rate: QuantityLike | None = None,
- unit: str = 'rad/s',
- normalize_gain: bool = False,
- inplace: bool = False,
- **kwargs,
Apply the given filter to this
Spectrogram.- Parameters:
- filt
numpy.ndarrayortuple 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.
- analog
bool, optional If
True, filter definition will be converted from Hertz to Z-domain digital representation, default:False.- sample_rate
float,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).- unit
str, optional For analogue ZPK filters, the units in which the zeros and poles are specified. Either
'Hz'or'rad/s'(default).- normalize_gain
bool, 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'ins2z()).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'ins2z()).
Only used for analogue filters in Hz (
analog=True, unit="Hz").- inplace
bool, optional If
True, this array will be overwritten with the filtered version, default:False.- kwargs
Additional keyword arguments passed to the filter function.
- filt
- Returns:
- result
Spectrogram The filtered version of the input
Spectrogram, ifinplace=Truewas given, this is just a reference to the modified input array.
- result
- Raises:
ValueErrorIf
filtarguments cannot be interpreted properly.
- classmethod from_spectra(
- *spectra: FrequencySeries,
- **kwargs,
Build a new
Spectrogramfrom a list of spectra.- Parameters:
- *spectra
FrequencySeries One or more frequency series to stack.
- dt
float,Quantity, optional Stride between given spectra.
- kwargs
Other keyword arguments to pass to the constructor.
- *spectra
- Returns:
SpectrogramA new
Spectrogramfrom a vertical stacking of the spectra The new object takes the metadata from the first givenFrequencySeriesif not given explicitly.
Notes
Each
FrequencySeriespassed to this constructor must be the same length.
- percentile(
- percentile: float,
Calculate a given spectral percentile for this
Spectrogram.- Parameters:
- percentile
float percentile (0 - 100) of the bins to compute
- percentile
- Returns:
- spectrum
FrequencySeries the given percentile
FrequencySeriescalculated from thisSpectralVaraicence
- spectrum
- plot(
- method: str = 'pcolormesh',
- figsize: tuple[float, float] = (12, 6),
- xscale: str = 'auto-gps',
- **kwargs,
Plot the data for this
Spectrogram.- Parameters:
- method
str, optional The plotting method to use to render this spectrogram, either
'pcolormesh'(default) or'imshow'.- figsize
tupleoffloat, optional (width, height)(inches) of the output figure.- xscale
str, optional The X-axis scale.
- kwargs
All keyword arguments are passed along to underlying functions, see below for references.
- method
- Returns:
- plot
Plot The
Plotcontaining the data.
- plot
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.
gwpy.plot.Axes.imshowgwpy.plot.Axes.pcolormeshFor documentation of keyword arguments used in rendering the
Spectrogramdata.
- ratio(
- operand: FrequencySeries | Quantity | Literal['mean', 'median'],
Calculate the ratio of this
Spectrogramagainst a reference.- Parameters:
- operand
str,FrequencySeries,Quantity A
FrequencySeriesorQuantityto 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.
- operand
- Returns:
- spectrogram
Spectrogram A new
Spectrogram.
- spectrogram
- Raises:
ValueErrorIf
operandis given as astrthat 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,
Calculate the
SpectralVarianceof thisSpectrogram.- Parameters:
- bins
ndarray, optional,defaultNone array of histogram bin edges, including the rightmost edge
- low
float, optional, default:None left edge of lowest amplitude bin, only read if
binsis not given- high
float, optional, default:None right edge of highest amplitude bin, only read if
binsis not given- nbins
int, optional, default:500 number of bins to generate, only read if
binsis not given- log
bool, optional, default:False calculate amplitude bins over a logarithmic scale, only read if
binsis not given- norm
bool, optional, default:False normalise bin counts to a unit sum
- density
bool, optional, default:False normalise bin counts to a unit integral
- bins
- Returns:
- specvar
SpectralVariance 2D-array of spectral frequency-amplitude counts
- specvar
See also
numpy.histogramfor 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,
Filter this
Spectrogramby applying a zero-pole-gain filter.- Parameters:
- zeros
array-like List of zero frequencies (in Hertz).
- poles
array-like List of pole frequencies (in Hertz).
- gain
float DC gain of filter.
- analog
bool, optional Type of ZPK being applied, if
analog=Trueall parameters will be converted in the Z-domain for digital filtering.- sample_rate
float,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).- unit
str, optional For analogue ZPK filters, the units in which the zeros and poles are specified. Either
'Hz'or'rad/s'(default).- normalize_gain
bool, 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'ins2z()).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'ins2z()).
Only used for analogue filters in Hz (
analog=True, unit="Hz").
- zeros
- Returns:
- specgram
Spectrogram The frequency-domain filtered version of the input data.
- specgram
See also
Spectrogram.filterFor 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)