pycbc.detector package

Submodules

pycbc.detector.ground module

This module provides utilities for calculating detector responses and timing between ground-based observatories.

class pycbc.detector.ground.Detector(detector_name, reference_time=1126259462.0)[source]

Bases: object

A gravitational wave detector

antenna_pattern(right_ascension, declination, polarization, t_gps, frequency=0, polarization_type='tensor')[source]

Return the detector response.

Parameters:
  • right_ascension (float or numpy.ndarray) – The right ascension of the source

  • declination (float or numpy.ndarray) – The declination of the source

  • polarization (float or numpy.ndarray) – The polarization angle of the source

  • polarization_type (string flag: Tensor, Vector or Scalar) – The gravitational wave polarizations. Default: ‘Tensor’

Returns:

  • fplus(default) or fx or fb (float or numpy.ndarray) – The plus or vector-x or breathing polarization factor for this sky location / orientation

  • fcross(default) or fy or fl (float or numpy.ndarray) – The cross or vector-y or longitudnal polarization factor for this sky location / orientation

effective_distance(distance, ra, dec, pol, time, inclination)[source]

Distance scaled to account for amplitude factors

The effective distance of the source. This scales the distance so that the amplitude is equal to a source which is optimally oriented with respect to the detector. For fixed detector-frame intrinsic parameters this is a measure of the expected signal strength.

Parameters:
  • distance (float) – Source luminosity distance in megaparsecs

  • ra (float) – The right ascension in radians

  • dec (float) – The declination in radians

  • pol (float) – Polarization angle of the gravitational wave in radians

  • time (float) – GPS time in seconds

  • inclination – The inclination of the binary’s orbital plane

Returns:

eff_dist – The effective distance of the source

Return type:

float

get_icrs_pos()[source]

Transforms GCRS frame to ICRS frame

Returns:

loc – ICRS coordinates in cartesian system

Return type:

numpy.ndarray shape (3,1) units: AU

gmst_estimate(gps_time)[source]
lal()[source]

Return lal data type detector instance

light_travel_time_to_detector(det)[source]

Return the light travel time from this detector :param det: The other detector to determine the light travel time to. :type det: Detector

Returns:

time – The light travel time in seconds

Return type:

float

optimal_orientation(t_gps)[source]
Return the optimal orientation in right ascension and declination

for a given GPS time.

Parameters:

t_gps (float) – Time in gps seconds

Returns:

  • ra (float) – Right ascension that is optimally oriented for the detector

  • dec (float) – Declination that is optimally oriented for the detector

project_wave(hp, hc, ra, dec, polarization, method='lal', reference_time=None)[source]

Return the strain of a waveform as measured by the detector. Apply the time shift for the given detector relative to the assumed geocentric frame and apply the antenna patterns to the plus and cross polarizations.

Parameters:
  • hp (pycbc.types.TimeSeries) – Plus polarization of the GW

  • hc (pycbc.types.TimeSeries) – Cross polarization of the GW

  • ra (float) – Right ascension of source location

  • dec (float) – Declination of source location

  • polarization (float) – Polarization angle of the source

  • method ({'lal', 'constant', 'vary_polarization'}) – The method to use for projecting the polarizations into the detector frame. Default is ‘lal’.

  • reference_time (float, Optional) – The time to use as, a reference for some methods of projection. Used by ‘constant’ and ‘vary_polarization’ methods. Uses average time if not provided.

set_gmst_reference()[source]
time_delay_from_detector(other_detector, right_ascension, declination, t_gps)[source]

Return the time delay from the given to detector for a signal with the given sky location; i.e. return t1 - t2 where t1 is the arrival time in this detector and t2 is the arrival time in the other detector. Note that this would return the same value as time_delay_from_earth_center if other_detector was geocentric. :param other_detector: A detector instance. :type other_detector: detector.Detector :param right_ascension: The right ascension (in rad) of the signal. :type right_ascension: float :param declination: The declination (in rad) of the signal. :type declination: float :param t_gps: The GPS time (in s) of the signal. :type t_gps: float

Returns:

The arrival time difference between the detectors.

Return type:

float

time_delay_from_earth_center(right_ascension, declination, t_gps)[source]

Return the time delay from the earth center

