Waveform models

Waveform models give the prediction for a GW signal emitted by a coalescing binary as a function of the parameters of the source.

gwfast uses Fourier domain waveform models. In particular, it provides a Python implementation for some selected state-of-the-art models, plus wrappers to use both the models present in the LIGO Algorithm Library, LAL, and the TEOBResumS models.

We separately release also WF4Py, a user-friendly package containing only the waveform models implemented in gwfast (avoiding some specific implementations needed for JAX derivatives). Our Python implementation is fully vectorised, thus giving the possibility to evaluate waveforms for multiple events at a time.

The WaveFormModel class

Waveforms in gwfast are built as classes, initialised as follows

class gwfast.waveforms.WaveFormModel(objType, fcutPar, is_newtonian=False, is_tidal=False, is_HigherModes=False, is_chi1chi2=True, is_Precessing=False, is_LAL=False, is_prec_ang=False, is_eccentric=False, is_holomorphic=False, apply_fcut=True)[source]

Abstract class to compute waveforms

Parameters
  • objType (str) – The kind of system the wf model is made for, can be 'BBH', 'BNS' or 'NSBH'.

  • fcutPar (float) – The cut frequency factor of the waveform. This can either be given in \(\rm Hz\), as for gwfast.waveforms.TaylorF2_RestrictedPN, or as an adimensional frequency (Mf), as for the IMR models.

  • is_newtonian (bool, optional) – Boolean specifying if the waveform is a simple Newtonian inspiral.

  • is_tidal (bool, optional) – Boolean specifying if the waveform includes tidal effects.

  • is_HigherModes (bool, optional) – Boolean specifying if the waveform includes the contribution of sub-dominant (higher-order) modes.

  • is_chi1chi2 (bool, optional) – Boolean specifying if, in the aligned spins only case, the individual spins are used in place of the 'chiS' and 'chiA' combinations.

  • is_Precessing (bool, optional) – Boolean specifying if the waveform includes spin-precession effects.

  • is_LAL (bool, optional) – Boolean specifying if the waveform comes from the LAL library.

  • is_prec_ang (bool, optional) – Boolean specifying if, in the precessing spin case, the angular variables of the spins are used, namely 'thetaJN', 'chi1', 'chi2', 'tilt1', 'tilt2', 'phiJL', 'phi12'.

  • is_eccentric (bool, optional) – Boolean specifying if the waveform includes orbital eccentricity.

  • is_holomorphic (bool, optional) – Boolean specifying if the waveform function is holomorphic (needed for derivatives handling).

  • apply_fcut (bool, optional) – Boolean specifying if the waveform has to be cut at the chosen maximum frequency specified by fcutPar (as in LAL) or not.

Each waveform includes four fundamental methods, to compute, given the parameter of the source(s), the phase of the signal, \(\Phi (f)\), the amplitude of the signal, \(A (f)\), the time to coalescence as a function of the frequency, \(\tau^* (f)\), and the cut frequency.

abstract WaveFormModel.Phi(f, **kwargs)[source]

Compute the phase of the GW as a function of frequency, given the events parameters.

We compute here only the GW phase, not the full phase of the signal, which also includes the reference phase and the time of coalescence.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the phase of, as in events.

Returns

GW phase for the chosen events evaluated on the frequency grid.

Return type

array

abstract WaveFormModel.Ampl(f, **kwargs)[source]

Compute the amplitude of the GW as a function of frequency, given the events parameters.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the amplitude of, as in events.

Returns

GW amplitude for the chosen events evaluated on the frequency grid.

Return type

array

WaveFormModel.tau_star(f, **kwargs)[source]

Compute the time to coalescence (in seconds) as a function of frequency (in \(\rm Hz\)), given the events parameters.

Parameters
  • f (array) – Frequency grid on which the time to coalescence will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the time to coalescence of, as in events.

Returns

time to coalescence for the chosen events evaluated on the frequency grid, in seconds.

Return type

array

Note

For gwfast.waveforms.NewtInspiral we use the \(\tau^* (f)\) expression in M. Maggiore – Gravitational Waves Vol. 1 eq. (4.21). For the other models instead we use the \(\tau^* (f)\) expression in arXiv:0907.0700 eq. (3.8b), valid up to 3.5 PN.

