Utilities needed for aligned-spin NR-hybrid surrogate models. More...
Prototypes | |
int | ReadHDF5DoubleDataset (REAL8 *param, LALH5File *sub, const char *name) |
Reads a REAL8 value from a H5 file/group. More... | |
int | ReadHDF5IntDataset (int *param, LALH5File *sub, const char *name) |
Reads an INT8 value from a H5 file/group. More... | |
int | NRHybSur_Init (NRHybSurData *data, LALH5File *file) |
Initialize a NRHybSurData structure from an open H5 file. More... | |
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. More... | |
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. More... | |
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. More... | |
void | NRHybSur_DestroyEvaluatedDataPieces (gsl_vector *phi_22, EvaluatedDataPieces **evaluated_mode_dps, const UINT4 num_modes_incl) |
Destroy phi_22 and an EvaluatedDataPieces structure. More... | |
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 frequencies, by ensuring Nyquist frequency is greater than the QNM frequency of the (max_ell,max_ell,0) mode, where max_ell is the maximum ell index among the included modes. More... | |
void | NRHybSur_set_default_modes (LALValue *ModeArray, const NRHybSurData *NR_hybsur_data) |
Activates all modes of an NRHybSur model. More... | |
int | NRHybSur_check_mode_array (UINT4 *num_modes_incl, UINT4 *max_ell, LALValue *ModeArray, const NRHybSurData *NR_hybsur_data) |
Sanity checks on ModeArray. More... | |
Utilities needed for aligned-spin NR-hybrid surrogate models.
Called from: LALSimIMRNRHybSur3dq8.h
Definition in file LALSimNRHybSurUtilities.h.
Go to the source code of this file.
Data Structures | |
struct | GPRHyperParams |
Data used in a single GPR fit. More... | |
struct | NRHybSurFitData |
Data used in a single NRHybSur fit. More... | |
struct | DataPiece |
NRHybSur data for a single waveform data piece. More... | |
struct | ModeDataPieces |
NRHybSur data pieces of a single mode. More... | |
struct | NRHybSurData |
NRHybSur surrogate data for all modes, to be loaded from a h5 file. More... | |
struct | EvaluatedDataPieces |
NRHybSur evaluated data for a single mode. More... | |
Reads a REAL8 value from a H5 file/group.
param | Output: REAL8 value from H5 file/group. |
sub | H5 file or group. |
name | Name of REAL8 dataset within file/group. |
Definition at line 68 of file LALSimNRHybSurUtilities.c.
Reads an INT8 value from a H5 file/group.
param | Output: int value from H5 file/group. |
sub | H5 file or group. |
name | Name of int dataset within file/group. |
Definition at line 92 of file LALSimNRHybSurUtilities.c.
int NRHybSur_Init | ( | NRHybSurData * | NR_hybsur_data, |
LALH5File * | file | ||
) |
Initialize a NRHybSurData structure from an open H5 file.
This will typically only be called once. For example from NRHybSur3dq8_Init_LALDATA.
The HDF5 file should have the following structure:
NR_hybsur_data | Output: Struct to save surrogate data. |
file | Opened HDF5 file. |
Definition at line 403 of file LALSimNRHybSurUtilities.c.
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.
fit_data | Data for fit. |
fit_params | Parameter space point to evaluate the fit at. size=D, the dimension of the model. |
x_train | Training set points. |
dummy_worker | Dummy worker array for computations. |
Definition at line 581 of file LALSimNRHybSurUtilities.c.
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.
The surrogate actually models the phase residual \( \phi^{res}_{22} \) defined in Eq.(44) of arxiv:1812.07865. Here we first evaluate that and then add the 0PN TaylorT3 phase to get the (2,2) mode phase.
Sets the orbital phase to phiRef at the reference frequency omegaM_22_ref. The orbital phase is obtained as phi_22/2, so this leaves a pi ambiguity. But the surrogate data is already aligned such that the heavier BH is on the +ve x-axis at t=-1000M. See Sec.VI.A.4 of arxiv:1812.07865, the resolves the pi ambiguity. This means that the after the realignment, the orbital phase at reference frequency omegaM_22_ref is phiRef, or the heavier BH is at azimuthal angle = phiRef from the +ve x-axis.
Only uses data at (2,2) mode frequencies >= omegaM_22_min. This determines the start time. The start time, along with the step size deltaTOverM, is used to determine the output_times. Uses cubic spline interpolation to interpolate from the surrogate's time array to output_times.
phi_22 | Output: (2,2) mode phase. |
output_times | Output: Time array. |
eta | Symmetric mass ratio. |
fit_params | Parameter space point to evaluate the fit at. size=D, the dimension of the model. |
omegaM_22_min | Start frequency of (2,2) mode in rad/M. |
deltaTOverM | Time step in M. |
phiRef | Orbital phase at reference frequency. |
omegaM_22_ref | Reference freq of (2,2) mode in rad/M. |
dummy_dp | Dummy vector to store phase evaluation. |
x_train | Training set points. |
dummy_worker | Dummy worker array for computations. |
NR_hybsur_data | Loaded surrogate data. |
Definition at line 887 of file LALSimNRHybSurUtilities.c.
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.
For (2,2) mode we model the amplitude and phase, but the phase is evaluated using NRHybSur_eval_phase_22, since it is required for all modes to transform from the coorbital frame to the inertial frame. For all other modes we evaluate the real and imaginary parts of the coorbital frame strain, defined in Eq.(39) of arxiv:1812.07865.
Only the required data pieces for a given mode will be evaluated.
this_mode_eval_dp | Output: evaluated waveform data pieces of a single mode. |
ell | \(\ell\) index of mode. |
m | m index of mode. |
data_pieces | Surrogate data pieces of this mode. |
output_times | Time vector to evaluate at. |
fit_params | Parameter space point to evaluate the fit at. size=D, the dimension of the model. |
dummy_dp | Dummy vector to store phase evaluation. |
x_train | Training set points. |
dummy_worker | Dummy worker array for computations. |
NR_hybsur_data | Loaded surrogate data. |
Definition at line 978 of file LALSimNRHybSurUtilities.c.
void NRHybSur_DestroyEvaluatedDataPieces | ( | gsl_vector * | phi_22, |
EvaluatedDataPieces ** | evaluated_mode_dps, | ||
const UINT4 | num_modes_incl | ||
) |
Destroy phi_22 and an EvaluatedDataPieces structure.
Free all associated memory.
phi_22 | \(\phi_{22}\) data piece. |
evaluated_mode_dps | All other data pieces. |
num_modes_incl | Number of models included. |
Definition at line 1056 of file LALSimNRHybSurUtilities.c.
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 frequencies, by ensuring Nyquist frequency is greater than the QNM frequency of the (max_ell,max_ell,0) mode, where max_ell is the maximum ell index among the included modes.
deltaT | Sampling interval (s). |
m1 | Mass of Bh1 (kg). |
m2 | Mass of Bh2 (kg). |
chi1z | Dimensionless spin of Bh1. |
chi2z | Dimensionless spin of Bh2. |
max_ell | Max ell index included. |
Definition at line 1101 of file LALSimNRHybSurUtilities.c.
void NRHybSur_set_default_modes | ( | LALValue * | ModeArray, |
const NRHybSurData * | NR_hybsur_data | ||
) |
Activates all modes of an NRHybSur model.
For NRHybSur3dq8 that is \( \ell \leq 4, m \geq 0 \), and (5,5), but not (4,1) or (4,0).
ModeArray | Output: Container for modes. |
NR_hybsur_data | Loaded surrogate data. |
Definition at line 1151 of file LALSimNRHybSurUtilities.c.
int NRHybSur_check_mode_array | ( | UINT4 * | num_modes_incl, |
UINT4 * | max_ell, | ||
LALValue * | ModeArray, | ||
const NRHybSurData * | NR_hybsur_data | ||
) |
Sanity checks on ModeArray.
Will raise an error if an unavailable mode is requested. Note that we only look for m>=0 modes in ModeArray, and will ignore m<0 modes even if present. The m<0 modes automatically get added when evaluting the waveform.
num_modes_incl | Output: Number of modes to include. |
max_ell | Output: Max ell index included. |
ModeArray | Container for modes. |
NR_hybsur_data | Loaded surrogate data. |
Definition at line 1175 of file LALSimNRHybSurUtilities.c.