pycbc.coordinates package

Submodules

pycbc.coordinates.base module

Base coordinate transformations, this module provides transformations between cartesian and spherical coordinates.

pycbc.coordinates.base.cartesian_to_spherical(x, y, z)[source]

Maps cartesian coordinates (x,y,z) to spherical coordinates (rho,phi,theta) where phi is in [0,2*pi] and theta is in [0,pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

  • rho ({numpy.array, float}) – The radial amplitude.

  • phi ({numpy.array, float}) – The azimuthal angle.

  • theta ({numpy.array, float}) – The polar angle.

pycbc.coordinates.base.cartesian_to_spherical_azimuthal(x, y)[source]

Calculates the azimuthal angle in spherical coordinates from Cartesian coordinates. The azimuthal angle is in [0,2*pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

Returns:

phi – The azimuthal angle.

Return type:

{numpy.array, float}

pycbc.coordinates.base.cartesian_to_spherical_polar(x, y, z)[source]

Calculates the polar angle in spherical coordinates from Cartesian coordinates. The polar angle is in [0,pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

theta – The polar angle.

Return type:

{numpy.array, float}

pycbc.coordinates.base.cartesian_to_spherical_rho(x, y, z)[source]

Calculates the magnitude in spherical coordinates from Cartesian coordinates.

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

rho – The radial amplitude.

Return type:

{numpy.array, float}

pycbc.coordinates.base.spherical_to_cartesian(rho, phi, theta)[source]

Maps spherical coordinates (rho,phi,theta) to cartesian coordinates (x,y,z) where phi is in [0,2*pi] and theta is in [0,pi].

Parameters:
  • rho ({numpy.array, float}) – The radial amplitude.

  • phi ({numpy.array, float}) – The azimuthal angle.

  • theta ({numpy.array, float}) – The polar angle.

Returns:

  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

pycbc.coordinates.space module

This module provides coordinate transformations related to space-borne detectors, such as coordinate transformations between space-borne detectors and ground-based detectors. Note that current LISA orbit used in this module is a circular orbit, need to be replaced by a more realistic and general orbit model in the near future.

pycbc.coordinates.space.earth_position_ssb(t_geo)[source]

Calculating the position vector and angular displacement of the Earth in the SSB frame, at a given time. By using Astropy.

Parameters:

t_geo (float) – The time when a GW signal arrives at the origin of geocentric frame, or any other time you want.

Returns:

  • (p, alpha) (tuple)

  • p (numpy.array) – The position vector of the Earth in the SSB frame. In the unit of ‘m’.

  • alpha (float) – The angular displacement of the Earth in the SSB frame. In the unit of ‘radian’.

pycbc.coordinates.space.geo_to_lisa(t_geo, longitude_geo, latitude_geo, polarization_geo, t0=7365189.431698299, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the geocentric frame to the LISA frame.

Parameters:
  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_lisa, longitude_lisa, latitude_lisa, polarization_lisa) (tuple)

  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

pycbc.coordinates.space.geo_to_ssb(t_geo, longitude_geo, latitude_geo, polarization_geo, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the geocentric frame to the SSB frame.

Parameters:
  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_ssb, longitude_ssb, latitude_ssb, polarization_ssb) (tuple)

  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

pycbc.coordinates.space.lisa_position_ssb(t_lisa, t0=7365189.431698299)[source]

Calculating the position vector and angular displacement of LISA in the SSB frame, at a given time. This function assumes LISA’s barycenter is orbiting around a circular orbit within the ecliptic behind the Earth. The period of it is one year.

Parameters:
  • t_lisa (float) – The time when a GW signal arrives at the origin of LISA frame, or any other time you want.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (p, alpha) (tuple)

  • p (numpy.array) – The position vector of LISA in the SSB frame. In the unit of ‘m’.

  • alpha (float) – The angular displacement of LISA in the SSB frame. In the unit of ‘radian’.

pycbc.coordinates.space.lisa_to_geo(t_lisa, longitude_lisa, latitude_lisa, polarization_lisa, t0=7365189.431698299, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the LISA frame to the geocentric frame.

Parameters:
  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_geo, longitude_geo, latitude_geo, polarization_geo) (tuple)

  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The ecliptic longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The ecliptic latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

pycbc.coordinates.space.lisa_to_ssb(t_lisa, longitude_lisa, latitude_lisa, polarization_lisa, t0=7365189.431698299)[source]

Converting the arrive time, the sky localization, and the polarization from the LISA frame to the SSB frame.

Parameters:
  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (t_ssb, longitude_ssb, latitude_ssb, polarization_ssb) (tuple)

  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

pycbc.coordinates.space.localization_to_propagation_vector(longitude, latitude, use_astropy=True, frame=None)[source]

Converting the sky localization to the corresponding propagation unit vector of a GW signal.

Parameters:
  • longitude (float) – The longitude, in the unit of ‘radian’.

  • latitude (float) – The latitude, in the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

Returns:

[[x], [y], [z]] – The propagation unit vector of that GW signal.

Return type:

numpy.array

pycbc.coordinates.space.polarization_newframe(polarization, k, rotation_matrix, use_astropy=True, old_frame=None, new_frame=None)[source]

Converting a polarization angle from a frame to a new frame by using rotation matrix method.

Parameters:
  • polarization (float) – The polarization angle in the old frame, in the unit of ‘radian’.

  • k (numpy.array) – The propagation unit vector of a GW signal in the old frame.

  • rotation_matrix (numpy.array) – The rotation matrix (of frame basis) from the old frame to the new frame.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • old_frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

  • new_frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None. The new frame for the new polarization angle.

Returns:

polarization_new_frame – The polarization angle in the new frame of that GW signal.

Return type:

float

pycbc.coordinates.space.propagation_vector_to_localization(k, use_astropy=True, frame=None)[source]

Converting the propagation unit vector to the corresponding sky localization of a GW signal.

Parameters:
  • k (numpy.array) – The propagation unit vector of a GW signal.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

Returns:

(longitude, latitude) – The sky localization of that GW signal.

Return type:

tuple

pycbc.coordinates.space.rotation_matrix_ssb_to_geo(epsilon=0.40909262775014904)[source]

The rotation matrix (of frame basis) from SSB frame to geocentric frame.

Parameters:

epsilon (float) – The Earth’s axial tilt (obliquity), in the unit of ‘radian’.

Returns:

r – A 3x3 rotation matrix from SSB frame to geocentric frame.

Return type:

numpy.array

pycbc.coordinates.space.rotation_matrix_ssb_to_lisa(alpha)[source]

The rotation matrix (of frame basis) from SSB frame to LISA frame. This function assumes the angle between LISA plane and the ecliptic is 60 degrees, and the period of LISA’s self-rotation and orbital revolution is both one year.

Parameters:

alpha (float) – The angular displacement of LISA in SSB frame. In the unit of ‘radian’.

Returns:

r_total – A 3x3 rotation matrix from SSB frame to LISA frame.

Return type:

numpy.array

pycbc.coordinates.space.ssb_to_geo(t_ssb, longitude_ssb, latitude_ssb, polarization_ssb, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the SSB frame to the geocentric frame.

Parameters:
  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_geo, longitude_geo, latitude_geo, polarization_geo) (tuple)

  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

pycbc.coordinates.space.ssb_to_lisa(t_ssb, longitude_ssb, latitude_ssb, polarization_ssb, t0=7365189.431698299)[source]

Converting the arrive time, the sky localization, and the polarization from the SSB frame to the LISA frame.

Parameters:
  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (t_lisa, longitude_lisa, latitude_lisa, polarization_lisa) (tuple)

  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

pycbc.coordinates.space.t_geo_from_ssb(t_ssb, longitude_ssb, latitude_ssb, use_astropy=True, frame=None)[source]

Calculating the time when a GW signal arrives at the barycenter of the Earth, by using the time and sky localization in SSB frame.

Parameters:
  • t_ssb (float) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

Returns:

t_geo – The time when a GW signal arrives at the origin of geocentric frame.

Return type:

float

pycbc.coordinates.space.t_lisa_from_ssb(t_ssb, longitude_ssb, latitude_ssb, t0=7365189.431698299)[source]

Calculating the time when a GW signal arrives at the barycenter of LISA, by using the time and sky localization in SSB frame.

Parameters:
  • t_ssb (float) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

t_lisa – The time when a GW signal arrives at the origin of LISA frame.

Return type:

float

pycbc.coordinates.space.t_ssb_from_t_geo(t_geo, longitude_ssb, latitude_ssb, use_astropy=True, frame=None)[source]

Calculating the time when a GW signal arrives at the barycenter of SSB, by using the time in geocentric frame and sky localization in SSB frame.

Parameters:
  • t_geo (float) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

Returns:

t_ssb – The time when a GW signal arrives at the origin of SSB frame.

Return type:

float

pycbc.coordinates.space.t_ssb_from_t_lisa(t_lisa, longitude_ssb, latitude_ssb, t0=7365189.431698299)[source]

Calculating the time when a GW signal arrives at the barycenter of SSB, by using the time in LISA frame and sky localization in SSB frame.

Parameters:
  • t_lisa (float) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

t_ssb – The time when a GW signal arrives at the origin of SSB frame.

Return type:

float

Module contents

This modules provides functions for base coordinate transformations, and more advanced transformations between ground-based detectors and space-borne detectors.

pycbc.coordinates.cartesian_to_spherical(x, y, z)[source]

Maps cartesian coordinates (x,y,z) to spherical coordinates (rho,phi,theta) where phi is in [0,2*pi] and theta is in [0,pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

  • rho ({numpy.array, float}) – The radial amplitude.

  • phi ({numpy.array, float}) – The azimuthal angle.

  • theta ({numpy.array, float}) – The polar angle.

pycbc.coordinates.cartesian_to_spherical_azimuthal(x, y)[source]

Calculates the azimuthal angle in spherical coordinates from Cartesian coordinates. The azimuthal angle is in [0,2*pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

Returns:

phi – The azimuthal angle.

Return type:

{numpy.array, float}

pycbc.coordinates.cartesian_to_spherical_polar(x, y, z)[source]

Calculates the polar angle in spherical coordinates from Cartesian coordinates. The polar angle is in [0,pi].

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

theta – The polar angle.

Return type:

{numpy.array, float}

pycbc.coordinates.cartesian_to_spherical_rho(x, y, z)[source]

Calculates the magnitude in spherical coordinates from Cartesian coordinates.

Parameters:
  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

Returns:

rho – The radial amplitude.

Return type:

{numpy.array, float}

pycbc.coordinates.earth_position_ssb(t_geo)[source]

Calculating the position vector and angular displacement of the Earth in the SSB frame, at a given time. By using Astropy.

Parameters:

t_geo (float) – The time when a GW signal arrives at the origin of geocentric frame, or any other time you want.

Returns:

  • (p, alpha) (tuple)

  • p (numpy.array) – The position vector of the Earth in the SSB frame. In the unit of ‘m’.

  • alpha (float) – The angular displacement of the Earth in the SSB frame. In the unit of ‘radian’.

pycbc.coordinates.geo_to_lisa(t_geo, longitude_geo, latitude_geo, polarization_geo, t0=7365189.431698299, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the geocentric frame to the LISA frame.

Parameters:
  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_lisa, longitude_lisa, latitude_lisa, polarization_lisa) (tuple)

  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

pycbc.coordinates.geo_to_ssb(t_geo, longitude_geo, latitude_geo, polarization_geo, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the geocentric frame to the SSB frame.

Parameters:
  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_ssb, longitude_ssb, latitude_ssb, polarization_ssb) (tuple)

  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

pycbc.coordinates.lisa_position_ssb(t_lisa, t0=7365189.431698299)[source]

Calculating the position vector and angular displacement of LISA in the SSB frame, at a given time. This function assumes LISA’s barycenter is orbiting around a circular orbit within the ecliptic behind the Earth. The period of it is one year.

Parameters:
  • t_lisa (float) – The time when a GW signal arrives at the origin of LISA frame, or any other time you want.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (p, alpha) (tuple)

  • p (numpy.array) – The position vector of LISA in the SSB frame. In the unit of ‘m’.

  • alpha (float) – The angular displacement of LISA in the SSB frame. In the unit of ‘radian’.

pycbc.coordinates.lisa_to_geo(t_lisa, longitude_lisa, latitude_lisa, polarization_lisa, t0=7365189.431698299, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the LISA frame to the geocentric frame.

Parameters:
  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_geo, longitude_geo, latitude_geo, polarization_geo) (tuple)

  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The ecliptic longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The ecliptic latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

pycbc.coordinates.lisa_to_ssb(t_lisa, longitude_lisa, latitude_lisa, polarization_lisa, t0=7365189.431698299)[source]

Converting the arrive time, the sky localization, and the polarization from the LISA frame to the SSB frame.

Parameters:
  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (t_ssb, longitude_ssb, latitude_ssb, polarization_ssb) (tuple)

  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

pycbc.coordinates.localization_to_propagation_vector(longitude, latitude, use_astropy=True, frame=None)[source]

Converting the sky localization to the corresponding propagation unit vector of a GW signal.

Parameters:
  • longitude (float) – The longitude, in the unit of ‘radian’.

  • latitude (float) – The latitude, in the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

Returns:

[[x], [y], [z]] – The propagation unit vector of that GW signal.

Return type:

numpy.array

pycbc.coordinates.polarization_newframe(polarization, k, rotation_matrix, use_astropy=True, old_frame=None, new_frame=None)[source]

Converting a polarization angle from a frame to a new frame by using rotation matrix method.

Parameters:
  • polarization (float) – The polarization angle in the old frame, in the unit of ‘radian’.

  • k (numpy.array) – The propagation unit vector of a GW signal in the old frame.

  • rotation_matrix (numpy.array) – The rotation matrix (of frame basis) from the old frame to the new frame.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • old_frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

  • new_frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None. The new frame for the new polarization angle.

Returns:

polarization_new_frame – The polarization angle in the new frame of that GW signal.

Return type:

float

pycbc.coordinates.propagation_vector_to_localization(k, use_astropy=True, frame=None)[source]

Converting the propagation unit vector to the corresponding sky localization of a GW signal.

Parameters:
  • k (numpy.array) – The propagation unit vector of a GW signal.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

  • frame (astropy.coordinates) – The frame from astropy.coordinates if use_astropy is True, the default is None.

Returns:

(longitude, latitude) – The sky localization of that GW signal.

Return type:

tuple

pycbc.coordinates.rotation_matrix_ssb_to_geo(epsilon=0.40909262775014904)[source]

The rotation matrix (of frame basis) from SSB frame to geocentric frame.

Parameters:

epsilon (float) – The Earth’s axial tilt (obliquity), in the unit of ‘radian’.

Returns:

r – A 3x3 rotation matrix from SSB frame to geocentric frame.

Return type:

numpy.array

pycbc.coordinates.rotation_matrix_ssb_to_lisa(alpha)[source]

The rotation matrix (of frame basis) from SSB frame to LISA frame. This function assumes the angle between LISA plane and the ecliptic is 60 degrees, and the period of LISA’s self-rotation and orbital revolution is both one year.

Parameters:

alpha (float) – The angular displacement of LISA in SSB frame. In the unit of ‘radian’.

Returns:

r_total – A 3x3 rotation matrix from SSB frame to LISA frame.

Return type:

numpy.array

pycbc.coordinates.spherical_to_cartesian(rho, phi, theta)[source]

Maps spherical coordinates (rho,phi,theta) to cartesian coordinates (x,y,z) where phi is in [0,2*pi] and theta is in [0,pi].

Parameters:
  • rho ({numpy.array, float}) – The radial amplitude.

  • phi ({numpy.array, float}) – The azimuthal angle.

  • theta ({numpy.array, float}) – The polar angle.

Returns:

  • x ({numpy.array, float}) – X-coordinate.

  • y ({numpy.array, float}) – Y-coordinate.

  • z ({numpy.array, float}) – Z-coordinate.

pycbc.coordinates.ssb_to_geo(t_ssb, longitude_ssb, latitude_ssb, polarization_ssb, use_astropy=True)[source]

Converting the arrive time, the sky localization, and the polarization from the SSB frame to the geocentric frame.

Parameters:
  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

  • use_astropy (bool) – Using Astropy to calculate the sky localization or not. Default is True.

Returns:

  • (t_geo, longitude_geo, latitude_geo, polarization_geo) (tuple)

  • t_geo (float or numpy.array) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_geo (float or numpy.array) – The longitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • latitude_geo (float or numpy.array) – The latitude of a GW signal in geocentric frame. In the unit of ‘radian’.

  • polarization_geo (float or numpy.array) – The polarization angle of a GW signal in geocentric frame. In the unit of ‘radian’.

pycbc.coordinates.ssb_to_lisa(t_ssb, longitude_ssb, latitude_ssb, polarization_ssb, t0=7365189.431698299)[source]

Converting the arrive time, the sky localization, and the polarization from the SSB frame to the LISA frame.

Parameters:
  • t_ssb (float or numpy.array) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float or numpy.array) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float or numpy.array) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • polarization_ssb (float or numpy.array) – The polarization angle of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

  • (t_lisa, longitude_lisa, latitude_lisa, polarization_lisa) (tuple)

  • t_lisa (float or numpy.array) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_lisa (float or numpy.array) – The longitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • latitude_lisa (float or numpy.array) – The latitude of a GW signal in LISA frame, in the unit of ‘radian’.

  • polarization_lisa (float or numpy.array) – The polarization angle of a GW signal in LISA frame. In the unit of ‘radian’.

pycbc.coordinates.t_geo_from_ssb(t_ssb, longitude_ssb, latitude_ssb, use_astropy=True, frame=None)[source]

Calculating the time when a GW signal arrives at the barycenter of the Earth, by using the time and sky localization in SSB frame.

Parameters:
  • t_ssb (float) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

Returns:

t_geo – The time when a GW signal arrives at the origin of geocentric frame.

Return type:

float

pycbc.coordinates.t_lisa_from_ssb(t_ssb, longitude_ssb, latitude_ssb, t0=7365189.431698299)[source]

Calculating the time when a GW signal arrives at the barycenter of LISA, by using the time and sky localization in SSB frame.

Parameters:
  • t_ssb (float) – The time when a GW signal arrives at the origin of SSB frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

t_lisa – The time when a GW signal arrives at the origin of LISA frame.

Return type:

float

pycbc.coordinates.t_ssb_from_t_geo(t_geo, longitude_ssb, latitude_ssb, use_astropy=True, frame=None)[source]

Calculating the time when a GW signal arrives at the barycenter of SSB, by using the time in geocentric frame and sky localization in SSB frame.

Parameters:
  • t_geo (float) – The time when a GW signal arrives at the origin of geocentric frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

Returns:

t_ssb – The time when a GW signal arrives at the origin of SSB frame.

Return type:

float

pycbc.coordinates.t_ssb_from_t_lisa(t_lisa, longitude_ssb, latitude_ssb, t0=7365189.431698299)[source]

Calculating the time when a GW signal arrives at the barycenter of SSB, by using the time in LISA frame and sky localization in SSB frame.

Parameters:
  • t_lisa (float) – The time when a GW signal arrives at the origin of LISA frame. In the unit of ‘s’.

  • longitude_ssb (float) – The ecliptic longitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • latitude_ssb (float) – The ecliptic latitude of a GW signal in SSB frame. In the unit of ‘radian’.

  • t0 (float) – The initial time offset of LISA, in the unit of ‘s’, default is 7365189.431698299. This makes sure LISA is behind the Earth by 19-23 degrees.

Returns:

t_ssb – The time when a GW signal arrives at the origin of SSB frame.

Return type:

float