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
Routines for IMRPhenomXAS | |
C code for IMRPhenomXAS phenomenological waveform model. 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.
Waveform flags: InsPhaseVersion: Determines the inspiral phase model.
IntPhaseVersion: Determines the intermediate phase model.
RDPhaseVersion: Determines the merger-ringdown phase model.
InsAmpVersion : Determines inspiral amplitude model.
IntAmpVersion : Determines intermediate amplitude model.
RDAmpVersion : Determines the merger-ringdown amplitude model.
| |
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. 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.
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:
PhenomXPExpansionOrder:
PhenomXPFinalSpinMod:
PhenomXPConvention (App. C and Table IV of arXiv:2004.06503):
| |
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. This is a C-code impementation of the PNR angle model outlined in arXiv:2107.08876. Available flags: PhenomXPNRUseTunedAngles:
| |
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. 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.
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 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.
AmpInterpol: Determines the gsl interpolation order for the amplitude.
| |
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. 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.
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.
UseModes: Determine how the polarizations hp, hc are computed.
ModesL0Frame: Determine in which frame the individual precessing modes are returned.
PrecModes: Determine which indiviual modes are returned, the non-precessing or the precessing.
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.
MBandPrecVersion: Determines the algorithm to build the non-uniform frequency grid for the Euler angles.
| |
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... | |
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.
[out] | htilde22 | FD waveform |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
distance | Luminosity distance (m) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = 0.3 | |
deltaF | Sampling frequency (Hz) | |
phi0 | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary |
Definition at line 175 of file LALSimIMRPhenomX.c.
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.
[out] | htilde22 | FD waveform |
[out] | freqs | Frequency series [Hz] |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
distance | Luminosity distance (m) | |
phi0 | Phase at reference frequency | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary |
Definition at line 319 of file LALSimIMRPhenomX.c.
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.
m1_SI | mass of companion 1 (kg) |
m2_SI | mass of companion 2 (kg) |
chi1L | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi2L | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
f_start | Initial frequency (Hz) |
Definition at line 408 of file LALSimIMRPhenomX.c.
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 | ||
) |
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
f_min | Starting GW frequency (Hz) | |
f_max | Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. | |
deltaF | Sampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0. | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 1038 of file LALSimIMRPhenomX.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
freqs | input Frequency series [Hz] | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 1255 of file LALSimIMRPhenomX.c.
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 | ||
) |
[out] | chi1L | Dimensionless aligned spin on companion 1 |
[out] | chi2L | Dimensionless aligned spin on companion 2 |
[out] | chi_p | Effective precession parameter: Schmidt, Ohme, Hannam, PRD, 91,024043 (2015) |
[out] | thetaJN | Angle between J0 and line of sight (z-direction) |
[out] | alpha0 | Initial value of alpha angle (azimuthal precession angle) |
[out] | phi_aligned | Initial phase to feed the underlying aligned-spin model |
[out] | zeta_polarization | Angle to rotate the polarizations |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
f_ref | Reference GW frequency (Hz) | |
phiRef | Reference phase (Hz) | |
incl | Inclination : angle between LN and the line of sight | |
chi1x | Initial value of chi1x: dimensionless spin of BH 1 in L frame | |
chi1y | Initial value of chi1y: dimensionless spin of BH 1 in L frame | |
chi1z | Initial value of chi1z: dimensionless spin of BH 1 in L frame | |
chi2x | Initial value of chi2x: dimensionless spin of BH 2 in L frame | |
chi2y | Initial value of chi2y: dimensionless spin of BH 2 in L frame | |
chi2z | Initial value of chi2z: dimensionless spin of BH 2 in L frame | |
lalParams | LAL Dictionary |
Definition at line 1417 of file LALSimIMRPhenomX.c.
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 | ||
) |
[out] | alpha_of_f | The azimuthal angle of L around J |
[out] | gamma_of_f | The third Euler angle describing L with respect to J. Fixed by minmal rotation condition. |
[out] | cosbeta_of_f | Cosine of polar angle between L and J |
freqs | Input Frequency series [Hz] | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
inclination | Inclination : angle between LN and the line of sight | |
fRef_In | Reference frequency (Hz) | |
mprime | Spherical harmonic order m | |
lalParams | LAL Dictionary struct |
Definition at line 1604 of file LALSimIMRPhenomX.c.
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 | ||
) |
[out] | alpha_of_f | Azimuthal angle of L w.r.t J |
[out] | gamma_of_f | Third Euler angle describing L w.r.t J, fixed by minimal rotation condition |
[out] | cosbeta_of_f | Cosine of polar angle between L and J |
freqs | Input Frequency series [Hz] | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
inclination | Inclination : angle between LN and the line of sight | |
fRef_In | Reference frequency (Hz) | |
mprime | Spherical harmonic order m | |
lalParams | LAL Dictionary struct |
Definition at line 1725 of file LALSimIMRPhenomX.c.
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.
[out] | alphaPNR | Alpha precession angle (rad) |
[out] | betaPNR | Beta precession angle (rad) |
[out] | gammaPNR | Gamma precession angle (rad) |
[out] | freqs | Frequency array (Hz) |
[out] | alphaPNR_ref | reference value of alpha (rad) |
[out] | gammaPNR_ref | reference value of gamma (rad) |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
inclination | Angle between orbital angular momentum and line-of-sight vector at reference frequency (rad) | |
deltaF | Frequency spacing (Hz) | |
f_min | Starting GW frequency (Hz) | |
f_max | Ending GW frequency (Hz) | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 73 of file LALSimIMRPhenomX_PNR.c.
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.
[out] | alphaPNR | Alpha precession angle (rad) |
[out] | betaPNR | Beta precession angle (rad) |
[out] | gammaPNR | Gamma precession angle (rad) |
[out] | freqs | Frequency array (Hz) |
[out] | alphaPNR_ref | Reference value of alpha (rad) |
[out] | gammaPNR_ref | Reference value of gamma (rad) |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
inclination | Angle between orbital angular momentum and line-of-sight vector at reference frequency (rad) | |
deltaF | Frequency spacing (Hz) | |
f_min | Starting GW frequency (Hz) | |
f_max | Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. | |
fRef_In | Reference frequency (Hz) | |
ell | Orbital index (int) | |
emmprime | Azimuthal index (int) | |
lalParams | LAL Dictionary struct |
Definition at line 296 of file LALSimIMRPhenomX_PNR.c.
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:
--> 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.
alphaPNR | alpha precession angle (rad) |
betaPNR | beta precession angle (rad) |
gammaPNR | gamma precession angle (rad) |
freqs | input frequency array (Hz) |
pWF | waveform struct |
pPrec | precession struct |
lalParams | LAL dictionary struct |
Definition at line 633 of file LALSimIMRPhenomX_PNR.c.
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.
--> 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
Definition at line 802 of file LALSimIMRPhenomX_PNR.c.
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.
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.
--> 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.
Definition at line 1025 of file LALSimIMRPhenomX_PNR.c.
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
m1_SI | mass of companion 1 (kg) |
m2_SI | mass of companion 2 (kg) |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
f_min | Starting GW frequency (Hz) |
f_max | Ending GW frequency (Hz) |
fRef | Reference frequency (Hz) |
ell | Orbital index |
emmprime | azimuthal index |
lalParams | LAL Dictionary struct |
Definition at line 1339 of file LALSimIMRPhenomX_PNR.c.
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.
m1_SI | mass of companion 1 (kg) |
m2_SI | mass of companion 2 (kg) |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
deltaF | Sampling Frequency (Hz) |
f_min | Starting GW frequency (Hz) |
f_max | Ending GW frequency (Hz) |
fRef_In | Reference frequency (Hz) |
lalParams | LAL Dictionary struct |
Definition at line 1411 of file LALSimIMRPhenomX_PNR.c.
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 | ||
) |
twospin | UNDOCUMENTED |
precversion | UNDOCUMENTED |
deltaF | UNDOCUMENTED |
m1_SI | mass of companion 1 (kg) |
m2_SI | mass of companion 2 (kg) |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) |
f_min | Starting GW frequency (Hz) |
f_max | Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. |
fRef | Reference frequency (Hz) |
lalParams | LAL Dictionary struct |
Definition at line 1489 of file LALSimIMRPhenomX_PNR.c.
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.
[out] | htildelm | FD waveform |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
ell | l index of the mode | |
emm | m index of the mode | |
distance | Luminosity distance (m) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = 0.3 | |
deltaF | Sampling frequency (Hz) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | Extra params |
Definition at line 287 of file LALSimIMRPhenomXHM.c.
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.
[out] | htildelm | FD waveform |
freqs | frequency array to evaluate model (positives) (Hz) | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
ell | l index of the mode | |
emm | m index of the mode | |
distance | Luminosity distance (m) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | UNDOCUMENTED |
Definition at line 751 of file LALSimIMRPhenomXHM.c.
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.
[out] | hlms | list with single modes h_lm |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
S1z | z-component of the dimensionless spin of object 1 | |
S2z | z-component of the dimensionless spin of object 2 | |
deltaF | frequency spacing (Hz) | |
f_min | starting GW frequency (Hz) | |
f_max | ending GW frequency (Hz) | |
f_ref | reference GW frequency (Hz) | |
phiRef | phase shift at reference frequency | |
distance | distance of source (m) | |
LALparams | LAL dictionary with extra options |
Definition at line 907 of file LALSimIMRPhenomXHM.c.
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 | ||
) |
[out] | phase | FD waveform |
freqs | frequency array to evaluate model (positives) | |
ell | l index of the mode | |
emm | m index of the mode | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
distance | Luminosity distance (m) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | UNDOCUMENTED |
Definition at line 1065 of file LALSimIMRPhenomXHM.c.
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 | ||
) |
[out] | amplitude | FD waveform |
freqs | frequency array to evaluate model (positives) | |
ell | l index of the mode | |
emm | m index of the mode | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
distance | Luminosity distance (m) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | UNDOCUMENTED |
Definition at line 1222 of file LALSimIMRPhenomXHM.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1L | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2L | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = 0.3 | |
deltaF | Sampling frequency (Hz) | |
distance | distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | reference orbital phase (rad) | |
fRef_In | Reference frequency | |
lalParams | linked list containing the extra parameters |
Definition at line 1398 of file LALSimIMRPhenomXHM.c.
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.
[out] | hptilde | Frequency domain h+ GW strain |
[out] | hctilde | Frequency domain hx GW strain |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = 0.3 | |
deltaF | Sampling frequency (Hz) | |
distance | Luminosity distance (m) | |
inclination | Inclination of the source | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary |
Definition at line 1485 of file LALSimIMRPhenomXHM.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
freqs | Input Frequency series [Hz] | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 1633 of file LALSimIMRPhenomXHM.c.
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.
[out] | amplitude | FD amp |
freqs | Input Frequency Array (Hz) | |
ell | l index of the mode | |
emm | m index of the mode | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
S1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
S2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
S2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Luminosity distance (m) | |
phiRef | Orbital amp at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | Extra params |
Definition at line 2396 of file LALSimIMRPhenomXHM.c.
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.
[out] | phase | FD amp |
freqs | Input Frequency Array (Hz) | |
ell | l index of the mode | |
emm | m index of the mode | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
S1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
S2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
S2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
S2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Luminosity distance (m) | |
phiRef | Orbital amp at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | Extra params |
Definition at line 2593 of file LALSimIMRPhenomXHM.c.
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.
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.
[out] | htildelm | FD waveform |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
ell | l index of the mode | |
emmIn | m index of the mode | |
distance | Luminosity distance (m) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = f_CUT | |
deltaF | Sampling frequency (Hz) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | Extra params |
Definition at line 448 of file LALSimIMRPhenomXHM_multiband.c.
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
[out] | htildelm | FD waveform |
htilde22 | Precomputed FD waveform of dominant mode | |
m1_SI | Mass of companion 1 (kg) | |
m2_SI | Mass of companion 2 (kg) | |
chi1L | Dimensionless aligned spin of companion 1 | |
chi2L | Dimensionless aligned spin of companion 2 | |
ell | l index of the mode | |
emmIn | m index of the mode | |
distance | Luminosity distance (m) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to Mf = f_CUT | |
deltaF | Sampling frequency (Hz) | |
phiRef | Orbital phase at fRef (rad) | |
fRef_In | Reference frequency (Hz) | |
lalParams | Extra params |
Definition at line 1242 of file LALSimIMRPhenomXHM_multiband.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
f_min | Starting GW frequency (Hz) | |
f_max | Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. | |
deltaF | Sampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0. | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 167 of file LALSimIMRPhenomXPHM.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
f_min | Starting GW frequency (Hz) | |
f_max | Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. | |
deltaF | Sampling frequency (Hz). To use non-uniform frequency grid, set deltaF <= 0. | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 382 of file LALSimIMRPhenomXPHM.c.
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.
[out] | hptilde | Frequency-domain waveform h+ |
[out] | hctilde | Frequency-domain waveform hx |
freqs | Input Frequency series (Hz) | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | Distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | Orbital phase (rad) at reference frequency | |
fRef_In | Reference frequency (Hz) | |
lalParams | LAL Dictionary struct |
Definition at line 582 of file LALSimIMRPhenomXPHM.c.
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.
[out] | hlmpos | Frequency-domain waveform hlm inertial frame positive frequencies |
[out] | hlmneg | Frequency-domain waveform hlm inertial frame negative frequencies |
l | First index of the (l,m) precessing mode | |
m | Second index of the (l,m) precessing mode | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | reference orbital phase (rad) | |
deltaF | Sampling frequency (Hz) | |
f_min | Starting GW frequency (Hz) | |
f_max | End frequency; 0 defaults to ringdown cutoff freq | |
fRef_In | Reference frequency | |
lalParams | LAL Dictionary |
Definition at line 2293 of file LALSimIMRPhenomXPHM.c.
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.
[out] | hlmpos | Frequency-domain waveform hlm inertial frame positive frequencies |
[out] | hlmneg | Frequency-domain waveform hlm inertial frame negative frequencies |
freqs | Input Frequency series [Hz] | |
l | First index of the (l,m) precessing mode | |
m | Second index of the (l,m) precessing mode | |
m1_SI | mass of companion 1 (kg) | |
m2_SI | mass of companion 2 (kg) | |
chi1x | x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1y | y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi1z | z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) | |
chi2x | x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2y | y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
chi2z | z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) | |
distance | distance of source (m) | |
inclination | inclination of source (rad) | |
phiRef | reference orbital phase (rad) | |
fRef_In | Reference frequency | |
lalParams | LAL Dictionary |
Definition at line 2559 of file LALSimIMRPhenomXPHM.c.