Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimIMRPhenomP.h
Go to the documentation of this file.
1#ifndef _LALSIM_IMR_PHENOMP_H
2#define _LALSIM_IMR_PHENOMP_H
3
4/*
5 * Copyright (C) 2013,2014,2015 Michael Puerrer, Alejandro Bohe
6 *
7 * This program is free software; you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation; either version 2 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with with program; see the file COPYING. If not, write to the
19 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
20 * MA 02110-1301 USA
21 */
22
23#include <lal/LALStdlib.h>
24#include <lal/LALSimIMR.h>
25#include <lal/LALConstants.h>
26
29
30#include <lal/FrequencySeries.h>
31#include <lal/LALSimInspiral.h>
32
33
34/* CONSTANTS */
35
36/**
37 * Tolerance used below which numbers are treated as zero for the calculation of atan2
38 */
39#define MAX_TOL_ATAN 1.0e-15
40
41
42/* ************************** PhenomP internal function prototypes *****************************/
43/* atan2 wrapper that returns 0 when both magnitudes of x and y are below tol, otherwise it returns
44 atan2(x, y) */
45static REAL8 atan2tol(REAL8 x, REAL8 y, REAL8 tol);
46
47/* PhenomC parameters for modified ringdown: Uses final spin formula of Barausse & Rezzolla, Astrophys.J.Lett.704:L40-L44, 2009 */
49 const REAL8 m1, /**< Mass of companion 1 (solar masses) */
50 const REAL8 m2, /**< Mass of companion 2 (solar masses) */
51 const REAL8 chi, /**< Reduced aligned spin of the binary chi = (m1*chi1 + m2*chi2)/M */
52 const REAL8 chip, /**< Dimensionless spin in the orbital plane */
53 LALDict *extraParams /**< linked list containing the extra testing GR parameters */
54);
55
56typedef struct tagNNLOanglecoeffs {
57 REAL8 alphacoeff1; /* Coefficient of omega^(-1) in alphaNNLO */
58 REAL8 alphacoeff2; /* Coefficient of omega^(-2/3) in alphaNNLO */
59 REAL8 alphacoeff3; /* Coefficient of omega^(-1/3) in alphaNNLO */
60 REAL8 alphacoeff4; /* Coefficient of log(omega) in alphaNNLO */
61 REAL8 alphacoeff5; /* Coefficient of omega^(1/3) in alphaNNLO */
62
63 REAL8 epsiloncoeff1; /* Coefficient of omega^(-1) in epsilonNNLO */
64 REAL8 epsiloncoeff2; /* Coefficient of omega^(-2/3) in epsilonNNLO */
65 REAL8 epsiloncoeff3; /* Coefficient of omega^(-1/3) in epsilonNNLO */
66 REAL8 epsiloncoeff4; /* Coefficient of log(omega) in epsilonNNLO */
67 REAL8 epsiloncoeff5; /* Coefficient of omega^(1/3) in epsilonNNLO */
69
71 NNLOanglecoeffs *angcoeffs, /**< Output: Structure to store results */
72 const REAL8 q, /**< Mass-ratio (convention q>1) */
73 const REAL8 chil, /**< Dimensionless aligned spin of the largest BH */
74 const REAL8 chip /**< Dimensionless spin component in the orbital plane */
75);
76
77typedef struct tagSpinWeightedSphericalHarmonic_l2 {
78 COMPLEX16 Y2m2, Y2m1, Y20, Y21, Y22;
80
81/* Internal core function to calculate PhenomP polarizations for a sequence of frequences. */
82static int PhenomPCore(
83 COMPLEX16FrequencySeries **hptilde, /**< Output: Frequency-domain waveform h+ */
84 COMPLEX16FrequencySeries **hctilde, /**< Output: Frequency-domain waveform hx */
85 const REAL8 chi1_l_in, /**< Dimensionless aligned spin on companion 1 */
86 const REAL8 chi2_l_in, /**< Dimensionless aligned spin on companion 2 */
87 const REAL8 chip, /**< Effective spin in the orbital plane */
88 const REAL8 thetaJ, /**< Angle between J0 and line of sight (z-direction) */
89 const REAL8 m1_SI_in, /**< Mass of companion 1 (kg) */
90 const REAL8 m2_SI_in, /**< Mass of companion 2 (kg) */
91 const REAL8 distance, /**< Distance of source (m) */
92 const REAL8 alpha0, /**< Initial value of alpha angle (azimuthal precession angle) */
93 const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
94 const REAL8 f_ref, /**< Reference frequency */
95 const REAL8Sequence *freqs, /**< Frequency points at which to evaluate the waveform (Hz) */
96 double deltaF, /**< Sampling frequency (Hz).
97 * If deltaF > 0, the frequency points given in freqs are uniformly spaced with
98 * spacing deltaF. Otherwise, the frequency points are spaced non-uniformly.
99 * Then we will use deltaF = 0 to create the frequency series we return. */
100 IMRPhenomP_version_type IMRPhenomP_version, /**< IMRPhenomPv1 uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal is a tidal version of IMRPhenomPv2 */
101 NRTidal_version_type NRTidal_version, /**< either NRTidal or NRTidalv2 for BNS waveform; NoNRT_V for BBH waveform */
102 LALDict *extraParams /**< linked list containing the extra testing GR parameters */
103);
104
105/* Internal core function to calculate PhenomP polarizations for a single frequency. */
107 const REAL8 fHz, /**< Frequency (Hz) */
108 const REAL8 eta, /**< Symmetric mass ratio */
109 const REAL8 distance, /**< Distance of source (m) */
110 const REAL8 M, /**< Total mass (Solar masses) */
111 const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
112 IMRPhenomDAmplitudeCoefficients *pAmp, /**< Internal IMRPhenomD amplitude coefficients */
113 IMRPhenomDPhaseCoefficients *pPhi, /**< Internal IMRPhenomD phase coefficients */
114 BBHPhenomCParams *PCparams, /**< Internal PhenomC parameters */
115 PNPhasingSeries *PNparams, /**< PN inspiral phase coefficients */
116 COMPLEX16 *hPhenom, /**< Output: IMRPhenom waveform (before precession) */
117 REAL8 *phasing, /**< Output: overall phasing */
118 IMRPhenomP_version_type IMRPhenomP_version, /**< Version number: 1 uses IMRPhenomC, 2 uses IMRPhenomD, NRTidal uses IMRPhenomPv2 with the NRTidal framework */
119 AmpInsPrefactors *amp_prefactors, /**< pre-calculated (cached for saving runtime) coefficients for amplitude. See LALSimIMRPhenomD_internals.c*/
120 PhiInsPrefactors *phi_prefactors /**< pre-calculated (cached for saving runtime) coefficients for phase. See LALSimIMRPhenomD_internals.*/
121);
122
124 const REAL8 fHz, /**< Frequency (Hz) */
125 COMPLEX16 hPhenom, /**< [in] IMRPhenom waveform (before precession) */
126 const REAL8 eta, /**< Symmetric mass ratio */
127 const REAL8 chi1_l, /**< Dimensionless aligned spin on companion 1 */
128 const REAL8 chi2_l, /**< Dimensionless aligned spin on companion 2 */
129 const REAL8 chip, /**< Dimensionless spin in the orbital plane */
130 const REAL8 M, /**< Total mass (Solar masses) */
131 NNLOanglecoeffs *angcoeffs, /**< Struct with PN coeffs for the NNLO angles */
132 SpinWeightedSphericalHarmonic_l2 *Y2m, /**< Struct of l=2 spherical harmonics of spin weight -2 */
133 const REAL8 alphaoffset, /**< f_ref dependent offset for alpha angle (azimuthal precession angle) */
134 const REAL8 epsilonoffset, /**< f_ref dependent offset for epsilon angle */
135 COMPLEX16 *hp, /**< [out] plus polarization \f$\tilde h_+\f$ */
136 COMPLEX16 *hc, /**< [out] cross polarization \f$\tilde h_x\f$ */
137 IMRPhenomP_version_type IMRPhenomP_version /**< IMRPhenomP(v1) uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
138);
139
140/* Simple 2PN version of L, without any spin terms expressed as a function of v */
142 const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
143 const REAL8 eta /**< Symmetric mass-ratio */
144);
145
147 const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
148 const REAL8 eta /**< Symmetric mass-ratio */
149);
150
152 REAL8 *cos_beta_half, /**< Output: cos(beta/2) */
153 REAL8 *sin_beta_half, /**< Output: sin(beta/2) */
154 const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
155 const REAL8 SL, /**< Dimensionfull aligned spin */
156 const REAL8 eta, /**< Symmetric mass-ratio */
157 const REAL8 Sp /**< Dimensionfull spin component in the orbital plane */
158);
159
161 REAL8 *cos_beta_half, /**< Output: cos(beta/2) */
162 REAL8 *sin_beta_half, /**< Output: sin(beta/2) */
163 const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
164 const REAL8 SL, /**< Dimensionfull aligned spin */
165 const REAL8 eta, /**< Symmetric mass-ratio */
166 const REAL8 Sp /**< Dimensionfull spin component in the orbital plane */
167);
168
170 const REAL8 m1, /**< Mass of companion 1 (solar masses) */
171 const REAL8 m2, /**< Mass of companion 2 (solar masses) */
172 const REAL8 chi1_l, /**< Aligned spin of BH 1 */
173 const REAL8 chi2_l, /**< Aligned spin of BH 2 */
174 const REAL8 chip /**< Dimensionless spin in the orbital plane */
175);
176
178 const REAL8 m1, /**< Mass of companion 1 (solar masses) */
179 const REAL8 m2, /**< Mass of companion 2 (solar masses) */
180 const REAL8 chi1_l, /**< Aligned spin of BH 1 */
181 const REAL8 chi2_l, /**< Aligned spin of BH 2 */
182 const REAL8 chip /**< Dimensionless spin in the orbital plane */
183);
184
186 const REAL8 nu, /**< Symmetric mass-ratio */
187 const REAL8 chi, /**< Effective aligned spin of the binary: chi = (m1*chi1 + m2*chi2)/M */
188 const REAL8 chip /**< Dimensionless spin in the orbital plane */
189);
190
191static REAL8 FinalSpinBarausse2009( /* Barausse & Rezzolla, Astrophys.J.Lett.704:L40-L44, 2009, arXiv:0904.2577 */
192 const REAL8 nu, /**< Symmetric mass-ratio */
193 const REAL8 a1, /**< |a_1| norm of dimensionless spin vector for BH 1 */
194 const REAL8 a2, /**< |a_2| norm of dimensionless spin vector for BH 2 */
195 const REAL8 cos_alpha, /**< cos(alpha) = \\hat a_1 . \\hat a_2 (Eq. 7) */
196 const REAL8 cos_beta_tilde, /**< cos(\\tilde beta) = \\hat a_1 . \\hat L (Eq. 9) */
197 const REAL8 cos_gamma_tilde /**< cos(\\tilde gamma) = \\hat a_2 . \\hat L (Eq. 9)*/
198);
199
200static bool approximately_equal(REAL8 x, REAL8 y, REAL8 epsilon);
201static void nudge(REAL8 *x, REAL8 X, REAL8 epsilon);
202
204 const REAL8 fHz, /**< Frequency (Hz) */
205 const REAL8 window, /**< planck window */
206 const REAL8 ampTidal, /**< tidal amplitude at a frequency sample */
207 const REAL8 phaseTidal, /**< tidal phasing at a frequency sample from NRTidal infrastructure*/
208 const REAL8 distance, /**< Distance of source (m) */
209 const REAL8 M, /**< Total mass (Solar masses) */
210 const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
211 IMRPhenomDAmplitudeCoefficients *pAmp, /**< Internal IMRPhenomD amplitude coefficients */
212 IMRPhenomDPhaseCoefficients *pPhi, /**< Internal IMRPhenomD phase coefficients */
213 PNPhasingSeries *PNparams, /**< PN inspiral phase coefficients */
214 COMPLEX16 *hPhenom, /**< [out] IMRPhenom waveform (before precession) */
215 REAL8 *phasing, /**< [out] overall phasing */
216 AmpInsPrefactors *amp_prefactors, /**< pre-calculated (cached for saving runtime) coefficients for amplitude. See LALSimIMRPhenomD_internals.c*/
217 PhiInsPrefactors *phi_prefactors /**< pre-calculated (cached for saving runtime) coefficients for phase. See LALSimIMRPhenomD_internals.*/
218);
219
220#endif // of #ifndef _LALSIM_IMR_PHENOMP_H
Internal function for IMRPhenomD phenomenological waveform model. See LALSimIMRPhenom....
static void ComputeNNLOanglecoeffs(NNLOanglecoeffs *angcoeffs, const REAL8 q, const REAL8 chil, const REAL8 chip)
static void nudge(REAL8 *x, REAL8 X, REAL8 epsilon)
static void WignerdCoefficients(REAL8 *cos_beta_half, REAL8 *sin_beta_half, const REAL8 v, const REAL8 SL, const REAL8 eta, const REAL8 Sp)
static REAL8 atan2tol(REAL8 x, REAL8 y, REAL8 tol)
static int PhenomPCore(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8 chi1_l_in, const REAL8 chi2_l_in, const REAL8 chip, const REAL8 thetaJ, const REAL8 m1_SI_in, const REAL8 m2_SI_in, const REAL8 distance, const REAL8 alpha0, const REAL8 phic, const REAL8 f_ref, const REAL8Sequence *freqs, double deltaF, IMRPhenomP_version_type IMRPhenomP_version, NRTidal_version_type NRTidal_version, LALDict *extraParams)
static void CheckMaxOpeningAngle(const REAL8 m1, const REAL8 m2, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip)
static int PhenomPCoreTwistUp(const REAL8 fHz, COMPLEX16 hPhenom, const REAL8 eta, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip, const REAL8 M, NNLOanglecoeffs *angcoeffs, SpinWeightedSphericalHarmonic_l2 *Y2m, const REAL8 alphaoffset, const REAL8 epsilonoffset, COMPLEX16 *hp, COMPLEX16 *hc, IMRPhenomP_version_type IMRPhenomP_version)
static BBHPhenomCParams * ComputeIMRPhenomCParamsRDmod(const REAL8 m1, const REAL8 m2, const REAL8 chi, const REAL8 chip, LALDict *extraParams)
static int PhenomPCoreOneFrequency(const REAL8 fHz, const REAL8 eta, const REAL8 distance, const REAL8 M, const REAL8 phic, IMRPhenomDAmplitudeCoefficients *pAmp, IMRPhenomDPhaseCoefficients *pPhi, BBHPhenomCParams *PCparams, PNPhasingSeries *PNparams, COMPLEX16 *hPhenom, REAL8 *phasing, IMRPhenomP_version_type IMRPhenomP_version, AmpInsPrefactors *amp_prefactors, PhiInsPrefactors *phi_prefactors)
static REAL8 L2PNR_v1(const REAL8 v, const REAL8 eta)
static REAL8 FinalSpinBarausse2009(const REAL8 nu, const REAL8 a1, const REAL8 a2, const REAL8 cos_alpha, const REAL8 cos_beta_tilde, const REAL8 cos_gamma_tilde)
static int PhenomPCoreOneFrequency_withTides(const REAL8 fHz, const REAL8 window, const REAL8 ampTidal, const REAL8 phaseTidal, const REAL8 distance, const REAL8 M, const REAL8 phic, IMRPhenomDAmplitudeCoefficients *pAmp, IMRPhenomDPhaseCoefficients *pPhi, PNPhasingSeries *PNparams, COMPLEX16 *hPhenom, REAL8 *phasing, AmpInsPrefactors *amp_prefactors, PhiInsPrefactors *phi_prefactors)
static REAL8 L2PNR(const REAL8 v, const REAL8 eta)
static REAL8 FinalSpinBarausse2009_all_spin_on_larger_BH(const REAL8 nu, const REAL8 chi, const REAL8 chip)
static REAL8 FinalSpinIMRPhenomD_all_in_plane_spin_on_larger_BH(const REAL8 m1, const REAL8 m2, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip)
static void WignerdCoefficients_SmallAngleApproximation(REAL8 *cos_beta_half, REAL8 *sin_beta_half, const REAL8 v, const REAL8 SL, const REAL8 eta, const REAL8 Sp)
static bool approximately_equal(REAL8 x, REAL8 y, REAL8 epsilon)
REAL8 M
Definition: bh_qnmode.c:133
const double a2
double complex COMPLEX16
double REAL8
IMRPhenomP_version_type
Definition: LALSimIMR.h:73
NRTidal_version_type
Definition: LALSimIMR.h:80
static const INT4 q
used to cache the recurring (frequency-independent) prefactors of AmpInsAnsatz.
Structure holding all coefficients for the amplitude.
Structure holding all coefficients for the phase.
used to cache the recurring (frequency-independent) prefactors of PhiInsAnsatzInt.