Source

class legwork.source.Source(m_1, m_2, ecc, dist, n_proc=1, f_orb=None, a=None, position=None, polarisation=None, inclination=None, weights=None, gw_lum_tol=0.05, stat_tol=0.01, interpolate_g=True, interpolate_sc=True, sc_params={})[source]

Bases: object

Class for generic GW sources

This class is for analysing a generic set of sources that may be stationary/evolving and circular/eccentric. If the type of sources are known, then a more specific subclass may be more useful

Parameters
m_1float/array

Primary mass. Must have astropy units of mass.

m_2float/array

Secondary mass. Must have astropy units of mass.

eccfloat/array

Initial eccentricity

distfloat/array

Luminosity distance to source. Must have astropy units of distance.

n_procint

Number of processors to split eccentric evolution over if needed

f_orbfloat/array

Orbital frequency (either a or f_orb must be supplied). This takes precedence over a. Must have astropy units of frequency.

afloat/array

Semi-major axis (either a or f_orb must be supplied). Must have astropy units of length.

positionSkyCoord/array, optional

Sky position of source. Must be specified using Astropy’s astropy.coordinates.SkyCoord class.

polarisationfloat/array, optional

GW polarisation angle of the source. Must have astropy angular units.

inclinationfloat/array, optional

Inclination of the source. Must have astropy angular units.

weightsfloat/array, optional

Statistical weights associated with each sample (used for plotted), default is equal weights

gw_lum_tolfloat

Allowed error on the GW luminosity when calculating SNRs. This is used to calculate maximum harmonics needed and transition between ‘eccentric’ and ‘circular’. This variable should be updated using the function legwork.source.Source.update_gw_lum_tol() (not Source._gw_lum_tol =) to ensure the cached calculations match the current tolerance.

stat_tolfloat

Fractional change in frequency over mission length above which a binary should be considered to be stationary

interpolate_gboolean

Whether to interpolate the g(n,e) function from Peters (1964). This results in a faster runtime for large collections of sources.

interpolate_scboolean

Whether to interpolate the LISA sensitivity curve

sc_paramsdict

Parameters for interpolated sensitivity curve. Include any of instrument, custom_psd, t_obs, L, approximate_R and confusion_noise. Default values are: “LISA”, None, “auto”, “auto”, False and “auto”.

Raises
ValueError

If both f_orb and a are missing. If only part of the position, inclination, and polarization are supplied. If array-like parameters don’t have the same length.

AssertionError

If a parameter is missing units

Attributes
m_cfloat/array

Chirp mass. Set using m_1 and m_2 in legwork.utils.chirp_mass()

ecc_tolfloat

Eccentricity above which a binary is considered eccentric. Set by legwork.source.Source.find_eccentric_transition()

snrfloat/array

Signal-to-noise ratio. Set by legwork.source.Source.get_snr()

max_snr_harmonicint/array

Harmonic with the maximum snr. Set by legwork.source.Source.get_snr()

n_sourcesint

Number of sources in class

Methods Summary

create_harmonics_functions()

Create two harmonics related functions as methods for the Source class

evolve_sources(t_evol[, create_new_class])

Evolve sources forward in time for t_evol amount of time.

find_eccentric_transition()

Find the eccentricity at which we must treat binaries at eccentric.

get_h_0_n(harmonics[, which_sources])

Computes the strain for binaries for the given harmonics.

get_h_c_n(harmonics[, which_sources])

Computes the characteristic strain for binaries for the given harmonics.

get_merger_time([save_in_class, ...])

Get the merger time for each source.

get_snr([t_obs, instrument, custom_psd, L, ...])

Computes the SNR for a generic binary.

get_snr_evolving([t_obs, instrument, ...])

Computes the SNR assuming an evolving binary

get_snr_stationary([t_obs, instrument, ...])

Computes the SNR assuming a stationary binary

get_source_mask([circular, stationary, t_obs])

Produce a mask of the sources.

plot_source_variables(xstr[, ystr, ...])

Plot distributions of Source variables.

