The Conversion module

For the case of GW specific result files, certain parameters can be derived given a certain set of samples. For instance, if you have samples for the mass of the bigger and smaller black holes, then the chirp mass, mass ratio, total mass and symmetric mass ratio can all be derived.

PESummary boasts a user friendly conversion module that handles all this for you. The user simply passes the parameters and samples either as a dictionary or a list and the conversion class will calculate all possible derived quantities.

The convert function

The convert function can be initalized as follows,

>>> from pesummary.gw.conversions import convert
>>> posterior = {"mass_1": [10, 2, 40], "mass_2": [5, 1, 20]}
>>> data = convert(posterior)
>>> print(data.keys())
dict_keys(['mass_1', 'mass_2', 'a_1', 'a_2', 'mass_ratio', 'inverted_mass_ratio', 'total_mass', 'chirp_mass', 'symmetric_mass_ratio'])
>>> print(data.parameters.added)
['a_1', 'a_2', 'mass_ratio', 'inverted_mass_ratio', 'reference_frequency', 'total_mass', 'chirp_mass', 'symmetric_mass_ratio']

Alternatively,

>>> parameters = ["mass_1", "mass_2"]
>>> samples = [[10, 5], [2, 1], [40, 20]]
>>> data = convert(parameters, samples)
>>> print(data.keys())
dict_keys(['mass_1', 'mass_2', 'a_1', 'a_2', 'mass_ratio', 'inverted_mass_ratio', 'total_mass', 'chirp_mass', 'symmetric_mass_ratio'])

The convert function takes a range of keyword arguments, some of which are required for certain conversions, and others allow the user to run often time-consuming conversions. A full list of the kwargs can be seen below,

pesummary.gw.conversions.convert(*args, restart_from_checkpoint=False, resume_file=None, **kwargs)[source]

Class to calculate all possible derived quantities

Parameters:
  • data (dict, list) – either a dictionary or samples or a list of parameters and a list of samples. See the examples below for details

  • extra_kwargs (dict, optional) – dictionary of kwargs associated with this set of posterior samples.

  • f_low (float, optional) – the low frequency cut-off to use when evolving the spins

  • f_ref (float, optional) – the reference frequency when spins are defined

  • f_final (float, optional) – the final frequency to use when integrating over frequencies

  • approximant (str, optional) – the approximant to use when evolving the spins

  • evolve_spins_forwards (float/str, optional) – the final velocity to evolve the spins up to.

  • evolve_spins_backwards (str, optional) – method to use when evolving the spins backwards to an infinite separation

  • return_kwargs (Bool, optional) – if True, return a modified dictionary of kwargs containing information about the conversion

  • NRSur_fits (float/str, optional) – the NRSurrogate model to use to calculate the remnant fits. If nothing passed, the average NR fits are used instead

  • multipole_snr (Bool, optional) – if True, the SNR in the (l, m) = [(2, 1), (3, 3), (4, 4)] multipoles are calculated from the posterior samples. samples.

  • precessing_snr (Bool, optional) – if True, the precessing SNR is calculated from the posterior samples.

  • psd (dict, optional) – dictionary containing a psd frequency series for each detector you wish to include in calculations

  • waveform_fits (Bool, optional) – if True, the approximant is used to calculate the remnant fits. Default is False which means that the average NR fits are used

  • multi_process (int, optional) – number of cores to use to parallelize the computationally expensive conversions

  • redshift_method (str, optional) – method you wish to use when calculating the redshift given luminosity distance samples. If redshift samples already exist, this method is not used. Default is ‘approx’ meaning that interpolation is used to calculate the redshift given N luminosity distance points.

  • cosmology (str, optional) – cosmology you wish to use when calculating the redshift given luminosity distance samples.

  • force_non_evolved (Bool, optional) – force non evolved remnant quantities to be calculated when evolved quantities already exist in the input. Default False

  • force_BBH_remnant_computation (Bool, optional) – force BBH remnant quantities to be calculated for systems that include tidal deformability parameters where BBH fits may not be applicable. Default False.

  • force_BH_spin_evolution (Bool, optional) – force BH spin evolution methods to be applied for systems that include tidal deformability parameters where these methods may not be applicable. Default False.

  • disable_remnant (Bool, optional) – disable all remnant quantities from being calculated. Default False.

  • add_zero_spin (Bool, optional) – if no spins are present in the posterior table, add spins with 0 value. Default False.

  • psd_default (str/pycbc.psd obj, optional) – Default PSD to use for conversions when no other PSD is provided.

  • regenerate (list, optional) – list of posterior distributions that you wish to regenerate

  • return_dict (Bool, optional) – if True, return a pesummary.utils.utils.SamplesDict object

  • resume_file (str, optional) – path to file to use for checkpointing. If not provided, checkpointing is not used. Default None

Examples

There are two ways of passing arguments to this conversion class, either a dictionary of samples or a list of parameters and a list of samples. See the examples below:

>>> samples = {"mass_1": 10, "mass_2": 5}
>>> converted_samples = convert(samples)
>>> parameters = ["mass_1", "mass_2"]
>>> samples = [10, 5]
>>> converted_samples = convert(parameters, samples)
>>> samples = {"mass_1": [10, 20], "mass_2": [5, 8]}
>>> converted_samples = convert(samples)
>>> parameters = ["mass_1", "mass_2"]
>>> samples = [[10, 5], [20, 8]]

Extra notes on specific conversions

Below are seperate pages which goes into extra detail about specific parameter conversions:

Core conversion functions

Of course, the core conversion functions can be also be used,

Cosmology

Below are conversion functions that are specific to a given cosmology,

pesummary.gw.conversions.cosmology.comoving_distance_from_z(redshift, cosmology='Planck15')[source]

Return the comoving distance given samples for the redshift

pesummary.gw.conversions.cosmology.dL_from_z(redshift, cosmology='Planck15')[source]

Return the luminosity distance given samples for the redshift

pesummary.gw.conversions.cosmology.m1_from_m1_source_z(mass_1_source, z)[source]

Return the detector-frame primary mass given samples for the source-frame primary mass and the redshift

