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
LALSimNRSur7dq4Remnant.c
Go to the documentation of this file.
1/* Copyright (C) 2019 Vijay Varma
2 * Evaluates NRSur7dq4Remnant model for remnant BH mass, spin and recoil kick
3 * for generically precessing 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 NRSur7dq4Remnant model for remnant BH mass, spin and recoil kick for
27 * generically precessing BBH.
28 *
29 * The binary data file (NRSur7dq4Remnant_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/1905.09300,
38 * https://journals.aps.org/prresearch/abstract/10.1103/PhysRevResearch.1.033015
39 *
40 * **Parameter ranges**:
41 *
42 * q = [1, 6.01]
43 *
44 * \f$\chi_{1}, \chi_{2}\f$ = [-1, 1]
45 *
46 * **Training parameter ranges**:
47 *
48 * q = [1, 4.01]
49 *
50 * \f$\chi_{1}, \chi_{2}\f$ = [-0.81, 0.81]
51 *
52 * But extrapolates reasonably to the above mass ratios and spins. However,
53 * if a guarantee of accuracy is required, this model should be used within
54 * the training parameter range.
55 */
56
57#ifdef __GNUC__
58#define UNUSED __attribute__ ((unused))
59#else
60#define UNUSED
61#endif
62
63
64#include <libgen.h>
65
66#include <gsl/gsl_bspline.h>
67#include <gsl/gsl_blas.h>
68#include <gsl/gsl_spline.h>
69#include <gsl/gsl_complex_math.h>
70
71#include <lal/SeqFactories.h>
72#include <lal/FileIO.h>
73#include <lal/XLALError.h>
74#include <lal/LALSimIMR.h>
75#ifdef LAL_HDF5_ENABLED
76#include <lal/H5FileIO.h>
77#endif
78
82
83
84#ifdef LAL_PTHREAD_LOCK
85#include <pthread.h>
86#endif
87
88#ifdef LAL_PTHREAD_LOCK
89static pthread_once_t NRSur7dq4Remnant_is_initialized = PTHREAD_ONCE_INIT;
90#endif
91
92
93/**
94 * Global surrogate data. This data will be loaded at most once. Any
95 * executable which calls NRSur7dq4Remnant_Init_LALDATA directly or by calling
96 * the XLALNRSur7dq4Remnant function will have a memory leak according to
97 * valgrind, because we never free this memory. This memory, however, gets freed
98 * after the executable terminates.
99 */
101
102//*************************************************************************/
103//************************* function definitions **************************/
104//*************************************************************************/
105
106/**
107 * Helper function to check if the NRSur7dq4Remnant model has been initialized.
108 */
109static bool NRSur7dq4Remnant_IsSetup(void) {
111 return true;
112 else
113 return false;
114}
115
116/**
117 * Surrogate initialization.
118 *
119 * This needs to be called once, before __lalsim_NRSur7dq4Remnant_data is used.
120 * It finds the H5 data file with the NRSur7dq4Remnant data and loads the
121 * surrogate. Can be called multiple times, will just return true if surrogate
122 * is already loaded.
123 */
125
126 if (NRSur7dq4Remnant_IsSetup()) return;
127
128 #if LAL_HDF5_ENABLED
129
130 LALH5File *file = NULL;
131 NRSurRemnant_LoadH5File(&file, NRSur7dq4Remnant_DATAFILE);
132
133 int ret = PrecessingNRSurRemnant_Init(&__lalsim_NRSur7dq4Remnant_data,
134 file);
135
136 ret |= ROM_check_canonical_file_basename(file,NRSur7dq4Remnant_DATAFILE,
137 "CANONICAL_FILE_BASENAME");
138
140 if (ret != XLAL_SUCCESS) {
141 XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s\n",
143 }
144
145 #else
146 XLAL_ERROR_VOID(XLAL_FAILURE, "Failure loading data from %s. HDF5 support is not enabled.\n",
148 #endif
149}
150
151/**
152 * Map from mass ratio and spins to surrogate fit parameters.
153 *
154 * The fit parameters are \f$[log_e(q), \chi_{1x}, \chi_{1y}, \hat{\chi},
155 * \chi_{2x}, \chi_{2y} \chi_a]\f$.
156 * \f$\hat{\chi}\f$ is defined in Eq.(8) of arxiv:1905.09300.
157 * \f$\chi_a = (\chi_{1z} - \chi_{2z})/2 \f$.
158 *
159 * The spins must be specified in the coorbital frame at t=-100M from
160 * the total amplitude peak, as described in Sec.IV.C of arxiv:1905.09300.
161 */
163 gsl_vector* fit_params, /**< Output: mapped fit parameters. */
164 const REAL8 q, /**< Mass ratio m1 / m2 >= 1. */
165 const REAL8 chiAx, /**< Dimless x-spin of heavier BH. */
166 const REAL8 chiAy, /**< Dimless y-spin of heavier BH. */
167 const REAL8 chiAz, /**< Dimless z-spin of heavier BH. */
168 const REAL8 chiBx, /**< Dimless x-spin of lighter BH. */
169 const REAL8 chiBy, /**< Dimless y-spin of lighter BH. */
170 const REAL8 chiBz, /**< Dimless z-spin of lighter BH. */
171 LALDict* LALparams /**< Dict with extra parameters */
172) {
173
174 // By default we do not allow unlimited_extrapolation
175 UINT4 unlim_extrap = 0;
176 if (LALparams != NULL &&
177 XLALDictContains(LALparams, "unlimited_extrapolation")) {
178 // Unless the user asks for it
179 unlim_extrap
180 = XLALDictLookupUINT4Value(LALparams, "unlimited_extrapolation");
181 }
182
183 //// Sanity check parameter ranges
184 REAL8 chiAmag = sqrt(chiAx*chiAx + chiAy*chiAy + chiAz*chiAz);
185 REAL8 chiBmag = sqrt(chiBx*chiBx + chiBy*chiBy + chiBz*chiBz);
186 // These are the limits of the training parameter space, so print a
187 // warning when extrapolated beyond.
188 REAL8 q_max_soft = 4.01;
189 REAL8 chi_max_soft = 0.81;
190 // These are the limits beyond which extrapolation is expected to be
191 // very bad, so raise an error.
192 REAL8 q_max_hard = 6.01;
193
194 if (q < 1) {
195 XLAL_ERROR(XLAL_FAILURE, "Invalid mass ratio q = %0.4f < 1\n", q);
196 }
197 if ((q > q_max_hard) && (unlim_extrap == 0)) {
199 "Too much extrapolation in mass ratio; q = %0.4f > %0.4f\n",
200 q, q_max_hard);
201 }
202 if (q > q_max_soft) {
204 "Extrapolating outside training range q = %0.4f > %0.4f\n",
205 q, q_max_soft);
206 }
207
208 if (chiAmag > 1) {
210 "Invalid spin magnitude |chiA| = %0.4f > 1\n", chiAmag);
211 }
212 if (chiBmag > 1) {
214 "Invalid spin magnitude |chiB| = %0.4f > 1\n", chiBmag);
215 }
216 if (chiAmag > chi_max_soft) {
218 "Extrapolating outside training range |chiA| = %0.4f > %0.4f\n",
219 chiAmag, chi_max_soft);
220 }
221 if (chiBmag > chi_max_soft) {
223 "Extrapolating outside training range |chiB| = %0.4f > %0.4f\n",
224 chiBmag, chi_max_soft);
225 }
226
227
228 const REAL8 eta = q/(1.+q)/(1.+q);
229 const REAL8 chi_wtAvg = (q*chiAz+chiBz)/(1.+q);
230 const REAL8 chi_hat
231 = (chi_wtAvg - 38.*eta/113.*(chiAz + chiBz))/(1. - 76.*eta/113.);
232 const REAL8 chi_a = (chiAz - chiBz)/2.;
233
234 XLAL_CHECK((fit_params != NULL) && (fit_params->size == 7), XLAL_EDIMS,
235 "Size of fit_params should be 7, not %zu.\n", fit_params->size);
236
237 gsl_vector_set(fit_params, 0, log(q));
238 gsl_vector_set(fit_params, 1, chiAx);
239 gsl_vector_set(fit_params, 2, chiAy);
240 gsl_vector_set(fit_params, 3, chi_hat);
241 gsl_vector_set(fit_params, 4, chiBx);
242 gsl_vector_set(fit_params, 5, chiBy);
243 gsl_vector_set(fit_params, 6, chi_a);
244
245 return XLAL_SUCCESS;
246}
247
248//************************ main functions ***************************/
249
250/**
251 * Returns the remnant BH's mass, spin, or kick according to NRSur7dq4Remnant
252 * model.
253 *
254 * Reference: arxiv:1905.09300.
255 *
256 * The input spins must be dimensionless, and given in the coorbital frame at
257 * t=-100M from the total amplitude peak, as described in the paper. The
258 * remnant mass is returned in units of total mass, the remnant spin is
259 * dimensionless, and the remnant kick is in units of c. The remnant spin and
260 * kick vectors are also returned in the same frame.
261 *
262 * remnant_property can be one of "mf", "chif" or "vf", respectively, for
263 * the remnant mass, spin or kick.
264 *
265 * Usage examples:
266 * mf = lalsim.NRSur7dq4Remnant(q, s1x, s1y, s1z, s2x, s2y, s2z, "mf").data[0]
267 * chif = lalsim.NRSur7dq4Remnant(q, s1x, s1y, s1z, s2x, s2y, s2z, "chif").data
268 */
270 gsl_vector **result, /**<Output: The requested remnant property. */
271 REAL8 q, /**< Mass ratio of Bh1/Bh2. q>=1. */
272 REAL8 s1x, /**< S1x in coorbital frame at t=-100M */
273 REAL8 s1y, /**< S1y in coorbital frame at t=-100M */
274 REAL8 s1z, /**< S1z in coorbital frame at t=-100M */
275 REAL8 s2x, /**< S2x in coorbital frame at t=-100M */
276 REAL8 s2y, /**< S2y in coorbital frame at t=-100M */
277 REAL8 s2z, /**< S2z in coorbital frame at t=-100M */
278 char *remnant_property, /**< One of "mf", "chif" or "vf" */
279 LALDict* LALparams /**< Dict with extra parameters */
280) {
281
282#ifdef LAL_PTHREAD_LOCK
283 (void) pthread_once(&NRSur7dq4Remnant_is_initialized,
285#else
287#endif
288
289 // Loaded surrogate data
291 if (!sur_data->setup) {
292 XLAL_ERROR(XLAL_EFAILED, "Error loading surrogate data.\n");
293 }
294
295 // assign size to dummy_worker
296 gsl_vector *dummy_worker = gsl_vector_alloc(sur_data->params_dim);
297
298 // Compute fit_params (initialize with dimension of the surrogate)
299 gsl_vector *fit_params = gsl_vector_alloc(sur_data->params_dim);
300 int ret = NRSur7dq4Remnant_fitParams(fit_params, q, s1x, s1y, s1z,
301 s2x, s2y, s2z, LALparams);
302 if(ret != XLAL_SUCCESS) {
303 XLAL_ERROR(XLAL_EFUNC, "Failed to evaluate fit_params.");
304 }
305
306 REAL8 tmp;
307 // For final mass, we have a single scalar fit
308 if (strcmp(remnant_property, "mf") == 0) {
309
310 ScalarFitData *fit_data = sur_data->mf_data;
311
312 *result = gsl_vector_alloc(1);
313 tmp = NRHybSur_eval_fit(fit_data, fit_params, sur_data->x_train,
314 dummy_worker);
315 gsl_vector_set(*result, 0, tmp);
316
317 return XLAL_SUCCESS;
318
319 // For final spin and kick, we have a vector fit
320 } else {
321
322 VectorFitData *vec_data;
323 if (strcmp(remnant_property, "chif") == 0) {
324 vec_data = sur_data->chif_data;
325 } else if (strcmp(remnant_property, "vf") == 0) {
326 vec_data = sur_data->vf_data;
327 } else {
328 XLAL_ERROR(XLAL_EINVAL, "Invalid remnant_property, should be one "
329 "of 'mf', 'chif' or 'vf'");
330 }
331
332 *result = gsl_vector_alloc(vec_data->vec_dim);
333 for (UINT4 i=0; i<vec_data->vec_dim; i++) {
334 tmp = NRHybSur_eval_fit(vec_data->fit_data[i], fit_params,
335 sur_data->x_train, dummy_worker);
336 gsl_vector_set(*result, i, tmp);
337 }
338
339 return XLAL_SUCCESS;
340 }
341}
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 bool NRSur7dq4Remnant_IsSetup(void)
Helper function to check if the NRSur7dq4Remnant model has been initialized.
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.
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.
static PrecessingRemnantFitData __lalsim_NRSur7dq4Remnant_data
Global surrogate data.
static void NRSur7dq4Remnant_Init_LALDATA(void)
Surrogate initialization.
NRSur7dq4Remnant model for remnant BH mass, spin and recoil kick for generically precessing BBH.
static const char NRSur7dq4Remnant_DATAFILE[]
Utils for NR surrogates for remnant BH mass, spin and recoil kick.
#define print_warning(...)
#define ScalarFitData
double i
Definition: bh_ringdown.c:118
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 generically precessing BBHs.
VectorFitData * vf_data
Fit data for recoil kick.
UINT4 setup
Indicates if NRSurRemnant has been initialized.
ScalarFitData * mf_data
Fit data for final mass.
UINT4 params_dim
Dimensions of the model.
VectorFitData * chif_data
Fit data for final spin.
gsl_matrix * x_train
Training set parameters, needed for GPR fits.
Data used in a single vector fit NOTE: basisFunctionOrders, coefs, componentIndices,...
FitData ** fit_data
Vector of FitData.
int vec_dim
Dimension of the vector.