C code for NRHybSur3dq8 waveform model, an NR-hybrid surrogate model for aligned-spin BBH. More...
Prototypes | |
static bool | NRHybSur3dq8_IsSetup (void) |
static void | NRHybSur3dq8_Init_LALDATA (void) |
int | NRHybSur3dq8_fitParams (gsl_vector *fit_params, const REAL8 q, const REAL8 chi1z, const REAL8 chi2z) |
Map from mass ratio and spins to surrogate fit parameters. More... | |
int | NRHybSur3dq8_core (gsl_vector **phi_22, EvaluatedDataPieces **evaluated_mode_dps, LIGOTimeGPS *epoch, const REAL8 deltaTOverM, const REAL8 fMin, const REAL8 fRef, REAL8 q, const REAL8 Mtot_sec, REAL8 chi1z, REAL8 chi2z, LALValue *ModeArray, LALDict *LALparams) |
This is the core function of the NRHybSur3dq8 model. More... | |
C code for NRHybSur3dq8 waveform model, an NR-hybrid surrogate model for aligned-spin BBH.
The binary data file is available at https://dcc.ligo.org/LIGO-T1900034. Get the lalsuite-extra repo or put the data into a location in your LAL_DATA_PATH.
Paper: https://arxiv.org/abs/1812.07865
Parameter ranges:
q = [1, 10.1] and \(\chi_{1z}, \chi_{2z}\) = [-0.81, 0.81] or q = [1, 9.1] and \(\chi_{1z}, \chi_{2z}\) = [-0.91, 0.91]
modes: \( \ell \leq 4, m \geq 0 \), and (5,5), but not (4,1) or (4,0). m<0 modes are determined from the m \(\geq0\) modes.
\(M \geq 2.25 M_{\odot} \), for fstart=20Hz, for all modes.
Training parameter ranges:
q = [1, 8]
\(\chi_{1z}, \chi_{2z}\) = [-0.8, 0.8]
But extrapolates reasonably to the above mass ratios and spins.
Definition in file LALSimIMRNRHybSur3dq8.h.
Go to the source code of this file.
Variables | |
static const char | NRHybSur3dq8_DATAFILE [] = "NRHybSur3dq8_lal.h5" |
|
static |
|
static |
int NRHybSur3dq8_fitParams | ( | gsl_vector * | fit_params, |
const REAL8 | q, | ||
const REAL8 | chi1z, | ||
const REAL8 | chi2z | ||
) |
Map from mass ratio and spins to surrogate fit parameters.
The fit parameters are \([log_e(q), \hat{\chi}, \chi_a]\). \(\hat{\chi}\) is defined in Eq.(46) of arxiv:1812.07865. \(\chi_a = (\chi_{1z} - \chi_{2z})/2 \). These are described in Sec.VI.C.3 of arxiv:1812.07865.
fit_params | Output: mapped fit parameters. |
q | Mass ratio m1 / m2 >= 1. |
chi1z | Dimless spin of heavier BH. |
chi2z | Dimless spin of lighter BH. |
Definition at line 155 of file LALSimIMRNRHybSur3dq8.c.
int NRHybSur3dq8_core | ( | gsl_vector ** | phi_22, |
EvaluatedDataPieces ** | evaluated_mode_dps, | ||
LIGOTimeGPS * | epoch, | ||
const REAL8 | deltaT, | ||
const REAL8 | fMin, | ||
const REAL8 | fRef, | ||
REAL8 | q, | ||
const REAL8 | Mtot_sec, | ||
REAL8 | chi1z, | ||
REAL8 | chi2z, | ||
LALValue * | ModeArray, | ||
LALDict * | LALparams | ||
) |
This is the core function of the NRHybSur3dq8 model.
It evaluates all required waveform modes. For each mode, it evaluates each waveform data piece. The different data pieces are described in Sec.VI.B of arxiv:1812.07865. For the (2,2) mode the data pieces are the amplitude and phase. Note that we model the phase residual but add back the 0PN term at evaluation time. For all other modes the data pieces are the real and imaginary parts of the strain in the coorbital frame.
The reference point is the time at which the (2,2) mode frequency equals fRef. If fRef=0, sets fRef=fMin. We set the orbital phase to 0 at the reference point. The orbital phase is obtained as \(\phi_{22}/2\), so this leaves a pi ambiguity. But the surrogate data is already aligned such that the heavier BH is on the +ve x-axis at t=-1000M. See Sec.VI.A.4 of arxiv:1812.07865. This resolves the pi ambiguity. This means that after the realignment, the orbital phase at reference frequency fRef is 0, or Bh1 is on the +ve x-axis. Note that this is alignment is done using only waveform quantities, so this doesn't necessarily correspond to the (gauge dependent) NR Bh trajectories. The modes are returned in this frame, which agrees with the LAL convention. When combining the modes to get the polarizations, the Ylms should be evaluated at (inclination, pi/2 - phiRef), following the LAL convention.
Only uses data at (2,2) mode frequencies >= fMin. This determines the start time. The start time, along with the step size deltaT, is used to determine the output_times. Uses cubic spline interpolation to interpolate from the surrogate's time array to output_times.
NOTE: If mass ratio q<1, the labels of the two BHs are swapped internally.
Returned values:
phi_22: contains the phase of the (2,2) mode. This is always evaluated as this is required for other modes as well to transform from coorbital frame to inertial frame.
evaluated_mode_dps: Contains all other data pieces. This is a list of size num_modes_incl, the number of modes to include. For each mode this contains the amplitude, and real and imaginary parts of the coorbital frame strain. For the (2,2) mode only the amplitude is required. For all other modes only the coorbital frame strain is required. So, we evaluate only the required data pieces of each mode. The amplitude and coorbital frame strain is in units of rh/M and needs to be rescaled to physical units.
epoch: Initial time value, w.r.t. the peak (t=0 at the peak) of the total waveform amplitude, as defined in Eq.38 of arxiv:1812.07865.
phi_22 | Output: phase of (2,2) mode. |
evaluated_mode_dps | Output: All other data pieces. |
epoch | Output: Initial time value, where t=0 is at the peak of the total waveform amplitude. |
deltaT | Sampling interval (s). |
fMin | Start GW frequency (Hz). |
fRef | Reference GW frequency (Hz). |
q | Mass ratio m1/m2. |
Mtot_sec | Total mass in geometric units (s). |
chi1z | Dimensionless spin of Bh1. |
chi2z | Dimensionless spin of Bh2. |
ModeArray | Container for (ell, m) modes to generate. |
LALparams | Dict with extra parameters |
Definition at line 230 of file LALSimIMRNRHybSur3dq8.c.
|
static |
Definition at line 73 of file LALSimIMRNRHybSur3dq8.h.