plot_sources_on_sc([snr_cutoff, fig, ax, ...])

Plot all sources in the class on the sensitivity curve

set_g(interpolate_g)

Set Source g function if user wants to interpolate g(n,e).

set_sc()

Set Source sensitivity curve function

update_gw_lum_tol(gw_lum_tol)

Update GW luminosity tolerance.

update_sc_params(sc_params)

Update sensitivity curve parameters

Methods Documentation

create_harmonics_functions()[source]

Create two harmonics related functions as methods for the Source class

The first function is stored at self.harmonics_required(ecc). This calculates the index of the highest harmonic required to calculate the SNR of a system with eccentricity ecc assuming the provided tolerance gw_lum_tol. This is equivalent to the total number of harmonics required since, when calculating SNR, harmonics in the range [1, harmonics_required(ecc)] are used. Note that the value returned by the function slightly conservative as we apply ceil to the interpolation result.

The second function is stored at self.max_strain_harmonic(ecc). This calculates the harmonic with the maximum strain for a system with eccentricity ecc.

evolve_sources(t_evol, create_new_class=False)[source]

Evolve sources forward in time for t_evol amount of time. If create_new_class is True then save the updated sources in a new Source class, otherwise, update the values in this class.

Parameters
t_evolfloat/array

Amount of time to evolve sources. Either a single value for all sources or an array of values corresponding to each source.

create_new_classbool, optional

Whether to save the evolved binaries in a new class or not. If not simply update the current class, by default False.

Returns
evolved_sourcesSource

The new class with evolved sources, only returned if create_new_class is True.

find_eccentric_transition()[source]

Find the eccentricity at which we must treat binaries at eccentric. We define this as the maximum eccentricity at which the n=2 harmonic is the total GW luminosity given the tolerance self._gw_lum_tol. Store the result in self.ecc_tol

get_h_0_n(harmonics, which_sources=None)[source]

Computes the strain for binaries for the given harmonics. Use which_sources to select a subset of the sources. Merged sources are set to have 0.0 strain.

Parameters
harmonicsint/array

Harmonic(s) at which to calculate the strain

which_sourcesboolean/array

Mask on which sources to compute values for (default is all)

Returns
h_0_nfloat/array

Dimensionless strain in the quadrupole approximation (unitless) shape of array is (number of sources, number of harmonics)

get_h_c_n(harmonics, which_sources=None)[source]

Computes the characteristic strain for binaries for the given harmonics. Use which_sources to select a subset of the sources. Merged sources are set to have 0.0 characteristic strain.

Parameters
harmonicsint/array

Harmonic(s) at which to calculate the strain

which_sources `boolean/array`

Mask on which sources to compute values for (default is all)

Returns
h_c_nfloat/array

Dimensionless characteristic strain in the quadrupole approximation shape of array is (number of sources, number of harmonics)

get_merger_time(save_in_class=True, which_sources=None, exact=True)[source]

Get the merger time for each source. Set save_in_class to true to save the values as an instance variable in the class. Use which_sources to select a subset of the sources in the class. Note that if save_in_class is set to True, which_sources will be ignored.

Parameters
save_in_classbool, optional

Whether the save the result into the class as an instance variable, by default True

which_sourcesbool/array, optional

A mask for the subset of sources for which to calculate the merger time, by default all sources (None)

exactboolean, optional

Whether to calculate the merger time exactly with numerical integration or to instead use a fit

Returns
t_mergefloat/array

Merger times

get_snr(t_obs=None, instrument=None, custom_psd=None, L=None, approximate_R=None, confusion_noise=None, n_step=100, verbose=False, re_interpolate_sc=True, which_sources=None)[source]

Computes the SNR for a generic binary. Also records the harmonic with maximum SNR for each binary in self.max_snr_harmonic.

Parameters
t_obsarray

Observation duration (default: value from sc_params)

instrument{{ ‘LISA’, ‘TianQin’, ‘custom’ }}

Instrument to observe with. If ‘custom’ then custom_psd must be supplied. (default: value from sc_params)

custom_psdfunction

