Detailed User Interface Information

Sensitivity

The main module in hasasia is the Sensitivity module which contains all of the needed methods for building a senstivity curve.

class hasasia.sensitivity.GWBSensitivityCurve(spectra, orf='hd', autocorr=False)[source]

Class to produce a sensitivity curve for a gravitational wave background, using Hellings-Downs spatial correlations.

Parameters:
orf : str, optional {‘hd’, ‘st’, ‘dipole’, ‘monopole’}

Overlap reduction function to be used in the sensitivity curve. Maybe be Hellings-Downs, Scalar-Tensor, Dipole or Monopole.

H_0

Hubble Constant. Assumed to be in units of km /(s Mpc) unless supplied as an astropy.quantity.

Omega_gw

Energy Density sensitivity

SNR(Sh)[source]

Calculate the signal-to-noise ratio of a given signal strain power spectral density, Sh. Must match frequency range and df of self.

S_eff

Strain power sensitivity.

S_effIJ

Strain power sensitivity.

h_c

Characteristic strain sensitivity

class hasasia.sensitivity.DeterSensitivityCurve(spectra, pulsar_term=True, include_corr=False, A_GWB=None)[source]
Parameters:
include_corr : bool

Whether to include cross correlations from the GWB as an additional noise source in full PTA correlation matrix. (Has little to no effect and adds a lot of computation time.)

A_GWB : float

Value of GWB amplitude for use in cross correlations.

H_0

Hubble Constant. Assumed to be in units of km /(s Mpc) unless supplied as an astropy.quantity.

NcalInvIJ

Inverse Noise Weighted Transmission Function that includes cross-correlation noise from GWB.

Omega_gw

Energy Density sensitivity

SNR(h0)[source]

Calculate the signal-to-noise ratio of a source given the strain amplitude. This is based on Equation (79) from Hazboun, et al., 2019 [1].

\[\rho(\hat{n})=h_0\sqrt{\frac{T_{\rm obs}}{S_{\rm eff}(f_0 ,\hat{k})}}\]
S_eff

Strain power sensitivity.

h_c

Characteristic strain sensitivity

class hasasia.sensitivity.Pulsar(toas, toaerrs, phi=None, theta=None, designmatrix=None, N=None, pdist=<Quantity 1. kpc>)[source]

Class to encode information about individual pulsars.

Parameters:
toas : array

Pulsar Times of Arrival [sec].

toaerrs : array

Pulsar TOA errors [sec].

phi : float

Ecliptic longitude of pulsar [rad].

theta : float

Ecliptic latitude of pulsar [rad].

designmatrix : array

Design matrix for pulsar’s timing model. N_TOA x N_param.

N : array

Covariance matrix for the pulsar. N_TOA x N_TOA. Made from toaerrs if not provided.

pdist : astropy.quantity, float

Earth-pulsar distance. Default units is kpc.

class hasasia.sensitivity.Spectrum(psr, nf=400, fmin=None, fmax=2e-07, freqs=None, tm_fit=True, **Tf_kwargs)[source]

Class to encode the spectral information for a single pulsar.

Parameters:
psr : hasasia.Pulsar

A hasasia.Pulsar instance.

nf : int, optional

Number of frequencies over which to build the various spectral densities.

fmin : float, optional [Hz]

Minimum frequency over which to build the various spectral densities. Defaults to the timespan/5 of the pulsar.

fmax : float, optional [Hz]

Minimum frequency over which to build the various spectral densities.

freqs : array, optional [Hz]

Optionally supply an array of frequencies over which to build the various spectral densities.

NcalInv

Inverse Noise Weighted Transmission Function.

Omega_gw

Energy Density sensitivity.

\[\Omega_{gw}=\frac{2\pi^2}{3\;H_0^2}f^3\;S_I\]
P_n

Inverse Noise Weighted Transmission Function.

S_I

Strain power sensitivity for this pulsar. Equation (74) in [1]

\[S_I=\frac{1}{\mathcal{N}^{-1}\;\mathcal{R}}\]
S_R

Residual power sensitivity for this pulsar.

\[S_R=\frac{1}{\mathcal{N}^{-1}}\]
add_noise_power(noise)[source]

Add any spectrum of noise. Must match length of frequency array.

