Reference on kuibit.gw_utils

The module gw_utils contains convenience functions and structures to analyze and work with gravitational waves.

First, the Detectors object is defined. Detectors is a named tuple with fields “hanford”, “livingston”, “virgo”. This is used every time we deal with specific detectors.

The functions provided are:

class kuibit.gw_utils.Detectors(hanford, livingston, virgo)

Create new instance of Detectors(hanford, livingston, virgo)

hanford

Alias for field number 0

livingston

Alias for field number 1

virgo

Alias for field number 2

kuibit.gw_utils.Schwarzschild_radius_to_tortoise(radii, mass)[source]

Transform radial coordinates radii to tortoise coordinates assuming mass mass.

Equation (26) in 1307.5307.

Parameters:
  • radii (float or 1D NumPy array) – Radius in Schwarzschild coordinates.

  • mass (float) – ADM mass.

Returns:

Tortoise radii.

Return type:

float or 1D NumPy array

kuibit.gw_utils.antenna_responses(theta, phi, polarization=0)[source]

Return the antenna response pattern of a detector on the z = 0 plane with the arms on the x and y directions for a given localization defined by the spherical angles theta and phi.

Parameters:
  • theta (float) – Meridional angle.

  • phi (float) – Azimuthal angle.

  • polarization (float) – Polarization angle of the wave.

Returns:

Antenna response for cross and plus polarizations (in this order).

Return type:

tuple of floats.

kuibit.gw_utils.antenna_responses_from_sky_localization(right_ascension, declination, time_utc, polarization=0)[source]

Return the antenna responses for Hanford, Livingston and Virgo for a given source.

See, http://research.physics.illinois.edu/cta/movies/bhbh_sim/wavestrain.html.

utc_time has to have the following formatting: %Y-%m-%d %H:%M, (eg 2015-09-14 09:50:45)

Parameters:
  • right_ascension (float) – Right ascension of the source in degrees

  • declination (str) – Declination of the source in degrees

  • time_utc – UTC time of the event

  • polarization (float) – Polarization of the wave

Returns antenna_pattern:

Cross and plus antenna pattern for the different interferometers.

Return type:

namedtuple with fields hanford, livingston, and virgo

kuibit.gw_utils.effective_amplitude_spectral_density(strain, *args, window_function=None, **kwargs)[source]

Return the effective amplitude spectral density for a given strain.

If window_function is not None, the window will be applied to the signal. All the unknown arguments are passed to the window function.

The effective amplitude spectral density is computed for a strain (h) as \(h_{\rm eff}(f) = f * \sqrt{(|h_{+}|^2 + |h_{\times}|^2) / 2.0}\) (see for instance 1604.00246)

Parameters:
  • strain (TimeSeries) – Gravitational-wave strain.

  • window_function (callable, str, or None) – If not None, apply window_function to the series before computing the strain.

  • kwargs (args,) – All the additional parameters are passed to the window function.

Returns:

Effective amplitude spectral density

Return type:

FrequencySeries

kuibit.gw_utils.luminosity_distance_to_redshift(luminosity_distance, Omega_m=0.309, Omega_L=0.691, initial_guess=0.1)[source]

Compute redshift from luminosity distance in Mpc assuming the ΛCDM cosmology.

This function is useful to correctly reproduce observed signals from cosmological sources (e.g., binary black holes far away).

The redshift is computed via root-finding, so an initial guess is needed.

Parameters:
  • luminosity_distance (float) – Luminosity distance in megaparsec.

  • Omega_m (float) – \(\Omega_m\) (matter) cosmological parameter.

  • Omega_L (float) – \(\Omega_m\) (dark energy) cosmological parameter.

  • initial_guess (float) – Initial guess to the redshift for the root-finding routine.

Returns z:

Redshift

Rtype z:

float

kuibit.gw_utils.ra_dec_to_theta_phi(right_ascension, declination, time_utc)[source]

Compute the spherical angles theta and phi for Hanford, Livingston and Virgo for a given source localization.

utc_time has to have the following formatting: %Y-%m-%d %H:%M, (eg 2015-09-14 09:50:45)

Parameters:
  • right_ascension (float) – Right ascension of the source in degrees.

  • declination (str) – Declination of the source in degrees.

  • time_utc – UTC time of the event.

Returns spherical coordinates:

Theta, phi for the different detectors.

Return type:

namedtuple with fields hanford, livingston, and virgo

kuibit.gw_utils.retarded_times_to_coordinate_times(retarded_times, radii, mass)[source]
Compute the coordinate times corresponding to the retarded times at the

coordinate radii.

First, the tortoise radius is computed from radii, then the coordinate times are computed with \(t = u + r_{tortoise}(radii)\).

This function is used to extrapolate gravitational waves to infinity.

Parameters:
  • retarded_times (float or 1D NumPy array) – Retarded times.

  • radii (float or 1D NumPy array) – Radii of evaluation.

  • mass (float) – ADM mass (needed to compute tortoise radius).

Returns:

Coordinate times corresponding to the given retarded times and evaluation radii.

Return type:

float or 1D NumPy array

kuibit.gw_utils.sYlm(ss, ll, mm, theta, phi)[source]

Compute spin-weighted spherical harmonics at the angles theta and phi.

When ss = 0, these are spherical harmonics.

Parameters:
  • ss (int) – Spin weight.

  • ll (int) – \(l\) multipolar number.

  • mm (int) – \(m\) multipolar number.

  • theta (float) – Meridional angle.

  • phi (float) – Azimuthal angle.

Returns sYlm:

Spin-weighted spherical harmonic evaluated at theta and phi.

Rtype sYlm:

float

kuibit.gw_utils.signal_to_noise_ratio_from_strain(h, *args, noise=None, fmin=0, fmax=inf, window_function=None, **kwargs)[source]

Return the signal to noise ratio given a strain and a power spectal density distribution for a detector.

If window_function is not None, the window will be applied to the signal. All the unknown arguments are passed to the window function.

The SNR is computed as \(sqrt of 4 \int_fmin^fmax |\tilde{h} f|^2 / Sn(f) d f\) (equation from 1408.0740)

Parameters:
  • h (TimeSeries) – Gravitational-wave strain.

  • noise (FrequencySeries) – Power spectral density of the noise of the detector.

  • fmin (float) – Minimum frequency over which to compute the SNR.

  • fmax (float) – Maximum frequency over which to compute the SNR.

  • window_function (callable, str, or None) – If not None, apply window_function to the series before computing the strain.

  • kwargs (args,) – All the additional parameters are passed to the window function.

Returns:

Signal-to-noise ratio.

Return type:

float