58#define UNUSED __attribute__ ((unused))
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>
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>
84#ifdef LAL_PTHREAD_LOCK
88#ifdef LAL_PTHREAD_LOCK
89static pthread_once_t NRSur7dq4Remnant_is_initialized = PTHREAD_ONCE_INIT;
137 "CANONICAL_FILE_BASENAME");
163 gsl_vector* fit_params,
175 UINT4 unlim_extrap = 0;
176 if (LALparams != NULL &&
184 REAL8 chiAmag = sqrt(chiAx*chiAx + chiAy*chiAy + chiAz*chiAz);
185 REAL8 chiBmag = sqrt(chiBx*chiBx + chiBy*chiBy + chiBz*chiBz);
188 REAL8 q_max_soft = 4.01;
189 REAL8 chi_max_soft = 0.81;
192 REAL8 q_max_hard = 6.01;
197 if ((
q > q_max_hard) && (unlim_extrap == 0)) {
199 "Too much extrapolation in mass ratio; q = %0.4f > %0.4f\n",
202 if (
q > q_max_soft) {
204 "Extrapolating outside training range q = %0.4f > %0.4f\n",
210 "Invalid spin magnitude |chiA| = %0.4f > 1\n", chiAmag);
214 "Invalid spin magnitude |chiB| = %0.4f > 1\n", chiBmag);
216 if (chiAmag > chi_max_soft) {
218 "Extrapolating outside training range |chiA| = %0.4f > %0.4f\n",
219 chiAmag, chi_max_soft);
221 if (chiBmag > chi_max_soft) {
223 "Extrapolating outside training range |chiB| = %0.4f > %0.4f\n",
224 chiBmag, chi_max_soft);
229 const REAL8 chi_wtAvg = (
q*chiAz+chiBz)/(1.+
q);
231 = (chi_wtAvg - 38.*eta/113.*(chiAz + chiBz))/(1. - 76.*eta/113.);
232 const REAL8 chi_a = (chiAz - chiBz)/2.;
235 "Size of fit_params should be 7, not %zu.\n", fit_params->size);
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);
278 char *remnant_property,
282#ifdef LAL_PTHREAD_LOCK
283 (void) pthread_once(&NRSur7dq4Remnant_is_initialized,
291 if (!sur_data->
setup) {
296 gsl_vector *dummy_worker = gsl_vector_alloc(sur_data->
params_dim);
299 gsl_vector *fit_params = gsl_vector_alloc(sur_data->
params_dim);
301 s2x, s2y, s2z, LALparams);
308 if (strcmp(remnant_property,
"mf") == 0) {
312 *result = gsl_vector_alloc(1);
315 gsl_vector_set(*result, 0, tmp);
323 if (strcmp(remnant_property,
"chif") == 0) {
325 }
else if (strcmp(remnant_property,
"vf") == 0) {
329 "of 'mf', 'chif' or 'vf'");
332 *result = gsl_vector_alloc(vec_data->
vec_dim);
335 sur_data->
x_train, dummy_worker);
336 gsl_vector_set(*result,
i, tmp);
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(...)
void XLALH5FileClose(LALH5File UNUSED *file)
#define XLAL_ERROR_VOID(...)
#define XLAL_CHECK(assertion,...)
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.