time_delay_from_location(other_location, right_ascension, declination, t_gps)[source]

Return the time delay from the given location to detector for a signal with the given sky location In other words return t1 - t2 where t1 is the arrival time in this detector and t2 is the arrival time in the other location.

Parameters:
  • other_location (numpy.ndarray of coordinates) – A detector instance.

  • right_ascension (float) – The right ascension (in rad) of the signal.

  • declination (float) – The declination (in rad) of the signal.

  • t_gps (float) – The GPS time (in s) of the signal.

Returns:

The arrival time difference between the detectors.

Return type:

float

pycbc.detector.ground.add_detector_on_earth(name, longitude, latitude, yangle=0, xangle=None, height=0, xlength=4000, ylength=4000, xaltitude=0, yaltitude=0)[source]

Add a new detector on the earth

Parameters:
  • name (str) – two-letter name to identify the detector

  • longitude (float) – Longitude in radians using geodetic coordinates of the detector

  • latitude (float) – Latitude in radians using geodetic coordinates of the detector

  • yangle (float) – Azimuthal angle of the y-arm (angle drawn from pointing north)

  • xangle (float) – Azimuthal angle of the x-arm (angle drawn from point north). If not set we assume a right angle detector following the right-hand rule.

  • xaltitude (float) – The altitude angle of the x-arm measured from the local horizon.

  • yaltitude (float) – The altitude angle of the y-arm measured from the local horizon.

  • height (float) – The height in meters of the detector above the standard reference ellipsoidal earth

pycbc.detector.ground.get_available_detectors()[source]

List the available detectors

pycbc.detector.ground.get_available_lal_detectors()[source]

Return list of detectors known in the currently sourced lalsuite. This function will query lalsuite about which detectors are known to lalsuite. Detectors are identified by a two character string e.g. ‘K1’, but also by a longer, and clearer name, e.g. KAGRA. This function returns both. As LAL doesn’t really expose this functionality we have to make some assumptions about how this information is stored in LAL. Therefore while we hope this function will work correctly, it’s possible it will need updating in the future. Better if lal would expose this information properly.

pycbc.detector.ground.gmst_accurate(gps_time)[source]
pycbc.detector.ground.load_detector_config(config_files)[source]

Add custom detectors from a configuration file

Parameters:

config_files (str or list of strs) – The config file(s) which specify new detectors

pycbc.detector.ground.overhead_antenna_pattern(right_ascension, declination, polarization)[source]

Return the antenna pattern factors F+ and Fx as a function of sky location and polarization angle for a hypothetical interferometer located at the north pole. Angles are in radians. Declinations of ±π/2 correspond to the normal to the detector plane (i.e. overhead and underneath) while the point with zero right ascension and declination is the direction of one of the interferometer arms. :param right_ascension: :type right_ascension: float :param declination: :type declination: float :param polarization: :type polarization: float

Returns:

  • f_plus (float)

  • f_cros (float)

pycbc.detector.ground.ppdets(ifos, separator=', ')[source]

Pretty-print a list (or set) of detectors: return a string listing the given detectors alphabetically and separated by the given string (comma by default).

pycbc.detector.ground.single_arm_frequency_response(f, n, arm_length)[source]

The relative amplitude factor of the arm response due to signal delay. This is relevant where the long-wavelength approximation no longer applies)

pycbc.detector.space module

This module provides utilities for simulating the GW response of space-based observatories.

class pycbc.detector.space.SpaceDetector(detector_name, reference_time=None, backend='LDC', **kwargs)[source]

Bases: AbsSpaceDet

Space-based detector.

Parameters:
  • detector_name (str) – The name of the detector. Accepts any output from get_available_space_detectors.

  • reference_time (float (optional)) – The reference time in seconds of the signal in the SSB frame. This is defined such that the detector mission start time corresponds to 0. Default None.

  • backend (str (optional)) – The backend architecture to use for generating TDI. Accepts ‘LDC’ or ‘FLR’. Default ‘LDC’.

project_wave(hp, hc, lamb, beta, *args, **kwargs)[source]

Placeholder for evaluating the TDI channels from the GW projections.

property sky_coords

List the sky coordinate names for the detector class.

pycbc.detector.space.get_available_space_detectors()[source]

List the available space detectors

Module contents