pesummary.gw.conversions.cosmology.m1_source_from_m1_z(mass_1, z)[source]

Return the source-frame primary mass given samples for the detector-frame primary mass and the redshift

pesummary.gw.conversions.cosmology.m2_from_m2_source_z(mass_2_source, z)[source]

Return the detector-frame secondary mass given samples for the source-frame secondary mass and the redshift

pesummary.gw.conversions.cosmology.m2_source_from_m2_z(mass_2, z)[source]

Return the source-frame secondary mass given samples for the detector-frame secondary mass and the redshift

pesummary.gw.conversions.cosmology.m_total_source_from_mtotal_z(total_mass, z)[source]

Return the source-frame total mass of the binary given samples for the detector-frame total mass and redshift

pesummary.gw.conversions.cosmology.mchirp_from_mchirp_source_z(mchirp_source, z)[source]

Return the detector-frame chirp mass of the binary given samples for the source-frame chirp mass and redshift

pesummary.gw.conversions.cosmology.mchirp_source_from_mchirp_z(mchirp, z)[source]

Return the source-frame chirp mass of the binary given samples for detector-frame chirp mass and redshift

pesummary.gw.conversions.cosmology.mtotal_from_mtotal_source_z(total_mass_source, z)[source]

Return the detector-frame total mass of the binary given samples for the source-frame total mass and redshift

pesummary.gw.conversions.cosmology.z_from_dL_approx(luminosity_distance, N=100, cosmology='Planck15', **kwargs)[source]

Return the approximate redshift given samples for the luminosity distance. This technique uses interpolation to estimate the redshift

pesummary.gw.conversions.cosmology.z_from_dL_exact(luminosity_distance, cosmology='Planck15', multi_process=1)[source]

Return the redshift given samples for the luminosity distance

Evolve

Below are functions that evolve BBH spins up to a specified value,

pesummary.gw.conversions.evolve.evolve_angles_backwards(mass_1, mass_2, a_1, a_2, tilt_1, tilt_2, phi_12, f_low, f_ref, approximant, method='precession_averaged', multi_process=1, return_fits_used=False, version='v2', evolution_approximant='SpinTaylorT5', **kwargs)[source]

Evolve BBH tilt angles backwards to infinite separation

Parameters:
  • mass_1 (float/np.ndarray) – float/array of primary mass samples of the binary

  • mass_2 (float/np.ndarray) – float/array of secondary mass samples of the binary

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • f_ref (float) – reference frequency where spins are defined

  • approximant (str) – Approximant used to generate the posterior samples

  • method (str) – Method to use when evolving tilts to infinity. Possible options are ‘precession_averaged’ and ‘hybrid_orbit_averaged’. ‘precession_averaged’ computes tilt angles at infinite separation assuming that precession averaged spin evolution from Gerosa et al. is valid starting from f_ref. ‘hybrid_orbit_averaged’ combines orbit-averaged evolution and ‘precession_averaged’ evolution as in Johnson-McDaniel et al. This is more accurate but slower than the ‘precession_averaged’ method.

  • multi_process (int, optional) – number of cores to run on when evolving the spins. Default: 1

  • return_fits_used (Bool, optional) – return a dictionary of fits used. Default False

  • version (str, optional) – version of the tilts_at_infinity.hybrid_spin_evolution.calc_tilts_at_infty_hybrid_evolve function to use within the lalsimulation library. Default ‘v2’. If an old version of lalsimulation is installed where ‘v2’ is not available, fall back to ‘v1’.

  • evolution_approximant (str, optional) – the approximant to use when evolving the spins. For allowed approximants see tilts_at_infinity.hybrid_spin_evolution.calc_tilts_at_infty_hybrid_evolve Default ‘SpinTaylorT5’

  • **kwargs (dict, optional) – all kwargs passed to the tilts_at_infinity.hybrid_spin_evolution.calc_tilts_at_infty_hybrid_evolve function in the lalsimulation library

pesummary.gw.conversions.evolve.evolve_angles_forwards(mass_1, mass_2, a_1, a_2, tilt_1, tilt_2, phi_12, f_low, f_ref, approximant, final_velocity='ISCO', tolerance=0.001, dt=0.1, multi_process=1, evolution_approximant='SpinTaylorT5')[source]

Evolve the BBH spin angles forwards to a specified value using lalsimulation.SimInspiralSpinTaylorPNEvolveOrbit. By default this is the Schwarzchild ISCO velocity.

Parameters:
  • mass_1 (float/np.ndarray) – float/array of primary mass samples of the binary

  • mass_2 (float/np.ndarray) – float/array of secondary mass samples of the binary

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • f_low (float) – low frequency cutoff used in the analysis

  • f_ref (float) – reference frequency where spins are defined

  • approximant (str) – Approximant used to generate the posterior samples

  • final_velocity (str, float) – final orbital velocity for the evolution. This can either be the Schwarzschild ISCO velocity 6**-0.5 ~= 0.408 (‘ISCO’) or a fraction of the speed of light

  • tolerance (float) – Only evolve spins if at least one spin’s magnitude is greater than tolerance

  • dt (float) – steps in time for the integration, in terms of the mass of the binary

  • multi_process (int, optional) – number of cores to run on when evolving the spins. Default: 1

  • evolution_approximant (str) – name of the approximant you wish to use to evolve the spins. Default is SpinTaylorT5. Other choices are SpinTaylorT1 or SpinTaylorT4

pesummary.gw.conversions.evolve.evolve_spins(*args, evolve_limit='ISCO', **kwargs)[source]

Evolve spins to a given limit.

Parameters:
  • *args (tuple) – all arguments passed to either evolve_angles_forwards or evolve_angles_backwards

  • evolve_limit (str/float, optional) – limit to evolve frequencies. If evolve_limit==’infinite_separation’ or evolve_limit==0, evolve spins to infinite separation. if evolve_limit==’ISCO’, evolve spins to ISCO frequency. If any other float, evolve spins to that frequency.

  • **kwargs (dict, optional) – all kwargs passed to either evolve_angles_forwards or evolve_angles_backwards

Mass

Below are conversion functions specific to the binaries component masses,

pesummary.gw.conversions.mass.component_mass_product(mass1, mass2)[source]

