30#define UNUSED __attribute__ ((unused))
39#include <gsl/gsl_bspline.h>
40#include <gsl/gsl_blas.h>
41#include <gsl/gsl_spline.h>
42#include <gsl/gsl_complex_math.h>
44#include <lal/SeqFactories.h>
45#include <lal/FileIO.h>
47#ifdef LAL_HDF5_ENABLED
48#include <lal/H5FileIO.h>
51#include <lal/XLALError.h>
52#include <lal/LALConstants.h>
53#include <lal/LALSimIMR.h>
69void NRSurRemnant_LoadH5File(
71 const char* NRSurRemnant_DATAFILE
77 "Unable to resolve data file '%s' in $LAL_DATA_PATH.\n"
78 "Note: LALSuite versions >= 7.25 require data files that are publicly available at:\n"
79 "https://git.ligo.org/waveforms/software/lalsuite-waveform-data\n"
80 "and on Zenodo at: https://zenodo.org/records/14999310.\n"
81 "For earlier LALSuite versions, use the files in lalsuite-extra, available at:\n"
82 "https://git.ligo.org/lscsoft/lalsuite-extra\n",
83 NRSurRemnant_DATAFILE);
86 char *dir = dirname(path);
87 const UINT4 size = strlen(dir) + strlen(NRSurRemnant_DATAFILE) + 2;
89 snprintf(file_path, size,
"%s/%s", dir, NRSurRemnant_DATAFILE);
94 "Unable to load data file %s in $LAL_DATA_PATH."
95 " File may be corrupted.\n", NRSurRemnant_DATAFILE);
105int NRSurRemnant_LoadScalarFit(
112 if (fit_data == NULL || *fit_data != NULL) {
131 sub,
"constant_value");
138 ret = ReadHDF5DoubleDataset(&(hyperparams->
y_train_mean),
139 sub,
"y_train_mean");
145 ret = ReadHDF5DoubleDataset(&((*fit_data)->data_mean),
152 ret = ReadHDF5DoubleDataset(&((*fit_data)->data_std),
159 ret = ReadHDF5DoubleDataset(&((*fit_data)->lin_intercept),
160 sub,
"lin_intercept");
171 ret = ReadHDF5RealVectorDataset(sub,
"length_scale",
181 hyperparams->
alpha = NULL;
182 ret = ReadHDF5RealVectorDataset(sub,
"alpha",
183 &(hyperparams->
alpha));
189 (*fit_data)->lin_coef = NULL;
190 ret = ReadHDF5RealVectorDataset(sub,
"lin_coef", &((*fit_data)->lin_coef));
195 (*fit_data)->hyperparams = hyperparams;
207int NRSurRemnant_LoadVectorFit(
215 if (vector_data == NULL || *vector_data != NULL) {
222 const size_t str_size = 20;
224 UNUSED
size_t nwritten;
230 nwritten = snprintf(sub_grp_name, str_size,
"%s/comp_%u", grp_name,
i);
234 ret = NRSurRemnant_LoadScalarFit(&fit_data,
file, sub_grp_name);
235 (*vector_data)->fit_data[
i] = fit_data;
237 (*vector_data)->vec_dim = vec_dim;
243#ifdef LAL_HDF5_ENABLED
250int PrecessingNRSurRemnant_Init(
255 if (sur_data == NULL) {
262 if (sur_data->
setup) {
264 "Model was already initialized. Exiting.");
268 gsl_matrix *x_train = NULL;
269 int ret = ReadHDF5RealMatrixDataset(
file,
"GPR_X_train", &x_train);
280 ret = NRSurRemnant_LoadScalarFit(&(sur_data->
mf_data),
file,
"mf");
284 NRSurRemnant_LoadVectorFit(&(sur_data->
chif_data), 3,
file,
"chif");
288 NRSurRemnant_LoadVectorFit(&(sur_data->
vf_data), 3,
file,
"vf");
307int AlignedSpinNRSurRemnant_Init(
312 if (sur_data == NULL) {
319 if (sur_data->
setup) {
321 "Model was already initialized. Exiting.");
325 gsl_matrix *x_train = NULL;
326 int ret = ReadHDF5RealMatrixDataset(
file,
"GPR_X_train", &x_train);
337 ret = NRSurRemnant_LoadScalarFit(&(sur_data->
mf_data),
file,
"mf");
341 ret = NRSurRemnant_LoadScalarFit(&(sur_data->
chifz_data),
file,
"chifz");
345 ret = NRSurRemnant_LoadScalarFit(&(sur_data->
vfx_data),
file,
"vfx");
349 ret = NRSurRemnant_LoadScalarFit(&(sur_data->
vfy_data),
file,
"vfy");
struct tagLALH5File LALH5File
Auxiliary functions for SEOBNRv1/v2 reduced order modeling codes in LALSimIMRSEOBNRv2ROMDoubleSpin....
Utils for NR surrogates for remnant BH mass, spin and recoil kick.
#define XLAL_FILE_RESOLVE_PATH(fname)
void XLALH5FileClose(LALH5File UNUSED *file)
LALH5File * XLALH5GroupOpen(LALH5File UNUSED *file, const char UNUSED *name)
LALH5File * XLALH5FileOpen(const char UNUSED *path, const char UNUSED *mode)
void * XLALMalloc(size_t n)
#define XLAL_ERROR_VOID(...)
#define XLAL_CHECK_ABORT(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.
Data used in a single GPR fit.
REAL8 constant_value
in kernel.
gsl_vector * alpha
Precomputed .
gsl_vector * length_scale
in kernel.
REAL8 y_train_mean
Mean value before GPR fit, usually zero.
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,...