LALSimulation  5.4.0.1-b72065a

Detailed Description

Routines to produce IMRPhenomX-family of phenomenological inspiral-merger-ringdown waveforms.

These are frequency-domain models for compact binaries at comparable and extreme mass ratios, tuned to numerical-relativity simulations.

The previous models are for binary black holes. There are also versions for binary neutron stars, using the extension of the binary black hole models given in https://arxiv.org/abs/1905.06011 ; DCC link: https://dcc.ligo.org/P1900148

Review status:
IMRPhenomXAS & IMRPhenomXHM reviewed by Maria Haney, Patricia Schmidt, Roberto Cotesta, Anuradha Samajdar, Jonathan Thompson, N.V. Krishnendu. IMRPhenomXP & IMRPhenomXPHM reviewed by Maria Haney, Jonathan Thompson, Marta Colleoni, David Keitel. Combined review wiki: https://git.ligo.org/waveforms/reviews/imrphenomx/-/wikis/home
Review status:
IMRPhenomXAS_NRTidalv2 & IMRPhenomXP_NRTidalv2 plus SpinTaylor precession option also applicable to IMRPhenomXP & IMRPhenomXPHM reviewed by Maria Haney, Sarp Akcay, N.V. Krishnendu, Shubhanshu Tiwari. Review wiki: https://git.ligo.org/waveforms/reviews/imrphenomxp_nrtidalv2/-/wikis/home

Routines for IMRPhenomXAS

C code for IMRPhenomXAS phenomenological waveform model.

Author
Geraint Pratten

This is an aligned-spin frequency domain model for the 22 mode. See G.Pratten et al arXiv:2001.11412 for details. Any studies that use this waveform model should include a reference to this paper.

Note
The model was calibrated to mass-ratios from 1 to 1000. The calibration points will be given in forthcoming papers.
Attention
The model is usable outside this parameter range, and in tests to date gives sensible physical results, but conclusive statements on the physical fidelity of the model for these parameters await comparisons against further numerical-relativity simulations. For more information, see the review wiki under https://git.ligo.org/waveforms/reviews/imrphenomx/wikis/home

Waveform flags: InsPhaseVersion: Determines the inspiral phase model.

  • 104 : Canonical TaylorF2 at 3.5PN including cubic-in-spin and quadratic-in-spin corrections. Uses 4 pseudo-PN coefficients. (RECOMMENDED).
  • 105 : Canonical TaylorF2 at 3.5PN including cubic-in-spin and quadratic-in-spin corrections. Uses 5 pseudo-PN coefficients.
  • 114 : Extended TaylorF2. Includes cubic-in-spin, quadratic-in-spin corrections, 4PN and 4.5PN orbital corrections. Uses 4 pseudo-PN coefficients.
  • 115 : Extended TaylorF2. Includes cubic-in-spin, quadratic-in-spin corrections, 4PN and 4.5PN orbital corrections. Uses 5 pseudo-PN coefficients.

IntPhaseVersion: Determines the intermediate phase model.

  • 104 : 4th order polynomial ansatz.
  • 105 : 5th order polynomial ansatz. (RECOMMENDED).

RDPhaseVersion: Determines the merger-ringdown phase model.

  • 105 : Deformed Lorentzian using 5 coefficients. (RECOMMENDED).

InsAmpVersion : Determines inspiral amplitude model.

  • 103 : Canonical PN re-expanded TaylorF2 amplitude with pseudo-PN corrections. (RECOMMENDED).

IntAmpVersion : Determines intermediate amplitude model.

  • 104 : Based on a 4th order polynomial ansatz. Less accurate but stable extrapolation. (RECOMMENDED).
  • 105 : Based on 5th order polynomial ansatz. More accurate in calibration domain, more unstable extrapolation.

RDAmpVersion : Determines the merger-ringdown amplitude model.

  • 103 : Deformed Lorentzian with 3 free coefficients. Uses 1 calibrated collocation point and 2 calibrated phenomenological coefficients. (RECOMMENDED).
int XLALSimIMRPhenomXASGenerateFD (COMPLEX16FrequencySeries **htilde22, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 distance, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 phi0, REAL8 fRef_In, LALDict *lalParams)
 Driver routine to calculate an IMRPhenomX aligned-spin, inspiral-merger-ringdown phenomenological waveform model in the frequency domain. More...
 
int XLALSimIMRPhenomXASFrequencySequence (COMPLEX16FrequencySeries **htilde22, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 distance, REAL8 phi0, REAL8 fRef_In, LALDict *lalParams)
 Compute waveform in LAL format at specified frequencies for the IMRPhenomX model. More...
 
double XLALSimIMRPhenomXASDuration (const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 chi1L, const REAL8 chi2L, const REAL8 f_start)
 Compute the duration of IMRPhenomXAS using the approximate SPA relation \(t_f \sim \frac{1}{2 \pi} \frac{d \varphi}{d f} \). More...
 

Routines for IMRPhenomXP

C code for IMRPhenomXP phenomenological waveform model.

Author
Geraint Pratten, Cecilio García Quirós

This is a precessing frequency domain model. See Pratten, García-Quirós, Colleoni et al arXiv:2004.06503 for details. Studies using this model are kindly asked to cite Pratten et al arXiv:2001.11412, García-Quirós et al arXiv:2001.10914 and Pratten, García-Quirós, Colleoni et al arXiv:2004.06503.

Note
The underlying aligned-spin model was calibrated for mass-ratios 1 to 1000.
Attention
The model can be called outside this parameter space and in all tests to date gives sensible physical results but conclusive statements on the physical fidelity of the model requires further detailed comparisons against numerical-relativity simulations. For more information, see the review wiki under https://git.ligo.org/waveforms/reviews/imrphenomx/wikis/home

IMRPhenomXP/HM is based on a modular framework. User can specify flags to control how Euler angles are calculated, the final spin parameterization and the conventions used in reconstructing the waveform in the LAL frame. A detailed discussion can be found in arXiv:2004.06503. The various flags are detailed below.

