LALSimulation  5.4.0.1-fe68b98
LALSimNRHybSurUtilities.h File Reference

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...
 

Detailed Description

Utilities needed for aligned-spin NR-hybrid surrogate models.

Author
Vijay Varma

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...
 

Function Documentation

◆ ReadHDF5DoubleDataset()

int ReadHDF5DoubleDataset ( REAL8 param,
LALH5File sub,
const char name 
)

Reads a REAL8 value from a H5 file/group.

Parameters
paramOutput: REAL8 value from H5 file/group.
subH5 file or group.
nameName of REAL8 dataset within file/group.

Definition at line 68 of file LALSimNRHybSurUtilities.c.

◆ ReadHDF5IntDataset()

int ReadHDF5IntDataset ( int *  param,
LALH5File sub,
const char name 
)

Reads an INT8 value from a H5 file/group.

Parameters
paramOutput: int value from H5 file/group.
subH5 file or group.
nameName of int dataset within file/group.

Definition at line 92 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_Init()

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:

  • Top level:
    'GPR_X_train': h5 dataset, shape=(N,dim) where N=size of training dataset.
    'domain': h5 dataset, size = number of time samples, can be sparse.
    'TaylorT3_t_ref': h5 dataset, INT8, used in TaylorT3 piece.
    'mode_list': h5 dataset, INT8, shape=(M,2), where M=num. of modes. Each element is (ell,m) of that mode. Only m>=0 modes are modeled. First mode should always be (2,2).
    'phaseAlignIdx: h5 dataset, INT8. Needed in TaylorT3 piece.
    'l2_m0': h5 group, surrogate data needed for this mode.
    'l2_m1':
    'l2_m2':
    • Subgroups for each mode:
      • 'l2_m2':
        'amp': h5 group, surrogate data needed for amplitude data piece.
        'phase' : h5 group, surrogate data needed for phase data piece.
      • 'l2_m1' and all other modes:
        're': h5 group, data needed for real part of coorbital frame waveform data piece.
        'im': h5 group, data needed for imag part of coorbital frame waveform data piece.
      • Subgroups for each data piece:
        'ei_basis': h5 dataset, shape=(n_nodes, L), where L=len(domain).
        'n_nodes': h5 dataset, INT8. Num. of empirical time nodes for this data piece.
        'node_num_0': h5 group, fit data for this node.
        'node_num_1':
        • Subgroups for each node:
          'alpha': h5 dataset, size=N, where N=size of training dataset.
          'constant_value': h5 dataset, REAL8, constant value in kernel.
          'data_mean': h5 dataset, REAL8. Mean of fit data.
          'data_std': h5 dataset, REAL8. Standard deviation of fit data.
          'length_scale': h5 dataset, size=dom, where dim=dimensions of the model. Length scale parameters in the kernel.
          'lin_coef': h5 dataset, size=dim. Coefs of of linear fit.
          'lin_intercept': h5 dataset, REAL8. Intercept of linear fit.
          'L': h5 dataset, shape=(N,N), where N=size of training dataset. Not used right now, needed to evaluate error estimate.
          'noise_level': h5 dataset, REAL8. Noise parameter, not needed for evaluating fit.
          'y_train_mean': h5 dataset, REAL8. Mean value before doing GPR fit, usually zero. We already remove the mean and store it in data_mean.
Parameters
NR_hybsur_dataOutput: Struct to save surrogate data.
fileOpened HDF5 file.

Definition at line 403 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_eval_fit()

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.

Parameters
fit_dataData for fit.
fit_paramsParameter space point to evaluate the fit at. size=D, the dimension of the model.
x_trainTraining set points.
dummy_workerDummy worker array for computations.

Definition at line 581 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_eval_phase_22()

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.

Parameters
phi_22Output: (2,2) mode phase.
output_timesOutput: Time array.
etaSymmetric mass ratio.
fit_paramsParameter space point to evaluate the fit at. size=D, the dimension of the model.
omegaM_22_minStart frequency of (2,2) mode in rad/M.
deltaTOverMTime step in M.
phiRefOrbital phase at reference frequency.
omegaM_22_refReference freq of (2,2) mode in rad/M.
dummy_dpDummy vector to store phase evaluation.
x_trainTraining set points.
dummy_workerDummy worker array for computations.
NR_hybsur_dataLoaded surrogate data.

Definition at line 887 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_eval_mode_data_pieces()

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.

Parameters
this_mode_eval_dpOutput: evaluated waveform data pieces of a single mode.
ell\(\ell\) index of mode.
mm index of mode.
data_piecesSurrogate data pieces of this mode.
output_timesTime vector to evaluate at.
fit_paramsParameter space point to evaluate the fit at. size=D, the dimension of the model.
dummy_dpDummy vector to store phase evaluation.
x_trainTraining set points.
dummy_workerDummy worker array for computations.
NR_hybsur_dataLoaded surrogate data.

Definition at line 978 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_DestroyEvaluatedDataPieces()

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.

Parameters
phi_22\(\phi_{22}\) data piece.
evaluated_mode_dpsAll other data pieces.
num_modes_inclNumber of models included.

Definition at line 1056 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_sanity_check_sample_rate()

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.

Parameters
deltaTSampling interval (s).
m1Mass of Bh1 (kg).
m2Mass of Bh2 (kg).
chi1zDimensionless spin of Bh1.
chi2zDimensionless spin of Bh2.
max_ellMax ell index included.

Definition at line 1101 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_set_default_modes()

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).

Parameters
ModeArrayOutput: Container for modes.
NR_hybsur_dataLoaded surrogate data.

Definition at line 1151 of file LALSimNRHybSurUtilities.c.

◆ NRHybSur_check_mode_array()

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.

Parameters
num_modes_inclOutput: Number of modes to include.
max_ellOutput: Max ell index included.
ModeArrayContainer for modes.
NR_hybsur_dataLoaded surrogate data.

Definition at line 1175 of file LALSimNRHybSurUtilities.c.