Return the product of mass1 and mass2

pesummary.gw.conversions.mass.component_masses_from_mchirp_q(mchirp, q)[source]

Return the primary and secondary masses given samples for the chirp mass and mass ratio

pesummary.gw.conversions.mass.component_masses_from_mtotal_q(mtotal, q)[source]

Return the primary and secondary masses given samples for the total mass and mass ratio

pesummary.gw.conversions.mass.eta_from_m1_m2(mass1, mass2)[source]

Return the symmetric mass ratio given the samples for mass1 and mass2

pesummary.gw.conversions.mass.eta_from_mtotal_q(total_mass, mass_ratio)[source]

Return the symmetric mass ratio given samples for the total mass and mass ratio

pesummary.gw.conversions.mass.invq_from_m1_m2(mass1, mass2)[source]

Return the inverted mass ratio (mass1/mass2 for mass1 > mass2) given the samples for mass1 and mass2

pesummary.gw.conversions.mass.invq_from_q(mass_ratio)[source]

Return the inverted mass ratio (mass1/mass2 for mass1 > mass2) given the samples for mass ratio (mass2/mass1)

pesummary.gw.conversions.mass.m1_from_mchirp_q(mchirp, q)[source]

Return the mass of the larger component given the chirp mass and mass ratio

pesummary.gw.conversions.mass.m1_from_mtotal_q(mtotal, q)[source]

Return the mass of the larger component given the total mass and mass ratio

pesummary.gw.conversions.mass.m2_from_mchirp_q(mchirp, q)[source]

Return the mass of the smaller component given the chirp mass and mass ratio

pesummary.gw.conversions.mass.m2_from_mtotal_q(mtotal, q)[source]

Return the mass of the smaller component given the total mass and mass ratio

pesummary.gw.conversions.mass.m_total_from_m1_m2(mass1, mass2)[source]

Return the total mass given the samples for mass1 and mass2

pesummary.gw.conversions.mass.mass2_from_m1_q(mass1, q)[source]

Return the secondary mass given samples for mass1 and mass ratio

pesummary.gw.conversions.mass.mchirp_from_m1_m2(mass1, mass2)[source]

Return the chirp mass given the samples for mass1 and mass2

pesummary.gw.conversions.mass.mchirp_from_mtotal_q(total_mass, mass_ratio)[source]

Return the chirp mass given samples for total mass and mass ratio

pesummary.gw.conversions.mass.q_from_eta(symmetric_mass_ratio)[source]

Return the mass ratio given samples for symmetric mass ratio

pesummary.gw.conversions.mass.q_from_m1_m2(mass1, mass2)[source]

Return the mass ratio given the samples for mass1 and mass2

NRUtils

Below are a series of helper functions that calculate the remnant properties of the binary,

pesummary.gw.conversions.nrutils.NRSur_fit(mass_1, mass_2, a_1, a_2, tilt_1, tilt_2, phi_12, phi_jl, theta_jn, phi_ref, f_low=20.0, f_ref=20.0, model='NRSur7dq4Remnant', fits=['final_mass', 'final_spin', 'final_kick'], return_fits_used=False, approximant=None, **kwargs)[source]

Return the NR fits based on a chosen NRSurrogate model

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object. In units of solar mass

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object. In units of solar mass

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • phi_jl (float/np.ndarray) – float/array of samples for the azimuthal angle of the orbital angular momentum around the total orbital angular momentum

  • theta_jn (float/np.ndarray) – float/array of samples for the angle between the total angular momentum and the line of sight

  • phi_ref (float/np.ndarray) – float/array of samples for the reference phase used in the analysis

  • f_low (float) – the low frequency cut-off used in the analysis

  • f_ref (float/np.ndarray, optional) – the reference frequency used in the analysis

  • model (str, optional) – NRSurrogate model that you wish to use for the calculation

  • fits (list, optional) – list of fits that you wish to evaluate

  • approximant (str, optional) – The approximant that was used to generate the posterior samples

  • kwargs (dict, optional) – optional kwargs that are passed directly to the lalsimulation.nrfits.eval_fits.eval_nrfit function

pesummary.gw.conversions.nrutils.bbh_final_mass_average(*args, fits=['UIB2016', 'Healyetal'], return_fits_used=False)[source]

Return the final mass averaged across multiple fits

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • fits (list, optional) – list of fits that you wish to use

  • return_fits_used (Bool, optional) – if True, return the fits that were used to calculate the average

pesummary.gw.conversions.nrutils.bbh_final_mass_non_precessing_Healyetal(*args, **kwargs)[source]

Return the final mass of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from Healy et al. If no version specified, the default fit is Healy and Lousto: arXiv:1610.09713

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • version (str, optional) – version of the fitting coefficients you wish to use. Default are the fits from 2016

  • final_spin (float/np.ndarray, optional) – float/array of precomputed final spins

pesummary.gw.conversions.nrutils.bbh_final_mass_non_precessing_Husaetal(*args)[source]

Return the final mass of the BH resulting from the merge of a BBH for an aligned-spin system using the fit from Husa et al: arXiv:1508.07250

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_final_mass_non_precessing_UIB2016(*args, **kwargs)[source]

Return the final mass of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from https://arxiv.org/abs/1611.00332

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • version (str, optional) – version of the fitting coefficients you wish to use

pesummary.gw.conversions.nrutils.bbh_final_mass_non_spinning_Panetal(*args)[source]

Return the final mass of the BH resulting from the merger of a non spinning BBH using the fit from Pan et al: Phys Rev D 84, 124052 (2011).

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

pesummary.gw.conversions.nrutils.bbh_final_spin_Panetal(*args)[source]

Return the final spin of the BH resulting from the merger of a BBH for a precessing system using the fit from Pan et al: Phys Rev D 84, 124052 (2011)

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_final_spin_average_non_precessing(*args, fits=['UIB2016', 'Healyetal', 'HBR2016'], return_fits_used=False)[source]

Return the final spin averaged across multiple non-precessing fits

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • fits (list, optional) – list of fits that you wish to use

  • return_fits_used (Bool, optional) – if True, return the fits that were used to calculate the average

