Localization package

Submodules

simple_pe.localization.event module

class simple_pe.localization.event.Event(dist, ra, dec, phi, psi, cosi, mchirp, t_gps)[source]

Bases: object

Class to hold the details of the event. This contains the sky location, orientation, mass of the event. The Event class can also have details of the active GW network, including the detectors, sensivitities and SNRs added. These can be used to calculate the localization of the event.

add_network(network)[source]

Calculate the sensitivities and SNRs for the various detectors in network

Parameters

network: network.Network

an object containing details of the network

alpha_net(mirror=False)[source]

Get the relative network sensitivity to the second polarization

Parameters

mirror: boolean

indicating whether we are considering the mirror location

Returns

float

value of alpha_network

calculate_mirror()[source]

Calculate the mirror location and detector sensitivity there

calculate_sensitivity()[source]

Calculate the network sensitivity to the event, given the sky location and masses.

combined_loc(method)[source]

Calculate the area from original and mirror locations for the given method

Parameters

method: str

localization method to use

classmethod from_params(params)[source]

Give a set of parameters, as used in the first 2 years paper, and use these to initialize an event

Parameters

params: dict

parameters in form used by first 2 years paper

classmethod from_snrs(net, snrs, times, mchirp, ra=None, dec=None)[source]

Give a network with SNR and time in each detector and use this to populate the event information. If ra and dec are provided, they are used. If not, then sky location is inferred from the time of arrival in different detectors. For fewer than 3 detectors, sky location is chosen arbitrarily among possible points.

Parameters

net: network.Network

a network with SNR and time for each detector

snrs: dict

the complex snr in each detector

times: dict

the time in each detector

mchirp: float

the chirp mass of the event

ra: float

the right ascenscion of the source (optional)

dec: float

the declination of the source (optional)

get_data(data)[source]

Get the relevant data for each detector and return it as an array

Parameters

data: str

string giving name of data field

Returns

np.array

with the data (for all ifos)

get_f(mirror=False)[source]

Get the network sensitivity to plus and cross in dominant polarization

Parameters

mirror: boolean

indicating whether we are considering the mirror location

Returns

np.array

length 2 array containing F_+, F_x response

get_fsig(mirror=False)[source]

Get F_plus/F_cross multiplied by sigma (sensitivity) for each detector

Parameters

mirror: boolean

indicating whether we are considering the mirror location

Returns

np.array

with the sensitivities of the detectors

get_snr(dt_i=None)[source]

Calculate the snr for each detector at the requested time offset

Parameters

dt_i: np.array

time shift to be applied in each detector

Returns

np.array

the complex snr for each detector

localization_factors(method, mirror=False)[source]

Calculate all the localization factors for a given source and network of detectors, given the complex snr, sensitivity, bandwidth, mean frequency, location of the detectors. Definition of terms given in equations (B.5) and (B.6) of “Localization of transient gravitational wave sources: beyond triangulation”, Class. Quantum Grav. 35 (2018) 105002.

Parameters

method: string

localization method to use

mirror: boolean

indicating whether we are considering the mirror location

Returns

b_i: np.array

the localization factor B_i (equation (B.5))

c_ij: np.array

the localization matrix C_ij (equation (B.6))

c_i: np.array

the localization factor c_i (equation (B.17))

c: float

the localization factor c (equation (B.17))

localize(method, mirror=False, p=0.9)[source]

Calculate localization of source at given probability with chosen method

Parameters

method: str

localization method to use

mirror: boolean

indicating whether we are considering the mirror location

p: float probability region to calculate (default 0.9)

localize_all(p=0.9, methods=None)[source]

Calculate localization with given set of methods.

Parameters

p: float

probability region to calculate

methods: list

A list of localization methods to use from ‘time’, ‘coh’, ‘left’, ‘right’, ‘marg’. Default is to calculate all.

marg_loc(mirror=False, p=0.9)[source]

Calculate the marginalized localization. Method described in “Localization of transient gravitational wave sources: beyond triangulation”, Class. Quantum Grav. 35 (2018).

