IMRPhenomD

Header file for utilities

Variables

const double lambda_num_params[19][11]

Numerically calibrated parameters from arXiv:1508.07253 see the table in the data directory for labeled version from lalsuite

template<class T>
struct lambda_parameters
#include <IMRPhenomD.h>

Structure to facilitate IMRPhenomD parameter transfers

Public Members

T rho[4]
T v2
T gamma[4]
T sigma[5]
T beta[5]
T alpha[7]
template<class T>
class IMRPhenomD
#include <IMRPhenomD.h>

Subclassed by IMRPhenomD_NRT< T >, IMRPhenomPv2< T >, gIMRPhenomD< T >, ppE_IMRPhenomD_Inspiral< T >

Public Functions

virtual ~IMRPhenomD() = default
virtual void fisher_calculation_sky_averaged(double *frequency, int length, gen_params *parameters, double **amplitude_deriv, double **phase_deriv, double *amplitude, int *amp_tapes, int *phase_tapes)
virtual void change_parameter_basis(T *old_param, T *new_param, bool sky_average)

Convience method to change parameter basis between common Fisher parameters and the intrinsic parameters of IMRPhenomD.

Takes input array of old parameters and ouputs array of transformed parameters

Parameters:
  • old_param – array of old params, order {A0, tc, phic, chirpmass, eta, spin1, spin2}

  • new_param – output new array: order {m1,m2,DL, spin1,spin2,phic,tc}

virtual void construct_amplitude_derivative(double *frequencies, int length, int dimension, double **amplitude_derivative, source_parameters<double> *input_params, int *tapes = NULL)

Construct the derivative of the amplitude for a given source evaluated by the given frequency.

Order of output: dh/d \theta : \theta \el {A0,tc, phic, chirp mass, eta, symmetric spin, antisymmetric spin}

Parameters:
  • frequencies – input array of frequency

  • length – length of the frequency array

  • amplitude_derivative – < dimension of the fisher output array for all the derivatives double[dimension][length]

  • input_params – Source parameters structure for the source

  • tapes – int array of tape ids, if NULL, these will be calculated

virtual void construct_phase_derivative(double *frequencies, int length, int dimension, double **phase_derivative, source_parameters<double> *input_params, int *tapes = NULL)

Construct the derivative of the phase for a given source evaluated by the given frequency.

Order of output: dh/d \theta : \theta \el {A0,tc, phic, chirp mass, eta, symmetric spin, antisymmetric spin}

Parameters:
  • frequencies – input array of frequency

  • length – length of the frequency array

  • phase_derivative – < dimension of the fisher output array for all the derivatives double[dimension][length]

  • input_params – Source parameters structure for the source

  • tapes – int array of tape ids, if NULL, these will be calculated

virtual void amplitude_tape(source_parameters<double> *input_params, int *tape)

Creates the tapes for derivatives of the amplitude.

For efficiency in long runs of large sets of fishers, the tapes can be precomputed and reused

Parameters:
  • input_params – source parameters structure of the desired source

  • tape – tape ids

virtual void phase_tape(source_parameters<double> *input_params, int *tape)

Creates the tapes for derivatives of phase.

For efficiency in long runs of large sets of fishers, the tapes can be precomputed and reused

Parameters:
  • input_params – source parameters structure of the desired source

  • tape – tape ids

virtual int construct_waveform(T *frequencies, int length, std::complex<T> *waveform, source_parameters<T> *params)

Constructs the waveform as outlined by.

arguments: array of frequencies, length of that array, a complex array for the output waveform, and a source_parameters structure

Parameters:
  • frequencies – T array of frequencies the waveform is to be evaluated at

  • length – integer length of the array of frequencies and the waveform

  • waveform – complex T array for the waveform to be output

virtual std::complex<T> construct_waveform(T frequency, source_parameters<T> *params)

overloaded method to evaluate the waveform for one frequency instead of an array

Parameters:

frequency – T array of frequencies the waveform is to be evaluated at

virtual int construct_amplitude(T *frequencies, int length, T *amplitude, source_parameters<T> *params)

Constructs the Amplitude as outlined by IMRPhenomD.

arguments: array of frequencies, length of that array, T array for the output amplitude, and a source_parameters structure

Parameters:
  • frequencies – T array of frequencies the waveform is to be evaulated at

  • length – integer length of the input array of frequencies and the output array

  • amplitude – output T array for the amplitude

  • params – Structure of source parameters to be initilized before computation

virtual int construct_phase(T *frequencies, int length, T *phase, source_parameters<T> *params)

Constructs the Phase as outlined by IMRPhenomD.

