LALSimulation  5.4.0.1-fe68b98
LALSimNRSur3dq8Remnant.c File Reference

NRSur3dq8Remnant model for remnant BH mass, spin and recoil kick for aligned-spin BBH. More...

Prototypes

static bool NRSur3dq8Remnant_IsSetup (void)
 Helper function to check if the NRSur3dq8Remnant model has been initialized. More...
 
static void NRSur3dq8Remnant_Init_LALDATA (void)
 Surrogate initialization. More...
 
static int NRSur3dq8Remnant_fitParams (gsl_vector *fit_params, const REAL8 q, const REAL8 chiAz, const REAL8 chiBz, LALDict *LALparams)
 Map from mass ratio and spins to surrogate fit parameters. More...
 
int XLALNRSur3dq8Remnant (REAL8 *result, REAL8 q, REAL8 s1z, REAL8 s2z, char *remnant_property, LALDict *LALparams)
 Returns the remnant BH's mass, spin, or kick according to NRSur3dq8Remnant model. More...
 

Detailed Description

NRSur3dq8Remnant model for remnant BH mass, spin and recoil kick for aligned-spin BBH.

Author
Vijay Varma

The binary data file is available at: https://dcc.ligo.org/LIGO-T1900034/public. Get the lalsuite-extra repo or put the data into a location in your LAL_DATA_PATH.

Paper: https://arxiv.org/abs/1809.09125, https://journals.aps.org/prl/abstract/10.1103/PhysRevLett.122.011101. The model is referred to as surfinBH3dq8 in the paper.

Parameter ranges:

q = [1, 9.1] \(\chi_{1z}, \chi_{2z}\) = [-0.91, 0.91]

OR

q = [1, 10.1] \(\chi_{1z}, \chi_{2z}\) = [-0.81, 0.81]

Training parameter ranges:

q = [1, 8] \(\chi_{1z}, \chi_{2z}\) = [-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 LALSimNRSur3dq8Remnant.c.

Go to the source code of this file.

Variables

static AlignedSpinRemnantFitData __lalsim_NRSur3dq8Remnant_data
 Global surrogate data. More...
 

Function Documentation

◆ NRSur3dq8Remnant_IsSetup()

static bool NRSur3dq8Remnant_IsSetup ( void  )
static

Helper function to check if the NRSur3dq8Remnant model has been initialized.

Definition at line 110 of file LALSimNRSur3dq8Remnant.c.

◆ NRSur3dq8Remnant_Init_LALDATA()

static void NRSur3dq8Remnant_Init_LALDATA ( void  )
static

Surrogate initialization.

This needs to be called once, before __lalsim_NRSur3dq8Remnant_data is used. It finds the H5 data file with the NRSur3dq8Remnant data and loads the surrogate. Can be called multiple times, will just return true if surrogate is already loaded.

Definition at line 125 of file LALSimNRSur3dq8Remnant.c.

◆ NRSur3dq8Remnant_fitParams()

static int NRSur3dq8Remnant_fitParams ( gsl_vector *  fit_params,
const REAL8  q,
const REAL8  chiAz,
const REAL8  chiBz,
LALDict *  LALparams 
)
static

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.(S5) of arxiv:1809.09125. \(\chi_a = (\chi_{1z} - \chi_{2z})/2 \).

Parameters
fit_paramsOutput: mapped fit parameters.
qMass ratio m1 / m2 >= 1.
chiAzDimless z-spin of heavier BH.
chiBzDimless z-spin of lighter BH.
LALparamsDict with extra parameters

Definition at line 149 of file LALSimNRSur3dq8Remnant.c.

◆ XLALNRSur3dq8Remnant()

int XLALNRSur3dq8Remnant ( REAL8 result,
REAL8  q,
REAL8  s1z,
REAL8  s2z,
char remnant_property,
LALDict *  LALparams 
)

Returns the remnant BH's mass, spin, or kick according to NRSur3dq8Remnant model.

Reference: arxiv:1809.09125. This model is referred to as surfinBH3dq8 in the paper.

The remnant mass is returned in units of total mass, the remnant spin is dimensionless, and the recoil kick is in units of c.

remnant_property can be one of "mf", "chifz", "vfx" or "vfy", respectively, for the remnant mass, the z-component of remnant spin, and the x and y components of the kick. The kick components are given in the coorbital frame at t=-100M from the total amplitude peak, as described in the paper. All other components of the spin and kick are zero due to the symmetries of aligned-spin systems and are not modelled by the surrogate.

Usage examples: mf = lalsim.NRSur3dq8Remnant(q, s1z, s2z, "mf") chifz = lalsim.NRSur3dq8Remnant(q, s1z, s2z, "chifz")

Parameters
resultOutput: The requested remnant property.
qMass ratio of Bh1/Bh2. q>=1.
s1zS1z z-spin of Bh1
s2zS2z z-spin of Bh2
remnant_propertyOne of "mf", "chifz", "vfx" or "vfy"
LALparamsDict with extra parameters

Definition at line 260 of file LALSimNRSur3dq8Remnant.c.

Variable Documentation

◆ __lalsim_NRSur3dq8Remnant_data

AlignedSpinRemnantFitData __lalsim_NRSur3dq8Remnant_data
static

Global surrogate data.

This data will be loaded at most once. Any executable which calls NRSur3dq8Remnant_Init_LALDATA directly or by calling the XLALNRSur3dq8Remnant 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 101 of file LALSimNRSur3dq8Remnant.c.