Parameters

mirror: boolean

indicating whether we are considering the mirror location

p: float

probability region to calculate (default=0.9)

projected_snr(method, mirror=False, dt_i=None)[source]

Calculate the projected SNR for a given method at either original or mirror sky location

Parameters

method: str

localization method to use

mirror: boolean

indicating whether we are considering the mirror location

dt_i: time shift to be applied in each detector

Returns

np.array

the complex snr for each detector

classmethod random_values(d_max=1000, mass=1.4, t_gps=1000000000)[source]

Generate an event with random distance, orientation at given time and mass

Parameters

d_max: float

maximum distance

mass: float

component mass (assumed equal mass)

t_gps: float

GPS time of event

simple_pe.localization.event.f(r, keys, localization, p)[source]

Function used in marginalization of localizations TODO: Fix documentation, and check if this function is used

Parameters

r: float

radius

keys: list

the different localization methods to consider

localization: loc.Localization

object storing localization information

p: float

probability

Returns

float

f

simple_pe.localization.event.snr_projection(f_sig, method)[source]

Function to calculate the SNR projection matrix p for a given set of detector responses, f_sig

Parameters

f_sig: array

a Nx2 array of detector responses [F+, Fx] x sigma

method: str

the way we project (one of “time”, “coh”, “left”, “right”)

Returns

p: array

SNR projection matrix

simple_pe.localization.loc module

class simple_pe.localization.loc.Localization(method, event, mirror=False, p=0.9, d_max=1000, area=0)[source]

Bases: object

class to hold the details and results of localization based on a given method

approx_like(d_max=1000)[source]

Calculate the approximate likelihood, based on equations A.18 and A.29 from “Localization of transient gravitational wave sources: beyond triangulation”, Class. Quantum Grav. 35 (2018) 105002.

Parameters

d_max: float

maximum distance, used for normalization

calc_area()[source]

Calculate the localization area using equations (12) and (15) of “Source localization with an advanced gravitational wave detector network”, DOI 10.1088/0264-9381/28/10/10502

calculate_dt_0(dt_i)[source]

Calculate the overall time offset at a point offset from the source point by dt_i that maximizes the overall SNR

Parameters

dt_i: np.array

The time offsets at which to evaluate the SNR

Returns

dt_0: float

The overall time offset that maximizes the SNR

calculate_m()[source]

Calculate the localization matrix as given in “Localization of transient gravitational wave sources: beyond triangulation”, Class. Quantum Grav. 35 (2018) 105002 this is a generalization of the timing based localization

calculate_max_snr()[source]

Calculate the maximum SNR nearby the given point. For this, calculate the localization factors and the SNR projection and see whether there is a higher network SNR at a nearby point. For timing, the peak will be at the initial point, but for coherent or left/right circular polarizations, the peak can be offset.

calculate_snr(dt_i)[source]

Calculate the SNR at a point offset from the source point dt_i is an array of time offsets for the detectors Note: we maximize over the overall time offset dt_0. If the time offsets are too large to trust the leading order approximation, then return the original SNR

Parameters

dt_i: np.array

The time offsets at which to evaluate the SNR

Returns

z: np.array

The individual SNRs at the offset time, consistent with a signal (either coherent or left/right circular)

snr: float

the network SNR consistent with a signal

generate_loc_grid(npts=10, scale=1.0)[source]

Generate a grid of points with extent governed by the localization

Parameters

npts: int

number of points in each dimension of the grid

scale: float

factor by which to scale grid, relative to default given by the localization eigenvectors

Returns

grid_points: np.array

containing a grid of points (theta and phi)

generate_samples(npts=100000, sky_weight=True)[source]

Generate a set of samples based on Gaussian distribution in localization eigendirections

Parameters

npts: int

number of points to generate

sky_weight: bool

weight points to be uniform on the sky, rather than uniform over the localization ellipse

Returns

samples: phi, theta samples

make_ellipse(npts=101, scale=1.0)[source]

Generate points that lie on the localization ellipse

