62#define UNUSED __attribute__ ((unused))
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>
75#include <lal/SeqFactories.h>
76#include <lal/FileIO.h>
78#include <lal/H5FileIO.h>
80#include <lal/XLALError.h>
81#include <lal/LALSimIMR.h>
90#ifdef LAL_PTHREAD_LOCK
94#ifdef LAL_PTHREAD_LOCK
95static pthread_once_t NRSur3dq8Remnant_is_initialized = PTHREAD_ONCE_INIT;
142 "CANONICAL_FILE_BASENAME");
164 gsl_vector* fit_params,
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";
176 UINT4 unlim_extrap = 0;
177 if (LALparams != NULL &&
187 if (fabs(chiAz) > 1) {
189 "Invalid spin magnitude |chiA| = %0.4f > 1\n", fabs(chiAz));
191 if (fabs(chiBz) > 1) {
193 "Invalid spin magnitude |chiB| = %0.4f > 1\n", fabs(chiBz));
195 if ((
q > 10.1) && (unlim_extrap == 0)) {
197 "Too much extrapolation in mass ratio; q=%0.4f > 10.1\n%s\n",
q,
202 "Extrapolating outside training range q=%0.4f > 8\n",
q);
204 if (fabs(chiAz) > 0.91 && (unlim_extrap == 0)) {
206 "Too much extrapolation; |chiAz|=%0.4f > 0.91\n%s\n", fabs(chiAz),
209 if (fabs(chiBz) > 0.91 && (unlim_extrap == 0)) {
211 "Too much extrapolation; |chiBz|=%0.4f > 0.91\n%s\n", fabs(chiBz),
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);
221 "Extrapolating outside training range |chiAz|=%0.4f > 0.81\n",
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);
231 "Extrapolating outside training range |chiBz|=%0.4f > 0.81\n",
236 const REAL8 chi_wtAvg = (
q*chiAz+chiBz)/(1.+
q);
238 = (chi_wtAvg - 38.*eta/113.*(chiAz + chiBz))/(1. - 76.*eta/113.);
239 const REAL8 chi_a = (chiAz - chiBz)/2.;
242 "Size of fit_params should be 3, not %zu.\n", fit_params->size);
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);
279 char *remnant_property,
283#ifdef LAL_PTHREAD_LOCK
284 (void) pthread_once(&NRSur3dq8Remnant_is_initialized,
292 if (!sur_data->
setup) {
297 gsl_vector *dummy_worker = gsl_vector_alloc(sur_data->
params_dim);
300 gsl_vector *fit_params = gsl_vector_alloc(sur_data->
params_dim);
308 if (strcmp(remnant_property,
"mf") == 0) {
311 }
else if (strcmp(remnant_property,
"chifz") == 0) {
314 }
else if (strcmp(remnant_property,
"vfx") == 0) {
317 }
else if (strcmp(remnant_property,
"vfy") == 0) {
322 "of 'mf', 'chifz', 'vfx' or 'vfy'");
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(...)
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 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.