LALSimulation  5.4.0.1-fe68b98
LALSimNRSurRemnantUtils.h
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 
30 
31 
32 
33 
34 /*
35  * Custom Warning function because XLAL_PRINT_WARNING is not enabled by
36  * default.
37  */
38 #define print_warning(...) \
39  if (lalDebugLevel & LALERRORBIT) \
40  printf("Warning - %s (%s:%d): ", __func__, __FILE__, __LINE__); \
41  printf(__VA_ARGS__);
42 
43 //*************************************************************************/
44 //************************* struct definitions ****************************/
45 //*************************************************************************/
46 
47 // These GPR fits need the same struct as NRHybSurFitData, but renaming it
48 // here to avoid confusion
49 #define ScalarFitData NRHybSurFitData
50 
51 /**
52  * Vector of fit data, where each component is a ScalarFitData
53  */
54 typedef struct tagVectorFitData {
55  UINT4 vec_dim; /**< Length of the vector. */
56  ScalarFitData **fit_data; /**< One ScalarFitData per component */
58 
59 /**
60  * NRSurRemnant GPR fit data for the mass, spin, and recoil kick for
61  * generically precessing BBHs.
62  *
63  * The final mass is included as a sclar, while the final spin and kick are
64  * included as 3-vectors.
65  */
66 typedef struct tagPrecessingRemnantFitData {
67  UINT4 setup; /**< Indicates if NRSurRemnant has been initialized */
68  UINT4 params_dim; /**< Dimensions of the model. */
69  ScalarFitData *mf_data; /**< Fit data for final mass. */
70  VectorFitData *chif_data; /**< Fit data for final spin. */
71  VectorFitData *vf_data; /**< Fit data for recoil kick. */
72  gsl_matrix *x_train; /**< Training set parameters, needed for GPR fits. */
74 
75 /**
76  * NRSurRemnant GPR fit data for the mass, spin, and recoil kick for
77  * aligned-spin BBHs.
78  *
79  * The final mass, z-component of the final spin, and x and y components of the
80  * recoil kick are included as scalars. The other components of the spin and
81  * kick are zero due to the symmetries of aligned-spin systems and are not
82  * modelled by the surrogate.
83  */
84 typedef struct tagAlignedSpinRemnantFitData {
85  UINT4 setup; /**< Indicates if NRSurRemnant has been initialized */
86  UINT4 params_dim; /**< Dimensions of the model. */
87  ScalarFitData *mf_data; /**< Fit data for final mass. */
88  ScalarFitData *chifz_data; /**< Fit data for final mass. */
89  ScalarFitData *vfx_data; /**< Fit data for final mass. */
90  ScalarFitData *vfy_data; /**< Fit data for final mass. */
91  gsl_matrix *x_train; /**< Training set parameters, needed for GPR fits. */
93 
94 //*************************************************************************/
95 //************************* function declarations *************************/
96 //*************************************************************************/
97 
99  LALH5File **file,
100  const char* NRSurRemnant_DATAFILE
101 );
102 
104  ScalarFitData **fit_data,
105  LALH5File *file,
106  const char *grp_name
107 );
108 
110  VectorFitData **vector_data,
111  UINT4 vec_dim,
112  LALH5File *file,
113  const char *grp_name
114 );
115 
117  PrecessingRemnantFitData *sur_data,
118  LALH5File *file
119 );
120 
122  AlignedSpinRemnantFitData *sur_data,
123  LALH5File *file
124 );
struct tagLALH5File LALH5File
Utilities needed for aligned-spin NR-hybrid surrogate models.
int NRSurRemnant_LoadScalarFit(ScalarFitData **fit_data, LALH5File *file, const char *grp_name)
Loads a single NRSurRemnant GPR fit, as described in the supplementary materials of arxiv:1809....
void NRSurRemnant_LoadH5File(LALH5File **file, const char *NRSurRemnant_DATAFILE)
Loads H5 file for a NRSurRemnant model.
int NRSurRemnant_LoadVectorFit(VectorFitData **vector_data, UINT4 vec_dim, LALH5File *file, const char *grp_name)
Loads a vector of NRSurRemnant GPR fits.
int AlignedSpinNRSurRemnant_Init(AlignedSpinRemnantFitData *sur_data, LALH5File *file)
Initializes fit data for an aligned-spin NRSurRemnant.
#define ScalarFitData
int PrecessingNRSurRemnant_Init(PrecessingRemnantFitData *sur_data, LALH5File *file)
Initializes fit data for a precessing NRSurRemnant.
uint32_t UINT4
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.
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,...
UINT4 vec_dim
Length of the vector.
ScalarFitData ** fit_data
One ScalarFitData per component.