pesummary.gw.conversions.nrutils.bbh_final_spin_average_precessing(*args, fits=['UIB2016', 'Healyetal', 'HBR2016'], return_fits_used=False)[source]

Return the final spin averaged across multiple fits

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • fits (list, optional) – list of fits that you wish to use

  • return_fits_used (Bool, optional) – if True, return the fits that were used to calculate the average

pesummary.gw.conversions.nrutils.bbh_final_spin_non_precessing_HBR2016(*args, **kwargs)[source]

Return the final spin of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from Hofmann, Barausse, and Rezzolla ApJL 825, L19 (2016)

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • version (str, optional) – version of the fitting coefficients you wish to use

pesummary.gw.conversions.nrutils.bbh_final_spin_non_precessing_Healyetal(*args, **kwargs)[source]

Return the final spin of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from Healy and Lousto: arXiv:1610.09713

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • version (str, optional) – version of the fitting coefficients you wish to use. Default are the fits from 2016

pesummary.gw.conversions.nrutils.bbh_final_spin_non_precessing_Husaetal(*args)[source]

Return the final spin of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from Husa et al: arXiv:1508.07250

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_final_spin_non_precessing_UIB2016(*args, **kwargs)[source]

Return the final spin of the BH resulting from the merger of a BBH for an aligned-spin system using the fit from https://arxiv.org/abs/1611.00332

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • version (str, optional) – version of the fitting coefficients you wish to use

pesummary.gw.conversions.nrutils.bbh_final_spin_non_spinning_Panetal(*args)[source]

Return the final spin of the BH resulting from the merger of a non spinning BBH using the fit from Pan et al: Phys Rev D 84, 124052 (2011).

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

pesummary.gw.conversions.nrutils.bbh_final_spin_precessing_HBR2016(*args, **kwargs)[source]

Return the final spin of the BH resulting from the merger of a BBH for a precessing system using the fit from Hofmann, Barausse, and Rezzolla ApJL 825, L19 (2016)

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • version (str, optional) – version of the fitting coefficients you wish to use

pesummary.gw.conversions.nrutils.bbh_final_spin_precessing_Healyetal(*args)[source]

Return the final spin of the BH resulting from the merger of a BBH for a precessing using the fit from Healy et al. If no version specified, the default fit is Healy and Lousto: arXiv:1610.09713

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_final_spin_precessing_UIB2016(*args)[source]

Return the final spin of the BH resulting from the merger of a BBH for a precessing using the fit from David Keitel, Xisco Jimenez Forteza, Sascha Husa, Lionel London et al. arxiv:1612.09566v1

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_final_spin_precessing_projected_Healyetal(*args)[source]

Return the final spin of the BH calculated from projected spins using the fit from Healy and Lousto: arXiv:1610.09713 augmenting it with the leading contribution from the in-plane spins

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

pesummary.gw.conversions.nrutils.bbh_final_spin_precessing_projected_UIB2016(*args)[source]

Return the final spin of the BH calculated from projected spins using the fit by David Keitel, Xisco Jimenez Forteza, Sascha Husa, Lionel London et al. arxiv:1612.09566v1 augmenting it with the leading contribution from the in-plane spins

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

pesummary.gw.conversions.nrutils.bbh_peak_luminosity_average(*args, fits=['UIB2016', 'Healyetal'], return_fits_used=False)[source]

Return the peak luminosity (in units of 10^56 ergs/s) averaged across multiple fits.

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • fits (list, optional) – list of fits that you wish to use

  • return_fits_used (Bool, optional) – if True, return the fits that were used to calculate the average

pesummary.gw.conversions.nrutils.bbh_peak_luminosity_non_precessing_Healyetal(*args)[source]

Return the peak luminosity (in units of 10^56 ergs/s) of an aligned-spin BBH using the fit from Healy and Lousto: arXiv:1610.09713

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_peak_luminosity_non_precessing_T1600018(*args)[source]

Return the peak luminosity (in units of 10^56 ergs/s) of an aligned-spin BBH using the fit by Sascha Husa, Xisco Jimenez Forteza, David Keitel [LIGO-T1500598] using 5th order in chieff

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

pesummary.gw.conversions.nrutils.bbh_peak_luminosity_non_precessing_UIB2016(*args)[source]

Return the peak luminosity (in units of 10^56 ergs/s) of an aligned-spin BBH using the fit by David Keitel, Xisco Jimenez Forteza, Sascha Husa, Lionel London et al. arxiv:1612.09566v1

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

Remnant

Below are conversion functions that calculate the remnant properties of the binary,

pesummary.gw.conversions.remnant.final_kick_of_merger(*args, method='NR', approximant='SEOBNRv4', NRfit='average', return_fits_used: False, model='NRSur7dq4Remnant')[source]

Return the final kick velocity of the remnant resulting from a BBH merger

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • method (str) – The method you wish to use to calculate the final kick of merger. Either NR, NRSurrogate or waveform

  • approximant (str) – Name of the approximant you wish to use if the chosen method is waveform or NRSurrogate

  • NRFit (str) – Name of the NR fit you wish to use if chosen method is NR

  • return_fits_used (Bool, optional) – if True, return the NR fits that were used. Only used when NRFit=’average’ or when method=’NRSurrogate’

  • model (str, optional) – The NRSurrogate model to use when evaluating the fits

pesummary.gw.conversions.remnant.final_kick_of_merger_from_NRSurrogate(*args, model='NRSur7dq4Remnant', return_fits_used=False, approximant='SEOBNRv4PHM')[source]

Return the final kick of the remnant resulting from a BBH merger using NRSurrogate fits

pesummary.gw.conversions.remnant.final_mass_of_merger(*args, method='NR', approximant='SEOBNRv4', NRfit='average', final_spin=None, return_fits_used=False, model='NRSur7dq4Remnant')[source]

Return the final mass resulting from a BBH merger

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • method (str) – The method you wish to use to calculate the final mass of merger. Either NR, NRSurrogate or waveform

  • approximant (str) – Name of the approximant you wish to use if the chosen method is waveform or NRSurrogate

  • NRFit (str) – Name of the NR fit you wish to use if chosen method is NR

  • return_fits_used (Bool, optional) – if True, return the NR fits that were used. Only used when NRFit=’average’ or when method=’NRSurrogate’

  • model (str, optional) – The NRSurrogate model to use when evaluating the fits