Parameters

npts: int

number of points

scale: float

factor by which to scale the ellipse. Default (scaling of 1) is to use the probability stored in the localization.

Returns

points: np.array

a set of points marking the boundary of the localization ellipse the points are projected onto the sky (so in theta/phi coordinates)

sky_project()[source]

Project localization matrix to zero out components in direction of source. This is implementing equations 10 and 11 from “Source localization with an advanced gravitational wave detector network”, DOI 10.1088/0264-9381/28/10/10502

simple_pe.localization.loc.evec_sigma(m)[source]

Calculate the eigenvalues and vectors of the localization matrix M. sigma is defined as the reciprocal of the square-root eigenvalue. Definitions from “Triangulation of gravitational wave sources with a network of detectors”, New J. Phys. 11 123006.

Parameters

m: np.array

square matrix for which we calculate the eigen-vectors and sigmas

Returns

evec: np.array

An array of eigenvectors of M, giving principal directions for localization

sigma: np.array

localization accuracy along eigen-directions.

simple_pe.localization.loc.project_to_sky(x, y, event_xyz, gmst, evec, ellipse=False, sky_weight=False)[source]

Project a set of points onto the sky.

Parameters

x: np.array

x coordinates of points (relative to sky location)

y: np.array

y coordinate of points (relative to sky location)

event_xyz: np.array

xyz location of event

gmst: float

gmst of event

evec: np.array

localization eigenvectors

ellipse: bool

is this an ellipse

sky_weight: bool

re-weight to uniform on sky

Returns

phi: np.array

phi coordinate of points

theta: np.array

theta coordinate of points

simple_pe.localization.sky module

simple_pe.localization.sky.CartesianToSpherical(x, system='equatorial')[source]

Convert 3-tuple Cartesian sky position x to SkyPosition object in the given system

simple_pe.localization.sky.CircularGrid(resolution, radius)[source]

Generates a SkyPositionTable of circular grids around the North Pole with given angular resolution and a maximal radius.

simple_pe.localization.sky.ConvertSkyPosition(skyPoint, system, zenith=None, gpstime=None)[source]

Convert the SkyPosition object skyPoint from it’s current system to the new system.

Valid systems are : ‘horizon’, ‘geographic’, ‘equatorial’, ‘ecliptic’, ‘galactic’

SkyPosition zenith and gpstime should be given as appropriate for ‘horizon’ and ‘geographic’ conversions respectively.

simple_pe.localization.sky.EclipticToEquatorial(input)[source]

Convert the SkyPosition object input from the inherited ‘eliptic’ system to to ‘equatorial’.

simple_pe.localization.sky.EquatorialToEcliptic(input)[source]

Convert the SkyPosition object input from the inherited ‘equatorial’ system to to ‘ecliptic’.

simple_pe.localization.sky.EquatorialToGalactic(input)[source]

Convert the SkyPosition object input from the inherited ‘equatorial’ system to to ‘galactic’.

simple_pe.localization.sky.EquatorialToGeographic(input, gpstime)[source]

Convert the SkyPosition object input from the inherited ‘equatorial’ system to to ‘geographic’.

simple_pe.localization.sky.GeographicToEquatorial(input, gpstime)[source]

Convert the SkyPosition object input from the inherited ‘geographical’ system to to ‘equatorial’.

simple_pe.localization.sky.HorizonToSystem(input, zenith)[source]

Convert the SkyPosition object input from ‘horizon’ to the inherited system using the SkyPosition zenith

simple_pe.localization.sky.ISOTimeDelayLine(ifos, ra, dec, gpstime=None, n=3)[source]

Returns the n-point SkyPositionTable describing a line of constant time delay through the given ra and dec. for the given 2-tuple ifos.

simple_pe.localization.sky.MaxTimeDelayLine(ifos, ra, dec, gpstime=None, n=3)[source]

Return the n-point SkyPositionTable describing the line perpendicular to the line of constant time delay through the given ra and dec for the 2-tuple ifos.