Custom function for computing the PSD. Must take the same arguments as legwork.psd.lisa_psd() even if it ignores some. (default: function from sc_params)

Lfloat

LISA arm length in metres

approximate_Rboolean

Whether to approximate the response function (default: no)

confusion_noisevarious

Galactic confusion noise. Acceptable inputs are either one of the values listed in legwork.psd.get_confusion_noise(), “auto” (automatically selects confusion noise based on instrument - ‘robson19’ if LISA and ‘huang20’ if TianQin), or a custom function that gives the confusion noise at each frequency for a given mission length where it would be called by running noise(f, t_obs) and return a value with units of inverse Hertz

n_stepint

Number of time steps during observation duration

verboseboolean

Whether to print additional information to user

re_interpolate_scboolean

Whether to re-interpolate the sensitivity curve if the observation time or instrument changes. If False, warning will instead be given

which_sourcesboolean/array

Mask of which sources to calculate the SNR for. If None then calculate SNR for all sources.

Returns
SNRarray

The signal-to-noise ratio

get_snr_evolving(t_obs=None, instrument=None, custom_psd=None, L=None, approximate_R=None, confusion_noise=None, re_interpolate_sc=True, n_step=100, which_sources=None, verbose=False)[source]

Computes the SNR assuming an evolving binary

Parameters
t_obsarray

Observation duration (default: follow sc_params)

instrument{{ ‘LISA’, ‘TianQin’, ‘custom’ }}

Instrument to observe with. If ‘custom’ then custom_psd must be supplied.

custom_psdfunction

Custom function for computing the PSD. Must take the same arguments as legwork.psd.lisa_psd() even if it ignores some.

Lfloat

LISA arm length in metres

approximate_Rboolean

Whether to approximate the response function (default: no)

confusion_noisevarious

Galactic confusion noise. Acceptable inputs are either one of the values listed in legwork.psd.get_confusion_noise(), “auto” (automatically selects confusion noise based on instrument - ‘robson19’ if LISA and ‘huang20’ if TianQin), or a custom function that gives the confusion noise at each frequency for a given mission length where it would be called by running noise(f, t_obs) and return a value with units of inverse Hertz

re_interpolate_scboolean

Whether to re-interpolate the sensitivity curve if the observation time or instrument changes. If False, warning will instead be given

n_stepint

Number of time steps during observation duration

which_sourcesbool/array

Mask on which sources to consider evolving and calculate (default is all sources in Class)

verboseboolean

Whether to print additional information to user

Returns
SNRarray

The signal-to-noise ratio

get_snr_stationary(t_obs=None, instrument=None, custom_psd=None, L=None, approximate_R=None, confusion_noise=None, re_interpolate_sc=True, which_sources=None, verbose=False)[source]

Computes the SNR assuming a stationary binary

Parameters
t_obsarray

Observation duration (default: follow sc_params)

instrument{{ ‘LISA’, ‘TianQin’, ‘custom’ }}

Instrument to observe with. If ‘custom’ then custom_psd must be supplied.

custom_psdfunction

Custom function for computing the PSD. Must take the same arguments as legwork.psd.lisa_psd() even if it ignores some.

Lfloat

LISA arm length in metres

approximate_Rboolean

Whether to approximate the response function (default: no)

confusion_noisevarious

Galactic confusion noise. Acceptable inputs are either one of the values listed in legwork.psd.get_confusion_noise(), “auto” (automatically selects confusion noise based on instrument - ‘robson19’ if LISA and ‘huang20’ if TianQin), or a custom function that gives the confusion noise at each frequency for a given mission length where it would be called by running noise(f, t_obs) and return a value with units of inverse Hertz

re_interpolate_scboolean

Whether to re-interpolate the sensitivity curve if the observation time or instrument changes. If False, warning will instead be given

which_sourcesbool/array

Mask on which sources to consider stationary and calculate (default is all sources in Class)

verboseboolean

Whether to print additional information to user

Returns
SNRarray

The signal-to-noise ratio

get_source_mask(circular=None, stationary=None, t_obs=None)[source]