Precession flags: PhenomXPrecVersion:

  • 101 : NNLO PN Euler angles and a 2PN non-spinning approximation to L
  • 102 : NNLO PN Euler angles and a 3PN spinning approximation to L
  • 103 : NNLO PN Euler angles and a 4PN spinning approximation to L
  • 104 : NNLO PN Euler angles and a 4PN spinning approximation to L augmeneted with leading PN order at all order in spin terms. See N. Siemonsen et al, PRD, 97, 124046, (2018), arXiv:1712.08603
  • 220 : MSA Euler angles and a 3PN spinning approximation to L, see K. Chatziioannou et al, PRD, 95, 104004, (2017), arXiv:1703.03967. Defaults to NNLO version 102 if MSA fails to initialize.
  • 221 : MSA Euler angles and a 3PN spinning approximation to L.
  • 222 : MSA Euler angles as implemented in LALSimInspiralFDPrecAngles.
  • 223 : MSA Euler angles as implemented in LALSimInspiralFDPrecAngles. Defaults to NNLO version 102 if MSA fails to initialize. [Default].
  • 224 : As version 220 but using the \(\phi_{z,0}\) and \(\zeta_{z,0}\) prescription from 223.
  • 310 : Numerical integration of SpinTaylor equations with constant angles in merger-ringdown
  • 311 : Numerical integration of SpinTaylor equations with constant angles in merger-ringdown, without tidal and non-black hole spin-induced quadrupole terms in the SpinTaylor equations when used in IMRPhenomXP_NRTidalv2, for comparison
  • 320 : Numerical integration of SpinTaylor equations, analytical continuation in merger-ringdown
  • 321 : Numerical integration of SpinTaylor equations, analytical continuation in merger-ringdown, without tidal and non-black hole spin-induced quadrupole terms in the SpinTaylor equations when used in IMRPhenomXP_NRTidalv2, for comparison

PhenomXPExpansionOrder:

  • -1, 0, 1, 2, 3, 4, 5. Controls the expansion order of the leading-order MSA terms for both \(\zeta\) and \(\phi_z\). [Default is 5].

PhenomXPFinalSpinMod:

  • 0 : Modify final spin based on \(\chi_p\). [Recommended default for NNLO angles].
  • 1 : Modify final spin using \(\chi_{1x}\). This is pathological. Do not use. Implemented to compare to PhenomPv3 pre-bug fix.
  • 2 : Modify final spin using norm of total in-plane spin vector.
  • 3 : Modify final spin using precession-averaged couplings from MSA analysis. Only works with MSA Euler angles (versions 220, 221, 222, 223 and 224). If MSA fails to initialize or called with NNLO angles, default to version 0. [Default]
  • 4: Modify final spin estimating the total in-plane spin from the PN spin-evolution equations.

PhenomXPConvention (App. C and Table IV of arXiv:2004.06503):

  • 0 : Conventions defined as following https://dcc.ligo.org/LIGO-T1500602
  • 1 : Convention defined following App. C, see Table II of arXiv:2004.06503 for specific details. [Default]
  • 5 : Conventions as used in PhenomPv3/HM
  • 6 : Conventions defined following App. C, see Table II of arXiv:2004.06503 for specific details.
  • 7 : Conventions defined following App. C, see Table II of arXiv:2004.06503 for specific details.
int XLALSimIMRPhenomXPGenerateFD (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, REAL8 f_min, REAL8 f_max, const REAL8 deltaF, REAL8 fRef_In, LALDict *lalParams)
 