WaveFormModel.fcut(**kwargs)[source]

Compute the cut frequency of the waveform as a function of the events parameters, in \(\rm Hz\).

Parameters

kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the cut frequency of, as in events.

Returns

Cut frequency of the waveform for the chosen events, in \(\rm Hz\).

Return type

array

Some models also have a method denoted as hphc to compute directly the \(\tilde{h}_+\) and \(\tilde{h}_{\times}\) polarisations of the gravitational wave, see below.

Each waveform further includes a dictionary storing the order of the parameters in the Fisher matrix

WaveFormModel.ParNums

Dictionary containing the number of the rows/columns in which the parameters will appear in the Fisher matrix.

Type

dict(int)

Python waveform models

We here report and briefly describe the waveform models implemented in pure Python in gwfast. We carefully checked that our Python implementation accurately reproduces the original LAL waveforms.

Note

Only when using these models it will be possible to compute the derivatives using Automatic Differentiation through JAX.

Leading order inspiral

This model includes only the leading order term in the inspiral.

class gwfast.waveforms.NewtInspiral(**kwargs)[source]

Leading order inspiral only waveform model.

Relevant references: M. Maggiore – Gravitational Waves Vol. 1, chapter 4.

Parameters

kwargs – Optional arguments to be passed to the parent class WaveFormModel.

TaylorF2

New in version 1.0.2: Added the possibility to terminate the waveform at the ISCO frequency of a remnant Kerr BH. Added the spin-induced quadrupole due to tidal effects.

This model includes contributions to the phase up to 3.5 order in the Post Newtonian, PN, expansion, and can thus be used to describe the inspiral. The amplitude is the same as in Newtonian approximation. Our implementation can include both the tidal terms at 5 and 6 PN (see arXiv:1402.5156) and a moderate eccentricity in the orbit \(e_0\lesssim 0.1\) up to 3 PN (see arXiv:1605.00304). There is no limitation in the parameters range, but, being an inspiral-only model, gwfast.waveforms.TaylorF2_RestrictedPN is better suited for BNS systems.

class gwfast.waveforms.TaylorF2_RestrictedPN(fHigh=None, is_tidal=False, use_3p5PN_SpinHO=False, phiref_vlso=False, is_eccentric=False, fRef_ecc=None, which_ISCO='Schw', use_QuadMonTid=False, **kwargs)[source]

TaylorF2 restricted PN waveform model, with coefficients up to 3.5 PN. The amplitude is thus the same as in Newtonian approximation, and the model is valid only in the inspiral.

This model can include both the contribution of tidal effects at 5 and 6 PN and the contribution of eccentricity up to 3 PN.

Relevant references:

[1] arXiv:0907.0700

[2] arXiv:1107.1267

[3] arXiv:1402.5156

[4] arXiv:1601.05588

[5] arXiv:1605.00304

Parameters
  • fHigh (float) –

    The cut frequency factor of the waveform, in \(\rm Hz\). By default this is set to two times the Innermost Stable Circular Orbit, ISCO, frequency of a remnant Schwarzschild BH having a mass equal to the total mass of the binary, see gwfast.gwfastGlobals.f_isco. Another useful value, whose coefficient is provided in gwfast.gwfastGlobals.f_qK, is the limit of the quasi-Keplerian approximation, defined as in arXiv:2108.05861 (see also arXiv:1605.00304), which is more conservative than two times the Schwarzschild ISCO.

  • is_tidal (bool, optional) – Boolean specifying if the waveform has to include tidal effects.

  • use_3p5PN_SpinHO (bool, optional) – Boolean specifying if the waveform has to include the quadratic- and cubic-in-spin contributions at 3.5 PN, which are not included in LAL.

  • phiref_vlso (bool, optional) – Boolean specifying if the reference frequency of the waveform has to be set to the Last Stable Orbit, LSO, frequency.

  • is_eccentric (bool, optional) – Boolean specifying if the waveform has to include orbital eccentricity.

  • fRef_ecc (float, optional) – The reference frequency for the provided eccentricity, \(f_{e_{0}}\).

  • which_ISCO (str, optional) –

    String specifying if the waveform has to be cut at two times the ISCO frequency of a remnant Schwarzschild (non-rotating) BH or of a Kerr BH, as in arXiv:2108.05861 (see in particular App. C), with the fits from arXiv:1605.01938. The Schwarzschild ISCO can be selected passing 'Schw', while the Kerr ISCO passing 'Kerr'. NOTE: the Kerr option pushes the validity of the model to the limit, and is not the default option.

  • use_QuadMonTid (bool, optional) – Boolean specifying if the waveform has to include the spin-induced quadrupole due to tidal effects, with the fits from arXiv:1608.02582.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