Note: All noise information is furnished by the covariance matrix in the hasasia.Pulsar object, this is simply useful for bookkeeping and plots.

add_red_noise_power(A=None, gamma=None, vals=False)[source]

Add power law red noise to the prefit residual power spectral density. As \(P=A^2(f/fyr)^{-\gamma}\).

Note: All noise information is furnished by the covariance matrix in the hasasia.Pulsar object, this is simply useful for bookkeeping and plots.

Parameters:
A : float

Amplitude of red noise.

gamma : float

Spectral index of red noise powerlaw.

vals : bool

Whether to return the psd values as an array. Otherwise just added to self.psd_prefit.

add_white_noise_power(sigma=None, dt=None, vals=False)[source]

Add power law red noise to the prefit residual power spectral density.

Note: All noise information is furnished by the covariance matrix in the hasasia.Pulsar object, this is simply useful for bookkeeping and plots.

Parameters:
sigma : float

TOA error.

dt : float

Time between observing epochs in [seconds].

vals : bool

Whether to return the psd values as an array. Otherwise just added to self.psd_prefit.

h_c

Characteristic strain sensitivity for this pulsar.

\[h_c=\sqrt{f\;S_I}\]
psd_postfit

Postfit Residual Power Spectral Density

psd_prefit

Prefit Residual Power Spectral Density

hasasia.sensitivity.R_matrix(designmatrix, N)[source]

Create R matrix as defined in Ellis et al (2013) and Demorest et al (2012)

Parameters:
designmatrix : array

Design matrix of timing model.

N : array

TOA uncertainties [s]

Returns:
R matrix
hasasia.sensitivity.G_matrix(designmatrix)[source]

Create G matrix as defined in van Haasteren 2013

Parameters:
designmatrix : array

Design matrix for a pulsar timing model.

Returns:
G matrix
hasasia.sensitivity.get_Tf(designmatrix, toas, N=None, nf=200, fmin=None, fmax=2e-07, freqs=None, exact_astro_freqs=False, from_G=True, twofreqs=False)[source]

Calculate the transmission function for a given pulsar design matrix, TOAs and TOA errors.

Parameters:
designmatrix : array

Design matrix for a pulsar timing model, N_TOA x N_param.

toas : array

Times-of-arrival for pulsar, N_TOA long.

N : array

Covariance matrix for pulsar time-of-arrivals, N_TOA x N_TOA. Often just a diagonal matrix of inverse TOA errors squared.

nf : int, optional

Number of frequencies at which to calculate transmission function.

fmin : float, optional

Minimum frequency at which to calculate transmission function.

fmax : float, optional

Maximum frequency at which to calculate transmission function.

exact_astro_freqs : bool, optional

Whether to use exact 1/year and 2/year frequency values in calculation.

from_G : bool, optional

Whether to use G matrix for transmission function calculate. If False R-matrix is used.

hasasia.sensitivity.get_NcalInv(psr, nf=200, fmin=None, fmax=2e-07, freqs=None, exact_yr_freqs=False, full_matrix=False, return_Gtilde_Ncal=False, tm_fit=True)[source]

