inspiral_range¶
-
gwpy.astro.inspiral_range(psd, snr=
8
, mass1=1.4
, mass2=1.4
, fmin=None
, fmax=None
, horizon=False
, **kwargs)[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
FrequencySeries
the instrumental power-spectral-density data
- snr
float
, optional the signal-to-noise ratio for which to calculate range, default:
8
- mass1
float
,Quantity
, optional the mass (
float
assumed in solar masses) of the first binary component, default:1.4
- mass2
float
,Quantity
, optional the mass (
float
assumed in solar masses) of the second binary component, default:1.4
- 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, default:False
- **kwargs
dict
, optional additional keyword arguments to
CBCWaveform
- psd
- Returns:¶
- range
Quantity
the calculated inspiral range [Mpc]
- range
See also
sensemon_range
for the method based on LIGO-T030276, also known as LIGO SenseMonitor
inspiral-range
the package which does heavy lifting for waveform simulation and cosmology calculations
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=30) >>> print(r) 70.4612102889 Mpc