Calculating Remnant Fits

PESummary has multiple methods for calculating the remnant properties of a binary. These include averaging multiple fits tuned to numerical relativity, and using waveform specific fitting functions. Below we go into more details about how each of these can be used with the summarypages executable and/or the convert function directly.

Average NR fits

By default, the average of the fits tuned to numerical relativity will always be calculated for you. We input the spins defined at a given reference frequency and return the final_mass_non_evolved, final_spin_non_evolved among other quantities. However, we also offer the ability to evolve the spins up to the ISCO frequency and input these to the fitting functions. The returned posterior distributions are then final_mass, final_spin among others (note the lack of non_evolved in the parameter name) and will consequently be more accurate than the corresponding non_evolved quantities. The evolved spin fits can be achieve with the following command summarypages executable:

$ summarypages --webdir ./evolved_spin --samples example.hdf5 \
               --evolve_spins --gw

or:

>>> from pesummary.gw.conversions import convert
>>> data = convert(data_table, evolve_spins='ISCO')

This then calls the following functions:

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_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_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_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

NRSurrogate fits

Alternatively, we may use an NRSurrogate remnant fit to evolve the spins and calculate the remnant properties. This can be done with the following summarypages executable:

$ summarypages --webdir ./NRSurrogate --samples example.hdf5 \
               --NRSur_fits --gw

By default, the NRSur7dq4Remnant is used. The user may choose their own NRSurrogate remnant model by passing it from the command line:

$ summarypages --webdir ./NRSurrogate --samples example.hdf5 \
               --NRSur_fits NRSur3dq8Remnant --gw

Using the convert function:

>>> data = convert(data_table, NRSur_fits='NRSur7dq4Remnant')

This then calls the following function:

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

A list of available NRSurrogate remnant models can be found here.

Note

This will require the LAL_DATA_PATH to be set correctly. This will require cloning https://git.ligo.org/lscsoft/lalsuite-extra and pointing to lalsuite-extra/data/lalsimulation:

export LAL_DATA_PATH=~/lalsuite-extra/data/lalsimulation

Waveform specific waveform fits

Finally specific waveform approximants can be used to evaluate the remnant fits. This means that the functions used when evaluating the waveform are directly used for the posterior samples. For example, we can use the SEOBNRv4PHM approximant to calculate the remnant properties with the following summarypages executable:

$ summarypages --webdir ./waveform --samples example.hdf5 \
               --waveform_fits --gw --approximant SEOBNRv4PHM

Sometimes, it can take a while to evaluate these fits for computationally expensive waveform models. Therefore we offer parallelisation to reduce wall time. This can be done with the following:

$ summarypages --webdir ./waveform --samples example.hdf5 \
               --waveform_fits --gw --approximant SEOBNRv4PHM \
               --multi_process 20

f_low is required for this conversion. If f_low cannot be extracted from the result file, the code will return a ValueError and ask for an f_low to be passed from the command line. This can be done with the following:

$ summarypages --webdir ./waveform --samples example.hdf5 \
               --waveform_fits --gw --approximant SEOBNRv4PHM \
               --multi_process 20 --f_low 3.0

Using the convert function:

>>> data = convert(data_table, waveform_fits=True, approximant="SEOBNRv4PHM")

This then calls the following function:

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_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