pesummary.gw.conversions.remnant.final_mass_of_merger_from_NR(*args, NRfit='average', final_spin=None, return_fits_used=False)[source]

Return the final mass resulting from a BBH merger using NR fits

Parameters:
  • NRfit (str) – Name of the fit you wish to use. If you wish to use a precessing fit please use the syntax ‘precessing_{}’.format(fit_name). If you wish to have an average NR fit, then pass ‘average’

  • final_spin (float/np.ndarray, optional) – precomputed final spin of the remnant.

  • return_fits_used (Bool, optional) – if True, return the fits that were used. Only used when NRfit=’average’

pesummary.gw.conversions.remnant.final_mass_of_merger_from_NRSurrogate(*args, model='NRSur7dq4Remnant', return_fits_used=False, approximant='SEOBNRv4PHM')[source]

Return the final mass resulting from a BBH merger using NRSurrogate fits

pesummary.gw.conversions.remnant.final_mass_of_merger_from_NSBH(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the final mass resulting from an NSBH merger using NSBH waveform models given samples for mass_1, mass_2, spin_1z and lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.remnant.final_mass_of_merger_from_waveform(*args, NSBH=False, **kwargs)[source]

Return the final mass resulting from a BBH/NSBH merger using a given approximant

Parameters:

NSBH (Bool, optional) – if True, use NSBH waveform fits. Default False

pesummary.gw.conversions.remnant.final_remnant_properties_from_NRSurrogate(*args, f_low=20.0, f_ref=20.0, model='NRSur7dq4Remnant', return_fits_used=False, properties=['final_mass', 'final_spin', 'final_kick'], approximant='SEOBNRv4PHM')[source]

Return the properties of the final remnant resulting from a BBH merger using NRSurrogate fits

Parameters:
  • f_low (float/np.ndarray) – The low frequency cut-off used in the analysis. Default is 20Hz

  • f_ref (float/np.ndarray) – The reference frequency used in the analysis. Default is 20Hz

  • model (str, optional) – The name of the NRSurrogate model you wish to use

  • return_fits_used (Bool, optional) – if True, return the approximant that was used.

  • properties (list, optional) – The list of properties you wish to calculate

  • approximant (str, optional) – The approximant that was used to generate the posterior samples

pesummary.gw.conversions.remnant.final_spin_of_merger(*args, method='NR', approximant='SEOBNRv4', NRfit='average', return_fits_used=False, model='NRSur7dq4Remnant')[source]

Return the final mass resulting from a BBH merger

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • a_1 (float/np.ndarray) – float/array of primary spin magnitudes

  • a_2 (float/np.ndarray) – float/array of secondary spin magnitudes

  • tilt_1 (float/np.ndarray) – float/array of primary spin tilt angle from the orbital angular momentum

  • tilt_2 (float/np.ndarray) – float/array of secondary spin tilt angle from the orbital angular momentum

  • phi_12 (float/np.ndarray) – float/array of samples for the angle between the in-plane spin components

  • method (str) – The method you wish to use to calculate the final mass of merger. Either NR, NRSurrogate or waveform

  • approximant (str) – Name of the approximant you wish to use if the chosen method is waveform or NRSurrogate

  • NRFit (str) – Name of the NR fit you wish to use if chosen method is NR

  • return_fits_used (Bool, optional) – if True, return the NR fits that were used. Only used when NRFit=’average’ or when method=’NRSurrogate’

  • model (str, optional) – The NRSurrogate model to use when evaluating the fits

pesummary.gw.conversions.remnant.final_spin_of_merger_from_NR(*args, NRfit='average', return_fits_used=False)[source]

Return the final spin resulting from a BBH merger using NR fits

Parameters:
  • NRfit (str) – Name of the fit you wish to use. If you wish to use a precessing fit please use the syntax ‘precessing_{}’.format(fit_name). If you wish to have an average NR fit, then pass ‘average’

  • return_fits_used (Bool, optional) – if True, return the fits that were used. Only used when NRfit=’average’

pesummary.gw.conversions.remnant.final_spin_of_merger_from_NRSurrogate(*args, model='NRSur7dq4Remnant', return_fits_used=False, approximant='SEOBNRv4PHM')[source]

Return the final spin resulting from a BBH merger using NRSurrogate fits

pesummary.gw.conversions.remnant.final_spin_of_merger_from_NSBH(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the final spin resulting from an NSBH merger using NSBH waveform models given samples for mass_1, mass_2, spin_1z and lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.remnant.final_spin_of_merger_from_waveform(*args, NSBH=False, **kwargs)[source]

Return the final spin resulting from a BBH/NSBH merger using a given approximant.

Parameters:

NSBH (Bool, optional) – if True, use NSBH waveform fits. Default False

pesummary.gw.conversions.remnant.peak_luminosity_of_merger(*args, NRfit='average', return_fits_used=False)[source]

Return the peak luminosity of an aligned-spin BBH using NR fits

Parameters:
  • mass_1 (float/np.ndarray) – float/array of masses for the primary object

  • mass_2 (float/np.ndarray) – float/array of masses for the secondary object

  • spin_1z (float/np.ndarray) – float/array of primary spin aligned with the orbital angular momentum

  • spin_2z (float/np.ndarray) – float/array of secondary spin aligned with the orbital angular momentum

  • NRFit (str) – Name of the NR fit you wish to use if chosen method is NR

  • return_fits_used (Bool, optional) – if True, return the NR fits that were used. Only used when NRFit=’average’

SNR

Below are functions to calculate conversions specific to the signal-to-noise ratio (SNR),

pesummary.gw.conversions.snr.multipole_snr(mass_1, mass_2, spin_1z, spin_2z, psi, iota, ra, dec, time, distance, phase, f_low=20.0, f_final=1024.0, psd={}, approx='IMRPhenomXHM', f_ref=None, return_data_used=False, multi_process=1, df=0.00390625, multipole=[21, 33, 44], psd_default='aLIGOZeroDetHighPower')[source]

Calculate the multipole snr as defined in Mills et al. arXiv:

Parameters:
  • mass_1 (float/np.ndarray) – Primary mass of the bianry

  • mass_2 (float/np.ndarray) – Secondary mass of the binary

  • spin_1z (float/np.ndarray) – The primary spin aligned with the total orbital angular momentum

  • spin_2z (float/np.ndarray) – The secondary spin aligned with the total orbital angular momentum

  • psi (float/np.ndarray) – The polarization angle of the binary

  • iota (float/np.ndarray) – The angle between the total orbital angular momentum and the line of sight

  • ra (float/np.ndarray) – The right ascension of the source

  • dec (float/np.ndarray) – The declination of the source

  • time (float/np.ndarray) – The merger time of the binary

  • distance (float/np.ndarray) – The luminosity distance of the source

  • phase (float/np.ndarray) – The phase of the source

  • f_low (float, optional) – Low frequency to use for integration. Default 20Hz

  • f_final (float, optional) – Final frequency to use for integration. Default 1024Hz

  • psd (dict, optional) – Dictionary of pycbc.types.frequencyseries.FrequencySeries objects, one for each detector. Default is to use the aLIGOZeroDetHighPower PSD

  • approx (str, optional) – The aligned spin higher order multipole approximant you wish to use. Default IMRPhenomXHM

  • f_ref (float, optional) – Reference frequency where the spins are defined. Default is f_low

  • return_data_used (Bool, optional) – if True, return a dictionary containing information about what data was used. Default False

  • multi_process (int, optional) – The number of cpus to use when computing the precessing_snr. Default 1

  • df (float, optional) – Frequency spacing between frequency samples. Default 1./256

  • multipole (list, optional) – List of multipoles to calculate the SNR for. Default [21, 33, 44]

pesummary.gw.conversions.snr.network_matched_filter_snr(IFO_matched_filter_snrs, IFO_optimal_snrs)[source]

Return the network matched filter SNR for a given detector network. Code adapted from Christopher Berry’s python notebook

Parameters:
  • IFO_matched_filter_snrs (list) – list of matched filter SNRs for each IFO in the network

  • IFO_optimal_snrs (list) – list of optimal SNRs

pesummary.gw.conversions.snr.network_snr(snrs)[source]

Return the network SNR for N IFOs

Parameters:

snrs (list) – list of numpy.array objects containing the snrs samples for a particular IFO

pesummary.gw.conversions.snr.precessing_snr(mass_1, mass_2, beta, psi_J, a_1, a_2, tilt_1, tilt_2, phi_12, theta_jn, ra, dec, time, phi_jl, distance, phase, f_low=20.0, psd={}, spin_1z=None, spin_2z=None, chi_eff=None, approx='IMRPhenomPv2', f_final=1024.0, f_ref=None, return_data_used=False, multi_process=1, duration=None, df=0.00390625, psd_default='aLIGOZeroDetHighPower', debug=True)[source]

Calculate the precessing snr as defined in Fairhurst et al. arXiv:1908.05707

Parameters:
  • mass_1 (np.ndarray) – Primary mass of the bianry

  • mass_2 (np.ndarray) – Secondary mass of the binary

  • beta (np.ndarray) – The angle between the total angular momentum and the total orbital angular momentum

  • psi_J (np.ndarray) – The polarization angle as defined with respect to the total angular momentum

  • a_1 (np.ndarray) – The spin magnitude on the larger object

  • a_2 (np.ndarray) – The spin magnitude on the secondary object

  • tilt_1 (np.ndarray) – The angle between the total orbital angular momentum and the primary spin

  • tilt_2 (np.ndarray) – The angle between the total orbital angular momentum and the secondary spin

  • phi_12 (np.ndarray) – The angle between the primary spin and the secondary spin

  • theta_jn (np.ndarray) – The angle between the total orbital angular momentum and the line of sight

  • ra (np.ndarray) – The right ascension of the source

  • dec (np.ndarray) – The declination of the source

  • time (np.ndarray) – The merger time of the binary

  • phi_jl (np.ndarray) – the precession phase

  • f_low (float, optional) – The low frequency cut-off to use for integration. Default is 20Hz

  • psd (dict, optional) – Dictionary of pycbc.types.frequencyseries.FrequencySeries objects, one for each detector. Default is to use the aLIGOZeroDetHighPower PSD

  • spin_1z (np.ndarray, optional) – The primary spin aligned with the total orbital angular momentum

  • spin_2z (np.ndarray, optional) – The secondary spin aligned with the total orbital angular momentum

  • chi_eff (np.ndarray, optional) – Effective spin of the binary

  • approx (str, optional) – The approximant you wish to use. Default IMRPhenomPv2

  • f_final (float, optional) – Final frequency to use for integration. Default 1024Hz

  • f_ref (float, optional) – Reference frequency where the spins are defined. Default is f_low

  • return_data_used (Bool, optional) – if True, return a dictionary containing information about what data was used. Default False

  • multi_process (int, optional) – The number of cpus to use when computing the precessing_snr. Default 1

  • duration (float, optional) – maximum IMR duration to use to estimate delta_f when PSD is not provided.

  • debug (Bool, optional) – if True, return posteriors for b_bar and the overlap between the 0th and 1st harmonics. These are useful for debugging.

Spins

Below are functions that deal with the binaries spins,

pesummary.gw.conversions.spins.chi_eff(mass1, mass2, spin1z, spin2z)[source]

Return chi_eff given samples for mass1, mass2, spin1z, spin2z

pesummary.gw.conversions.spins.chi_p(mass1, mass2, spin1x, spin1y, spin2x, spin2y)[source]

Return chi_p given samples for mass1, mass2, spin1x, spin1y, spin2x, spin2y

pesummary.gw.conversions.spins.chi_p_2spin(mass1, mass2, spin1x, spin1y, spin2x, spin2y)[source]

Return the magnitude of a modified chi_p which takes into account precessing spin information from both compact objects given samples for mass1, mass2, spin1x, spin1y, spin2x, spin2y. See Eq.9 of arXiv:2012.02209

pesummary.gw.conversions.spins.chi_p_from_tilts(mass1, mass2, a_1, tilt_1, a_2, tilt_2)[source]

Return chi_p given samples for mass1, mass2, a_1, tilt_2, a_2, tilt_2

pesummary.gw.conversions.spins.component_spins(theta_jn, phi_jl, tilt_1, tilt_2, phi_12, a_1, a_2, mass_1, mass_2, f_ref, phase)[source]

Return the component spins given samples for theta_jn, phi_jl, tilt_1, tilt_2, phi_12, a_1, a_2, mass_1, mass_2, f_ref, phase

pesummary.gw.conversions.spins.opening_angle(mass_1, mass_2, phi_jl, tilt_1, tilt_2, phi_12, a_1, a_2, f_ref, phase)[source]

Return the opening angle of the system given samples for mass_1, mass_2, cartesian spins and a reference frequency

pesummary.gw.conversions.spins.phi1_from_spins(spin_1x, spin_1y)[source]

Return phi_1 given samples for spin_1x and spin_1y

pesummary.gw.conversions.spins.phi2_from_spins(spin_2x, spin_2y)[source]

Return phi_2 given samples for spin_2x and spin_2y

pesummary.gw.conversions.spins.phi_12_from_phi1_phi2(phi1, phi2)[source]

Return the difference in azimuthal angle between S1 and S2 given samples for phi1 and phi2

pesummary.gw.conversions.spins.spin_angles(mass_1, mass_2, inc, spin1x, spin1y, spin1z, spin2x, spin2y, spin2z, f_ref, phase)[source]

Return the spin angles given samples for mass_1, mass_2, inc, spin1x, spin1y, spin1z, spin2x, spin2y, spin2z, f_ref, phase

pesummary.gw.conversions.spins.spin_angles_from_azimuthal_and_polar_angles(a_1, a_2, a_1_azimuthal, a_1_polar, a_2_azimuthal, a_2_polar)[source]

Return the spin angles given samples for a_1, a_2, a_1_azimuthal, a_1_polar, a_2_azimuthal, a_2_polar

pesummary.gw.conversions.spins.tilt_angles_and_phi_12_from_spin_vectors_and_L(a_1, a_2, Ln)[source]

Return the tilt angles and phi_12 given samples for the spin vectors and the orbital angular momentum

Parameters:
  • a_1 (np.ndarray) – Spin vector for the larger object

  • a_2 (np.ndarray) – Spin vector for the smaller object

  • Ln (np.ndarray) – Orbital angular momentum of the binary

pesummary.gw.conversions.spins.viewing_angle_from_inclination(inclination)[source]

Return the viewing angle of the binary given samples for the source inclination angle. For a precessing system, the source inclination angle is theta_jn: the angle between the total angular momentum J and the line of sight N. For a non-precessing system, the source inclination angle is iota: the angle between the total orbital angular momentum L and the line of sight N

Tidal

Below are functions that are specific to objects that experience tidal forces,

pesummary.gw.conversions.tidal.NSBH_baryonic_torus_mass(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the baryonic torus mass given samples for mass_1, mass_2, spin_1z, lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.tidal.NSBH_merger_type(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH', percentages=True, percentage_round=2, _ringdown=None, _disruption=None, _torus=None)[source]

Determine the merger type based on the disruption frequency, ringdown frequency and torus mass. If percentages = True, a dictionary is returned showing the number of samples which fall in each category. If percentages = False, an array of length mass_1 is returned with elements indicating the merger type for each sample

pesummary.gw.conversions.tidal.NSBH_ringdown_frequency(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the ringdown frequency given samples for mass_1, mass_2, spin_1z, lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.tidal.NSBH_tidal_disruption_frequency(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the tidal disruption frequency given samples for mass_1, mass_2, spin_1z, lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.tidal.NS_baryonic_mass(compactness, NS_mass)[source]

Calculate the neutron star baryonic mass from its compactness and gravitational mass in solar masses

pesummary.gw.conversions.tidal.NS_compactness_from_NSBH(mass_1, mass_2, spin_1z, lambda_2, approximant='IMRPhenomNSBH')[source]

Calculate the neutron star compactness given samples for mass_1, mass_2, spin_1z, lambda_2. mass_1 and mass_2 should be in solar mass units.

pesummary.gw.conversions.tidal.NS_compactness_from_lambda(lambda_x)[source]

Calculate neutron star compactness from its tidal deformability

pesummary.gw.conversions.tidal.delta_lambda_from_lambda1_lambda2(lambda1, lambda2, mass1, mass2)[source]

Return the second dominant tidal term given samples for lambda1 and lambda 2

pesummary.gw.conversions.tidal.lambda1_from_4_parameter_piecewise_polytrope_equation_of_state(log_pressure, gamma_1, gamma_2, gamma_3, mass_1, mass_2)[source]

Convert 4 parameter piecewise polytrope EOS parameters to the tidal deformability parameters lambda_1

pesummary.gw.conversions.tidal.lambda1_from_lambda_tilde(lambda_tilde, mass1, mass2)[source]

Return the individual tidal parameter given samples for lambda_tilde

pesummary.gw.conversions.tidal.lambda1_lambda2_from_4_parameter_piecewise_polytrope_equation_of_state(log_pressure, gamma_1, gamma_2, gamma_3, mass_1, mass_2, multi_process=1)[source]

Convert 4 parameter piecewise polytrope EOS parameters to the tidal deformability parameters lambda_1, lambda_2

pesummary.gw.conversions.tidal.lambda1_lambda2_from_spectral_decomposition(spectral_decomposition_gamma_0, spectral_decomposition_gamma_1, spectral_decomposition_gamma_2, spectral_decomposition_gamma_3, mass_1, mass_2, multi_process=1)[source]

Convert spectral decomposition parameters to the tidal deformability parameters lambda_1, lambda_2

pesummary.gw.conversions.tidal.lambda1_minus_lambda2(lambda1, lambda2)[source]

Return the primary objects tidal deformability minus the secondary objests tidal deformability

pesummary.gw.conversions.tidal.lambda1_plus_lambda2(lambda1, lambda2)[source]

Return the sum of the primary objects tidal deformability and the secondary objects tidal deformability

pesummary.gw.conversions.tidal.lambda2_from_4_parameter_piecewise_polytrope_equation_of_state(log_pressure, gamma_1, gamma_2, gamma_3, mass_1, mass_2)[source]

Convert 4 parameter piecewise polytrope EOS parameters to the tidal deformability parameters lambda_2

pesummary.gw.conversions.tidal.lambda2_from_lambda1(lambda1, mass1, mass2)[source]

Return the individual tidal parameter given samples for lambda1

pesummary.gw.conversions.tidal.lambda_tilde_from_lambda1_lambda2(lambda1, lambda2, mass1, mass2)[source]

Return the dominant tidal term given samples for lambda1 and lambda2

pesummary.gw.conversions.tidal.wrapper_for_lambda1_lambda2_from_spectral_decomposition(args)[source]

Wrapper function to calculate the tidal deformability parameters from the spectral decomposition parameters for a pool of workers

pesummary.gw.conversions.tidal.wrapper_for_lambda1_lambda2_polytrope_EOS(args)[source]

Wrapper function to calculate the tidal deformability parameters from the 4_parameter_piecewise_polytrope_equation_of_state parameters for a pool of workers

Time

Below are functions that deal with the event time,

pesummary.gw.conversions.time.time_in_each_ifo(detector, ra, dec, time_gps)[source]

Return the event time in a given detector, given samples for ra, dec, time

TGR

Below are functions that generate parameters required for testing General Relativity

pesummary.gw.conversions.tgr.generate_imrct_deviation_parameters(samples, evolve_spins_forward=True, inspiral_string='inspiral', postinspiral_string='postinspiral', approximant=None, f_low=None, return_samples_used=False, **kwargs)[source]

Generate deviation parameter pdfs for the IMR Consistency Test

Parameters:
  • samples (MultiAnalysisSamplesDict) – Dictionary containing inspiral and postinspiral samples

  • evolve_spins_forward (bool) – If True, evolve spins to the ISCO frequency. Default: True.

  • inspiral_string (string) – Identifier for the inspiral samples

  • postinspiral_string (string) – Identifier for the post-inspiral samples

  • approximant (dict, optional) – The approximant used for the inspiral and postinspiral analyses. Keys of the dictionary must be the same as the inspiral_string and postinspiral_string. Default None

  • f_low (dict, optional) – The low frequency cut-off used for the inspiral and postinspiral analyses. Keys of the dictionary must be the same as the inspiral_string and postinspiral_string. Default None

  • return_samples_used (Bool, optional) – if True, return the samples which were used to generate the IMRCT deviation parameters. These samples will match the input but may include remnant samples if they were not previously present

  • kwargs (dict, optional) – Keywords to be passed to imrct_deviation_parameters_from_final_mass_final_spin

Returns:

  • imrct_deviations (ProbabilityDict2d) – 2d pdf of the IMRCT deviation parameters

  • data (dict) – Metadata

pesummary.gw.conversions.tgr.imrct_deviation_parameters_from_final_mass_final_spin(final_mass_inspiral, final_spin_inspiral, final_mass_postinspiral, final_spin_postinspiral, N_bins=101, final_mass_deviation_lim=1, final_spin_deviation_lim=1, multi_process=1, use_kde=False, kde=<class 'scipy.stats._kde.gaussian_kde'>, kde_kwargs={}, interp_method=<class 'scipy.interpolate._interpolate.interp2d'>, interp_kwargs={'bounds_error': False, 'fill_value': 0.0}, vectorize=False)[source]

Compute the IMR Consistency Test deviation parameters. Code borrows from the implementation in lalsuite: https://git.ligo.org/lscsoft/lalsuite/-/blob/master/lalinference/python/lalinference/imrtgr/imrtgrutils.py and https://git.ligo.org/lscsoft/lalsuite/-/blob/master/lalinference/bin/imrtgr_imr_consistency_test.py

Parameters:
  • final_mass_inspiral (np.ndarray) – values of final mass calculated from the inspiral part

  • final_spin_inspiral (np.ndarray) – values of final spin calculated from the inspiral part

  • final_mass_postinspiral (np.ndarray) – values of final mass calculated from the post-inspiral part

  • final_spin_postinspiral (np.ndarray) – values of final spin calculated from the post-inspiral part

  • final_mass_deviation_lim (float) – Maximum absolute value of the final mass deviation parameter. Default 1.

  • final_spin_deviation_lim (float) – Maximum absolute value of the final spin deviation parameter. Default 1.

  • N_bins (int, optional) – Number of equally spaced bins between [-final_mass_deviation_lim, final_mass_deviation_lim] and [-final_spin_deviation_lim, final_spin_deviation_lim]. Default is 101.

  • multi_process (int) – Number of parallel processes. Default is 1.

  • use_kde (bool) – If True, uses kde instead of interpolation. Default is False.

  • kde (method) – KDE method to use. Default is scipy.stats.gaussian_kde

  • kde_kwargs (dict) – Arguments to be passed to the KDE method

  • interp_method (method) – Interpolation method to use. Default is scipy.interpolate.interp2d

  • interp_kwargs (dict, optional) – Arguments to be passed to the interpolation method Default is dict(fill_value=0.0, bounds_error=False)

  • vectorize (bool) – if True, use vectorized imrct_deviation_parameters_integrand function. This is quicker but does consume more memory. Default: False

Returns:

imrct_deviations – Contains the 2d pdf of the IMRCT deviation parameters

Return type:

ProbabilityDict2d

pesummary.gw.conversions.tgr.imrct_deviation_parameters_integrand(*args, vectorize=False, **kwargs)[source]

Compute the final mass and final spin deviation parameters

Parameters:
  • *args (tuple) – all args passed to either _imrct_deviation_parameters_integrand_vectorized or _imrct_deviation_parameters_integrand_series

  • vectorize (bool) – Vectorize the calculation. Note that vectorize=True uses up a lot of memory

  • kwargs (dict, optional) – all kwargs passed to either _imrct_deviation_parameters_integrand_vectorized or _imrct_deviation_parameters_integrand_series