Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimNRSur3dq8Remnant.c
Go to the documentation of this file.
1/* Copyright (C) 2019 Vijay Varma
2 * Evaluates NRSur3dq8Remnant model for remnant BH mass, spin and recoil kick
3 * for aligned-spin BBH.
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with with program; see the file COPYING. If not, write to the
17 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
18 * MA 02110-1301 USA
19 */
20
21/**
22 * \author Vijay Varma
23 *
24 * \file
25 *
26 * \brief NRSur3dq8Remnant model for remnant BH mass, spin and recoil kick for
27 * aligned-spin BBH.
28 *
29 * The binary data file (NRSur3dq8Remnant_v1.0.h5) is available at:
30 * https://git.ligo.org/waveforms/software/lalsuite-waveform-data
31 * and on Zenodo at: https://zenodo.org/records/14999310.
32 * Get the lalsuite-waveform-data repo or put the data into a location in your
33 * LAL_DATA_PATH.
34 * The data is also available on CIT at /home/lalsimulation_data and via CVMFS
35 * at /cvmfs/shared.storage.igwn.org/igwn/shared/auxiliary/obs_sci/cbc/waveform/lalsimulation_data
36 *
37 * **Paper**: https://arxiv.org/abs/1809.09125,
38 * https://journals.aps.org/prl/abstract/10.1103/PhysRevLett.122.011101.
39 * The model is referred to as surfinBH3dq8 in the paper.
40 *
41 * **Parameter ranges**:
42 *
43 * q = [1, 9.1]
44 * \f$\chi_{1z}, \chi_{2z}\f$ = [-0.91, 0.91]
45 *
46 * OR
47 *
48 * q = [1, 10.1]
49 * \f$\chi_{1z}, \chi_{2z}\f$ = [-0.81, 0.81]
50 *
51 * **Training parameter ranges**:
52 *
53 * q = [1, 8]
54 * \f$\chi_{1z}, \chi_{2z}\f$ = [-0.81, 0.81]
55 *
56 * But extrapolates reasonably to the above mass ratios and spins. However,
57 * if a guarantee of accuracy is required, this model should be used within
58 * the training parameter range.
59 */
60
61#ifdef __GNUC__
62#define UNUSED __attribute__ ((unused))
63#else
64#define UNUSED
65#endif
66
67
68#include <libgen.h>
69
70#include <gsl/gsl_bspline.h>
71#include <gsl/gsl_blas.h>
72#include <gsl/gsl_spline.h>
73#include <gsl/gsl_complex_math.h>
74
75#include <lal/SeqFactories.h>
76#include <lal/FileIO.h>
77#if LAL_HDF5_ENABLED
78#include <lal/H5FileIO.h>
79#endif
80#include <lal/XLALError.h>
81#include <lal/LALSimIMR.h>
82
85
87
88
89
90#ifdef LAL_PTHREAD_LOCK
91#include <pthread.h>
92#endif
93
94#ifdef LAL_PTHREAD_LOCK
95static pthread_once_t NRSur3dq8Remnant_is_initialized = PTHREAD_ONCE_INIT;
96#endif
97
98
99/**
100 * Global surrogate data. This data will be loaded at most once. Any
101 * executable which calls NRSur3dq8Remnant_Init_LALDATA directly or by calling
102 * the XLALNRSur3dq8Remnant function will have a memory leak according to
103 * valgrind, because we never free this memory. This memory, however, gets freed
104 * after the executable terminates.
105 */
107
108//*************************************************************************/
109//************************* function definitions **************************/
110//*************************************************************************/
111
112/**
113 * Helper function to check if the NRSur3dq8Remnant model has been initialized.
114 */
115static bool NRSur3dq8Remnant_IsSetup(void) {
117 return true;
118 else
119 return false;
120}
121
122/**
123 * Surrogate initialization.
124 *
125 * This needs to be called once, before __lalsim_NRSur3dq8Remnant_data is used.
126 * It finds the H5 data file with the NRSur3dq8Remnant data and loads the
127 * surrogate. Can be called multiple times, will just return true if surrogate
128 * is already loaded.
129 */
131
132 if (NRSur3dq8Remnant_IsSetup()) return;
133
134 #if LAL_HDF5_ENABLED
135 LALH5File *file = NULL;
136 NRSurRemnant_LoadH5File(&file, NRSur3dq8Remnant_DATAFILE);
137
138 int ret = AlignedSpinNRSurRemnant_Init(&__lalsim_NRSur3dq8Remnant_data,
139 file);
140
141 ret |= ROM_check_canonical_file_basename(file,NRSur3dq8Remnant_DATAFILE,
142 "CANONICAL_FILE_BASENAME");
143
145 if (ret != XLAL_SUCCESS) {
146 XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s\n",
148 }
149 #else
150 XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s. HDF5 support is not enabled.\n",
152 #endif
153
154}
155
156/**
157 * Map from mass ratio and spins to surrogate fit parameters.
158 *
159 * The fit parameters are \f$[log_e(q), \hat{\chi}, \chi_a]\f$.
160 * \f$\hat{\chi}\f$ is defined in Eq.(S5) of arxiv:1809.09125.
161 * \f$\chi_a = (\chi_{1z} - \chi_{2z})/2 \f$.
162 */
164 gsl_vector* fit_params, /**< Output: mapped fit parameters. */
165 const REAL8 q, /**< Mass ratio m1 / m2 >= 1. */
166 const REAL8 chiAz, /**< Dimless z-spin of heavier BH. */
167 const REAL8 chiBz, /**< Dimless z-spin of lighter BH. */
168 LALDict* LALparams /**< Dict with extra parameters */
169) {
170
171 //// Sanity check parameter ranges
172 const char *param_validity = "This model is valid for q <= 9.1 & "
173 "|chiAz,chiBz| <= 0.91, or q <= 10.1 & |chiAz,chiBz| <= 0.81";
174
175 // By default we do not allow unlimited_extrapolation
176 UINT4 unlim_extrap = 0;
177 if (LALparams != NULL &&
178 XLALDictContains(LALparams, "unlimited_extrapolation")) {
179 // Unless the user asks for it
180 unlim_extrap
181 = XLALDictLookupUINT4Value(LALparams, "unlimited_extrapolation");
182 }
183
184 if (q < 1) {
185 XLAL_ERROR(XLAL_FAILURE, "Invalid mass ratio q = %0.4f < 1\n", q);
186 }
187 if (fabs(chiAz) > 1) {
189 "Invalid spin magnitude |chiA| = %0.4f > 1\n", fabs(chiAz));
190 }
191 if (fabs(chiBz) > 1) {
193 "Invalid spin magnitude |chiB| = %0.4f > 1\n", fabs(chiBz));
194 }
195 if ((q > 10.1) && (unlim_extrap == 0)) {
197 "Too much extrapolation in mass ratio; q=%0.4f > 10.1\n%s\n", q,
198 param_validity);
199 }
200 if (q > 8) {
202 "Extrapolating outside training range q=%0.4f > 8\n", q);
203 }
204 if (fabs(chiAz) > 0.91 && (unlim_extrap == 0)) {
206 "Too much extrapolation; |chiAz|=%0.4f > 0.91\n%s\n", fabs(chiAz),
207 param_validity);
208 }
209 if (fabs(chiBz) > 0.91 && (unlim_extrap == 0)) {
211 "Too much extrapolation; |chiBz|=%0.4f > 0.91\n%s\n", fabs(chiBz),
212 param_validity);
213 }
214 if (fabs(chiAz) > 0.81) {
215 if ((q > 9.1) && (unlim_extrap == 0)) {
217 "Too much extrapolation; q=%0.4f > 9.1 & |chiAz|=%.04f"
218 " >0.81\n%s\n", q, fabs(chiAz), param_validity);
219 }
221 "Extrapolating outside training range |chiAz|=%0.4f > 0.81\n",
222 fabs(chiAz));
223 }
224 if (fabs(chiBz) > 0.81) {
225 if ((q > 9.1) && (unlim_extrap == 0)) {
227 "Too much extrapolation; q=%0.4f > 9.1 & |chiBz|=%.04f"
228 " >0.81\n%s\n", q, fabs(chiBz), param_validity);
229 }
231 "Extrapolating outside training range |chiBz|=%0.4f > 0.81\n",
232 fabs(chiBz));
233 }
234
235 const REAL8 eta = q/(1.+q)/(1.+q);
236 const REAL8 chi_wtAvg = (q*chiAz+chiBz)/(1.+q);
237 const REAL8 chi_hat
238 = (chi_wtAvg - 38.*eta/113.*(chiAz + chiBz))/(1. - 76.*eta/113.);
239 const REAL8 chi_a = (chiAz - chiBz)/2.;
240
241 XLAL_CHECK((fit_params != NULL) && (fit_params->size == 3), XLAL_EDIMS,
242 "Size of fit_params should be 3, not %zu.\n", fit_params->size);
243
244 gsl_vector_set(fit_params, 0, log(q));
245 gsl_vector_set(fit_params, 1, chi_hat);
246 gsl_vector_set(fit_params, 2, chi_a);
247
248 return XLAL_SUCCESS;
249}
250
251//************************ main functions ***************************/
252
253/**
254 * Returns the remnant BH's mass, spin, or kick according to NRSur3dq8Remnant
255 * model.
256 *
257 * Reference: arxiv:1809.09125. This model is referred to as surfinBH3dq8 in
258 * the paper.
259 *
260 * The remnant mass is returned in units of total mass, the remnant spin is
261 * dimensionless, and the recoil kick is in units of c.
262 *
263 * remnant_property can be one of "mf", "chifz", "vfx" or "vfy", respectively,
264 * for the remnant mass, the z-component of remnant spin, and the x and y
265 * components of the kick. The kick components are given in the coorbital frame
266 * at t=-100M from the total amplitude peak, as described in the paper. All
267 * other components of the spin and kick are zero due to the symmetries of
268 * aligned-spin systems and are not modelled by the surrogate.
269 *
270 * Usage examples:
271 * mf = lalsim.NRSur3dq8Remnant(q, s1z, s2z, "mf")
272 * chifz = lalsim.NRSur3dq8Remnant(q, s1z, s2z, "chifz")
273 */
275 REAL8 *result, /**<Output: The requested remnant property. */
276 REAL8 q, /**< Mass ratio of Bh1/Bh2. q>=1. */
277 REAL8 s1z, /**< S1z z-spin of Bh1 */
278 REAL8 s2z, /**< S2z z-spin of Bh2 */
279 char *remnant_property, /**< One of "mf", "chifz", "vfx" or "vfy" */
280 LALDict* LALparams /**< Dict with extra parameters */
281) {
282
283#ifdef LAL_PTHREAD_LOCK
284 (void) pthread_once(&NRSur3dq8Remnant_is_initialized,
286#else
288#endif
289
290 // Loaded surrogate data
292 if (!sur_data->setup) {
293 XLAL_ERROR(XLAL_EFAILED, "Error loading surrogate data.\n");
294 }
295
296 // assign size to dummy_worker
297 gsl_vector *dummy_worker = gsl_vector_alloc(sur_data->params_dim);
298
299 // Compute fit_params (initialize with dimension of the surrogate)
300 gsl_vector *fit_params = gsl_vector_alloc(sur_data->params_dim);
301 int ret = NRSur3dq8Remnant_fitParams(fit_params, q, s1z, s2z, LALparams);
302 if(ret != XLAL_SUCCESS) {
303 XLAL_ERROR(XLAL_EFUNC, "Failed to evaluate fit_params.");
304 }
305
306 // All quantities are treated as scalars
307 ScalarFitData *fit_data = NULL;
308 if (strcmp(remnant_property, "mf") == 0) {
309 // final mass
310 fit_data = sur_data->mf_data;
311 } else if (strcmp(remnant_property, "chifz") == 0) {
312 // z-component of final spin (other components are zero)
313 fit_data = sur_data->chifz_data;
314 } else if (strcmp(remnant_property, "vfx") == 0) {
315 // x-component of recoil kick
316 fit_data = sur_data->vfx_data;
317 } else if (strcmp(remnant_property, "vfy") == 0) {
318 // y-component of recoil kick, z-component is zero
319 fit_data = sur_data->vfy_data;
320 } else {
321 XLAL_ERROR(XLAL_EINVAL, "Invalid remnant_property, should be one "
322 "of 'mf', 'chifz', 'vfx' or 'vfy'");
323 }
324
325 *result = NRHybSur_eval_fit(fit_data, fit_params, sur_data->x_train,
326 dummy_worker);
327 return XLAL_SUCCESS;
328
329}
struct tagLALH5File LALH5File
int XLALDictContains(const LALDict *dict, const char *key)
UINT4 XLALDictLookupUINT4Value(const LALDict *dict, const char *key)
Auxiliary functions related to HDF5 waveform data files.
REAL8 NRHybSur_eval_fit(const NRHybSurFitData *fit_data, const gsl_vector *fit_params, const gsl_matrix *x_train, gsl_vector *dummy_worker)
Evaluate a NRHybSur fit.
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.
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.
static void NRSur3dq8Remnant_Init_LALDATA(void)
Surrogate initialization.
static bool NRSur3dq8Remnant_IsSetup(void)
Helper function to check if the NRSur3dq8Remnant model has been initialized.
static AlignedSpinRemnantFitData __lalsim_NRSur3dq8Remnant_data
Global surrogate data.
NRSur3dq8Remnant model for remnant BH mass, spin and recoil kick for aligned-spin BBH.
static const char NRSur3dq8Remnant_DATAFILE[]
Utils for NR surrogates for remnant BH mass, spin and recoil kick.
#define print_warning(...)
#define ScalarFitData
void XLALH5FileClose(LALH5File UNUSED *file)
double REAL8
uint32_t UINT4
static const INT4 q
#define XLAL_ERROR_VOID(...)
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
XLAL_SUCCESS
XLAL_EFUNC
XLAL_EDIMS
XLAL_EINVAL
XLAL_EFAILED
XLAL_FAILURE
NRSurRemnant GPR fit data for the mass, spin, and recoil kick for aligned-spin BBHs.
ScalarFitData * mf_data
Fit data for final mass.
UINT4 params_dim
Dimensions of the model.
gsl_matrix * x_train
Training set parameters, needed for GPR fits.
ScalarFitData * vfy_data
Fit data for final mass.
ScalarFitData * vfx_data
Fit data for final mass.
ScalarFitData * chifz_data
Fit data for final mass.
UINT4 setup
Indicates if NRSurRemnant has been initialized.