simple_pe.localization.sky.MaxTimeDelayLine3(ifo1, ifo2, ra, dec, gpstime=None, n=3)[source]

Alternative implementation of MaxTimeDelayLine.

simple_pe.localization.sky.SkyPatch(ifos, ra, dec, radius, gpstime, dt=0.0005, sigma=1.65, grid='circular')[source]

Returns a SkyPositionTable of circular rings emanating from a given central ra and dec. out to the maximal radius.

class simple_pe.localization.sky.SkyPosition[source]

Bases: object

get_time_delay(ifo1, ifo2, gpstime)[source]

Return the time delay for this SkyPosition between arrival at ifo1 relative to ifo2, for the given gpstime.

latitude
longitude
normalize()[source]

Normalize this SkyPosition to be within standard limits [0 <= longitude < 2pi] and [-pi < latitude <= pi]

opening_angle(other)[source]

Calcalate the opening angle between this SkyPosition and the other SkyPosition

probability
process_id
rotate(R)[source]

Apply the 3x3 rotation matrix R to this SkyPosition.

solid_angle
system
class simple_pe.localization.sky.SkyPositionTable(attrs=None)[source]

Bases: Table

ligo.lw.table.Table holding SkyPosition objects.

normalize()[source]
parseTimeDelayDegeneracy(ifos, gpstime=LIGOTimeGPS(1395517690, 0), dt=0.0005)[source]
rotate(R)[source]
tableName = 'sky_positions:table'
valid_columns = {'latitude': 'real_4', 'longitude': 'real_4', 'probability': 'real_4', 'process_id': 'ilwd:char', 'solid_angle': 'real_4', 'system': 'lstring'}
simple_pe.localization.sky.SphericalToCartesian(skyPoint)[source]

Convert SkyPosition object into Cartesian 3-vector

simple_pe.localization.sky.SystemToHorizon(input, zenith)[source]

Convert the SkyPosition object input from the inherited system to ‘horizon’ using the SkyPosition zenith

simple_pe.localization.sky.ThreeSiteSkyGrid(ifos, gpstime, dt=0.0005, tiling='square', sky='half')[source]

Generates a SkyPositionTable for a three-site all-sky grid, covering either a hemisphere (sky=’half’) or full sphere (sky=’full’), using either ‘square’ or ‘hexagonal tiling.

simple_pe.localization.sky.TwoSiteSkyGrid(ifos, gpstime, dt=0.0005, sky='half', point=None)[source]

Generates a SkyPositionTable for a two-site all-sky grid, covering either a hemisphere (sky=’half’) or full sphere (sky=’full’).

The grid can be forced to pass through the given SkyPosition point if required.

simple_pe.localization.sky.fromfile(fobj, loncol=0, latcol=1, probcol=None, angcol=None, system=None, degrees=False)[source]

Returns a SkyPositionTable as read from the file object fobj.

simple_pe.localization.sky.parse_sites(ifos)[source]

Returns a new list of interferometers containing one only per site. I.e. this removes H2 if included.

simple_pe.localization.sky.tofile(fobj, grid, delimiter=' ', degrees=False)[source]

Writes the SkyPositionTable object grid to the file object fobj

simple_pe.localization.sky_loc module

simple_pe.localization.sky_loc.chisq_loc(t_phi_theta, t_i, d_i, f_band_i)[source]

Calculate expression (A.3) from 2nd localization paper

Param:

t_phi_theta: the time of arrival and (phi, theta) for source location

Param:

t_i: array with times of arrival in each detector

Param:

d_i: matrix with locations of detectors

Param:

f_band_i: array with bandwidths in each detector

simple_pe.localization.sky_loc.localization_from_timing(ifos, arrival_times, bandwidths)[source]

Calculate RA and dec based upon time of arrival in a network of ifos

Param:

ifos: list of ifos

Param:

arrival_times: dictionary of arrival times in different ifos

Param:

bandwidths: dictionary of signal bandwidth in each ifo

Return ra:

the right ascension of the signal

Return dec:

the declination of the signal