IMRPhenomD

This is a full inspiral–merger-ringdown model tuned with NR simulations, which can be used to simulate signals coming from BBH mergers, with non–precessing spins up to \(|\chi_z|\sim 0.85\) and mass ratios up to \(q = m_1/m_2 \sim 18\).

class gwfast.waveforms.IMRPhenomD(fRef=None, **kwargs)[source]

IMRPhenomD waveform model.

Relevant references:

[1] arXiv:1508.07250

[2] arXiv:1508.07253

Parameters
  • fRef (float, optional) – Reference frequency of the waveform, in \(\rm Hz\). If not provided, the minimum of the frequency grid will be used.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

IMRPhenomD_NRTidalv2

This is a full inspiral-merger-ringdown model tuned with NR simulations, which extends the IMRPhenomD model to include tidal effects, and can thus be used to accurately describe signals coming from BNS mergers. The validity has been assessed for masses between \(1\,{\rm M}_{\odot}\) and \(3\,{\rm M}_{\odot}\), spins up to \(|\chi_z|\sim 0.6\) and tidal deformabilities up to \(\Lambda_i\sim 5000\).

class gwfast.waveforms.IMRPhenomD_NRTidalv2(fRef=None, **kwargs)[source]

IMRPhenomD_NRTidalv2 waveform model.

Relevant references:

[1] arXiv:1508.07250

[2] arXiv:1508.07253

[3] arXiv:1905.06011

Parameters
  • fRef (float, optional) – Reference frequency of the waveform, in \(\rm Hz\). If not provided, the minimum of the frequency grid will be used.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

The model includes a Planck taper filter to terminate the waveform after merger, we thus the cut frequency slightly before the end of the filter, for numerical stability.

IMRPhenomHM

This is a full inspiral–merger-ringdown model tuned with NR simulations, which takes into account not only the \((2,2)\) quadrupole of the signal, but also the sub–dominant multipoles \((l,m) = (2,1),\, (3,2),\, (3,3),\, (4,3)\), and \((4,4)\), that can be particularly relevant to better describe the signal coming from BBH systems. The calibration range is the same of the IMRPhenomD model.

class gwfast.waveforms.IMRPhenomHM(fRef=None, **kwargs)[source]

IMRPhenomHM waveform model.

Relevant references:

[1] arXiv:1508.07250

[2] arXiv:1508.07253

[3] arXiv:1708.00404

[4] arXiv:1909.10010

Parameters
  • fRef (float, optional) – Reference frequency of the waveform, in \(\rm Hz\). If not provided, the minimum of the frequency grid will be used.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

To avoid for loops and recomputing coefficients (given that in this case the amplitude cannot be computed separately from the phase), this model features a gwfast.waveforms.IMRPhenomHM.hphc method, to compute directly \(\tilde{h}_+\) and \(\tilde{h}_{\times}\)

IMRPhenomHM.hphc(f, **kwargs)[source]

Compute the plus and cross polarisations of the GW as a function of frequency, given the events parameters, avoiding for loops over the modes.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the phase of, as in events.

Returns

Plus and cross polarisations of the GW for the chosen events evaluated on the frequency grid.

Return type

tuple(array, array)

Also, in this case the gwfast.waveforms.IMRPhenomHM.Phi and gwfast.waveforms.IMRPhenomHM.Ampl methods return dictionaries containing the phases and amplitudes of the various modes, respectively. The dictionaries have keys '21', '22', '32', '33', '43' and '44'.

IMRPhenomHM.Phi(f, **kwargs)[source]

Compute the phase of the GW as a function of frequency, given the events parameters.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the phase of, as in events.

Returns

GW phases of the various modes for the chosen events evaluated on the frequency grid.

Return type