int XLALSimIMRPhenomXPFrequencySequence (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Compute waveform in LAL format at specified frequencies for the IMRPhenomXP model. More...
 
int XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame (REAL8 *chi1L, REAL8 *chi2L, REAL8 *chi_p, REAL8 *thetaJN, REAL8 *alpha0, REAL8 *phi_aligned, REAL8 *zeta_polarization, REAL8 m1_SI, REAL8 m2_SI, REAL8 f_ref, REAL8 phiRef, REAL8 incl, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, LALDict *lalParams)
 
int XLALSimIMRPhenomXPMSAAngles (REAL8Sequence **alpha_of_f, REAL8Sequence **gamma_of_f, REAL8Sequence **cosbeta_of_f, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 fRef_In, INT4 mprime, LALDict *lalParams)
 
int XLALSimIMRPhenomXPPNAngles (REAL8Sequence **alpha_of_f, REAL8Sequence **gamma_of_f, REAL8Sequence **cosbeta_of_f, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 fRef_In, INT4 mprime, LALDict *lalParams)
 

Routines for PNR angles

C code for routines to implement PNR angles in IMRPhenomXP/IMRPhenomXPHM.

Author
Eleanor Hamilton, Sebastian Khan, Jonathan E. Thompson

This is a C-code impementation of the PNR angle model outlined in arXiv:2107.08876.

Available flags: PhenomXPNRUseTunedAngles:

  • 0 : Disable PNR angles. Use the default NNLO/MSA precession angles in IMRPhenomXP/HM without tuning.
  • 1 : Enable PNR angles.
int XLALSimIMRPhenomX_PNR_GeneratePNRAngles (REAL8Sequence **alphaPNR, REAL8Sequence **betaPNR, REAL8Sequence **gammaPNR, REAL8Sequence **freqs, REAL8 *alphaPNR_ref, REAL8 *gammaPNR_ref, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, LALDict *lalParams)
 This is an external wrapper to generate the (2,2) PNR angles, following the prescription outlined in arXiv:2107.08876, given the standard inputs given to generate FD waveforms. More...
 
int XLALSimIMRPhenomX_PNR_GeneratePNRAnglesHM (REAL8Sequence **alphaPNR, REAL8Sequence **betaPNR, REAL8Sequence **gammaPNR, REAL8Sequence **freqs, REAL8 *alphaPNR_ref, REAL8 *gammaPNR_ref, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, INT4 ell, INT4 emmprime, LALDict *lalParams)
 This is an external wrapper to generate the (l,m) PNR angles, following the prescriptions outlined in arXiv:2107.08876 and arXiv:##### FIXME: add reference, given the standard inputs given to generate FD waveforms. More...
 
int IMRPhenomX_PNR_GeneratePNRAngles (REAL8Sequence *alphaPNR, REAL8Sequence *betaPNR, REAL8Sequence *gammaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
 Generate the tuned precession angles outlined in arXiv:2107.08876. More...
 
int IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies (REAL8Sequence *alphaPNR, REAL8Sequence *betaPNR, REAL8Sequence *gammaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF_SingleSpin, IMRPhenomXPrecessionStruct *pPrec_SingleSpin, IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomX_PNR_beta_parameters *betaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
 Generate the tuned precession angles outlined in arXiv:2107.08876 specifically on a uniform frequency grid. More...
 
int IMRPhenomX_PNR_GeneratePNRAlphaForAntisymmetry (REAL8Sequence *alphaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
 Internal helper function to generate PNR alpha for the antisymmetric waveform. More...
 
int IMRPhenomX_PNR_GeneratePNRAngleInterpolants (IMRPhenomX_PNR_angle_spline *angle_spline, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalparams)
 Generate the tuned precession angles outlined in arXiv:2107.08876 specifically populate the angle interpolant struct. More...
 
REAL8 XLALSimIMRPhenomX_PNR_GenerateEffectiveRingdownFreq (REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 f_min, REAL8 f_max, REAL8 fRef, UINT4 ell, UINT4 emmprime, LALDict *lalParams)
 External wrapper for IMRPhenomX_PNR_GenerateEffectiveRingdownFreq. More...
 
REAL8 XLALSimIMRPhenomX_PNR_GenerateRingdownPNRBeta (REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, LALDict *lalParams)
 External wrapper for IMRPhenomX_PNR_GenerateRingdownPNRBeta. More...
 
int XLALSimIMRPhenomX_PNR_InterpolationDeltaF (INT4 *twospin, INT4 *precversion, REAL8 *deltaF, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 f_min, REAL8 f_max, REAL8 fRef, LALDict *lalParams)
 

Routines for IMRPhenomXHM

C code for IMRPhenomXHM phenomenological waveform model.

Author
Cecilio García Quirós, Marta Colleoni, Sascha Husa

This is an aligned-spin frequency domain model with subdomimant modes: 2|2|, 2|1|, 3|3|, 3|2|, 4|4| . See C.García-Quirós et al arXiv:2001.10914 for details. Any studies that use this waveform model should include a reference to this paper.

Note
The model was calibrated to mass-ratios from 1 to 1000. The calibration points will be given in forthcoming papers.
Attention
The model is usable outside this parameter range, and in tests to date gives sensible physical results, but conclusive statements on the physical fidelity of the model for these parameters await comparisons against further numerical-relativity simulations. For more information, see the review wiki under https://git.ligo.org/waveforms/reviews/imrphenomx/wikis/home. DCC link to the paper and supplementary material: https://dcc.ligo.org/P2000011-v2

Waveform flags:

PhenomXHMReleaseVersion: this flag was introduced after the recalibration of the higher modes amplitudes in 122022. The default value is 122022, pointing to the new release. The old release can be recovered setting this option to 122019. For future releases this flag will be updated.

NOTE: The following flags are only intended for further model development and not of general interest to the user. Changes from the implemented default values (https://git.ligo.org/waveforms/reviews/imrphenomxhm-amplitude-recalibration/-/wikis/home) are generally not advised. IMRPhenomXHMInspiral(Intermediate)(Ringdown)Amp(Phase)FreqsVersion: controls the cutting frequencies between regions and the collocation points frequencies IMRPhenomXHMInspiral(Intermediate)(Ringdown)Amp(Phase)FitsVersion: controls the version of the parameter space fits for collocation points values and coefficients IMRPhenomXHMInspiral(Intermediate)(Ringdown)Amp(Phase)Version: controls the version of the ansatz. In the 122022 release the new format allows for more flexibility to choose the collocation points to be used The 122022 release does not modify the phases, so all the Phase flags point to the old value of PhaseVersion. Notice that the PhaseFreqs(Fits)Version flags cannot be modified through LALDict options.

int XLALSimIMRPhenomXHMGenerateFDOneMode (COMPLEX16FrequencySeries **htildelm, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, UINT4 ell, INT4 emm, REAL8 distance, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns the Fourier domain strain of just one mode: h_lm. More...
 
int XLALSimIMRPhenomXHMFrequencySequenceOneMode (COMPLEX16FrequencySeries **htildelm, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, UINT4 ell, INT4 emm, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Get the model evaluated in a prebuilt frequency array. More...
 
int XLALSimIMRPhenomXHMModes (SphHarmFrequencySeries **hlms, REAL8 m1_SI, REAL8 m2_SI, REAL8 S1z, REAL8 S2z, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 f_ref, REAL8 phiRef, REAL8 distance, LALDict *LALparams)
 Function to obtain a SphHarmFrequencySeries with the individual modes h_lm. More...
 
int XLALSimIMRPhenomXHM_SpheroidalPhase (REAL8Sequence **phase, const REAL8Sequence *freqs, UINT4 ell, INT4 emm, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 
int XLALSimIMRPhenomXHM_SpheroidalAmplitude (REAL8Sequence **amplitude, const REAL8Sequence *freqs, UINT4 ell, INT4 emm, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 
int XLALSimIMRPhenomXHM (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns the hptilde and hctilde of the multimode waveform for positive frequencies. More...
 
int XLALSimIMRPhenomXHM2 (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns the hptilde and hctilde of the multimode waveform for positive frequencies. More...
 
int XLALSimIMRPhenomXHMFrequencySequence (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1z, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies specified in the REAL8Sequence *freqs (which can be unequally spaced). More...
 
int XLALSimIMRPhenomXHMAmplitude (REAL8Sequence **amplitude, const REAL8Sequence *freqs, UINT4 ell, INT4 emm, REAL8 m1_SI, REAL8 m2_SI, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns amplitude of one single mode in a custom frequency array. More...
 
int XLALSimIMRPhenomXHMPhase (REAL8Sequence **phase, const REAL8Sequence *freqs, UINT4 ell, INT4 emm, REAL8 m1_SI, REAL8 m2_SI, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns phase of one single mode in a custom frequency array. More...
 
INT4 IMRPhenomXHM_Phase_for_Initialization (double *phase, double *dphase, double Mf, INT4 ell, INT4 emm, IMRPhenomXWaveformStruct *pWF, LALDict *lalParams)
 

Routines for IMRPhenomXHM Multibanding

C code for applying Multibanding to IMRPhenomXHM_Multimode

Author
Cecilio García Quirós, Sascha Husa

This is a technique to make the evaluation of waveform models faster by evaluating the model in a coarser non-uniform grid and interpolate this to the final fine uniform grid. We apply this technique to the fourier domain model IMRPhenomXHM as described in this paper: https://arxiv.org/abs/2001.10897

Multibanding flags: ThresholdMband: Determines the strength of the Multibanding algorithm. The lower this value is, the slower is the evaluation but more accurate is the final waveform compared to the one without multibanding.

  • 0.001: DEFAULT value
  • 0: switch off the multibanding

AmpInterpol: Determines the gsl interpolation order for the amplitude.

  • 1: linear interpolation (DEFAULT)
  • 3: cubic interpolation
int XLALSimIMRPhenomXHMMultiBandOneMode (COMPLEX16FrequencySeries **htildelm, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, UINT4 ell, INT4 emmIn, REAL8 distance, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Return htildelm, the waveform of one mode without mode-mixing. More...
 
int XLALSimIMRPhenomXHMMultiBandOneModeMixing (COMPLEX16FrequencySeries **htildelm, COMPLEX16FrequencySeries *htilde22, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, UINT4 ell, INT4 emmIn, REAL8 distance, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns htildelm the waveform of one mode that present mode-mixing. More...
 

Routines for IMRPhenomXPHM

C code for IMRPhenomXPHM phenomenological waveform model.

Author
Cecilio García Quirós, Geraint Pratten

This is a frequency domain precessing model based on the twisting-up of the aligned spin model with higher modes IMRPhenomXHM. See G.Pratten et al arXiv:2004.06503 for details. Any studies that use this waveform model should include a reference to this paper.

Note
DCC link to the paper: https://dcc.ligo.org/LIGO-P2000039. This paper will be refered in the code as the "Precessing paper".

Waveform flags:

All the flags for IMRPhenomXP apply here plus the following ones:

TwistPhenomHM: option to twist-up the AS model PhenomHM instead of PhenomXHM. It is only available for the polarizations, not for individual modes.

  • 0: (DEFAULT) twist-up PhenomXHM
  • 1: twist-up PhenomHM

UseModes: Determine how the polarizations hp, hc are computed.

  • 0: (DEFAULT) Compute the non-precessing modes once and do the twistin up as in eq. 3.5-3.7 in the Precessing paper.
  • 1: Compute first the individual precessing modes in the inertial J-frame and sum them to get the polarizations.

ModesL0Frame: Determine in which frame the individual precessing modes are returned.

  • 0: inertial J-frame (DEFAULT).
  • 1: inertial L0-frame (only working near the aligned spin limit).

PrecModes: Determine which indiviual modes are returned, the non-precessing or the precessing.

  • 0: (DEFAULT) Return the precessing individual mode in the J-frame.
  • 1: Return the non-precessing individual mode before the twisting-up with the modified final spin.

Multibanding flags:

PrecThresholdMband: Determines the accuracy and speed of the Multibanding algorithm for the Euler angles. The higher the threshold the faster is the algorithm but also less accurate.

  • 0.001 (DEFAULT)
  • 0: Switch off the multibanding.

MBandPrecVersion: Determines the algorithm to build the non-uniform frequency grid for the Euler angles.

  • 0: (DEFAULT) Not use multibanding. Activated to 1 when PrecThresholdMband is non-zero.
  • 1: Use the same grid that for the non-precessing modes. Activated when PrecThresholdMband is non-zero.
int XLALSimIMRPhenomXPHM (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 fRef_In, LALDict *lalParams)
 Returns hptilde and hctilde of the multimode precessing waveform for positive frequencies in an equally spaced grid. More...
 
int XLALSimIMRPhenomXPHMFromModes (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 f_min, REAL8 f_max, REAL8 deltaF, REAL8 fRef_In, LALDict *lalParams)
 Returns hptilde and hctilde of the multimode precessing waveform for positive frequencies in an equally space grid. More...
 
int XLALSimIMRPhenomXPHMFrequencySequence (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
 Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies specified in the REAL8Sequence *freqs (which can be unequally spaced). More...
 
int XLALSimIMRPhenomXPHMOneMode (COMPLEX16FrequencySeries **hlmpos, COMPLEX16FrequencySeries **hlmneg, const UINT4 l, const INT4 m, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, const REAL8 deltaF, const REAL8 f_min, const REAL8 f_max, const REAL8 fRef_In, LALDict *lalParams)
 Function to compute one hlm precessing mode in an uniform frequency grid. More...
 
int XLALSimIMRPhenomXPHMFrequencySequenceOneMode (COMPLEX16FrequencySeries **hlmpos, COMPLEX16FrequencySeries **hlmneg, const REAL8Sequence *freqs, const UINT4 l, const INT4 m, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, const REAL8 fRef_In, LALDict *lalParams)
 Function to compute one hlm precessing mode on a custom frequency grid. More...
 

Function Documentation

◆ XLALSimIMRPhenomXASGenerateFD()

int XLALSimIMRPhenomXASGenerateFD ( COMPLEX16FrequencySeries **  htilde22,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  distance,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  phi0,
REAL8  fRef_In,
LALDict *  lalParams 
)

Driver routine to calculate an IMRPhenomX aligned-spin, inspiral-merger-ringdown phenomenological waveform model in the frequency domain.

arXiv:2001.11412, https://arxiv.org/abs/2001.11412

All input parameters should be in SI units. Angles should be in radians.

XLALSimIMRPhenomXASGenerateFD() returns the strain of the 2-2 mode as a complex frequency series with equal spacing deltaF and contains zeros from zero frequency to the starting frequency and zeros beyond the cutoff frequency in the ringdown.

Parameters
[out]htilde22FD waveform
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
distanceLuminosity distance (m)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = 0.3
deltaFSampling frequency (Hz)
phi0Orbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary

Definition at line 175 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXASFrequencySequence()

int XLALSimIMRPhenomXASFrequencySequence ( COMPLEX16FrequencySeries **  htilde22,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  distance,
REAL8  phi0,
REAL8  fRef_In,
LALDict *  lalParams 
)

Compute waveform in LAL format at specified frequencies for the IMRPhenomX model.

All input parameters should be in SI units. Angles should be in radians.

XLALSimIMRPhenomXASFrequencySequence() returns the strain of the 2-2 mode as a complex frequency series with entries exactly at the frequencies specified in the sequence freqs (which can be unequally spaced). No zeros are added. Assumes positive frequencies.

Parameters
[out]htilde22FD waveform
[out]freqsFrequency series [Hz]
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
distanceLuminosity distance (m)
phi0Phase at reference frequency
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary

Definition at line 319 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXASDuration()

double XLALSimIMRPhenomXASDuration ( const REAL8  m1_SI,
const REAL8  m2_SI,
const REAL8  chi1L,
const REAL8  chi2L,
const REAL8  f_start 
)

Compute the duration of IMRPhenomXAS using the approximate SPA relation \(t_f \sim \frac{1}{2 \pi} \frac{d \varphi}{d f} \).

All input parameters should be in SI units. Angles should be in radians.

XLALSimIMRPhenomXASDuration() returns the duration in s of IMRPhenomXAS from the specified starting frequency in Hz up to the peak ringdown frequency as defined in Eq. 5.14 of https://arxiv.org/abs/2001.11412.

Parameters
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1Lz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2Lz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
f_startInitial frequency (Hz)

Definition at line 408 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXPGenerateFD()

int XLALSimIMRPhenomXPGenerateFD ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
const REAL8  distance,
const REAL8  inclination,
const REAL8  phiRef,
REAL8  f_min,
REAL8  f_max,
const REAL8  deltaF,
REAL8  fRef_In,
LALDict *  lalParams 
)
Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified.
deltaFSampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0.
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 1038 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXPFrequencySequence()

int XLALSimIMRPhenomXPFrequencySequence ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
const REAL8  distance,
const REAL8  inclination,
const REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Compute waveform in LAL format at specified frequencies for the IMRPhenomXP model.

XLALSimIMRPhenomXPGenerateFD() returns the plus and cross polarizations as a complex frequency series with equal spacing deltaF and contains zeros from zero frequency to the starting frequency and zeros beyond the cutoff frequency in the ringdown.

In contrast, XLALSimIMRPhenomXPFrequencySequence() returns a complex frequency series with entries exactly at the frequencies specified in the sequence freqs (which can be unequally spaced). No zeros are added.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
freqsinput Frequency series [Hz]
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 1255 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame()

int XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame ( REAL8 chi1L,
REAL8 chi2L,
REAL8 chi_p,
REAL8 thetaJN,
REAL8 alpha0,
REAL8 phi_aligned,
REAL8 zeta_polarization,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  f_ref,
REAL8  phiRef,
REAL8  incl,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
LALDict *  lalParams 
)
Parameters
[out]chi1LDimensionless aligned spin on companion 1
[out]chi2LDimensionless aligned spin on companion 2
[out]chi_pEffective precession parameter: Schmidt, Ohme, Hannam, PRD, 91,024043 (2015)
[out]thetaJNAngle between J0 and line of sight (z-direction)
[out]alpha0Initial value of alpha angle (azimuthal precession angle)
[out]phi_alignedInitial phase to feed the underlying aligned-spin model
[out]zeta_polarizationAngle to rotate the polarizations
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
f_refReference GW frequency (Hz)
phiRefReference phase (Hz)
inclInclination : angle between LN and the line of sight
chi1xInitial value of chi1x: dimensionless spin of BH 1 in L frame
chi1yInitial value of chi1y: dimensionless spin of BH 1 in L frame
chi1zInitial value of chi1z: dimensionless spin of BH 1 in L frame
chi2xInitial value of chi2x: dimensionless spin of BH 2 in L frame
chi2yInitial value of chi2y: dimensionless spin of BH 2 in L frame
chi2zInitial value of chi2z: dimensionless spin of BH 2 in L frame
lalParamsLAL Dictionary

Definition at line 1417 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXPMSAAngles()

int XLALSimIMRPhenomXPMSAAngles ( REAL8Sequence **  alpha_of_f,
REAL8Sequence **  gamma_of_f,
REAL8Sequence **  cosbeta_of_f,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  inclination,
REAL8  fRef_In,
INT4  mprime,
LALDict *  lalParams 
)
Parameters
[out]alpha_of_fThe azimuthal angle of L around J
[out]gamma_of_fThe third Euler angle describing L with respect to J. Fixed by minmal rotation condition.
[out]cosbeta_of_fCosine of polar angle between L and J
freqsInput Frequency series [Hz]
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
inclinationInclination : angle between LN and the line of sight
fRef_InReference frequency (Hz)
mprimeSpherical harmonic order m
lalParamsLAL Dictionary struct

Definition at line 1604 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomXPPNAngles()

int XLALSimIMRPhenomXPPNAngles ( REAL8Sequence **  alpha_of_f,
REAL8Sequence **  gamma_of_f,
REAL8Sequence **  cosbeta_of_f,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  inclination,
REAL8  fRef_In,
INT4  mprime,
LALDict *  lalParams 
)
Parameters
[out]alpha_of_fAzimuthal angle of L w.r.t J
[out]gamma_of_fThird Euler angle describing L w.r.t J, fixed by minimal rotation condition
[out]cosbeta_of_fCosine of polar angle between L and J
freqsInput Frequency series [Hz]
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
inclinationInclination : angle between LN and the line of sight
fRef_InReference frequency (Hz)
mprimeSpherical harmonic order m
lalParamsLAL Dictionary struct

Definition at line 1725 of file LALSimIMRPhenomX.c.

◆ XLALSimIMRPhenomX_PNR_GeneratePNRAngles()

int XLALSimIMRPhenomX_PNR_GeneratePNRAngles ( REAL8Sequence **  alphaPNR,
REAL8Sequence **  betaPNR,
REAL8Sequence **  gammaPNR,
REAL8Sequence **  freqs,
REAL8 alphaPNR_ref,
REAL8 gammaPNR_ref,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  inclination,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  fRef_In,
LALDict *  lalParams 
)

This is an external wrapper to generate the (2,2) PNR angles, following the prescription outlined in arXiv:2107.08876, given the standard inputs given to generate FD waveforms.

Parameters
[out]alphaPNRAlpha precession angle (rad)
[out]betaPNRBeta precession angle (rad)
[out]gammaPNRGamma precession angle (rad)
[out]freqsFrequency array (Hz)
[out]alphaPNR_refreference value of alpha (rad)
[out]gammaPNR_refreference value of gamma (rad)
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
inclinationAngle between orbital angular momentum and line-of-sight vector at reference frequency (rad)
deltaFFrequency spacing (Hz)
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz)
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 73 of file LALSimIMRPhenomX_PNR.c.

◆ XLALSimIMRPhenomX_PNR_GeneratePNRAnglesHM()

int XLALSimIMRPhenomX_PNR_GeneratePNRAnglesHM ( REAL8Sequence **  alphaPNR,
REAL8Sequence **  betaPNR,
REAL8Sequence **  gammaPNR,
REAL8Sequence **  freqs,
REAL8 alphaPNR_ref,
REAL8 gammaPNR_ref,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  inclination,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  fRef_In,
INT4  ell,
INT4  emmprime,
LALDict *  lalParams 
)

This is an external wrapper to generate the (l,m) PNR angles, following the prescriptions outlined in arXiv:2107.08876 and arXiv:##### FIXME: add reference, given the standard inputs given to generate FD waveforms.

Parameters
[out]alphaPNRAlpha precession angle (rad)
[out]betaPNRBeta precession angle (rad)
[out]gammaPNRGamma precession angle (rad)
[out]freqsFrequency array (Hz)
[out]alphaPNR_refReference value of alpha (rad)
[out]gammaPNR_refReference value of gamma (rad)
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
inclinationAngle between orbital angular momentum and line-of-sight vector at reference frequency (rad)
deltaFFrequency spacing (Hz)
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified.
fRef_InReference frequency (Hz)
ellOrbital index (int)
emmprimeAzimuthal index (int)
lalParamsLAL Dictionary struct

Definition at line 296 of file LALSimIMRPhenomX_PNR.c.

◆ IMRPhenomX_PNR_GeneratePNRAngles()

int IMRPhenomX_PNR_GeneratePNRAngles ( REAL8Sequence alphaPNR,
REAL8Sequence betaPNR,
REAL8Sequence gammaPNR,
const REAL8Sequence freqs,
IMRPhenomXWaveformStruct pWF,
IMRPhenomXPrecessionStruct pPrec,
LALDict *  lalParams 
)

Generate the tuned precession angles outlined in arXiv:2107.08876.

The general flow is as follows:

  • First check if deltaF > 0:

--> if it is, then we generate angles sampled on a uniform frequency grid and then compute the values of the angles at fRef using linear interpolation.

--> otherwise we expect non-uniform frequencies, so we generate interpolants and then evaluate them at fRef.

  • The reference values of alpha and gamma are stored in the pPrec struct for use later, and then we re-map the J-frame sky position using the reference value of beta.
  • NOTE: alpha0 is recomputed in IMRPhenomX_PNR_RemapThetaJSF, hence it is not used here unlike epsilon0. Instead it is applied inside IMRPhenomX_PNR_RemapThetaJSF.
Parameters
alphaPNRalpha precession angle (rad)
betaPNRbeta precession angle (rad)
gammaPNRgamma precession angle (rad)
freqsinput frequency array (Hz)
pWFwaveform struct
pPrecprecession struct
lalParamsLAL dictionary struct

Definition at line 633 of file LALSimIMRPhenomX_PNR.c.

◆ IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies()

int IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies ( REAL8Sequence alphaPNR,
REAL8Sequence betaPNR,
REAL8Sequence gammaPNR,
const REAL8Sequence freqs,
IMRPhenomXWaveformStruct pWF_SingleSpin,
IMRPhenomXPrecessionStruct pPrec_SingleSpin,
IMRPhenomX_PNR_alpha_parameters alphaParams,
IMRPhenomX_PNR_beta_parameters betaParams,
IMRPhenomXWaveformStruct pWF,
IMRPhenomXPrecessionStruct pPrec,
LALDict *  lalParams 
)

Generate the tuned precession angles outlined in arXiv:2107.08876 specifically on a uniform frequency grid.

  • First, we populate the required alpha and beta parameter structs, along with the two-spin structs if needed.
  • Next we store the two connection frequencies possibly needed for the HM angle frequency mapping.
  • Check if we should be attaching the MR contributions to beta: this is determined by calling the function IMRPhenomX_PNR_AttachMRBeta.

--> If yes, loop through the frequencies and generate full alpha and beta

--> If no, loop through the frequencies and generate full alpha but only the inspiral beta

  • Finally, if an initialized struct is passed for gamma, compute it.

Definition at line 802 of file LALSimIMRPhenomX_PNR.c.

◆ IMRPhenomX_PNR_GeneratePNRAlphaForAntisymmetry()

int IMRPhenomX_PNR_GeneratePNRAlphaForAntisymmetry ( REAL8Sequence alphaPNR,
const REAL8Sequence freqs,
IMRPhenomXWaveformStruct pWF,
IMRPhenomXPrecessionStruct pPrec,
LALDict *  lalParams 
)

Internal helper function to generate PNR alpha for the antisymmetric waveform.

Definition at line 930 of file LALSimIMRPhenomX_PNR.c.

◆ IMRPhenomX_PNR_GeneratePNRAngleInterpolants()

int IMRPhenomX_PNR_GeneratePNRAngleInterpolants ( IMRPhenomX_PNR_angle_spline angle_spline,
IMRPhenomXWaveformStruct pWF,
IMRPhenomXPrecessionStruct pPrec,
LALDict *  lalparams 
)

Generate the tuned precession angles outlined in arXiv:2107.08876 specifically populate the angle interpolant struct.

  • First, we check for an activated mode array.

--> If it exists, then we're likely computing HM angles with these interpolants and we need to extend the frequency range to capture the HM frequency map. Use the activated modes to determine this.

--> Otherwise we use the default f_min and f_max.

  • Get the required deltaF from IMRPhenomX_PNR_HMInterpolationDeltaF.
  • Compute uniform alpha and beta
  • Populate interpolants for alpha and beta, then compute gamma using these interpolants

Definition at line 1025 of file LALSimIMRPhenomX_PNR.c.

◆ XLALSimIMRPhenomX_PNR_GenerateEffectiveRingdownFreq()

REAL8 XLALSimIMRPhenomX_PNR_GenerateEffectiveRingdownFreq ( REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  f_min,
REAL8  f_max,
REAL8  fRef,
UINT4  ell,
UINT4  emmprime,
LALDict *  lalParams 
)

External wrapper for IMRPhenomX_PNR_GenerateEffectiveRingdownFreq.

Computes the effective ringdown frequency for a given (l,m) multipole following the prescription given in arXiv:##### FIXME: add citation

Parameters
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz)
fRefReference frequency (Hz)
ellOrbital index
emmprimeazimuthal index
lalParamsLAL Dictionary struct

Definition at line 1339 of file LALSimIMRPhenomX_PNR.c.

◆ XLALSimIMRPhenomX_PNR_GenerateRingdownPNRBeta()

REAL8 XLALSimIMRPhenomX_PNR_GenerateRingdownPNRBeta ( REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  fRef_In,
LALDict *  lalParams 
)

External wrapper for IMRPhenomX_PNR_GenerateRingdownPNRBeta.

Generate PNR beta as described in Eq. 60 or arXiv:2107.08876 and evaluate at the upper connection frequency f_f described in Sec. 8C.

Parameters
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
deltaFSampling Frequency (Hz)
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz)
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 1411 of file LALSimIMRPhenomX_PNR.c.

◆ XLALSimIMRPhenomX_PNR_InterpolationDeltaF()

int XLALSimIMRPhenomX_PNR_InterpolationDeltaF ( INT4 twospin,
INT4 precversion,
REAL8 deltaF,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  f_min,
REAL8  f_max,
REAL8  fRef,
LALDict *  lalParams 
)
Parameters
twospinUNDOCUMENTED
precversionUNDOCUMENTED
deltaFUNDOCUMENTED
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified.
fRefReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 1489 of file LALSimIMRPhenomX_PNR.c.

◆ XLALSimIMRPhenomXHMGenerateFDOneMode()

int XLALSimIMRPhenomXHMGenerateFDOneMode ( COMPLEX16FrequencySeries **  htildelm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
UINT4  ell,
INT4  emm,
REAL8  distance,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns the Fourier domain strain of just one mode: h_lm.

Supports positive and negative m. Notice than even when the specified frequencies are positives, the m>0 only lives in the negative frequencies regime. With m>0 the mode h_lm is zero for positive frequencies and for the negative frequencies is equal to (-1)^l h*_l-m(-f). In the contrary, h_l-m is zero for negative frequencies and only lives for positive frequencies. This is a wrapper function that uses XLALSimIMRPhenomXASGenerateFD for the 22 mode and IMRPhenomXHMGenerateFDOneMode for the higher modes.

Parameters
[out]htildelmFD waveform
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
elll index of the mode
emmm index of the mode
distanceLuminosity distance (m)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = 0.3
deltaFSampling frequency (Hz)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsExtra params

Definition at line 287 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMFrequencySequenceOneMode()

int XLALSimIMRPhenomXHMFrequencySequenceOneMode ( COMPLEX16FrequencySeries **  htildelm,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
UINT4  ell,
INT4  emm,
REAL8  distance,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Get the model evaluated in a prebuilt frequency array.

It does not use Multibanding.

Parameters
[out]htildelmFD waveform
freqsfrequency array to evaluate model (positives) (Hz)
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
elll index of the mode
emmm index of the mode
distanceLuminosity distance (m)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsUNDOCUMENTED

Definition at line 751 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMModes()

int XLALSimIMRPhenomXHMModes ( SphHarmFrequencySeries **  hlms,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  S1z,
REAL8  S2z,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  f_ref,
REAL8  phiRef,
REAL8  distance,
LALDict *  LALparams 
)

Function to obtain a SphHarmFrequencySeries with the individual modes h_lm.

By default it returns all the modes available in the model, both positive and negatives. With the mode array option in the LAL dictionary, the user can specify a custom mode array. This function is to be used by ChooseFDModes.

Parameters
[out]hlmslist with single modes h_lm
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
S1zz-component of the dimensionless spin of object 1
S2zz-component of the dimensionless spin of object 2
deltaFfrequency spacing (Hz)
f_minstarting GW frequency (Hz)
f_maxending GW frequency (Hz)
f_refreference GW frequency (Hz)
phiRefphase shift at reference frequency
distancedistance of source (m)
LALparamsLAL dictionary with extra options

Definition at line 907 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHM_SpheroidalPhase()

int XLALSimIMRPhenomXHM_SpheroidalPhase ( REAL8Sequence **  phase,
const REAL8Sequence freqs,
UINT4  ell,
INT4  emm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  distance,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)
Parameters
[out]phaseFD waveform
freqsfrequency array to evaluate model (positives)
elll index of the mode
emmm index of the mode
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
distanceLuminosity distance (m)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsUNDOCUMENTED

Definition at line 1065 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHM_SpheroidalAmplitude()

int XLALSimIMRPhenomXHM_SpheroidalAmplitude ( REAL8Sequence **  amplitude,
const REAL8Sequence freqs,
UINT4  ell,
INT4  emm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  distance,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)
Parameters
[out]amplitudeFD waveform
freqsfrequency array to evaluate model (positives)
elll index of the mode
emmm index of the mode
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
distanceLuminosity distance (m)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsUNDOCUMENTED

Definition at line 1222 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHM()

int XLALSimIMRPhenomXHM ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns the hptilde and hctilde of the multimode waveform for positive frequencies.

XLALSimIMRPhenomXHM calls the function for a single mode that can be XLALSimIMRPhenomXHMGenerateFDOneMode or XLALSimIMRPhenomXHMMultiBandOneMode, depending on if the Multibanding is active or not.

By default XLALSimIMRPhenomXHM is only used when the Multibanding is activated, since each mode has a different coarse frequency array and we can not recycle the array.

This is just a wrapper of the function that actually carry out the calculations: IMRPhenomXHM_MultiMode2.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1Lz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2Lz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = 0.3
deltaFSampling frequency (Hz)
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
fRef_InReference frequency
lalParamslinked list containing the extra parameters

Definition at line 1398 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHM2()

int XLALSimIMRPhenomXHM2 ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns the hptilde and hctilde of the multimode waveform for positive frequencies.

XLALSimIMRPhenomXHM2 builds each mode explicitly in the loop over modes, recycling some common quantities between modes like the frequency array, the powers of frequencies, structures, etc so it has less overhead than XLALSimIMRPhenomXHM.

By default XLALSimIMRPhenomXHM2 is only used for the version without Multibanding.

This is just a wrapper of the function that actually carry out the calculations: IMRPhenomXHM_MultiMode2.

Parameters
[out]hptildeFrequency domain h+ GW strain
[out]hctildeFrequency domain hx GW strain
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = 0.3
deltaFSampling frequency (Hz)
distanceLuminosity distance (m)
inclinationInclination of the source
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary

Definition at line 1485 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMFrequencySequence()

int XLALSimIMRPhenomXHMFrequencySequence ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1z,
REAL8  chi2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies specified in the REAL8Sequence *freqs (which can be unequally spaced).

No zeros are added. Assumes positive frequencies.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
freqsInput Frequency series [Hz]
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 1633 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMAmplitude()

int XLALSimIMRPhenomXHMAmplitude ( REAL8Sequence **  amplitude,
const REAL8Sequence freqs,
UINT4  ell,
INT4  emm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns amplitude of one single mode in a custom frequency array.

If the in-plane spin components are non-zero, this is the amplitude of the modified co-precessing mode, where the ringdown/damping frequencies corresponds to those of the precessing final spin.

Parameters
[out]amplitudeFD amp
freqsInput Frequency Array (Hz)
elll index of the mode
emmm index of the mode
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
S2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
S2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceLuminosity distance (m)
phiRefOrbital amp at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsExtra params

Definition at line 2396 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMPhase()

int XLALSimIMRPhenomXHMPhase ( REAL8Sequence **  phase,
const REAL8Sequence freqs,
UINT4  ell,
INT4  emm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns phase of one single mode in a custom frequency array.

If the in-plane spin components are non-zero, this is the phase of the modified co-precessing mode, where the ringdown/damping frequencies corresponds to those of the precessing final spin.

Parameters
[out]phaseFD amp
freqsInput Frequency Array (Hz)
elll index of the mode
emmm index of the mode
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
S2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
S2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
S2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceLuminosity distance (m)
phiRefOrbital amp at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsExtra params

Definition at line 2593 of file LALSimIMRPhenomXHM.c.

◆ IMRPhenomXHM_Phase_for_Initialization()

INT4 IMRPhenomXHM_Phase_for_Initialization ( double *  phase,
double *  dphase,
double  Mf,
INT4  ell,
INT4  emm,
IMRPhenomXWaveformStruct pWF,
LALDict *  lalParams 
)

Definition at line 2813 of file LALSimIMRPhenomXHM.c.

◆ XLALSimIMRPhenomXHMMultiBandOneMode()

int XLALSimIMRPhenomXHMMultiBandOneMode ( COMPLEX16FrequencySeries **  htildelm,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
UINT4  ell,
INT4  emmIn,
REAL8  distance,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Return htildelm, the waveform of one mode without mode-mixing.

e^(I phi) is interpolated with linear order and an iterative procedure. Amplitude uses standard 1st (by default) or 3rd order interpolation with gsl.

Parameters
[out]htildelmFD waveform
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
elll index of the mode
emmInm index of the mode
distanceLuminosity distance (m)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = f_CUT
deltaFSampling frequency (Hz)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsExtra params

Definition at line 448 of file LALSimIMRPhenomXHM_multiband.c.

◆ XLALSimIMRPhenomXHMMultiBandOneModeMixing()

int XLALSimIMRPhenomXHMMultiBandOneModeMixing ( COMPLEX16FrequencySeries **  htildelm,
COMPLEX16FrequencySeries htilde22,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1L,
REAL8  chi2L,
UINT4  ell,
INT4  emmIn,
REAL8  distance,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns htildelm the waveform of one mode that present mode-mixing.

The multibanding is applied to the spherical part (inspiral and intermediate) and to the spheroidal part (ringdown). Both are interpolated and evaluated in the fine frequency grid and the ringdown part is rotated back to spherical

Parameters
[out]htildelmFD waveform
htilde22Precomputed FD waveform of dominant mode
m1_SIMass of companion 1 (kg)
m2_SIMass of companion 2 (kg)
chi1LDimensionless aligned spin of companion 1
chi2LDimensionless aligned spin of companion 2
elll index of the mode
emmInm index of the mode
distanceLuminosity distance (m)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to Mf = f_CUT
deltaFSampling frequency (Hz)
phiRefOrbital phase at fRef (rad)
fRef_InReference frequency (Hz)
lalParamsExtra params

Definition at line 1242 of file LALSimIMRPhenomXHM_multiband.c.

◆ XLALSimIMRPhenomXPHM()

int XLALSimIMRPhenomXPHM ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns hptilde and hctilde of the multimode precessing waveform for positive frequencies in an equally spaced grid.

This is the default function used when calling ChooseFDWaveform.

It computes the non-precessing modes just once and do the twisting up according to eqs. 3.5-3.7 in the Precessing paper.

It is just a wrapper of the internal function that actually carries out the calculation: IMRPhenomXPHM_hplushcross.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified.
deltaFSampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0.
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 167 of file LALSimIMRPhenomXPHM.c.

◆ XLALSimIMRPhenomXPHMFromModes()

int XLALSimIMRPhenomXPHMFromModes ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  f_min,
REAL8  f_max,
REAL8  deltaF,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns hptilde and hctilde of the multimode precessing waveform for positive frequencies in an equally space grid.

This function is equivalent to XLALSimIMRPhenomXPHM, gives the same result but it is computed by first calling all the precessing indiviudal modes in the J-frame and then sum them all with Ylm(thetaJN, 0) since the J-frame is aligned such that the line of sight N is in the x-z plane of the J-frame. See appendix C and in particular eq. C8 in Precessing paper.

This function is slower since it has to compute the non-precessing modes again for every precessing mode.

It is just a wrapper of the internal function that actually carries out the calculation: IMRPhenomXPHM_hplushcross_from_modes.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
f_minStarting GW frequency (Hz)
f_maxEnding GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified.
deltaFSampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0.
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 382 of file LALSimIMRPhenomXPHM.c.

◆ XLALSimIMRPhenomXPHMFrequencySequence()

int XLALSimIMRPhenomXPHMFrequencySequence ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
const REAL8Sequence freqs,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  fRef_In,
LALDict *  lalParams 
)

Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies specified in the REAL8Sequence *freqs (which can be unequally spaced).

No zeros are added. Assumes positive frequencies. This is the function used when calling ChooseFDWaveformSequence. It is a wrapper that calls the function IMRPhenomXPHM_hplushcross or IMRPhenomXPHM_hplushcross_from_modes and we can change which one is used by modifying the option 'UseModes' in the LAL dictionary. No multibanding is used since this technique is only for equal-spaced frequency grids.

Parameters
[out]hptildeFrequency-domain waveform h+
[out]hctildeFrequency-domain waveform hx
freqsInput Frequency series (Hz)
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distanceDistance of source (m)
inclinationinclination of source (rad)
phiRefOrbital phase (rad) at reference frequency
fRef_InReference frequency (Hz)
lalParamsLAL Dictionary struct

Definition at line 582 of file LALSimIMRPhenomXPHM.c.

◆ XLALSimIMRPhenomXPHMOneMode()

int XLALSimIMRPhenomXPHMOneMode ( COMPLEX16FrequencySeries **  hlmpos,
COMPLEX16FrequencySeries **  hlmneg,
const UINT4  l,
const INT4  m,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
const REAL8  distance,
const REAL8  inclination,
const REAL8  phiRef,
const REAL8  deltaF,
const REAL8  f_min,
const REAL8  f_max,
const REAL8  fRef_In,
LALDict *  lalParams 
)

Function to compute one hlm precessing mode in an uniform frequency grid.

By default the mode is given in the inertial J-frame. It can be transformed to the L0-frame with the option "PhenomXPHMModesL0Frame" which currently only works for cases near AS limit. It can return the co-precessing mode with the option "PhenomXPHMPrecModes": this corresponds to the XHM mode with the modified ringdown/damping frequencies of the precessing final spin. In this case only m<0 are supported, so hlmneg will be filled with zeros. Returns two frequency series, one for the positive frequencies and other for the negative frequencies since, as opposite to the aligned spin case, in the precessing case all the modes have support in the whole frequency regime. This is a wrapper of the internal core function that actually does the calculation IMRPhenomXPHM_OneMode.

Parameters
[out]hlmposFrequency-domain waveform hlm inertial frame positive frequencies
[out]hlmnegFrequency-domain waveform hlm inertial frame negative frequencies
lFirst index of the (l,m) precessing mode
mSecond index of the (l,m) precessing mode
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
deltaFSampling frequency (Hz)
f_minStarting GW frequency (Hz)
f_maxEnd frequency; 0 defaults to ringdown cutoff freq
fRef_InReference frequency
lalParamsLAL Dictionary

Definition at line 2293 of file LALSimIMRPhenomXPHM.c.

◆ XLALSimIMRPhenomXPHMFrequencySequenceOneMode()

int XLALSimIMRPhenomXPHMFrequencySequenceOneMode ( COMPLEX16FrequencySeries **  hlmpos,
COMPLEX16FrequencySeries **  hlmneg,
const REAL8Sequence freqs,
const UINT4  l,
const INT4  m,
REAL8  m1_SI,
REAL8  m2_SI,
REAL8  chi1x,
REAL8  chi1y,
REAL8  chi1z,
REAL8  chi2x,
REAL8  chi2y,
REAL8  chi2z,
const REAL8  distance,
const REAL8  inclination,
const REAL8  phiRef,
const REAL8  fRef_In,
LALDict *  lalParams 
)

Function to compute one hlm precessing mode on a custom frequency grid.

Equivalent options and behaviour to that of XLALSimIMRPhenomXPHMOneMode.

Parameters
[out]hlmposFrequency-domain waveform hlm inertial frame positive frequencies
[out]hlmnegFrequency-domain waveform hlm inertial frame negative frequencies
freqsInput Frequency series [Hz]
lFirst index of the (l,m) precessing mode
mSecond index of the (l,m) precessing mode
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
chi1xx-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1yy-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi1zz-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1)
chi2xx-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2yy-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
chi2zz-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1)
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
fRef_InReference frequency
lalParamsLAL Dictionary

Definition at line 2559 of file LALSimIMRPhenomXPHM.c.