NRSur7dq4Remnant model for remnant BH mass, spin and recoil kick for generically precessing BBH. More...
Prototypes | |
static bool | NRSur7dq4Remnant_IsSetup (void) |
Helper function to check if the NRSur7dq4Remnant model has been initialized. More... | |
static void | NRSur7dq4Remnant_Init_LALDATA (void) |
Surrogate initialization. More... | |
static int | NRSur7dq4Remnant_fitParams (gsl_vector *fit_params, const REAL8 q, const REAL8 chiAx, const REAL8 chiAy, const REAL8 chiAz, const REAL8 chiBx, const REAL8 chiBy, const REAL8 chiBz, LALDict *LALparams) |
Map from mass ratio and spins to surrogate fit parameters. More... | |
int | XLALNRSur7dq4Remnant (gsl_vector **result, REAL8 q, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, char *remnant_property, LALDict *LALparams) |
Returns the remnant BH's mass, spin, or kick according to NRSur7dq4Remnant model. More... | |
NRSur7dq4Remnant model for remnant BH mass, spin and recoil kick for generically precessing BBH.
The binary data file is available at: https://dcc.ligo.org/LIGO-T1900393/public. Get the lalsuite-extra repo or put the data into a location in your LAL_DATA_PATH.
Paper: https://arxiv.org/abs/1905.09300, https://journals.aps.org/prresearch/abstract/10.1103/PhysRevResearch.1.033015
Parameter ranges:
q = [1, 6.01]
\(\chi_{1}, \chi_{2}\) = [-1, 1]
Training parameter ranges:
q = [1, 4.01]
\(\chi_{1}, \chi_{2}\) = [-0.81, 0.81]
But extrapolates reasonably to the above mass ratios and spins. However, if a guarantee of accuracy is required, this model should be used within the training parameter range.
Definition in file LALSimNRSur7dq4Remnant.c.
Go to the source code of this file.
Variables | |
static PrecessingRemnantFitData | __lalsim_NRSur7dq4Remnant_data |
Global surrogate data. More... | |
|
static |
Helper function to check if the NRSur7dq4Remnant model has been initialized.
Definition at line 104 of file LALSimNRSur7dq4Remnant.c.
|
static |
Surrogate initialization.
This needs to be called once, before __lalsim_NRSur7dq4Remnant_data is used. It finds the H5 data file with the NRSur7dq4Remnant data and loads the surrogate. Can be called multiple times, will just return true if surrogate is already loaded.
Definition at line 119 of file LALSimNRSur7dq4Remnant.c.
|
static |
Map from mass ratio and spins to surrogate fit parameters.
The fit parameters are \([log_e(q), \chi_{1x}, \chi_{1y}, \hat{\chi}, \chi_{2x}, \chi_{2y} \chi_a]\). \(\hat{\chi}\) is defined in Eq.(8) of arxiv:1905.09300. \(\chi_a = (\chi_{1z} - \chi_{2z})/2 \).
The spins must be specified in the coorbital frame at t=-100M from the total amplitude peak, as described in Sec.IV.C of arxiv:1905.09300.
fit_params | Output: mapped fit parameters. |
q | Mass ratio m1 / m2 >= 1. |
chiAx | Dimless x-spin of heavier BH. |
chiAy | Dimless y-spin of heavier BH. |
chiAz | Dimless z-spin of heavier BH. |
chiBx | Dimless x-spin of lighter BH. |
chiBy | Dimless y-spin of lighter BH. |
chiBz | Dimless z-spin of lighter BH. |
LALparams | Dict with extra parameters |
Definition at line 147 of file LALSimNRSur7dq4Remnant.c.
int XLALNRSur7dq4Remnant | ( | gsl_vector ** | result, |
REAL8 | q, | ||
REAL8 | s1x, | ||
REAL8 | s1y, | ||
REAL8 | s1z, | ||
REAL8 | s2x, | ||
REAL8 | s2y, | ||
REAL8 | s2z, | ||
char * | remnant_property, | ||
LALDict * | LALparams | ||
) |
Returns the remnant BH's mass, spin, or kick according to NRSur7dq4Remnant model.
Reference: arxiv:1905.09300.
The input spins must be dimensionless, and given in the coorbital frame at t=-100M from the total amplitude peak, as described in the paper. The remnant mass is returned in units of total mass, the remnant spin is dimensionless, and the remnant kick is in units of c. The remnant spin and kick vectors are also returned in the same frame.
remnant_property can be one of "mf", "chif" or "vf", respectively, for the remnant mass, spin or kick.
Usage examples: mf = lalsim.NRSur7dq4Remnant(q, s1x, s1y, s1z, s2x, s2y, s2z, "mf").data[0] chif = lalsim.NRSur7dq4Remnant(q, s1x, s1y, s1z, s2x, s2y, s2z, "chif").data
result | Output: The requested remnant property. |
q | Mass ratio of Bh1/Bh2. q>=1. |
s1x | S1x in coorbital frame at t=-100M |
s1y | S1y in coorbital frame at t=-100M |
s1z | S1z in coorbital frame at t=-100M |
s2x | S2x in coorbital frame at t=-100M |
s2y | S2y in coorbital frame at t=-100M |
s2z | S2z in coorbital frame at t=-100M |
remnant_property | One of "mf", "chif" or "vf" |
LALparams | Dict with extra parameters |
Definition at line 254 of file LALSimNRSur7dq4Remnant.c.
|
static |
Global surrogate data.
This data will be loaded at most once. Any executable which calls NRSur7dq4Remnant_Init_LALDATA directly or by calling the XLALNRSur7dq4Remnant function will have a memory leak according to valgrind, because we never free this memory. This memory, however, gets freed after the executable terminates.
Definition at line 95 of file LALSimNRSur7dq4Remnant.c.