arguments: array of frequencies, length of that array, T array for the output phase, and a source_parameters structure

Parameters:
  • frequencies – T array of frequencies the waveform is to be evaluated at

  • length – integer length of the input and output arrays

  • phase – output T array for the phasee

  • params – structure of source parameters to be calculated before computation

virtual T build_amp(T f, lambda_parameters<T> *lambda, source_parameters<T> *params, useful_powers<T> *pows, T *amp_coeff, T *deltas)

constructs the IMRPhenomD amplitude for frequency f

arguments: numerical parameters from Khan et al lambda_parameters structure, source_parameters structure, useful_powers<T> structure, PN parameters for the inspiral portions of the waveform, and the delta parameters for the intermediate region, numerically solved for using the amp_connection_coeffs function

virtual T build_phase(T f, lambda_parameters<T> *lambda, source_parameters<T> *params, useful_powers<T> *pows, T *phase_coeff)

constructs the IMRPhenomD phase for frequency f

arguments: numerical parameters from Khan et al lambda_parameters structure, source_parameters structure, useful_powers structure, PN parameters for the inspiral portions of the waveform

virtual T assign_lambda_param_element(source_parameters<T> *source_param, int i)

Calculate the lambda parameters from Khan et al for element i.

virtual void assign_lambda_param(source_parameters<T> *source_param, lambda_parameters<T> *lambda)

Wrapper for the Lambda parameter assignment that handles the looping.

virtual void precalc_powers_ins(T f, T M, useful_powers<T> *Mf_pows)