Calculate the inverse-noise-wieghted transmission function for a given pulsar. This calculates \(\mathcal{N}^{-1}(f,f') , \; \mathcal{N}^{-1}(f)\) in [1], see Equations (19-20).

Parameters:
psr : array

Pulsar object.

nf : int, optional

Number of frequencies at which to calculate transmission function.

fmin : float, optional

Minimum frequency at which to calculate transmission function.

fmax : float, optional

Maximum frequency at which to calculate transmission function.

exact_yr_freqs : bool, optional

Whether to use exact 1/year and 2/year frequency values in calculation.

Returns:
inverse-noise-weighted transmission function
hasasia.sensitivity.resid_response(freqs)[source]

Returns the timing residual response function for a pulsar across as set of frequencies. See Equation (53) in [1].

\[\mathcal{R}(f)=\frac{1}{12\pi^2\;f^2}\]
hasasia.sensitivity.HellingsDownsCoeff(phi, theta, autocorr=False)[source]

Calculate Hellings and Downs coefficients from two lists of sky positions.

Parameters:
phi : array, list

Pulsar axial coordinate.

theta : array, list

Pulsar azimuthal coordinate.

Returns:
ThetaIJ : array

An Npair-long array of angles between pairs of pulsars.

chiIJ : array

An Npair-long array of Hellings and Downs relation coefficients.

pairs : array

A 2xNpair array of pair indices corresponding to input order of sky coordinates.

chiRSS : float

Root-sum-squared value of all Hellings-Downs coefficients.

hasasia.sensitivity.get_Tspan(psrs)[source]

Returns the total timespan from a list or arry of Pulsar objects, psrs.

hasasia.sensitivity.get_TspanIJ(psr1, psr2)[source]

Returns the overlapping timespan of two Pulsar objects, psr1/psr2.

hasasia.sensitivity.corr_from_psd(freqs, psd, toas, fast=True)[source]

Calculates the correlation matrix over a set of TOAs for a given power spectral density.

Parameters:
freqs : array

Array of freqs over which the psd is given.

psd : array

Power spectral density to use in calculation of correlation matrix.

toas : array

Pulsar times-of-arrival to use in correlation matrix.

fast : bool, optional

Fast mode uses a matix inner product, while the slower mode uses the numpy.trapz function which is slower, but more accurate.

Returns:
corr : array

A 2-dimensional array which represents the correlation matrix for the given set of TOAs.

hasasia.sensitivity.quantize_fast(toas, toaerrs, flags=None, dt=0.1)[source]

Function to quantize and average TOAs by observation epoch. Used especially for NANOGrav multiband data.

Pulled from [3].

Parameters:
times : array

TOAs for a pulsar.

flags : array, optional

Flags for TOAs.

dt : float

Coarse graining time [days].

hasasia.sensitivity.red_noise_powerlaw(A, freqs, gamma=None, alpha=None)[source]

Add power law red noise to the prefit residual power spectral density. As \(P=A^2(f/fyr)^{-\gamma}\)

Parameters:
A : float

Amplitude of red noise.

gamma : float

Spectral index of red noise powerlaw.

freqs : array

Frequencies at which to calculate the red noise power law.

hasasia.sensitivity.Agwb_from_Seff_plaw(freqs, Tspan, SNR, S_eff, gamma=4.333333333333333, alpha=None)[source]

Must supply numpy.ndarrays.

hasasia.sensitivity.PI_hc(freqs, Tspan, SNR, S_eff, N=200)[source]

Power law-integrated characteristic strain.

hasasia.sensitivity.nanograv_11yr_stoch()[source]

Returns a GWBSensitivityCurve object built using with the NANOGrav 11-year data set.

hasasia.sensitivity.nanograv_11yr_deter()[source]

Returns a DeterSensitivityCurve object built using with the NANOGrav 11-year data set.

Skymap

The hasasia.Skymap module contains the methods necessary for making sensitivity maps for pulsar timing array.

class hasasia.skymap.SkySensitivity(spectra, theta_gw, phi_gw, pulsar_term=False, pol='gr', iota=None, psi=None)[source]

Class to make sky maps for deterministic PTA gravitational wave signals. Calculated in terms of \(\hat{n}=-\hat{k}\).

Parameters:
theta_gw : list, array

Gravitational wave source sky location colatitude at which to calculate sky map.

phi_gw : list, array

Gravitational wave source sky location longitude at which to calculate sky map.

pulsar_term : bool, str, optional [True, False, ‘explicit’]

Flag for including the pulsar term in sky map sensitivity. True includes an idealized factor of two from Equation (36) of [1]. The ‘explicit’ flag turns on an explicit calculation of pulsar terms using pulsar distances. (This option takes considerably more computational resources.)

pol: str, optional [‘gr’,’scalar-trans’,’scalar-long’,’vector-long’]

Polarization of gravitational waves to be used in pulsar antenna patterns. Only one can be used at a time.

A_gwb(h_div_A, SNR=1)[source]

Method to return a skymap of amplitudes needed to see signal the specified signal, given the specified SNR.

Parameters:
h_div_A : array

An array that represents the frequency dependence of a signal that has been divided by the amplitude. Must cover the same frequencies covered by the Skymap.freqs.

SNR : float, optional

Desired signal-to-noise ratio.

Returns:
An array representing the skymap of amplitudes needed to see the
given signal.
H_0

Hubble Constant. Assumed to be in units of km /(s Mpc) unless supplied as an astropy.quantity.

NcalInvIJ

Inverse Noise Weighted Transmission Function that includes cross-correlation noise from GWB.

Omega_gw

Energy Density sensitivity

SNR(h0, iota=None, psi=None)[source]

Calculate the signal-to-noise ratio of a source given the strain amplitude. This is based on Equation (79) from Hazboun, et al., 2019 [1].

\[\rho(\hat{n})=h_0\sqrt{\frac{T_{\rm obs}}{S_{\rm eff}(f_0 ,\hat{k})}}\]
S_SkyI

Per Pulsar Strain power sensitivity.

S_SkyI_full(iota, psi)[source]

Per Pulsar Strain power sensitivity.

S_eff

Strain power sensitivity.

S_eff_full(iota, psi)[source]

Strain power sensitivity.

S_eff_mean

Strain power sensitivity.

h_c

Characteristic strain sensitivity

h_thresh(SNR=1, iota=None, psi=None)[source]

Method to return a skymap of amplitudes needed to see a circular binary, given the specified SNR. This is based on Equation (80) from Hazboun, et al., 2019 [1].

\[h_0=\rho(\hat{n})\sqrt{\frac{S_{\rm eff}(f_0 ,\hat{k})}{T_{\rm obs}}}\]
Parameters:
SNR : float, optional

Desired signal-to-noise ratio.

Returns:
An array representing the skymap of amplitudes needed to see the
given signal with the SNR threshold specified.
sky_response_full(iota, psi)[source]

Calculate the signal-to-noise ratio of a source given the strain amplitude. This is based on Equation (79) from Hazboun, et al., 2019 [1].

\[\rho(\hat{n})=h_0\sqrt{\frac{T_{\rm obs}}{S_{\rm eff}(f_0 ,\hat{k})}}\]
hasasia.skymap.h_circ(M_c, D_L, f0, Tspan, f)[source]

Convenience function that returns the Fourier domain representation of a single circular super-massive binary black hole.

Parameters:
M_c : float [M_sun]

Chirp mass of a SMBHB.

D_L : float [Mpc]

Luminosity distance to a SMBHB.

f0 : float [Hz]

Frequency of the SMBHB.

Tspan : float [sec]

Timespan that the binary has been observed. Usually taken as the timespan of the data set.

f : array [Hz]

Array of frequencies over which to model the Fourier domain signal.

Returns:
hcw : array [strain]

Array of strain values across the frequency range provided for a circular SMBHB.

Utils

hasasia.utils.create_design_matrix(toas, RADEC=False, PROPER=False, PX=False)[source]

Return designmatrix for quadratic spindown model + optional astrometric parameters

Parameters:
toas : array

TOA measurements [s]

RADEC : bool, optional

Includes RA/DEC fitting.

PROPER : bool, optional

Includes proper motion fitting.

PX : bool, optional

Includes parallax fitting.

Returns:
M : array

Design matrix for quadratic spin down + optional astrometry fit.

hasasia.utils.fap(F, Npsrs=None)[source]

False alarm probability of the F-statistic Use None for the Fe statistic and the number of pulsars for the Fp stat.

Simulations

hasasia.sim.sim_pta(timespan, cad, sigma, phi, theta, Npsrs=None, A_rn=None, alpha=None, freqs=None, uneven=False, A_gwb=None, kwastro={'PROPER': True, 'PX': True, 'RADEC': True})[source]

Make a simulated pulsar timing array. Using the available parameters, the function returns a list of pulsar objects encoding them.

Parameters:
timespan : float, array, list

Timespan of observations in [years].

cad : float, array, list

Cadence of observations [number/yr].

sigma : float, array, list

TOA RMS Error [sec]. Single float, Npsrs long array, or Npsrs x NTOA array excepted.

phi : array, list

Pulsar’s longitude in ecliptic coordinates.

theta : array, list

Pulsar’s colatitude in ecliptic coordinates.

Npsrs : int, optional

Number of pulsars. Only needed if all pulsars have the same noise characteristics.

A_rn : float, optional

Red noise amplitude to be injected for each pulsar.

alpha : float, optional

Red noise spectral index to be injected for each pulsar.

freqs : array, optional

Array of frequencies at which to calculate the red noise. Same array used for all pulsars.

uneven : bool, optional

Option to have the toas be unevenly sampled.

Returns:
psrs : list

List of hasasia.Pulsar() objects.