dict(array, array, array, array, array, array)

IMRPhenomHM.Ampl(f, **kwargs)[source]

Compute the amplitude of the GW as a function of frequency, given the events parameters.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the amplitude of, as in events.

Returns

GW amplitudes of the various modes for the chosen events evaluated on the frequency grid.

Return type

dict(array, array, array, array, array, array)

To combine the various modes and obtain the full amplitude and phase from these outputs, it is possible to use the function

gwfast.gwfastUtils.Add_Higher_Modes(Ampl, Phi, iota, phi=0.0)[source]

Compute the total signal from a collection of different modes.

Parameters
  • Ampl (dict(array, array, ...)) – Dictionary containing the amplitudes for each mode computed on a grid of frequencies. The keys are expected to be stings made up of \(l\) and \(m\), e.g. for \((2,2)\) –> key= '22'.

  • Phi (dict(array, array, ...)) – Dictionary containing the phases for each mode computed on a grid of frequencies.

  • iota (array or float) – The inclination angle(s) of the system(s) with respect to orbital angular momentum, \(\iota\), in \(\rm rad\).

  • phi (array or float) – The second angular direction of the spherical coordinate system.

Returns

Plus and cross polarisations of the GW for the chosen events evaluated on the frequency grid.

Return type

tuple(array, array)

IMRPhenomNSBH

This is a full inspiral–merger-ringdown model tuned with NR simulations, which can describe the signal coming from the merger of a NS and a BH, with mass ratios up to \(q = m_1/m_2 \sim 100\), also taking into account tidal effects and the impact of the possible tidal disruption of the NS.

class gwfast.waveforms.IMRPhenomNSBH(fRef=None, verbose=True, **kwargs)[source]

IMRPhenomNSBH waveform model

The inputs labelled as 1 refer to the BH (e.g. 'chi1z') and with 2 to the NS (e.g. 'Lambda2')

Relevant references:

[1] arXiv:1508.07250

[2] arXiv:1508.07253

[3] arXiv:1509.00512

[4] arXiv:1905.06011

Parameters
  • fRef (float, optional) – Reference frequency of the waveform, in \(\rm Hz\). If not provided, the minimum of the frequency grid will be used.

  • verbose (bool, optional) – Boolean specifying if the code has to print additional details during execution.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

Note

In LAL, to compute the parameter \(\xi_{\rm tide}\) in arXiv:1509.00512 eq. (8), the roots are extracted. In Python this would break the possibility to vectorise so, to circumvent the issue, we compute a grid of \(\xi_{\rm tide}\) as a function of the compactness, mass ratio and BH spin, and then use a 3D interpolator. The first time the code runs, if this interpolator is not already present, it will be computed. The base resolution of the grid is 200 pts per parameter, that we find sufficient to reproduce the LAL implementation with good precision, given the smooth behaviour of the function, but this can be raised if needed. In this case, it is necessary to change the name of the file assigned to self.path_xiTide_tab and the res input passed to the function that loads the grid.

IMRPhenomNSBH._make_xiTide_interpolator(res=200)[source]

Load the table of the parameter \(\xi_{\rm tide}\) if present or computes it if not, and builds the needed 3-D interpolator.

Parameters

res (int, optional) – Resolution of the grid in compactness, mass ratio and spin.

Wrappers to external waveform models

gwfast features wrappers to use the waveform models implemented in other libraries.

Note

When using these models the derivatives are computed using Numerical Differentiation (finite differences) through the numdifftools package.

LAL wrapper

This is a wrapper to use the waveform models present in the LIGO Algorithm Library

class gwfast.waveforms.LAL_WF(approximant, fcutPar=0.3, is_tidal=False, is_HigherModes=False, is_Precessing=False, is_eccentric=False, compute_sequence=True, fRef_ecc=None, fRef=None, **kwargs)[source]

Wrapper for using LAL waveforms.

