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
LALSimNRSurRemnantUtils.c
Go to the documentation of this file.
1/* Copyright (C) 2019 Vijay Varma
2 * Utils for evaluates NR surrogate remnant fits for the mass, spin and recoil
3 * kick of the final BH left behind in a BBH merger.
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 Utils for NR surrogates for remnant BH mass, spin and recoil kick.
27 */
28
29#ifdef __GNUC__
30#define UNUSED __attribute__ ((unused))
31#else
32#define UNUSED
33#endif
34
36
37#include <libgen.h>
38
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>
43
44#include <lal/SeqFactories.h>
45#include <lal/FileIO.h>
46
47#ifdef LAL_HDF5_ENABLED
48#include <lal/H5FileIO.h>
49#endif
50
51#include <lal/XLALError.h>
52#include <lal/LALConstants.h>
53#include <lal/LALSimIMR.h>
54
57
58//*************************************************************************/
59//************************* function definitions **************************/
60//*************************************************************************/
61
62//********************* Surrogate loading functions ************************/
63// These should only be called once, when initializing the surrogate
64
65#if LAL_HDF5_ENABLED
66/**
67 * Loads H5 file for a NRSurRemnant model.
68 */
69void NRSurRemnant_LoadH5File(
70 LALH5File **file, /**< Output: Returns the opened H5 file. */
71 const char* NRSurRemnant_DATAFILE /**< H5 file name to load. */
72) {
73
74 char *path = XLAL_FILE_RESOLVE_PATH(NRSurRemnant_DATAFILE);
75 if (path==NULL){
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);
84 }
85
86 char *dir = dirname(path);
87 const UINT4 size = strlen(dir) + strlen(NRSurRemnant_DATAFILE) + 2;
88 char *file_path = XLALMalloc(size);
89 snprintf(file_path, size, "%s/%s", dir, NRSurRemnant_DATAFILE);
90
91 *file = XLALH5FileOpen(file_path, "r");
92 if (*file==NULL) {
94 "Unable to load data file %s in $LAL_DATA_PATH."
95 " File may be corrupted.\n", NRSurRemnant_DATAFILE);
96 }
97 XLALFree(path);
98 XLALFree(file_path);
99}
100
101/**
102 * Loads a single NRSurRemnant GPR fit, as described in the supplementary
103 * materials of arxiv:1809.09125.
104 */
105int NRSurRemnant_LoadScalarFit(
106 ScalarFitData **fit_data, /**< Output: Fit data. *fit_data
107 should be NULL. Space will be allocated. */
108 LALH5File *file, /**< Opened H5 file. */
109 const char *grp_name /**< H5 group name. */
110) {
111
112 if (fit_data == NULL || *fit_data != NULL) {
113 XLAL_ERROR(XLAL_EFAULT, "*fit_data should be NULL");
114 }
115 if (file == NULL) {
116 XLAL_ERROR(XLAL_EFAULT, "file should not be NULL");
117 }
118
119 // Open h5 group
120 LALH5File *sub;
121 sub = XLALH5GroupOpen(file, grp_name);
122 *fit_data = XLALMalloc(sizeof(ScalarFitData));
123
124 // Load fit data
125 GPRHyperParams *hyperparams = XLALMalloc(sizeof(GPRHyperParams));
126
127 //// Load scalars needed for fit
128
129 // constant_value = \sigma_k^2 in Eq. S3 of arxiv:1809.09125
130 int ret = ReadHDF5DoubleDataset(&(hyperparams->constant_value),
131 sub, "constant_value");
132 if (ret != XLAL_SUCCESS) {
133 XLAL_ERROR(XLAL_EFUNC, "Failed to load constant_value.");
134 }
135
136 // Mean value before doing GPR fit, usually zero. We already remove the
137 // mean and store it in data_mean.
138 ret = ReadHDF5DoubleDataset(&(hyperparams->y_train_mean),
139 sub, "y_train_mean");
140 if (ret != XLAL_SUCCESS) {
141 XLAL_ERROR(XLAL_EFUNC, "Failed to load y_train_mean.");
142 }
143
144 // Mean of fit data.
145 ret = ReadHDF5DoubleDataset(&((*fit_data)->data_mean),
146 sub, "data_mean");
147 if (ret != XLAL_SUCCESS) {
148 XLAL_ERROR(XLAL_EFUNC, "Failed to load data_mean.");
149 }
150
151 // Standard deviation of fit data.
152 ret = ReadHDF5DoubleDataset(&((*fit_data)->data_std),
153 sub, "data_std");
154 if (ret != XLAL_SUCCESS) {
155 XLAL_ERROR(XLAL_EFUNC, "Failed to load data_std.");
156 }
157
158 // Intercept of linear fit.
159 ret = ReadHDF5DoubleDataset(&((*fit_data)->lin_intercept),
160 sub, "lin_intercept");
161 if (ret != XLAL_SUCCESS) {
162 XLAL_ERROR(XLAL_EFUNC, "Failed to load lin_intercept.");
163 }
164
165
166 //// Load arrays needed for fit
167
168 // Length scale parameters in the kernel.
169 // Array of size=dom, where dim=dimensions of the model.
170 hyperparams->length_scale = NULL;
171 ret = ReadHDF5RealVectorDataset(sub, "length_scale",
172 &(hyperparams->length_scale));
173 if (ret != XLAL_SUCCESS) {
174 XLAL_ERROR(XLAL_EFUNC, "Failed to load length_scale.");
175 }
176
177 // alpha = \f$ K_{x x}^{-1} {\bf f}\f$, which is used in the mean value
178 // in Eq. S2 of arxiv:1809.09125. This can be precomputed as it does not
179 // depend on the x_{*} values.
180 // Array of size=N, where N=size of training dataset.
181 hyperparams->alpha = NULL;
182 ret = ReadHDF5RealVectorDataset(sub, "alpha",
183 &(hyperparams->alpha));
184 if (ret != XLAL_SUCCESS) {
185 XLAL_ERROR(XLAL_EFUNC, "Failed to load alpha.");
186 }
187
188 // Coefs of of linear fit.
189 (*fit_data)->lin_coef = NULL;
190 ret = ReadHDF5RealVectorDataset(sub, "lin_coef", &((*fit_data)->lin_coef));
191 if (ret != XLAL_SUCCESS) {
192 XLAL_ERROR(XLAL_EFUNC, "Failed to load lin_coef.");
193 }
194
195 (*fit_data)->hyperparams = hyperparams;
196
197 XLALH5FileClose(sub);
198 return ret;
199}
200#endif
201
202
203#if LAL_HDF5_ENABLED
204/**
205 * Loads a vector of NRSurRemnant GPR fits
206 */
207int NRSurRemnant_LoadVectorFit(
208 VectorFitData **vector_data, /**< Output: Vector of fit data. *vector_data
209 should be NULL. Space will be allocated. */
210 UINT4 vec_dim, /**< Length of the vector */
211 LALH5File *file, /**< Opened H5 file. */
212 const char *grp_name /**< H5 group name. */
213) {
214
215 if (vector_data == NULL || *vector_data != NULL) {
216 XLAL_ERROR(XLAL_EFAULT, "*vector_data should be NULL");
217 }
218
219 *vector_data = XLALMalloc(sizeof(VectorFitData));
220 (*vector_data)->fit_data = XLALMalloc( vec_dim * sizeof(ScalarFitData *) );
221
222 const size_t str_size = 20;
223 char *sub_grp_name = XLALMalloc(str_size);
224 UNUSED size_t nwritten;
225
226 // For each component load a Scalar fit
227 int ret = XLAL_SUCCESS;
228 for (UINT4 i=0; i<vec_dim; i++){
229
230 nwritten = snprintf(sub_grp_name, str_size, "%s/comp_%u", grp_name, i);
231 XLAL_CHECK_ABORT(nwritten < str_size);
232
233 ScalarFitData *fit_data = NULL;
234 ret = NRSurRemnant_LoadScalarFit(&fit_data, file, sub_grp_name);
235 (*vector_data)->fit_data[i] = fit_data;
236 }
237 (*vector_data)->vec_dim = vec_dim;
238
239 return ret;
240}
241#endif
242
243#ifdef LAL_HDF5_ENABLED
244/**
245 * Initializes fit data for a precessing NRSurRemnant.
246 *
247 * The data includes a ScalarFitData for the final mass, and a VectorFitData
248 * for the final spin and kick 3-vectors.
249 */
250int PrecessingNRSurRemnant_Init(
251 PrecessingRemnantFitData *sur_data, /**< Output: Loaded surrogate data. */
252 LALH5File *file /**< Opened H5 file. */
253) {
254
255 if (sur_data == NULL) {
256 XLAL_ERROR(XLAL_EFAULT, "sur_data should not be NULL");
257 }
258 if (file == NULL) {
259 XLAL_ERROR(XLAL_EFAULT, "file should not be NULL");
260 }
261
262 if (sur_data->setup) {
264 "Model was already initialized. Exiting.");
265 }
266
267 // x_train contains the training parameters and are common for all fits
268 gsl_matrix *x_train = NULL;
269 int ret = ReadHDF5RealMatrixDataset(file, "GPR_X_train", &x_train);
270 if (ret != XLAL_SUCCESS) {
271 XLAL_ERROR(XLAL_EFUNC, "Failed to load GPR_X_train.");
272 }
273 sur_data->x_train = x_train;
274
275 // dimension of the surrogate
276 sur_data->params_dim = x_train->size2;
277
278 // final mass
279 sur_data->mf_data = NULL;
280 ret = NRSurRemnant_LoadScalarFit(&(sur_data->mf_data), file, "mf");
281
282 // final spin vector
283 sur_data->chif_data = NULL;
284 NRSurRemnant_LoadVectorFit(&(sur_data->chif_data), 3, file, "chif");
285
286 // recoil kick vector
287 sur_data->vf_data = NULL;
288 NRSurRemnant_LoadVectorFit(&(sur_data->vf_data), 3, file, "vf");
289
290 if (ret == XLAL_SUCCESS){
291 sur_data->setup = 1;
292 }
293
294 return ret;
295}
296#endif
297
298#if LAL_HDF5_ENABLED
299/**
300 * Initializes fit data for an aligned-spin NRSurRemnant.
301 *
302 * The data includes a ScalarFitData for the final mass, z-component of the
303 * final spin, and x and y components of the recoil kick. The other components
304 * of the spin and kick are zero due to the symmetries of aligned-spin systems
305 * and are not modelled by the surrogate.
306 */
307int AlignedSpinNRSurRemnant_Init(
308 AlignedSpinRemnantFitData *sur_data, /**< Output: Loaded surrogate data. */
309 LALH5File *file /**< Opened H5 file. */
310) {
311
312 if (sur_data == NULL) {
313 XLAL_ERROR(XLAL_EFAULT, "sur_data should not be NULL");
314 }
315 if (file == NULL) {
316 XLAL_ERROR(XLAL_EFAULT, "file should not be NULL");
317 }
318
319 if (sur_data->setup) {
321 "Model was already initialized. Exiting.");
322 }
323
324 // x_train contains the training parameters and are common for all fits
325 gsl_matrix *x_train = NULL;
326 int ret = ReadHDF5RealMatrixDataset(file, "GPR_X_train", &x_train);
327 if (ret != XLAL_SUCCESS) {
328 XLAL_ERROR(XLAL_EFUNC, "Failed to load GPR_X_train.");
329 }
330 sur_data->x_train = x_train;
331
332 // dimension of the surrogate
333 sur_data->params_dim = x_train->size2;
334
335 // final mass
336 sur_data->mf_data = NULL;
337 ret = NRSurRemnant_LoadScalarFit(&(sur_data->mf_data), file, "mf");
338
339 // z-component of final spin (other components are zero)
340 sur_data->chifz_data = NULL;
341 ret = NRSurRemnant_LoadScalarFit(&(sur_data->chifz_data), file, "chifz");
342
343 // x-component of recoil kick
344 sur_data->vfx_data = NULL;
345 ret = NRSurRemnant_LoadScalarFit(&(sur_data->vfx_data), file, "vfx");
346
347 // y-component of recoil kick, z-component is zero
348 sur_data->vfy_data = NULL;
349 ret = NRSurRemnant_LoadScalarFit(&(sur_data->vfy_data), file, "vfy");
350
351 if (ret == XLAL_SUCCESS){
352 sur_data->setup = 1;
353 }
354
355 return ret;
356}
357#endif
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 ScalarFitData
double i
Definition: bh_ringdown.c:118
#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)
uint32_t UINT4
void * XLALMalloc(size_t n)
void XLALFree(void *p)
#define XLAL_ERROR_VOID(...)
#define XLAL_ERROR(...)
#define XLAL_CHECK_ABORT(assertion)
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EIO
XLAL_FAILURE
path
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,...