Produce a mask of the sources.

Create a mask based on whether binaries are circular or eccentric and stationary or evolving. Tolerance levels are defined in the class.

Parameters
circularbool

None means either, True means only circular binaries and False means only eccentric

stationarybool

None means either, True means only stationary binaries and False means only evolving

t_obsfloat

Observation time, default is value from self._sc_params

Returns
maskbool/array

Mask for the sources

plot_source_variables(xstr, ystr=None, which_sources=None, exclude_merged_sources=True, log_scale=False, **kwargs)[source]

Plot distributions of Source variables. If two variables are specified then produce a 2D distribution, otherwise a 1D distribution.

Parameters
xstr{ ‘m_1’, ‘m_2’, ‘m_c’, ‘ecc’, ‘dist’, ‘f_orb’, ‘f_GW’, ‘a’, ‘snr’ }

Which variable to plot on the x axis

ystr{ ‘m_1’, ‘m_2’, ‘m_c’, ‘ecc’, ‘dist’, ‘f_orb’, ‘f_GW’, ‘a’, snr’ }

Which variable to plot on the y axis (if None then a 1D distribution is made using xstr)

which_sourcesboolean array

Mask for which sources should be plotted (default is all sources)

exclude_merged_sourcesboolean

Whether to exclude merged sources in distributions (default is True)

log_scalebool or tuple of bools

Whether to use a log scale for the axes. For a 1D plot, only a bool can be supplied and it applies to the x-axis. For a 2D plot, a single bool is applied to both axes, a tuple is applied to the x- and y-axis respectively.

**kwargsvarious

When only xstr is provided, the kwargs are the same as legwork.visualisation.plot_1D_dist(). When both xstr and ystr are provided, the kwargs are the same as legwork.visualisation.plot_2D_dist(). Note that if xlabel or ylabel is not passed then this function automatically creates one using a default string and (if applicable) the Astropy units of the variable.

Returns
figmatplotlib Figure

The figure on which the distribution is plotted

axmatplotlib Axis

The axis on which the distribution is plotted

plot_sources_on_sc(snr_cutoff=0, fig=None, ax=None, show=True, label='Stationary sources', sc_vis_settings={}, **kwargs)[source]

Plot all sources in the class on the sensitivity curve

Parameters
snr_cutofffloat

SNR below which sources will not be plotted (default is to plot all sources)

fig: `matplotlib Figure`

A figure on which to plot the distribution. Both ax and fig must be supplied for either to be used

ax: `matplotlib Axis`

An axis on which to plot the distribution. Both ax and fig must be supplied for either to be used

showboolean

Whether to immediately show the plot

labelstr

Label to use for the plotted points

sc_vis_settingsdict

Dictionary of parameters to pass to plot_sensitivity_curve(), e.g. {“y_quantity”: “h_c”} will plot characteristic strain instead of ASD

**kwargsvarious

Keyword arguments to be passed to plotting functions

Returns
figmatplotlib Figure

The figure on which the sources are plotted

axmatplotlib Axis

The axis on which the sources are plotted

Notes

Warning

Note that this function is not yet implemented for evolving sources. Evolving sources will not be plotted and a warning will be shown instead. We are working on implementing this soon!

set_g(interpolate_g)[source]

Set Source g function if user wants to interpolate g(n,e). Otherwise just leave the function as None.

Parameters
interpolate_gboolean

Whether to interpolate the g(n,e) function from Peters (1964)

set_sc()[source]

Set Source sensitivity curve function

If user wants to interpolate then perform interpolation of LISA sensitivity curve using sc_params. Otherwise just leave the function as None.

update_gw_lum_tol(gw_lum_tol)[source]

Update GW luminosity tolerance. Use the updated value to recalculate harmonics_required function and transition to eccentric

Parameters
gw_lum_tolfloat

Allowed error on the GW luminosity when calculating SNRs

update_sc_params(sc_params)[source]

Update sensitivity curve parameters

Update the parameters used to interpolate sensitivity curve and perform interpolation again to match new params