Parameters
  • approximant (str) – Name of the waveform model to use, as reported in LAL. If an invalid name is provided, the code will list the names of all the available Fourier domain approximants.

  • fcutPar (float, optional) – The cut frequency factor of the waveform as an adimensional frequency (Mf).

  • is_tidal (bool, optional) – Boolean specifying if the waveform includes tidal effects.

  • is_HigherModes (bool, optional) – Boolean specifying if the waveform includes the contribution of sub-dominant (higher-order) modes.

  • is_Precessing (bool, optional) – Boolean specifying if the waveform includes spin-precession effects.

  • is_eccentric (bool, optional) – Boolean specifying if the waveform includes orbital eccentricity.

  • fRef_ecc (float, optional) – The reference frequency for the provided eccentricity, \(f_{e_{0}}\).

  • fRef (float, optional) – Reference frequency of the waveform, in \(\rm Hz\). If not provided, the minimum of the frequency grid will be used.

  • compute_sequence (bool, optional) – Boolean to specify which LAL function to use among SimInspiralChooseFDWaveformSequence (True) and SimInspiralChooseFDWaveform (False).

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

Note

As a defualt, we use the LAL function SimInspiralChooseFDWaveformSequence, which computes the waveform on a given frequency grid. However, this shows numerical issues with some waveform models (e.g. IMRPhenomXHM), we thus also give the possibility to use the function SimInspiralChooseFDWaveform which appears more stable. This performs the computation on a LAL defined grid, which then has to be interpolated, resulting in less accurate evaluation at low frequencies and slower execution time. This can be chosen with the boolean compute_sequence: setting it to True means that the function will perform the computation directly on the user grid, setting it to False it will let LAL choose the grid and then extrapolate.

Warning

In this case it is necessary to specify if the chosen waveform includes tidal effects, higher order modes, precessing spins or orbital eccentricity through the is_tidal, is_HigherModes, is_Precessing and is_eccentric arguments, respectively.

Given that LAL outputs directly \(\tilde{h}_+\) and \(\tilde{h}_{\times}\), this class has a gwfast.waveforms.LAL_WF.hphc method

LAL_WF.hphc(f, **kwargs)[source]

Compute the plus and cross polarisations of the GW as a function of frequency, given the events parameters.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the phase of, as in events.

Returns

Plus and cross polarisations of the GW for the chosen events evaluated on the frequency grid.

Return type

tuple(array, array)

TEOBResumS wrapper

New in version 1.0.2.

This is a wrapper to use the TEOBResumS waveform models, available here, in their Fourier domain version TEOBResumSPA

class gwfast.waveforms.TEOBResumSPA_WF(modes=[[2, 1], [2, 2], [3, 1], [3, 2], [3, 3], [4, 1], [4, 2], [4, 3], [4, 4]], fcutPar=0.3, is_tidal=False, is_Precessing=False, **kwargs)[source]

Wrapper for using TEOBResumSPA waveforms.

Relevant references:

[1] arXiv:1406.6913

[2] arXiv:1506.08457

[3] arXiv:1806.01772

[4] arXiv:1904.09550

[5] arXiv:2001.09082

[6] arXiv:2012.00027

[7] arXiv:2104.07533

Parameters
  • modes (list(list)) – List of modes to use in the waveform, each one being a list of two elements, representing \(l\) and \(m\), respectively.

  • fcutPar (float, optional) – The cut frequency factor of the waveform as an adimensional frequency (Mf).

  • is_tidal (bool, optional) – Boolean specifying if the waveform includes tidal effects.

  • is_Precessing (bool, optional) – Boolean specifying if the waveform includes spin-precession effects.

  • kwargs – Optional arguments to be passed to the parent class WaveFormModel, such as is_chi1chi2.

Warning

In this case it is necessary to specify if the waveform has to include tidal effects or precessing spins through the is_tidal and is_Precessing arguments, respectively.

Given that TEOBResumS can output directly \(\tilde{h}_+\) and \(\tilde{h}_{\times}\), this class has a gwfast.waveforms.TEOBResumSPA_WF.hphc method

TEOBResumSPA_WF.hphc(f, **kwargs)[source]

Compute the plus and cross polarisations of the GW as a function of frequency, given the events parameters.

Parameters
  • f (array) – Frequency grid on which the phase will be computed, in \(\rm Hz\).

  • kwargs (dict(array, array, ...)) – Dictionary with arrays containing the parameters of the events to compute the phase of, as in events.

Returns

Plus and cross polarisations of the GW for the chosen events evaluated on the frequency grid.

Return type

tuple(array, array)