Astrophysics Utilitiesο
Stability: Stable
Astrophysics-related tools including binary inspiral range calculations.
- gwexpy.astro.burst_range(psd: FrequencySeries, snr: float = 8, energy: float = 0.01, fmin: float = 100, fmax: float = 500) units.Quantity[source]
Calculate the integrated GRB-like GW burst range from a strain PSD.
- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data.
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
energy (float, optional) β The relative energy output of the GW burst.
fmin (float, optional) β The lower frequency cutoff of the burst range integral.
fmax (float, optional) β The upper frequency cutoff of the burst range integral.
- Returns:
range β The GRB-like-burst sensitive range [Mpc (default)].
- Return type:
~astropy.units.Quantity
Examples
Grab some data for LIGO-Livingston around GW150914 and generate a PSD
>>> from gwpy.timeseries import TimeSeries >>> hoft = TimeSeries.fetch_open_data("H1", 1126259446, 1126259478) >>> hoff = hoft.psd(fftlength=4)
Now we can calculate the
burst_range():>>> from gwpy.astro import burst_range >>> r = burst_range(hoff, fmin: float = 30) >>> print(r) 42.5055584195 Mpc
- gwexpy.astro.burst_range_spectrum(psd: FrequencySeries, snr: float = 8, energy: float = 0.01) FrequencySeries[source]
Calculate the frequency-dependent GW burst range from a strain PSD.
- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data.
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
energy (float, optional) β The relative energy output of the GW burst.
- Returns:
rangespec β The burst range FrequencySeries [Mpc (default)].
- Return type:
~gwpy.frequencyseries.FrequencySeries
- gwexpy.astro.inspiral_range(psd: FrequencySeries, snr: float = 8, mass1: float = 1.4, mass2: float = 1.4, fmin: float | None = None, fmax: float | None = None, *, horizon: bool = False, **kwargs) units.Quantity[source]
Calculate the cosmology-corrected inspiral sensitive distance.
This method returns the distance (in megaparsecs) to which a compact binary inspiral with the given component masses would be detectable given the instrumental PSD. The calculation is defined in Belczynski et. al (2014):
https://dx.doi.org/10.1088/0004-637x/789/2/120
- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data.
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
mass1 (float, ~astropy.units.Quantity, optional) β The mass of the first binary component (float assumed in solar masses).
mass2 (float, ~astropy.units.Quantity, optional) β The mass of the second binary component (float assumed in solar masses).
fmin (float, optional) β The lower frequency cut-off of the integral, default: psd.df.
fmax (float, optional) β The maximum frequency limit of the integral, defaults to the rest-frame innermost stable circular orbit (ISCO) frequency.
horizon (bool, optional) β If True, return the maximal βhorizonβ luminosity distance, otherwise return the angle-averaged comoving distance.
**kwargs β Additional keyword arguments to ~inspiral_range.waveform.CBCWaveform.
- Returns:
range β The calculated inspiral range [Mpc].
- Return type:
~astropy.units.Quantity
Examples
Grab some data for LIGO-Livingston around GW150914 and generate a PSD:
>>> from gwpy.timeseries import TimeSeries >>> hoft = TimeSeries.fetch_open_data("H1", 1126259446, 1126259478) >>> hoff = hoft.psd(fftlength=4)
Now, we can calculate the
inspiral_range():>>> from gwpy.astro import inspiral_range >>> r = inspiral_range(hoff, fmin: float = 30) >>> print(r) 70.4612102889 Mpc
See also
sensemon_rangeFor the method based on
LIGO-T030276, also known as LIGO SenseMonitor.inspiral-rangeThe package which does heavy lifting for waveform simulation and cosmology calculations.
- gwexpy.astro.inspiral_range_psd(psd: FrequencySeries, snr: float = 8, mass1: float = 1.4, mass2: float = 1.4, *, horizon: bool = False, **kwargs) FrequencySeries[source]
Calculate the cosmology-corrected inspiral sensitive distance PSD.
This method returns the power spectral density (in
Mpc**2 / Hz) to which a compact binary inspiral with the given component masses would be detectable given the instrumental PSD. The calculation is defined in Belczynski et. al (2014):https://dx.doi.org/10.1088/0004-637x/789/2/120
- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data.
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
mass1 (float, ~astropy.units.Quantity, optional) β The mass of the first binary component (float assumed in solar masses).
mass2 (float, ~astropy.units.Quantity, optional) β The mass of the second binary component (float assumed in solar masses).
horizon (bool, optional) β If True, return the maximal βhorizonβ luminosity distance, otherwise return the angle-averaged comoving distance.
**kwargs β Additional keyword arguments to ~inspiral_range.waveform.CBCWaveform.
- Returns:
rspec β The calculated inspiral sensitivity PSD [Mpc^2 / Hz].
- Return type:
~gwpy.frequencyseries.FrequencySeries
See also
sensemon_range_psdFor the method based on
LIGO-T030276, also known as LIGO SenseMonitor.inspiral-rangeThe package which does heavy lifting for waveform simulation and cosmology calculations.
- gwexpy.astro.range_spectrogram(hoft: TimeSeries | Spectrogram, stride: float | None = None, fftlength: float | None = None, overlap: float | None = None, window: str | numpy.ndarray = 'hann', method: str = 'median', nproc: int = 1, range_func: Callable | None = None, **rangekwargs) Spectrogram[source]
Calculate the average range spectrogram (Mpc or Mpc^2 / Hz) from strain.
- Parameters:
hoft (~gwpy.timeseries.TimeSeries or ~gwpy.spectrogram.Spectrogram) β record of gravitational-wave strain output from a detector
stride (float, optional) β number of seconds in a single PSD (i.e., step size of spectrogram), required if hoft is an instance of TimeSeries
fftlength (float, optional) β number of seconds in a single FFT
overlap (float, optional) β number of seconds of overlap between FFTs, defaults to the recommended overlap for the given window (if given), or 0
window (str, numpy.ndarray, optional) β window function to apply to timeseries prior to FFT, see
scipy.signal.get_window()for details on acceptable formatsmethod (str, optional) β FFT-averaging method, defaults to median averaging, see
spectrogram()for more detailsnproc (int, optional) β number of CPUs to use in parallel processing of FFTs, default: 1
fmin (float, optional) β low frequency cut-off (Hz), defaults to 1/fftlength
fmax (float, optional) β high frequency cut-off (Hz), defaults to Nyquist frequency of hoft
range_func (callable, optional) β the function to call to generate the range for each stride, defaults to
inspiral_rangeunlessenergyis given as a keyword argument to the range function**rangekwargs (dict, optional) β additional keyword arguments to
burst_range_spectrum()orinspiral_range_psd()(see βNotesβ below), defaults to inspiral range with mass1 = mass2 = 1.4 solar masses
- Returns:
out β time-frequency spectrogram of astrophysical range
- Return type:
~gwpy.spectrogram.Spectrogram
Notes
This method is designed to show the contribution to a gravitational-wave detectorβs sensitive range across frequency bins as a function of time. It supports the range to compact binary inspirals and to unmodelled GW bursts, each a class of transient event.
If inspiral range is requested and fmax exceeds the frequency of the innermost stable circular orbit (ISCO), the output will extend only up to the latter.
See also
gwpy.timeseries.TimeSeries.spectrogramfor the underlying power spectral density estimator
inspiral_range_psdfor the function that computes inspiral range integrand
burst_range_spectrumfor the function that computes burst range integrand
range_timeseriesfor TimeSeries trends of the astrophysical range
- gwexpy.astro.range_timeseries(hoft: TimeSeries | Spectrogram, stride: float | None = None, fftlength: float | None = None, overlap: float | None = None, window: str | numpy.ndarray = 'hann', method: str = 'median', nproc: int = 1, range_func: Callable | None = None, **rangekwargs) TimeSeries[source]
Estimate timeseries of astrophysical range (Mpc).
- Parameters:
hoft (~gwpy.timeseries.TimeSeries, ~gwpy.spectrogram.Spectrogram) β Detector (strain) data from which to estimate range.
stride (float, optional) β Desired step size (seconds) of range timeseries, required if hoft is an instance of TimeSeries.
fftlength (float, optional) β Number of seconds in a single FFT.
overlap (float, optional) β Number of seconds of overlap between FFTs, defaults to the recommended overlap for the given window (if given), or
0.window (str, numpy.ndarray, optional) β Window function to apply to timeseries prior to FFT, see
scipy.signal.get_window()for details on acceptable formats.method (str, optional) β FFT-averaging method, defaults to median averaging, see
spectrogram()for more detailsnproc (int, optional) β number of CPUs to use in parallel processing of FFTs, default: 1
range_func (callable, optional) β the function to call to generate the range for each stride, defaults to
inspiral_rangeunlessenergyis given as a keyword argument to the range function**rangekwargs (dict, optional) β additional keyword arguments to
burst_range()orinspiral_range()(see βNotesβ below), defaults to inspiral range with mass1 = mass2 = 1.4 solar masses
- Returns:
out β timeseries trends of astrophysical range
- Return type:
~gwpy.timeseries.TimeSeries
Notes
This method is designed to quantify a gravitational-wave detectorβs sensitive range as a function of time. It supports the range to compact binary inspirals and to unmodelled GW bursts, each a class of transient event.
See also
gwpy.timeseries.TimeSeries.spectrogramfor the underlying power spectral density estimator
inspiral_rangefor the function that computes inspiral range
burst_rangefor the function that computes burst range
range_spectrogramfor a ~gwpy.spectrogram.Spectrogram of the range integrand
- gwexpy.astro.sensemon_range(psd: FrequencySeries, snr: float = 8, mass1: float = 1.4, mass2: float = 1.4, fmin: float | None = None, fmax: float | None = None, *, horizon: bool = False) units.Quantity[source]
Approximate the inspiral sensitive distance from a GW strain PSD.
This method returns the distance (in megaparsecs) to which a compact binary inspiral with the given component masses would be detectable given the instrumental PSD. The calculation is as defined in:
LIGO-T030276.- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
mass1 (float, ~astropy.units.Quantity, optional) β The mass of the first binary component (float assumed in solar masses).
mass2 (float, ~astropy.units.Quantity, optional) β The mass of the second binary component (float assumed in solar masses).
fmin (float, optional) β The lower frequency cut-off of the integral, default: psd.df.
fmax (float, optional) β The maximum frequency limit of the integral, defaults to innermost stable circular orbit (ISCO) frequency
horizon (bool, optional) β If True, return the maximal βhorizonβ sensitive distance, otherwise return the angle-averaged range.
- Returns:
range β the calculated inspiral range [Mpc]
- Return type:
~astropy.units.Quantity
Examples
Grab some data for LIGO-Livingston around GW150914 and generate a PSD:
>>> from gwpy.timeseries import TimeSeries >>> hoft = TimeSeries.fetch_open_data("H1", 1126259446, 1126259478) >>> hoff = hoft.psd(fftlength=4)
Now we can calculate the
sensemon_range():>>> from gwpy.astro import sensemon_range >>> r = sensemon_range(hoff, fmin: float = 30) >>> print(r) 70.4612102889 Mpc
- gwexpy.astro.sensemon_range_psd(psd: FrequencySeries, snr: float = 8, mass1: float = 1.4, mass2: float = 1.4, *, horizon: bool = False) FrequencySeries[source]
Approximate the inspiral sensitive distance PSD from a GW strain PSD.
This method returns the power spectral density (in
Mpc**2 / Hz) to which a compact binary inspiral with the given component masses would be detectable given the instrumental PSD. The calculation is defined in:LIGO-T030276.- Parameters:
psd (~gwpy.frequencyseries.FrequencySeries) β The instrumental power-spectral-density data.
snr (float, optional) β The signal-to-noise ratio for which to calculate range.
mass1 (float, ~astropy.units.Quantity, optional) β The mass of the first binary component (float assumed in solar masses).
mass2 (float, ~astropy.units.Quantity, optional) β The mass of the second binary component (float assumed in solar masses).
horizon (bool, optional) β If True, return the maximal βhorizonβ sensitive distance, otherwise return the angle-averaged range.
- Returns:
rspec β The calculated inspiral sensitivity PSD [Mpc^2 / Hz].
- Return type:
~gwpy.frequencyseries.FrequencySeries