Pre-calculate powers of Mf, to speed up calculations for the inspiral waveform (both amplitude and phase.

It seems the pow() function is very slow, so to speed things up, powers of Mf will be precomputed and passed to the functions within the frequency loops

virtual void precalc_powers_PI(useful_powers<T> *PI_pows)

Pre-calculate powers of pi, to speed up calculations for the inspiral phase.

It seems the pow() function is very slow, so to speed things up, powers of PI will be precomputed and passed to the functions within the frequency loops

virtual void precalc_powers_ins_phase(T f, T M, useful_powers<T> *Mf_pows)

Pre-calculate powers of Mf, to speed up calculations for the inspiral phase.

It seems the pow() function is very slow, so to speed things up, powers of Mf will be precomputed and passed to the functions within the frequency loops

virtual void precalc_powers_ins_amp(T f, T M, useful_powers<T> *Mf_pows)

Pre-calculate powers of Mf, to speed up calculations for the inspiral amplitude.

It seems the pow() function is very slow, so to speed things up, powers of Mf will be precomputed and passed to the functions within the frequency loops

virtual void assign_pn_amplitude_coeff(source_parameters<T> *source_param, T *coeff)

Calculates the static PN coeffecients for the amplitude.

virtual void assign_static_pn_phase_coeff(source_parameters<T> *source_param, T *coeff)

Calculates the static PN coeffecients for the phase - coeffecients 0,1,2,3,4,7.

virtual void assign_nonstatic_pn_phase_coeff(source_parameters<T> *source_param, T *coeff, T f)

Calculates the dynamic PN phase coefficients 5,6.

f is in Hz

virtual void assign_nonstatic_pn_phase_coeff_deriv(source_parameters<T> *source_param, T *Dcoeff, T f)

Calculates the derivative of the dynamic PN phase coefficients 5,6.

f is in Hz

virtual void post_merger_variables(source_parameters<T> *source_param)

Calculates the post-merger ringdown frequency and dampening frequency.

Returns in Hz - assigns fRD to var[0] and fdamp to var[1]

virtual void calc_fring(source_parameters<T> *source_params)
virtual void calc_fdamp(source_parameters<T> *source_params)
virtual void _calc_fring(source_parameters<T> *source_params)

Deprecated version, but ADOL-C version &#8212; can be forced to use if dep_postmerger flag set to TRUE.

virtual void _calc_fdamp(source_parameters<T> *source_params)

Deprecated version, but ADOL-C version &#8212; can be forced to use if dep_postmerger flag set to TRUE.

virtual T final_spin(source_parameters<T> *params)
virtual T FinalSpin0815_s(T eta, T s)

Formula to predict the final spin. Equation 3.6 arXiv:1508.07250 s defined around Equation 3.6.

virtual T FinalSpin0815(T eta, T chi1, T chi2)

Wrapper function for FinalSpin0815_s.

virtual T EradRational0815_s(T eta, T s)

Formula to predict the total radiated energy. Equation 3.7 and 3.8 arXiv:1508.07250 Input parameter s defined around Equation 3.7 and 3.8.

virtual T EradRational0815(T eta, T chi1, T chi2)

Wrapper function for EradRational0815_s.

virtual T fpeak(source_parameters<T> *params, lambda_parameters<T> *lambda)

Solves for the peak frequency, where the waveform transitions from intermediate to merger-ringdown.

returns Hz

virtual T amp_ins(T f, source_parameters<T> *param, T *pn_coeff, lambda_parameters<T> *lambda, useful_powers<T> *pow)

Calculates the scaled inspiral amplitude A/A0 for frequency f with precomputed powers of MF and PI.

return a T

additional argument contains useful powers of MF and PI in structure userful_powers

virtual T Damp_ins(T f, source_parameters<T> *param, T *pn_coeff, lambda_parameters<T> *lambda)

Calculates the derivative wrt frequency for the scaled inspiral amplitude A/A0 for frequency f.

This is an analytic derivative for the smoothness condition on the amplitude connection

return a T

virtual T phase_ins(T f, source_parameters<T> *param, T *pn_coeff, lambda_parameters<T> *lambda, useful_powers<T> *pow)

Calculates the inspiral phase for frequency f with precomputed powers of MF and PI for speed.

return a T

extra argument of precomputed powers of MF and pi, contained in the structure useful_powers<T>

virtual T Dphase_ins(T f, source_parameters<T> *param, T *pn_coeff, lambda_parameters<T> *lambda)

Calculates the derivative of the inspiral phase for frequency f.

For phase continuity and smoothness return a T

virtual T amp_mr(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the scaled merger-ringdown amplitude A/A0 for frequency f.

return a T

virtual T phase_mr(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the merger-ringdown phase for frequency f.

return a T

virtual T Damp_mr(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the derivative wrt frequency for the scaled merger-ringdown amplitude A/A0 for frequency f.

This is an analytic derivative for the smoothness condition on the amplitude connection

The analytic expression was obtained from Mathematica - See the mathematica folder for code

return a T

virtual T Dphase_mr(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the derivative of the merger-ringdown phase for frequency f.

For phase continuity and smoothness return a T

virtual T amp_int(T f, source_parameters<T> *param, lambda_parameters<T> *lambda, T *deltas)

Calculates the scaled intermediate range amplitude A/A0 for frequency f.

return a T

virtual T phase_int(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the intermediate phase for frequency f.

return a T

virtual T Dphase_int(T f, source_parameters<T> *param, lambda_parameters<T> *lambda)

Calculates the derivative of the intermediate phase for frequency f.

For phase continuity and smoothness return a T

virtual void phase_connection_coefficients(source_parameters<T> *param, lambda_parameters<T> *lambda, T *pn_coeffs)

Calculates the phase connection coefficients alpha{0,1} and beta{0,1}.

Note: these coefficients are stored in the lambda parameter structure, not a separate array

virtual T calculate_beta1(source_parameters<T> *param, lambda_parameters<T> *lambda, T *pn_coeffs)
virtual T calculate_beta0(source_parameters<T> *param, lambda_parameters<T> *lambda, T *pn_coeffs)
virtual T calculate_alpha1(source_parameters<T> *param, lambda_parameters<T> *lambda)
virtual T calculate_alpha0(source_parameters<T> *param, lambda_parameters<T> *lambda)
virtual void amp_connection_coeffs(source_parameters<T> *param, lambda_parameters<T> *lambda, T *pn_coeffs, T *coeffs)

Solves for the connection coefficients to ensure the transition from inspiral to merger ringdown is continuous and smooth.

virtual T calculate_delta_parameter_0(T f1, T f2, T f3, T v1, T v2, T v3, T dd1, T dd3, T M)

Calculates the delta_0 component.

Solved in Mathematica and imported to C

virtual T calculate_delta_parameter_1(T f1, T f2, T f3, T v1, T v2, T v3, T dd1, T dd3, T M)

Calculates the delta_1 component.

Solved in Mathematica and imported to C

virtual T calculate_delta_parameter_2(T f1, T f2, T f3, T v1, T v2, T v3, T dd1, T dd3, T M)

Calculates the delta_2 component.

Solved in Mathematica and imported to C

virtual T calculate_delta_parameter_3(T f1, T f2, T f3, T v1, T v2, T v3, T dd1, T dd3, T M)

Calculates the delta_3 component.

Solved in Mathematica and imported to C

virtual T calculate_delta_parameter_4(T f1, T f2, T f3, T v1, T v2, T v3, T dd1, T dd3, T M)

Calculates the delta_4 component.

Solved in Mathematica and imported to C