LALSimulation  5.4.0.1-fe68b98
LALSimNRHybSurUtilities.h
Go to the documentation of this file.
1 /* Copyright (C) 2018 Vijay Varma
2  * Utility functions that are useful for evaluating NRHybrid surrogates.
3  *
4  * This program is free software; you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation; either version 2 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with with program; see the file COPYING. If not, write to the
16  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17  * MA 02110-1301 USA
18  */
19 
20 /**
21  * \author Vijay Varma
22  *
23  * \file
24  *
25  * \brief Utilities needed for aligned-spin NR-hybrid surrogate models.
26  *
27  * Called from:
28  * LALSimIMRNRHybSur3dq8.h
29  */
30 
31 #include <gsl/gsl_vector.h>
32 #include <gsl/gsl_matrix.h>
33 
34 #include <lal/LALSimIMR.h>
35 
36 
37 //*************************************************************************/
38 //************************* struct definitions ****************************/
39 //*************************************************************************/
40 
41 /**
42  * Data used in a single GPR fit.
43  *
44  * We follow sklearn closely. We use the kernel given in
45  * Eq.(S3) of arxiv:1809.09125 :
46  * \f[
47  * K(x, x') = \sigma_k^2 exp(-1/2 \sum_i^D (x^{i} - x'^{i})^2/\sigma_i^2)
48  * \f]
49  * where D is the dimension of the model.
50  *
51  * Note that Eq.(S3) also includes the WhiteKernel, which has a noise
52  * parameter, but we don't need that here since we only need to evaluate \f$
53  * K_{x* x} \f$ when evaluating the fits (See the mean value in Eq.(S2) of same
54  * paper). The other term we need is alpha = \f$ K_{x x}^{-1} {\bf f}\f$, which
55  * involves the WhiteKernel, but is precomputed offline. alpha is a vector of
56  * size N, where N is the number of cases in the training data set.
57  */
58 typedef struct tagGPRHyperParams {
59  REAL8 constant_value; /**< \f$ \sigma_k^2 \f$ in kernel. */
60  REAL8 y_train_mean; /**< Mean value before GPR fit, usually zero. */
61  gsl_vector *length_scale; /**< \f$ \sigma_i \f$ in kernel. */
62  gsl_vector *alpha; /**< Precomputed \f$ K_{x x}^{-1} {\bf f}\f$. */
64 
65 /**
66  * Data used in a single NRHybSur fit.
67  *
68  * Described in supplemental materials of arxiv:1809.09125.
69  * We first subtract a linear fit to the data. Then we subtract the mean of
70  * the resulting values. Then we normalize by dividing by the standard
71  * deviation of the resulting values. Finally, we construct a GPR fit for the
72  * resulting values. Below, D = dimension of the model.
73  */
74 typedef struct tagNRHybSurFitData {
75  REAL8 data_mean; /**< Mean of data after linear fit. */
76  REAL8 data_std; /**< Standard deviation after removing mean. */
77  REAL8 lin_intercept; /**< Constant intercept from linear fit. */
78  gsl_vector *lin_coef; /**< Linear coefs from linear fit, size=D. */
79  GPRHyperParams *hyperparams; /**< Hyperparameters from GPR fit. */
81 
82 /**
83  * NRHybSur data for a single waveform data piece.
84  */
85 typedef struct tagDataPiece {
86  int n_nodes; /**< Number of empirical nodes. */
87  NRHybSurFitData **fit_data; /**< NRHybSurFitData at each empirical node. */
88  gsl_matrix *ei_basis; /**< Empirical interpolation matrix. */
89 } DataPiece;
90 
91 /**
92  * NRHybSur data pieces of a single mode.
93  *
94  * The different data pieces are described in Sec.VI.B of arxiv:1812.07865.
95  * For (2,2) mode we model the amplitude \f$ A_{22} \f$ and phase residual
96  * \f$ \phi^{res}_{22} \f$ defined in Eq.(44) of arxiv:1812.07865.
97  * For all other modes we model the real and imaginary parts of the coorbital
98  * frame strain, defined in Eq.(39) of arxiv:1812.07865.
99  * Only the required data pieces for a given mode will be loaded.
100  */
101 typedef struct tagModeDataPieces {
102  UINT4 ell; /**< Mode \f$ \ell \f$ value. */
103  UINT4 m; /**< Mode m value. */
104  DataPiece *ampl_data_piece; /**< Amplitude data piece. */
105  DataPiece *phase_res_data_piece; /**< Phase residual data piece. */
106  DataPiece *coorb_re_data_piece;/**< Coorbital frame real part data piece. */
107  DataPiece *coorb_im_data_piece;/**< Coorbital frame imag part data piece. */
109 
110 /**
111  * NRHybSur surrogate data for all modes, to be loaded from a h5 file.
112  */
113 typedef struct tagNRHybSurData {
114  UINT4 setup; /**< Indicates if NRHybSur has been initialized */
115  UINT4 num_modes_modeled; /**< Number of modeled modes. */
116  int phaseAlignIdx; /**< Alignment index for phase data piece. */
117  REAL8 params_dim; /**< Dimensions of the model. */
118  gsl_vector *domain; /**< Common time samples for the surrogate. */
119  gsl_vector *TaylorT3_factor_without_eta; /**< Precomputed quantity for use
120  in evaluating the 0PN TaylorT3 phase contribution. */
121  gsl_matrix_long *mode_list; /**< List of modeled modes, first element is
122  always (2,2). */
123  gsl_matrix *x_train; /**< Training set parameters, needed for GPR fits. */
124  ModeDataPieces **mode_data_pieces; /**< Data pieces of all modes, same
125  order as mode_list. */
126 } NRHybSurData;
127 
128 /**
129  * NRHybSur evaluated data for a single mode
130  *
131  * For the (2,2) mode only the amplitude is stored in this struct. The phase
132  * is computed separately as it is needed for all modes to transform from the
133  * coorbital frame to the inertial frame.
134  * For all other modes the real and imaginary parts of the strain in the
135  * coorbital frame are stored.
136  * Only the required data pieces for a given mode will be evaluated.
137  */
138 typedef struct tagEvaluatedDataPieces {
139  UINT4 ell; /**< Mode \f$ \ell \f$ value. */
140  UINT4 m; /**< Mode m value. */
141  gsl_vector *ampl_eval; /**< Amplitude data piece evaluation. */
142  gsl_vector *coorb_re_eval; /**< Coorbital frame real part evaluation. */
143  gsl_vector *coorb_im_eval; /**< Coorbital frame imag part evaluation. */
145 
146 
147 //*************************************************************************/
148 //************************* function declarations *************************/
149 //*************************************************************************/
150 
152  REAL8 *param,
153  LALH5File *sub,
154  const char *name
155 );
156 
158  int *param,
159  LALH5File *sub,
160  const char *name
161 );
162 
163 int NRHybSur_Init(
165  LALH5File *file
166  );
167 
169  const NRHybSurFitData *fit_data,
170  const gsl_vector *fit_params,
171  const gsl_matrix *x_train,
172  gsl_vector *dummy_worker
173 );
174 
176  gsl_vector **phi_22,
177  gsl_vector **output_times,
178  const REAL8 eta,
179  const gsl_vector *fit_params,
180  const REAL8 omegaM_22_min,
181  const REAL8 deltaTOverM,
182  const REAL8 phiRef,
183  const REAL8 omegaM_22_ref,
184  gsl_vector *dummy_dp,
185  const gsl_matrix *x_train,
186  gsl_vector *dummy_worker,
187  const NRHybSurData *NR_hybsur_data
188 );
189 
191  EvaluatedDataPieces **this_mode_eval_dp,
192  const UINT4 ell,
193  const UINT4 m,
194  const ModeDataPieces *data_pieces,
195  const gsl_vector *output_times,
196  const gsl_vector *fit_params,
197  gsl_vector *dummy_dp,
198  const gsl_matrix *x_train,
199  gsl_vector *dummy_worker,
200  const NRHybSurData *NR_hybsur_data
201 );
202 
204  gsl_vector *phi_22,
205  EvaluatedDataPieces **evaluated_mode_dps,
206  const UINT4 num_modes_incl
207 );
208 
210  REAL8 deltaT,
211  REAL8 m1,
212  REAL8 m2,
213  REAL8 chi1z,
214  REAL8 chi2z,
215  UINT4 max_ell
216 );
217 
219  LALValue *ModeArray,
220  const NRHybSurData *NR_hybsur_data
221  );
222 
224  UINT4 *num_modes_incl,
225  UINT4 *max_ell,
226  LALValue *ModeArray,
227  const NRHybSurData *NR_hybsur_data
228  );
struct tagLALH5File LALH5File
const char * name
int NRHybSur_check_mode_array(UINT4 *num_modes_incl, UINT4 *max_ell, LALValue *ModeArray, const NRHybSurData *NR_hybsur_data)
Sanity checks on ModeArray.
int ReadHDF5DoubleDataset(REAL8 *param, LALH5File *sub, const char *name)
Reads a REAL8 value from a H5 file/group.
void NRHybSur_set_default_modes(LALValue *ModeArray, const NRHybSurData *NR_hybsur_data)
Activates all modes of an NRHybSur model.
int ReadHDF5IntDataset(int *param, LALH5File *sub, const char *name)
Reads an INT8 value from a H5 file/group.
int NRHybSur_Init(NRHybSurData *data, LALH5File *file)
Initialize a NRHybSurData structure from an open H5 file.
int NRHybSur_eval_phase_22(gsl_vector **phi_22, gsl_vector **output_times, const REAL8 eta, const gsl_vector *fit_params, const REAL8 omegaM_22_min, const REAL8 deltaTOverM, const REAL8 phiRef, const REAL8 omegaM_22_ref, gsl_vector *dummy_dp, const gsl_matrix *x_train, gsl_vector *dummy_worker, const NRHybSurData *NR_hybsur_data)
Evaluate the phase of the (2,2) mode.
int NRHybSur_eval_mode_data_pieces(EvaluatedDataPieces **this_mode_eval_dp, const UINT4 ell, const UINT4 m, const ModeDataPieces *data_pieces, const gsl_vector *output_times, const gsl_vector *fit_params, gsl_vector *dummy_dp, const gsl_matrix *x_train, gsl_vector *dummy_worker, const NRHybSurData *NR_hybsur_data)
Evaluate waveform data pieces of a single mode.
int NRHybSur_sanity_check_sample_rate(REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 chi1z, REAL8 chi2z, UINT4 max_ell)
Sanity check (warning only, not error) that the sample rate is high enough to capture the ringdown fr...
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.
void NRHybSur_DestroyEvaluatedDataPieces(gsl_vector *phi_22, EvaluatedDataPieces **evaluated_mode_dps, const UINT4 num_modes_incl)
Destroy phi_22 and an EvaluatedDataPieces structure.
sigmaKerr data[0]
double REAL8
uint32_t UINT4
static const INT4 m
NRHybSur data for a single waveform data piece.
int n_nodes
Number of empirical nodes.
NRHybSurFitData ** fit_data
NRHybSurFitData at each empirical node.
gsl_matrix * ei_basis
Empirical interpolation matrix.
NRHybSur evaluated data for a single mode.
gsl_vector * coorb_re_eval
Coorbital frame real part evaluation.
gsl_vector * ampl_eval
Amplitude data piece evaluation.
gsl_vector * coorb_im_eval
Coorbital frame imag part evaluation.
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.
NRHybSur data pieces of a single mode.
DataPiece * phase_res_data_piece
Phase residual data piece.
UINT4 m
Mode m value.
DataPiece * ampl_data_piece
Amplitude data piece.
UINT4 ell
Mode value.
DataPiece * coorb_im_data_piece
Coorbital frame imag part data piece.
DataPiece * coorb_re_data_piece
Coorbital frame real part data piece.
NRHybSur surrogate data for all modes, to be loaded from a h5 file.
ModeDataPieces ** mode_data_pieces
Data pieces of all modes, same order as mode_list.
gsl_vector * domain
Common time samples for the surrogate.
gsl_vector * TaylorT3_factor_without_eta
Precomputed quantity for use in evaluating the 0PN TaylorT3 phase contribution.
REAL8 params_dim
Dimensions of the model.
UINT4 setup
Indicates if NRHybSur has been initialized.
gsl_matrix * x_train
Training set parameters, needed for GPR fits.
int phaseAlignIdx
Alignment index for phase data piece.
UINT4 num_modes_modeled
Number of modeled modes.
gsl_matrix_long * mode_list
List of modeled modes, first element is always (2,2).
Data used in a single NRHybSur fit.
REAL8 data_std
Standard deviation after removing mean.
REAL8 data_mean
Mean of data after linear fit.
GPRHyperParams * hyperparams
Hyperparameters from GPR fit.
gsl_vector * lin_coef
Linear coefs from linear fit, size=D.
REAL8 lin_intercept
Constant intercept from linear fit.
double deltaT
Definition: unicorn.c:24