LALSimulation  5.4.0.1-fe68b98
LALSimIMRPhenomP.c
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2013,2014,2015 Michael Puerrer, Alejandro Bohe
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 #include <stdlib.h>
21 #include <math.h>
22 #include <float.h>
23 #include <gsl/gsl_errno.h>
24 #include <gsl/gsl_spline.h>
25 #include <gsl/gsl_math.h>
26 #include <gsl/gsl_sf_trig.h>
27 
28 #include <lal/Date.h>
29 #include <lal/FrequencySeries.h>
30 #include <lal/LALAtomicDatatypes.h>
31 #include <lal/LALConstants.h>
32 #include <lal/LALDatatypes.h>
33 #include <lal/LALSimInspiral.h>
34 #include <lal/Units.h>
35 #include <lal/XLALError.h>
36 #include <lal/SphericalHarmonics.h>
37 #include <lal/Sequence.h>
38 #include <lal/LALStdlib.h>
39 #include <lal/LALStddef.h>
40 
41 #include "LALSimIMR.h"
43 /* This is ugly, but allows us to reuse internal IMRPhenomC and IMRPhenomD functions without making those functions XLAL */
46 
47 #include "LALSimIMRPhenomP.h"
48 
49 #ifndef _OPENMP
50 #define omp ignore
51 #endif
52 
53 /* Macro functions to rotate the components of a vector about an axis */
54 #define ROTATEZ(angle, vx, vy, vz)\
55 tmp1 = vx*cos(angle) - vy*sin(angle);\
56 tmp2 = vx*sin(angle) + vy*cos(angle);\
57 vx = tmp1;\
58 vy = tmp2
59 
60 #define ROTATEY(angle, vx, vy, vz)\
61 tmp1 = vx*cos(angle) + vz*sin(angle);\
62 tmp2 = - vx*sin(angle) + vz*cos(angle);\
63 vx = tmp1;\
64 vz = tmp2
65 
66 const double sqrt_6 = 2.44948974278317788;
67 
68 /* ************************************* Implementation ****************************************/
69 
70 /**
71  * @addtogroup LALSimIMRPhenom_c
72  * @{
73  *
74  * @name Routines for IMR Phenomenological Model "P"
75  * @{
76  *
77  * @author Michael Puerrer, Alejandro Bohe
78  *
79  * @brief Functions for producing IMRPhenomP waveforms for precessing binaries,
80  * as described in Hannam et al., arXiv:1308.3271 [gr-qc].
81  *
82  * @note Three versions of IMRPhenomP are available (selected by IMRPhenomP_version):
83  * * version 1 ("IMRPhenomP"): based on IMRPhenomC
84  * (outdated, not reviewed!)
85  * * version 2 ("IMRPhenomPv2"): based on IMRPhenomD
86  * (to be used, currently under review as of Dec 2015)
87  * * version NRTidal ("IMRPhenomPv2_NRTidal" and "IMRPhenomPv2_NRTidalv2"): based on IMRPhenomPv2
88  * (framework for NR-tuned tidal effects added to PhenomD aligned phasing and then twisted up).
89  * Two flavors of NRTidal models are available:
90  * original ("IMRPhenomPv2_NRTidal", based on https://arxiv.org/pdf/1706.02969.pdf)
91  * and an improved version 2 ("IMRPhenomPv2_NRTidalv2", based on https://arxiv.org/pdf/1905.06011.pdf).
92  * The different NRTidal versions employ different internal switches (selected by NRTidal_version).
93  *
94  * Each IMRPhenomP version inherits its range of validity
95  * over the parameter space from the respective aligned-spin waveform.
96  *
97  * @attention A time-domain implementation of IMRPhenomPv2 is available in XLALChooseTDWaveform().
98  * This is based on a straight-forward inverse Fourier transformation via XLALSimInspiralTDfromFD(),
99  * but it was not included in the IMRPhenomPv2 review. Use it at your own risk.
100  * IMRPhenomPv2_NRTidal is also available in the time domain through the same transformation.
101  * Visual checks have been performed during the review, and unphysical features may arise for
102  * mass ratios smaller than 1.5 and when both tidal parameters are greater than 2000. In this
103  * case, a warning is issued, both for the time and frequency domain version.
104  */
105 
106 static REAL8 atan2tol(REAL8 a, REAL8 b, REAL8 tol)
107 {
108  REAL8 c;
109  if (fabs(a) < tol && fabs(b) < tol)
110  c = 0.;
111  else
112  c = atan2(a, b);
113  return c;
114 }
115 
116 /**
117  * Deprecated : used the old convention (view frame for the spins)
118  * Function to map LAL parameters
119  * (masses, 6 spin components and Lhat at f_ref)
120  * into IMRPhenomP intrinsic parameters
121  * (chi1_l, chi2_l, chip, thetaJ, alpha0).
122  *
123  * All input masses and frequencies should be in SI units.
124  *
125  * See Fig. 1. in arxiv:1408.1810 for a diagram of the angles.
126  */
128  REAL8 *chi1_l, /**< [out] Dimensionless aligned spin on companion 1 */
129  REAL8 *chi2_l, /**< [out] Dimensionless aligned spin on companion 2 */
130  REAL8 *chip, /**< [out] Effective spin in the orbital plane */
131  REAL8 *thetaJ, /**< [out] Angle between J0 and line of sight (z-direction) */
132  REAL8 *alpha0, /**< [out] Initial value of alpha angle (azimuthal precession angle) */
133  const REAL8 m1_SI, /**< Mass of companion 1 (kg) */
134  const REAL8 m2_SI, /**< Mass of companion 2 (kg) */
135  const REAL8 f_ref, /**< Reference GW frequency (Hz) */
136  const REAL8 lnhatx, /**< Initial value of LNhatx: orbital angular momentum unit vector */
137  const REAL8 lnhaty, /**< Initial value of LNhaty */
138  const REAL8 lnhatz, /**< Initial value of LNhatz */
139  const REAL8 s1x, /**< Initial value of s1x: dimensionless spin of BH 1 */
140  const REAL8 s1y, /**< Initial value of s1y: dimensionless spin of BH 1 */
141  const REAL8 s1z, /**< Initial value of s1z: dimensionless spin of BH 1 */
142  const REAL8 s2x, /**< Initial value of s2x: dimensionless spin of BH 2 */
143  const REAL8 s2y, /**< Initial value of s2y: dimensionless spin of BH 2 */
144  const REAL8 s2z, /**< Initial value of s2z: dimensionless spin of BH 2 */
145  IMRPhenomP_version_type IMRPhenomP_version /**< IMRPhenomP(v1) uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD */
146 )
147 {
148  // Note that the angle phiJ defined below and alpha0 are degenerate. Therefore we do not output phiJ.
149 
150  /* Check arguments for sanity */
151  XLAL_CHECK(chi1_l != NULL, XLAL_EFAULT);
152  XLAL_CHECK(chi2_l != NULL, XLAL_EFAULT);
153  XLAL_CHECK(chip != NULL, XLAL_EFAULT);
154  XLAL_CHECK(thetaJ != NULL, XLAL_EFAULT);
155  XLAL_CHECK(alpha0 != NULL, XLAL_EFAULT);
156 
157  XLAL_CHECK(f_ref > 0, XLAL_EDOM, "Reference frequency must be positive.\n");
158  XLAL_CHECK(m1_SI > 0, XLAL_EDOM, "m1 must be positive.\n");
159  XLAL_CHECK(m2_SI > 0, XLAL_EDOM, "m2 must be positive.\n");
160  XLAL_CHECK(fabs(lnhatx*lnhatx + lnhaty*lnhaty + lnhatz*lnhatz - 1) < 1e-10, XLAL_EDOM, "Lnhat must be a unit vector.\n");
161  XLAL_CHECK(fabs(s1x*s1x + s1y*s1y + s1z*s1z) <= 1.0, XLAL_EDOM, "|S1/m1^2| must be <= 1.\n");
162  XLAL_CHECK(fabs(s2x*s2x + s2y*s2y + s2z*s2z) <= 1.0, XLAL_EDOM, "|S2/m2^2| must be <= 1.\n");
163 
164  const REAL8 m1 = m1_SI / LAL_MSUN_SI; /* Masses in solar masses */
165  const REAL8 m2 = m2_SI / LAL_MSUN_SI;
166  const REAL8 M = m1+m2;
167  const REAL8 m1_2 = m1*m1;
168  const REAL8 m2_2 = m2*m2;
169  const REAL8 eta = m1 * m2 / (M*M); /* Symmetric mass-ratio */
170 
171  /* Aligned spins */
172  *chi1_l = lnhatx*s1x + lnhaty*s1y + lnhatz*s1z; /* Dimensionless aligned spin on BH 1 */
173  *chi2_l = lnhatx*s2x + lnhaty*s2y + lnhatz*s2z; /* Dimensionless aligned spin on BH 2 */
174 
175  /* Spin components orthogonal to lnhat */
176  const REAL8 S1_perp_x = (s1x - *chi1_l*lnhatx) * m1_2;
177  const REAL8 S1_perp_y = (s1y - *chi1_l*lnhaty) * m1_2;
178  const REAL8 S1_perp_z = (s1z - *chi1_l*lnhatz) * m1_2;
179  const REAL8 S2_perp_x = (s2x - *chi2_l*lnhatx) * m2_2;
180  const REAL8 S2_perp_y = (s2y - *chi2_l*lnhaty) * m2_2;
181  const REAL8 S2_perp_z = (s2z - *chi2_l*lnhatz) * m2_2;
182  const REAL8 S1_perp = sqrt(S1_perp_x*S1_perp_x + S1_perp_y*S1_perp_y + S1_perp_z*S1_perp_z);
183  const REAL8 S2_perp = sqrt(S2_perp_x*S2_perp_x + S2_perp_y*S2_perp_y + S2_perp_z*S2_perp_z);
184 
185  const REAL8 A1 = 2 + (3*m2) / (2*m1);
186  const REAL8 A2 = 2 + (3*m1) / (2*m2);
187  const REAL8 ASp1 = A1*S1_perp;
188  const REAL8 ASp2 = A2*S2_perp;
189  const REAL8 num = (ASp2 > ASp1) ? ASp2 : ASp1;
190  const REAL8 den = (m2 > m1) ? A2*m2_2 : A1*m1_2;
191  *chip = num / den; /* chip = max(A1 Sp1, A2 Sp2) / (A_i m_i^2) for i index of larger BH (See Eqn. 32 in technical document) */
192 
193  /* Compute L, J0 and orientation angles */
194  const REAL8 m_sec = M * LAL_MTSUN_SI; /* Total mass in seconds */
195  const REAL8 piM = LAL_PI * m_sec;
196  const REAL8 v_ref = cbrt(piM * f_ref);
197 
198  REAL8 L0 = 0.0;
199  switch (IMRPhenomP_version) {
200  case IMRPhenomPv1_V:
201  L0 = M*M * L2PNR_v1(v_ref, eta); /* Use 2PN approximation for L. */
202  break;
203  case IMRPhenomPv2_V:
204  L0 = M*M * L2PNR(v_ref, eta); /* Use 2PN approximation for L. */
205  break;
206  default:
207  XLAL_ERROR( XLAL_EINVAL, "Unknown IMRPhenomP version!\nAt present only v1 and v2 are available." );
208  break;
209  }
210 
211  const REAL8 Jx0 = L0 * lnhatx + m1_2*s1x + m2_2*s2x;
212  const REAL8 Jy0 = L0 * lnhaty + m1_2*s1y + m2_2*s2y;
213  const REAL8 Jz0 = L0 * lnhatz + m1_2*s1z + m2_2*s2z;
214  const REAL8 J0 = sqrt(Jx0*Jx0 + Jy0*Jy0 + Jz0*Jz0);
215 
216  /* Compute thetaJ, the angle between J0 and line of sight (z-direction) */
217  if (J0 < 1e-10) {
218  XLAL_PRINT_WARNING("Warning: |J0| < 1e-10. Setting thetaJ = 0.\n");
219  *thetaJ = 0;
220  } else {
221  *thetaJ = acos(Jz0 / J0);
222  }
223 
224  REAL8 phiJ; // We only use this angle internally since it is degenerate with alpha0.
225  phiJ = atan2tol(Jy0, Jx0, MAX_TOL_ATAN); /* Angle of J0 in the plane of the sky */
226  /* Note: Compared to the similar code in SpinTaylorF2 we have defined phiJ as the angle between the positive
227  (rather than the negative) x-axis and the projection of J0, since this is a more natural definition of the angle.
228  We have also renamed the angle from psiJ to phiJ. */
229 
230  /* Rotate Lnhat back to frame where J is along z and the line of sight in the Oxz plane with >0 projection in x, to figure out initial alpha */
231  /* The rotation matrix is
232  {
233  {-cos(thetaJ)*cos(phiJ), -cos(thetaJ)*sin(phiJ), sin(thetaJ)},
234  {sin(phiJ), -cos(phiJ), 0},
235  {cos(phiJ)*sin(thetaJ), sin(thetaJ)*sin(phiJ),cos(thetaJ)}
236  }
237  */
238  const REAL8 rotLx = -lnhatx*cos(*thetaJ)*cos(phiJ) - lnhaty*cos(*thetaJ)*sin(phiJ) + lnhatz*sin(*thetaJ);
239  const REAL8 rotLy = lnhatx*sin(phiJ) - lnhaty*cos(phiJ);
240  *alpha0 = atan2tol(rotLy, rotLx, MAX_TOL_ATAN);
241 
242  return XLAL_SUCCESS;
243 }
244 
245 
246 
247 /**
248  * Function to map LAL parameters
249  * (masses, 6 spin components, phiRef and inclination at f_ref)
250  * (assumed to be in the source frame
251  * where LN points in the z direction
252  * i.e. lnhat = (0,0,1)
253  * and the separation vector n is in the x direction
254  * and the spherical angles of the line of sight N are (incl,Pi/2-phiRef))
255  * into IMRPhenomP intrinsic parameters
256  * (chi1_l, chi2_l, chip, thetaJN, alpha0 and phi_aligned).
257  *
258  * All input masses and frequencies should be in SI units.
259  *
260  * See Fig. 1. in arxiv:1408.1810 for a diagram of the angles.
261  */
263  REAL8 *chi1_l, /**< [out] Dimensionless aligned spin on companion 1 */
264  REAL8 *chi2_l, /**< [out] Dimensionless aligned spin on companion 2 */
265  REAL8 *chip, /**< [out] Effective spin in the orbital plane */
266  REAL8 *thetaJN, /**< [out] Angle between J0 and line of sight (z-direction) */
267  REAL8 *alpha0, /**< [out] Initial value of alpha angle (azimuthal precession angle) */
268  REAL8 *phi_aligned, /**< [out] Initial phase to feed the underlying aligned-spin model */
269  REAL8 *zeta_polariz, /**< [out] Angle to rotate the polarizations */
270  const REAL8 m1_SI, /**< Mass of companion 1 (kg) */
271  const REAL8 m2_SI, /**< Mass of companion 2 (kg) */
272  const REAL8 f_ref, /**< Reference GW frequency (Hz) */
273  const REAL8 phiRef, /**< Reference phase */
274  const REAL8 incl, /**< Inclination : angle between LN and the line of sight */
275  const REAL8 s1x, /**< Initial value of s1x: dimensionless spin of BH 1 */
276  const REAL8 s1y, /**< Initial value of s1y: dimensionless spin of BH 1 */
277  const REAL8 s1z, /**< Initial value of s1z: dimensionless spin of BH 1 */
278  const REAL8 s2x, /**< Initial value of s2x: dimensionless spin of BH 2 */
279  const REAL8 s2y, /**< Initial value of s2y: dimensionless spin of BH 2 */
280  const REAL8 s2z, /**< Initial value of s2z: dimensionless spin of BH 2 */
281  IMRPhenomP_version_type IMRPhenomP_version /**< IMRPhenomP(v1) uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
282 )
283 {
284  // Note that the angle phiJ defined below and alpha0 are degenerate. Therefore we do not output phiJ.
285 
286  /* Check arguments for sanity */
287  XLAL_CHECK(chi1_l != NULL, XLAL_EFAULT);
288  XLAL_CHECK(chi2_l != NULL, XLAL_EFAULT);
289  XLAL_CHECK(chip != NULL, XLAL_EFAULT);
290  XLAL_CHECK(thetaJN != NULL, XLAL_EFAULT);
291  XLAL_CHECK(alpha0 != NULL, XLAL_EFAULT);
292  XLAL_CHECK(phi_aligned != NULL, XLAL_EFAULT);
293 
294  XLAL_CHECK(f_ref > 0, XLAL_EDOM, "Reference frequency must be positive.\n");
295  XLAL_CHECK(m1_SI > 0, XLAL_EDOM, "m1 must be positive.\n");
296  XLAL_CHECK(m2_SI > 0, XLAL_EDOM, "m2 must be positive.\n");
297  XLAL_CHECK(fabs(s1x*s1x + s1y*s1y + s1z*s1z) <= 1.0, XLAL_EDOM, "|S1/m1^2| must be <= 1.\n");
298  XLAL_CHECK(fabs(s2x*s2x + s2y*s2y + s2z*s2z) <= 1.0, XLAL_EDOM, "|S2/m2^2| must be <= 1.\n");
299 
300  const REAL8 m1 = m1_SI / LAL_MSUN_SI; /* Masses in solar masses */
301  const REAL8 m2 = m2_SI / LAL_MSUN_SI;
302  const REAL8 M = m1+m2;
303  const REAL8 m1_2 = m1*m1;
304  const REAL8 m2_2 = m2*m2;
305  const REAL8 eta = m1 * m2 / (M*M); /* Symmetric mass-ratio */
306 
307  /* From the components in the source frame, we can easily determine
308  chi1_l, chi2_l, chip and phi_aligned, which we need to return.
309  We also compute the spherical angles of J,
310  which we need to transform to the J frame*/
311 
312  /* Aligned spins */
313  *chi1_l = s1z; /* Dimensionless aligned spin on BH 1 */
314  *chi2_l = s2z; /* Dimensionless aligned spin on BH 2 */
315 
316  /* Magnitude of the spin projections in the orbital plane */
317  const REAL8 S1_perp = m1_2*sqrt(s1x*s1x + s1y*s1y);
318  const REAL8 S2_perp = m2_2*sqrt(s2x*s2x + s2y*s2y);
319  /* From this we can compute chip*/
320  const REAL8 A1 = 2 + (3*m2) / (2*m1);
321  const REAL8 A2 = 2 + (3*m1) / (2*m2);
322  const REAL8 ASp1 = A1*S1_perp;
323  const REAL8 ASp2 = A2*S2_perp;
324  const REAL8 num = (ASp2 > ASp1) ? ASp2 : ASp1;
325  const REAL8 den = (m2 > m1) ? A2*m2_2 : A1*m1_2;
326  *chip = num / den; /* chip = max(A1 Sp1, A2 Sp2) / (A_i m_i^2) for i index of larger BH (See Eqn. 32 in technical document) */
327 
328  /* Compute L, J0 and orientation angles */
329  const REAL8 m_sec = M * LAL_MTSUN_SI; /* Total mass in seconds */
330  const REAL8 piM = LAL_PI * m_sec;
331  const REAL8 v_ref = cbrt(piM * f_ref);
332 
333  const int ExpansionOrder = 5; // Used in PhenomPv3 only
334 
335  REAL8 L0 = 0.0;
336  switch (IMRPhenomP_version) {
337  case IMRPhenomPv1_V:
338  L0 = M*M * L2PNR_v1(v_ref, eta); /* Use 2PN approximation for L. */
339  break;
340  case IMRPhenomPv2_V:
342  L0 = M*M * L2PNR(v_ref, eta); /* Use 2PN approximation for L. */
343  break;
344  case IMRPhenomPv3_V: /*Pv3 uses 3PN spinning for L but in non-precessing limit uses the simpler L2PNR function */
345  if ((s1x == 0. && s1y == 0. && s2x == 0. && s2y == 0.))
346  { // non-precessing case
347  L0 = M * M * L2PNR(v_ref, eta); /* Use 2PN approximation for L. */
348  } else { // precessing case
349  L0 = M * M * PhenomInternal_OrbAngMom3PN(f_ref / 2., m1_SI, m2_SI, s1x, s1y, s1z, s2x, s2y, s2z, f_ref, ExpansionOrder); /* Use 3PN spinning approximation for L. */
350  }
351  break;
352  default:
353  XLAL_ERROR( XLAL_EINVAL, "Unknown IMRPhenomP version!\nAt present only v1 and v2 are available." );
354  break;
355  }
356  // Below, _sf indicates source frame components. We will also use _Jf for J frame components
357  const REAL8 J0x_sf = m1_2*s1x + m2_2*s2x;
358  const REAL8 J0y_sf = m1_2*s1y + m2_2*s2y;
359  const REAL8 J0z_sf = L0 + m1_2*s1z + m2_2*s2z;
360  const REAL8 J0 = sqrt(J0x_sf*J0x_sf + J0y_sf*J0y_sf + J0z_sf*J0z_sf);
361 
362  /* Compute thetaJ, the angle between J0 and LN (z-direction) */
363  REAL8 thetaJ_sf;
364  if (J0 < 1e-10) {
365  XLAL_PRINT_WARNING("Warning: |J0| < 1e-10. Setting thetaJ = 0.\n");
366  thetaJ_sf = 0;
367  } else {
368  thetaJ_sf = acos(J0z_sf / J0);
369  }
370 
371  REAL8 phiJ_sf;
372  if (fabs(J0x_sf) < MAX_TOL_ATAN && fabs(J0y_sf) < MAX_TOL_ATAN)
373  phiJ_sf = LAL_PI/2. - phiRef; // aligned spin limit
374  else
375  phiJ_sf = atan2(J0y_sf, J0x_sf); /* azimuthal angle of J0 in the source frame */
376 
377  *phi_aligned = - phiJ_sf;
378 
379  /* We now have to rotate to the "J frame" where we can easily
380  compute alpha0, the azimuthal angle of LN,
381  as well as thetaJ, the angle between J and N.
382  The J frame is defined imposing that J points in the z direction
383  and the line of sight N is in the xz plane (with positive projection along x).
384  The components of any vector in the (new) J frame are obtained from those
385  in the (old) source frame by multiplying by RZ[kappa].RY[-thetaJ].RZ[-phiJ]
386  where kappa will be determined by rotating N with RY[-thetaJ].RZ[-phiJ]
387  (which brings J to the z axis) and taking the opposite of azimuthal angle of the rotated N.
388  */
389  REAL8 tmp1,tmp2;
390  // First we determine kappa
391  // in the source frame, the components of N are given in Eq (35c) of T1500606-v6
392  REAL8 Nx_sf = sin(incl)*cos(LAL_PI/2. - phiRef);
393  REAL8 Ny_sf = sin(incl)*sin(LAL_PI/2. - phiRef);
394  REAL8 Nz_sf = cos(incl);
395  REAL8 tmp_x = Nx_sf;
396  REAL8 tmp_y = Ny_sf;
397  REAL8 tmp_z = Nz_sf;
398  ROTATEZ(-phiJ_sf, tmp_x, tmp_y, tmp_z);
399  ROTATEY(-thetaJ_sf, tmp_x, tmp_y, tmp_z);
400  REAL8 kappa;
401  kappa = - atan2tol(tmp_y,tmp_x, MAX_TOL_ATAN);
402 
403  // Then we determine alpha0, by rotating LN
404  tmp_x = 0.;
405  tmp_y = 0.;
406  tmp_z = 1.; // in the source frame, LN=(0,0,1)
407  ROTATEZ(-phiJ_sf, tmp_x, tmp_y, tmp_z);
408  ROTATEY(-thetaJ_sf, tmp_x, tmp_y, tmp_z);
409  ROTATEZ(kappa, tmp_x, tmp_y, tmp_z);
410  if (fabs(tmp_x) < MAX_TOL_ATAN && fabs(tmp_y) < MAX_TOL_ATAN)
411  *alpha0 = LAL_PI; //this is the aligned spin case
412  else
413  *alpha0 = atan2(tmp_y,tmp_x);
414 
415  // Finally we determine thetaJ, by rotating N
416  tmp_x = Nx_sf;
417  tmp_y = Ny_sf;
418  tmp_z = Nz_sf;
419  ROTATEZ(-phiJ_sf, tmp_x, tmp_y, tmp_z);
420  ROTATEY(-thetaJ_sf, tmp_x, tmp_y, tmp_z);
421  ROTATEZ(kappa, tmp_x, tmp_y, tmp_z);
422  REAL8 Nx_Jf = tmp_x; // let's store those two since we will reuse them later (we don't need the y component)
423  REAL8 Nz_Jf = tmp_z;
424  *thetaJN = acos(Nz_Jf); // No normalization needed, we are dealing with a unit vector
425 
426  /* Finally, we need to redefine the polarizations :
427  PhenomP's polarizations are defined following Arun et al (arXiv:0810.5336)
428  i.e. projecting the metric onto the P,Q,N triad defined with P=NxJ/|NxJ| (see (2.6) in there).
429  By contrast, the triad X,Y,N used in LAL
430  ("waveframe" in the nomenclature of T1500606-v6)
431  is defined in e.g. eq (35) of this document
432  (via its components in the source frame; note we use the defautl Omega=Pi/2).
433  Both triads differ from each other by a rotation around N by an angle \zeta
434  and we need to rotate the polarizations accordingly by 2\zeta
435  */
436  REAL8 Xx_sf = -cos(incl)*sin(phiRef);
437  REAL8 Xy_sf = -cos(incl)*cos(phiRef);
438  REAL8 Xz_sf = sin(incl);
439  tmp_x = Xx_sf;
440  tmp_y = Xy_sf;
441  tmp_z = Xz_sf;
442  ROTATEZ(-phiJ_sf, tmp_x, tmp_y, tmp_z);
443  ROTATEY(-thetaJ_sf, tmp_x, tmp_y, tmp_z);
444  ROTATEZ(kappa, tmp_x, tmp_y, tmp_z);
445  //now the tmp_a are the components of X in the J frame
446  //we need the polar angle of that vector in the P,Q basis of Arun et al
447  // P=NxJ/|NxJ| and since we put N in the (pos x)z half plane of the J frame
448  REAL8 PArunx_Jf = 0.;
449  REAL8 PAruny_Jf = -1.;
450  REAL8 PArunz_Jf = 0.;
451  // Q=NxP
452  REAL8 QArunx_Jf = Nz_Jf;
453  REAL8 QAruny_Jf = 0.;
454  REAL8 QArunz_Jf = -Nx_Jf;
455  REAL8 XdotPArun = tmp_x*PArunx_Jf+tmp_y*PAruny_Jf+tmp_z*PArunz_Jf;
456  REAL8 XdotQArun = tmp_x*QArunx_Jf+tmp_y*QAruny_Jf+tmp_z*QArunz_Jf;
457  *zeta_polariz = atan2(XdotQArun , XdotPArun);
458 
459  return XLAL_SUCCESS;
460 }
461 
462 
463 /**
464  * Driver routine to compute the precessing inspiral-merger-ringdown
465  * phenomenological waveform IMRPhenomP in the frequency domain.
466  *
467  * Reference:
468  * - Hannam et al., arXiv:1308.3271 [gr-qc]
469  *
470  * \ref XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame should be called first
471  * to map LAL parameters into IMRPhenomP intrinsic parameters
472  * (chi1_l, chi2_l, chip, thetaJ, alpha0).
473  *
474  * This function can be used for equally-spaced frequency series.
475  * For unequal spacing, use \ref XLALSimIMRPhenomPFrequencySequence instead.
476  */
478  COMPLEX16FrequencySeries **hptilde, /**< [out] Frequency-domain waveform h+ */
479  COMPLEX16FrequencySeries **hctilde, /**< [out] Frequency-domain waveform hx */
480  const REAL8 chi1_l, /**< Dimensionless aligned spin on companion 1 */
481  const REAL8 chi2_l, /**< Dimensionless aligned spin on companion 2 */
482  const REAL8 chip, /**< Effective spin in the orbital plane */
483  const REAL8 thetaJ, /**< Angle between J0 and line of sight (z-direction) */
484  const REAL8 m1_SI, /**< Mass of companion 1 (kg) */
485  const REAL8 m2_SI, /**< Mass of companion 2 (kg) */
486  const REAL8 distance, /**< Distance of source (m) */
487  const REAL8 alpha0, /**< Initial value of alpha angle (azimuthal precession angle) */
488  const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
489  const REAL8 deltaF, /**< Sampling frequency (Hz) */
490  const REAL8 f_min, /**< Starting GW frequency (Hz) */
491  const REAL8 f_max, /**< End frequency; 0 defaults to ringdown cutoff freq */
492  const REAL8 f_ref, /**< Reference frequency */
493  IMRPhenomP_version_type IMRPhenomP_version, /**< IMRPhenomPv1 uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
494  NRTidal_version_type NRTidal_version, /**< either NRTidal or NRTidalv2 for BNS waveform; NoNRT_V for BBH waveform */
495  LALDict *extraParams) /**<linked list that may contain the extra testing GR parameters and/or tidal parameters */
496 {
497  // See Fig. 1. in arxiv:1408.1810 for diagram of the angles.
498  // Note that the angles phiJ which is calculated internally in XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame
499  // and alpha0 are degenerate. Therefore phiJ is not passed to this function.
500 
501  // Use f_min, f_max, deltaF to compute freqs sequence
502  // Instead of building a full sequency we only transfer the boundaries and let
503  // the internal core function do the rest (and properly take care of corner cases).
504  XLAL_CHECK (f_min > 0, XLAL_EDOM, "Minimum frequency must be positive.");
505  XLAL_CHECK (f_max >= 0, XLAL_EDOM, "Maximum frequency must be non-negative.");
506  XLAL_CHECK ( ( f_max == 0 ) || ( f_max > f_min ), XLAL_EDOM, "f_max <= f_min");
508  XLAL_CHECK(freqs != NULL, XLAL_EFAULT);
509  freqs->data[0] = f_min;
510  freqs->data[1] = f_max;
511 
512  int retcode = PhenomPCore(hptilde, hctilde,
513  chi1_l, chi2_l, chip, thetaJ, m1_SI, m2_SI, distance, alpha0, phic, f_ref, freqs, deltaF, IMRPhenomP_version, NRTidal_version, extraParams);
514  XLAL_CHECK(retcode == XLAL_SUCCESS, XLAL_EFUNC, "Failed to generate IMRPhenomP waveform.");
516  return (retcode);
517 }
518 
519 /**
520  * Driver routine to compute the precessing inspiral-merger-ringdown
521  * phenomenological waveform IMRPhenomP in the frequency domain.
522  *
523  * Reference:
524  * - Hannam et al., arXiv:1308.3271 [gr-qc]
525  *
526  * \ref XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame should be called first
527  * to map LAL parameters into IMRPhenomP intrinsic parameters
528  * (chi1_l, chi2_l, chip, thetaJ, alpha0).
529  *
530  * This function can be used for user-specified,
531  * potentially unequally-spaced frequency series.
532  * For equal spacing with a given deltaF, use \ref XLALSimIMRPhenomP instead.
533  */
535  COMPLEX16FrequencySeries **hptilde, /**< [out] Frequency-domain waveform h+ */
536  COMPLEX16FrequencySeries **hctilde, /**< [out] Frequency-domain waveform hx */
537  const REAL8Sequence *freqs, /**< Frequency points at which to evaluate the waveform (Hz) */
538  const REAL8 chi1_l, /**< Dimensionless aligned spin on companion 1 */
539  const REAL8 chi2_l, /**< Dimensionless aligned spin on companion 2 */
540  const REAL8 chip, /**< Effective spin in the orbital plane */
541  const REAL8 thetaJ, /**< Angle between J0 and line of sight (z-direction) */
542  const REAL8 m1_SI, /**< Mass of companion 1 (kg) */
543  const REAL8 m2_SI, /**< Mass of companion 2 (kg) */
544  const REAL8 distance, /**< Distance of source (m) */
545  const REAL8 alpha0, /**< Initial value of alpha angle (azimuthal precession angle) */
546  const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
547  const REAL8 f_ref, /**< Reference frequency */
548  IMRPhenomP_version_type IMRPhenomP_version, /**< IMRPhenomPv1 uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
549  NRTidal_version_type NRTidal_version, /**< either NRTidal or NRTidalv2 for BNS waveform; NoNRT_V for BBH waveform */
550  LALDict *extraParams) /**<linked list that may contain the extra testing GR parameters and/or tidal parameters */
551 {
552  // See Fig. 1. in arxiv:1408.1810 for diagram of the angles.
553  // Note that the angles phiJ which is calculated internally in XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame
554  // and alpha0 are degenerate. Therefore phiJ is not passed to this function.
555 
556  // Call the internal core function with deltaF = 0 to indicate that freqs is non-uniformly
557  // spaced and we want the strain only at these frequencies
558  int retcode = PhenomPCore(hptilde, hctilde,
559  chi1_l, chi2_l, chip, thetaJ, m1_SI, m2_SI, distance, alpha0, phic, f_ref, freqs, 0, IMRPhenomP_version, NRTidal_version, extraParams);
560  XLAL_CHECK(retcode == XLAL_SUCCESS, XLAL_EFUNC, "Failed to generate IMRPhenomP waveform.");
561  return(retcode);
562 }
563 
564 /** @} */
565 /** @} */
566 
567 
568 /**
569  * Internal core function to calculate
570  * plus and cross polarizations of the PhenomP model
571  * for a set of frequencies.
572  * This can handle either user-specified frequency points
573  * or create an equally-spaced frequency series.
574  */
575 static int PhenomPCore(
576  COMPLEX16FrequencySeries **hptilde, /**< [out] Frequency-domain waveform h+ */
577  COMPLEX16FrequencySeries **hctilde, /**< [out] Frequency-domain waveform hx */
578  const REAL8 chi1_l_in, /**< Dimensionless aligned spin on companion 1 */
579  const REAL8 chi2_l_in, /**< Dimensionless aligned spin on companion 2 */
580  const REAL8 chip, /**< Effective spin in the orbital plane */
581  const REAL8 thetaJ, /**< Angle between J0 and line of sight (z-direction) */
582  const REAL8 m1_SI_in, /**< Mass of companion 1 (kg) */
583  const REAL8 m2_SI_in, /**< Mass of companion 2 (kg) */
584  const REAL8 distance, /**< Distance of source (m) */
585  const REAL8 alpha0, /**< Initial value of alpha angle (azimuthal precession angle) */
586  const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
587  const REAL8 f_ref, /**< Reference frequency */
588  const REAL8Sequence *freqs_in, /**< Frequency points at which to evaluate the waveform (Hz) */
589  double deltaF, /**< Sampling frequency (Hz).
590  * If deltaF > 0, the frequency points given in freqs are uniformly spaced with
591  * spacing deltaF. Otherwise, the frequency points are spaced non-uniformly.
592  * Then we will use deltaF = 0 to create the frequency series we return. */
593  IMRPhenomP_version_type IMRPhenomP_version, /**< IMRPhenomPv1 uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
594  NRTidal_version_type NRTidal_version, /**< either NRTidal or NRTidalv2 for BNS waveform; NoNRT_V for BBH waveform */
595  LALDict *extraParams /**<linked list that may contain the extra testing GR parameters and/or tidal parameters */
596  )
597 {
598  /* Check inputs for sanity */
599  XLAL_CHECK(NULL != hptilde, XLAL_EFAULT);
600  XLAL_CHECK(NULL != hctilde, XLAL_EFAULT);
601  XLAL_CHECK(*hptilde == NULL, XLAL_EFAULT);
602  XLAL_CHECK(*hctilde == NULL, XLAL_EFAULT);
603  XLAL_CHECK(deltaF >= 0, XLAL_EDOM, "deltaF must be non-negative.\n");
604  XLAL_CHECK(m1_SI_in > 0, XLAL_EDOM, "m1 must be positive.\n");
605  XLAL_CHECK(m2_SI_in > 0, XLAL_EDOM, "m2 must be positive.\n");
606  XLAL_CHECK(f_ref > 0, XLAL_EDOM, "Reference frequency must be non-negative.\n");
607  XLAL_CHECK(distance > 0, XLAL_EDOM, "distance must be positive.\n");
608  XLAL_CHECK(fabs(chi1_l_in) <= 1.0, XLAL_EDOM, "Aligned spin chi1_l=%g must be <= 1 in magnitude!\n", chi1_l_in);
609  XLAL_CHECK(fabs(chi2_l_in) <= 1.0, XLAL_EDOM, "Aligned spin chi2_l=%g must be <= 1 in magnitude!\n", chi2_l_in);
610  XLAL_CHECK(fabs(chip) <= 1.0, XLAL_EDOM, "In-plane spin chip =%g must be <= 1 in magnitude!\n", chip);
611 
612  // See Fig. 1. in arxiv:1408.1810 for diagram of the angles.
613  // Note that the angles phiJ which is calculated internally in XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame
614  // and alpha0 are degenerate. Therefore phiJ is not passed to this function.
615  /* Phenomenological parameters */
616  IMRPhenomDAmplitudeCoefficients *pAmp = NULL;
617  IMRPhenomDPhaseCoefficients *pPhi = NULL;
618  BBHPhenomCParams *PCparams = NULL;
619  PNPhasingSeries *pn = NULL;
620  // Spline
621  gsl_interp_accel *acc_fixed = NULL;
622  gsl_spline *phiI_fixed = NULL;
623  REAL8Sequence *freqs_fixed = NULL;
624  REAL8Sequence *phase_fixed = NULL;
625  REAL8Sequence *freqs = NULL;
626  int errcode = XLAL_SUCCESS;
627  LALDict *extraParams_in = extraParams;
628  // Tidal corrections
629  REAL8Sequence *phi_tidal = NULL;
630  REAL8Sequence *amp_tidal = NULL;
631  REAL8Sequence *planck_taper = NULL;
632  REAL8Sequence *phi_tidal_fixed = NULL;
633  REAL8Sequence *amp_tidal_fixed = NULL;
634  REAL8Sequence *planck_taper_fixed = NULL;
635  int ret = 0;
636 
637  // Enforce convention m2 >= m1
638  REAL8 chi1_l, chi2_l;
639  REAL8 m1_SI, m2_SI;
640  REAL8 lambda1_in = 0.0;
641  REAL8 lambda2_in = 0.0;
642  REAL8 quadparam1_in = 1.0;
643  REAL8 quadparam2_in = 1.0;
644 
645  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
646  int retcode;
647  retcode = XLALSimInspiralSetQuadMonParamsFromLambdas(extraParams);
648  XLAL_CHECK(retcode == XLAL_SUCCESS, XLAL_EFUNC, "Failed to set quadparams from Universal relation.\n");
649  lambda1_in = XLALSimInspiralWaveformParamsLookupTidalLambda1(extraParams);
650  lambda2_in = XLALSimInspiralWaveformParamsLookupTidalLambda2(extraParams);
651  quadparam1_in = 1. + XLALSimInspiralWaveformParamsLookupdQuadMon1(extraParams);
652  quadparam2_in = 1. + XLALSimInspiralWaveformParamsLookupdQuadMon2(extraParams);
653  }
654 
655  REAL8 lambda1, lambda2;
656  REAL8 quadparam1, quadparam2;
657  /* declare HO 3.5PN spin-spin and spin-cubed terms added later to Pv2_NRTidalv2 */
658  REAL8 SS_3p5PN = 0., SSS_3p5PN = 0.;
659  REAL8 SS_3p5PN_n = 0., SSS_3p5PN_n = 0.;
660 
661  if (m2_SI_in >= m1_SI_in) {
662  m1_SI = m1_SI_in;
663  m2_SI = m2_SI_in;
664  chi1_l = chi1_l_in;
665  chi2_l = chi2_l_in;
666  lambda1 = lambda1_in;
667  lambda2 = lambda2_in;
668  quadparam1 = quadparam1_in;
669  quadparam2 = quadparam2_in;
670  }
671  else { // swap bodies 1 <-> 2
672  m1_SI = m2_SI_in;
673  m2_SI = m1_SI_in;
674  chi1_l = chi2_l_in;
675  chi2_l = chi1_l_in;
676  lambda1 = lambda2_in;
677  lambda2 = lambda1_in;
678  quadparam1 = quadparam2_in;
679  quadparam2 = quadparam1_in;
680  }
681 
683  XLAL_CHECK(XLAL_SUCCESS == errcode, errcode, "init_useful_powers() failed.");
684 
685  /* Find frequency bounds */
686  if (!freqs_in || !freqs_in->data) XLAL_ERROR(XLAL_EFAULT);
687  double f_min = freqs_in->data[0];
688  double f_max = freqs_in->data[freqs_in->length - 1];
689  XLAL_CHECK(f_min > 0, XLAL_EDOM, "Minimum frequency must be positive.\n");
690  XLAL_CHECK(f_max >= 0, XLAL_EDOM, "Maximum frequency must be non-negative.\n");
691 
692  /* External units: SI; internal units: solar masses */
693  const REAL8 m1 = m1_SI / LAL_MSUN_SI;
694  const REAL8 m2 = m2_SI / LAL_MSUN_SI;
695  const REAL8 M = m1 + m2;
696  const REAL8 m_sec = M * LAL_MTSUN_SI; /* Total mass in seconds */
697  REAL8 q = m2 / m1; /* q >= 1 */
698  REAL8 eta = m1 * m2 / (M*M); /* Symmetric mass-ratio */
699  const REAL8 piM = LAL_PI * m_sec;
700  /* New variables needed for the NRTidalv2 model */
701  REAL8 X_A = m1/M;
702  REAL8 X_B = m2/M;
703  REAL8 pn_fac = 3.*pow(piM,2./3.)/(128.*eta);
704 
705  LIGOTimeGPS ligotimegps_zero = LIGOTIMEGPSZERO; // = {0, 0}
706 
707  // Note:
708  // * IMRPhenomP uses chi_eff both in the aligned part and the twisting
709  // * IMRPhenomPv2 uses chi1_l, chi2_l in the aligned part and chi_eff in the twisting
710  const REAL8 chi_eff = (m1*chi1_l + m2*chi2_l) / M; /* Effective aligned spin */
711  const REAL8 chil = (1.0+q)/q * chi_eff; /* dimensionless aligned spin of the largest BH */
712 
713  switch (IMRPhenomP_version) {
714  case IMRPhenomPv1_V:
715  XLAL_PRINT_WARNING("Warning: IMRPhenomP(v1) is unreviewed.\n");
716  if (eta < 0.0453515) /* q = 20 */
717  XLAL_ERROR(XLAL_EDOM, "IMRPhenomP(v1): Mass ratio is way outside the calibration range. m1/m2 should be <= 20.\n");
718  else if (eta < 0.16) /* q = 4 */
719  XLAL_PRINT_WARNING("IMRPhenomP(v1): Warning: The model is only calibrated for m1/m2 <= 4.\n");
720  /* If spins are above 0.9 or below -0.9, throw an error. */
721  /* The rationale behind this is given at this page: https://www.lsc-group.phys.uwm.edu/ligovirgo/cbcnote/WaveformsReview/IMRPhenomCdevel-SanityCheck01 */
722  if (chi_eff > 0.9 || chi_eff < -0.9)
723  XLAL_ERROR(XLAL_EDOM, "IMRPhenomP(v1): Effective spin chi_eff = %g outside the range [-0.9,0.9] is not supported!\n", chi_eff);
724  break;
725  case IMRPhenomPv2_V:
727  if (q > 18.0)
728  XLAL_PRINT_WARNING("IMRPhenomPv2: Warning: The underlying non-precessing model is calibrated up to m1/m2 <= 18.\n");
729  else if (q > 100.0)
730  XLAL_ERROR(XLAL_EDOM, "IMRPhenomPv2: Mass ratio q > 100 which is way outside the calibration range q <= 18.\n");
731  if ((q < 1.5) && (lambda1 > 2000.0) && (lambda2 > 2000.0))
732  XLAL_PRINT_WARNING("NRTidal: Warning: Entering region of parameter space where waveform might not be reliable; q=%g,lambda1=%g, lambda2=%g\n",q, lambda1, lambda2);
733  CheckMaxOpeningAngle(m1, m2, chi1_l, chi2_l, chip);
734  break;
735  default:
736  XLAL_ERROR( XLAL_EINVAL, "Unknown IMRPhenomP version!\nAt present only v1 and v2 are available." );
737  break;
738  }
739 
740  if (eta > 0.25 || q < 1.0) {
741  nudge(&eta, 0.25, 1e-6);
742  nudge(&q, 1.0, 1e-6);
743  }
744 
745  NNLOanglecoeffs angcoeffs; /* Next-to-next-to leading order PN coefficients for Euler angles alpha and epsilon */
746  ComputeNNLOanglecoeffs(&angcoeffs,q,chil,chip);
747 
748  /* Compute the offsets due to the choice of integration constant in alpha and epsilon PN formula */
749  const REAL8 omega_ref = piM * f_ref;
750  const REAL8 logomega_ref = log(omega_ref);
751  const REAL8 omega_ref_cbrt = cbrt(piM * f_ref); // == v0
752  const REAL8 omega_ref_cbrt2 = omega_ref_cbrt*omega_ref_cbrt;
753  const REAL8 alphaNNLOoffset = (angcoeffs.alphacoeff1/omega_ref
754  + angcoeffs.alphacoeff2/omega_ref_cbrt2
755  + angcoeffs.alphacoeff3/omega_ref_cbrt
756  + angcoeffs.alphacoeff4*logomega_ref
757  + angcoeffs.alphacoeff5*omega_ref_cbrt);
758 
759  const REAL8 epsilonNNLOoffset = (angcoeffs.epsiloncoeff1/omega_ref
760  + angcoeffs.epsiloncoeff2/omega_ref_cbrt2
761  + angcoeffs.epsiloncoeff3/omega_ref_cbrt
762  + angcoeffs.epsiloncoeff4*logomega_ref
763  + angcoeffs.epsiloncoeff5*omega_ref_cbrt);
764 
765  /* Compute Ylm's only once and pass them to PhenomPCoreOneFrequency() below. */
767  const REAL8 ytheta = thetaJ;
768  const REAL8 yphi = 0;
769  Y2m.Y2m2 = XLALSpinWeightedSphericalHarmonic(ytheta, yphi, -2, 2, -2);
770  Y2m.Y2m1 = XLALSpinWeightedSphericalHarmonic(ytheta, yphi, -2, 2, -1);
771  Y2m.Y20 = XLALSpinWeightedSphericalHarmonic(ytheta, yphi, -2, 2, 0);
772  Y2m.Y21 = XLALSpinWeightedSphericalHarmonic(ytheta, yphi, -2, 2, 1);
773  Y2m.Y22 = XLALSpinWeightedSphericalHarmonic(ytheta, yphi, -2, 2, 2);
774 
775 
776  REAL8 fCut = 0.0;
777  REAL8 finspin = 0.0;
778  REAL8 f_final = 0.0;
779 
780  switch (IMRPhenomP_version) {
781  case IMRPhenomPv1_V:
782  XLAL_PRINT_INFO("*** IMRPhenomP version 1: based on IMRPhenomC ***");
783  // PhenomC with ringdown using Barausse 2009 formula for final spin
784  PCparams = ComputeIMRPhenomCParamsRDmod(m1, m2, chi_eff, chip, extraParams);
785  if (!PCparams) {
786  errcode = XLAL_EFUNC;
787  goto cleanup;
788  }
789  fCut = PCparams->fCut;
790  f_final = PCparams->fRingDown;
791  break;
792  case IMRPhenomPv2_V:
794  XLAL_PRINT_INFO("*** IMRPhenomP version 2: based on IMRPhenomD ***");
795  // PhenomD uses FinalSpin0815() to calculate the final spin if the spins are aligned.
796  // We use a generalized version of FinalSpin0815() that includes the in-plane spin chip.
797  finspin = FinalSpinIMRPhenomD_all_in_plane_spin_on_larger_BH(m1, m2, chi1_l, chi2_l, chip);
798  if( fabs(finspin) > 1.0 ) {
799  XLAL_PRINT_WARNING("Warning: final spin magnitude %g > 1. Setting final spin magnitude = 1.", finspin);
800  finspin = copysign(1.0, finspin);
801  }
802  // IMRPhenomD assumes that m1 >= m2.
804  ComputeIMRPhenomDAmplitudeCoefficients(pAmp, eta, chi2_l, chi1_l, finspin);
805  pPhi = XLALMalloc(sizeof(IMRPhenomDPhaseCoefficients));
806  ComputeIMRPhenomDPhaseCoefficients(pPhi, eta, chi2_l, chi1_l, finspin, extraParams);
807  if (extraParams==NULL)
808  {
809  extraParams=XLALCreateDict();
810  }
812  // Start making changes here: use XLALSimInspiralWaveformParamsInsertdQuadMon1() function
813  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
816  XLALSimInspiralWaveformParamsInsertdQuadMon1(extraParams, quadparam1-1.);
817  XLALSimInspiralWaveformParamsInsertdQuadMon2(extraParams, quadparam2-1.);
818  XLALSimInspiralTaylorF2AlignedPhasing(&pn, m1, m2, chi1_l, chi2_l, extraParams);
819  } else {
820  XLALSimInspiralTaylorF2AlignedPhasing(&pn, m1, m2, chi1_l, chi2_l, extraParams);
821  }
822 
823  if (!pAmp || !pPhi || !pn) {
824  errcode = XLAL_EFUNC;
825  goto cleanup;
826  }
827 
828  // Subtract 3PN spin-spin term below as this is in LAL's TaylorF2 implementation
829  // (LALSimInspiralPNCoefficients.c -> XLALSimInspiralPNPhasing_F2), but
830  // was not available when PhenomD was tuned.
831  pn->v[6] -= (Subtract3PNSS(m1, m2, M, eta, chi1_l, chi2_l) * pn->v[0]);
832 
833  PhiInsPrefactors phi_prefactors;
834  errcode = init_phi_ins_prefactors(&phi_prefactors, pPhi, pn);
835  XLAL_CHECK(XLAL_SUCCESS == errcode, errcode, "init_phi_ins_prefactors failed");
836 
837  ComputeIMRPhenDPhaseConnectionCoefficients(pPhi, pn, &phi_prefactors, 1.0, 1.0);
838  // This should be the same as the ending frequency in PhenomD
839  fCut = f_CUT / m_sec;
840  f_final = pAmp->fRD / m_sec;
841  break;
842  default:
843  XLALPrintError( "XLAL Error - %s: Unknown IMRPhenomP version!\nAt present only v1 and v2 are available.\n", __func__);
844  errcode = XLAL_EINVAL;
845  goto cleanup;
846  break;
847  }
848 
849  XLAL_CHECK ( fCut > f_min, XLAL_EDOM, "fCut = %.2g/M <= f_min", fCut );
850 
851  /* Default f_max to params->fCut */
852  REAL8 f_max_prime = f_max ? f_max : fCut;
853  f_max_prime = (f_max_prime > fCut) ? fCut : f_max_prime;
854  if (f_max_prime <= f_min){
855  XLALPrintError("XLAL Error - %s: f_max <= f_min\n", __func__);
856  errcode = XLAL_EDOM;
857  goto cleanup;
858  }
859 
860  /* Allocate hp, hc */
861 
862  UINT4 L_fCut = 0; // number of frequency points before we hit fCut
863  size_t n = 0;
864  UINT4 offset = 0; // Index shift between freqs and the frequency series
865  if (deltaF > 0) { // freqs contains uniform frequency grid with spacing deltaF; we start at frequency 0
866  /* Set up output array with size closest power of 2 */
867  if (f_max_prime < f_max) /* Resize waveform if user wants f_max larger than cutoff frequency */
868  n = NextPow2(f_max / deltaF) + 1;
869  else
870  n = NextPow2(f_max_prime / deltaF) + 1;
871 
872  /* coalesce at t=0 */
873  XLAL_CHECK(XLALGPSAdd(&ligotimegps_zero, -1. / deltaF), XLAL_EFUNC,
874  "Failed to shift coalescence time by -1.0/deltaF with deltaF=%g.", deltaF); // shift by overall length in time
875 
876  *hptilde = XLALCreateCOMPLEX16FrequencySeries("hptilde: FD waveform", &ligotimegps_zero, 0.0, deltaF, &lalStrainUnit, n);
877  if(!*hptilde) {
878  errcode = XLAL_ENOMEM;
879  goto cleanup;
880  }
881  *hctilde = XLALCreateCOMPLEX16FrequencySeries("hctilde: FD waveform", &ligotimegps_zero, 0.0, deltaF, &lalStrainUnit, n);
882  if(!*hctilde) {
883  errcode = XLAL_ENOMEM;
884  goto cleanup;
885  }
886 
887  // Recreate freqs using only the lower and upper bounds
888  size_t i_min = (size_t) (f_min / deltaF);
889  size_t i_max = (size_t) (f_max_prime / deltaF);
890  freqs = XLALCreateREAL8Sequence(i_max - i_min);
891  if (!freqs) {
892  errcode = XLAL_EFUNC;
893  XLALPrintError("XLAL Error - %s: Frequency array allocation failed.", __func__);
894  goto cleanup;
895  }
896  for (UINT4 i=i_min; i<i_max; i++)
897  freqs->data[i-i_min] = i*deltaF;
898  L_fCut = freqs->length;
899  offset = i_min;
900  } else { // freqs contains frequencies with non-uniform spacing; we start at lowest given frequency
901  n = freqs_in->length;
902 
903  *hptilde = XLALCreateCOMPLEX16FrequencySeries("hptilde: FD waveform", &ligotimegps_zero, f_min, 0, &lalStrainUnit, n);
904  if(!*hptilde) {
905  XLALPrintError("XLAL Error - %s: Failed to allocate frequency series for hptilde polarization with f_min=%f and %tu bins.", __func__, f_min, n);
906  errcode = XLAL_ENOMEM;
907  goto cleanup;
908  }
909  *hctilde = XLALCreateCOMPLEX16FrequencySeries("hctilde: FD waveform", &ligotimegps_zero, f_min, 0, &lalStrainUnit, n);
910  if(!*hctilde) {
911  XLALPrintError("XLAL Error - %s: Failed to allocate frequency series for hctilde polarization with f_min=%f and %tu bins.", __func__, f_min, n);
912  errcode = XLAL_ENOMEM;
913  goto cleanup;
914  }
915  offset = 0;
916 
917  // Enforce that FS is strictly increasing
918  // (This is needed for phase correction towards the end of this function.)
919  for (UINT4 i=1; i<n; i++)
920  {
921  if (!(freqs_in->data[i] > freqs_in->data[i-1]))
922  {
923  XLALPrintError("XLAL Error - %s: Frequency sequence must be strictly increasing!\n", __func__);
924  errcode = XLAL_EDOM;
925  goto cleanup;
926  }
927  }
928 
929  freqs = XLALCreateREAL8Sequence(n);
930  if (!freqs) {
931  XLALPrintError("XLAL Error - %s: Frequency array allocation failed.", __func__);
932  errcode = XLAL_ENOMEM;
933  goto cleanup;
934  }
935  // Restrict sequence to frequencies <= fCut
936  for (UINT4 i=0; i<n; i++)
937  if (freqs_in->data[i] <= fCut) {
938  freqs->data[i] = freqs_in->data[i];
939  L_fCut++;
940  }
941  }
942 
943 
944  memset((*hptilde)->data->data, 0, n * sizeof(COMPLEX16));
945  memset((*hctilde)->data->data, 0, n * sizeof(COMPLEX16));
946  XLALUnitMultiply(&((*hptilde)->sampleUnits), &((*hptilde)->sampleUnits), &lalSecondUnit);
947  XLALUnitMultiply(&((*hctilde)->sampleUnits), &((*hctilde)->sampleUnits), &lalSecondUnit);
948 
949  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
950  /* Generating the NR tidal amplitude and phase */
951  /* Get FD tidal phase correction and amplitude factor from arXiv:1706.02969 */
952  if (NRTidal_version == NRTidal_V) {
953  phi_tidal = XLALCreateREAL8Sequence(L_fCut);
954  planck_taper = XLALCreateREAL8Sequence(L_fCut);
955  ret = XLALSimNRTunedTidesFDTidalPhaseFrequencySeries(phi_tidal, amp_tidal, planck_taper, freqs, m1_SI, m2_SI, lambda1, lambda2, NRTidal_V);
956  XLAL_CHECK(XLAL_SUCCESS == ret, ret, "XLALSimNRTunedTidesFDTidalPhaseFrequencySeries Failed.");
957  }
958  else if (NRTidal_version == NRTidalv2_V) {
959  phi_tidal = XLALCreateREAL8Sequence(L_fCut);
960  amp_tidal = XLALCreateREAL8Sequence(L_fCut);
961  planck_taper = XLALCreateREAL8Sequence(L_fCut);
962  ret = XLALSimNRTunedTidesFDTidalPhaseFrequencySeries(phi_tidal, amp_tidal, planck_taper, freqs, m1_SI, m2_SI, lambda1, lambda2, NRTidalv2_V);
963  XLAL_CHECK(XLAL_SUCCESS == ret, ret, "XLALSimNRTunedTidesFDTidalPhaseFrequencySeries Failed.");
964  /* Get the PN SS-tail and SSS terms */
965  XLALSimInspiralGetHOSpinTerms(&SS_3p5PN, &SSS_3p5PN, X_A, X_B, chi1_l, chi2_l, quadparam1, quadparam2);
966  }
967  }
968 
969 // phis = XLALMalloc(L_fCut*sizeof(REAL8)); // array for waveform phase
970 // if(!phis) {
971 // errcode = XLAL_ENOMEM;
972 // goto cleanup;
973 // }
974 
975  AmpInsPrefactors amp_prefactors;
976  PhiInsPrefactors phi_prefactors;
977 
978  if (IMRPhenomP_version == IMRPhenomPv2_V || IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
979  errcode = init_amp_ins_prefactors(&amp_prefactors, pAmp);
980  XLAL_CHECK(XLAL_SUCCESS == errcode, errcode, "init_amp_ins_prefactors() failed.");
981  errcode = init_phi_ins_prefactors(&phi_prefactors, pPhi, pn);
982  XLAL_CHECK(XLAL_SUCCESS == errcode, errcode, "init_phi_ins_prefactors() failed.");
983  }
984 
985  /*
986  We can't call XLAL_ERROR() directly with OpenMP on.
987  Keep track of return codes for each thread and in addition use flush to get out of
988  the parallel for loop as soon as possible if something went wrong in any thread.
989  */
990  #pragma omp parallel for
991  for (UINT4 i=0; i<L_fCut; i++) { // loop over frequency points in sequence
992  COMPLEX16 hPhenom = 0.0; // IMRPhenom waveform (before precession) at a given frequency point
993  COMPLEX16 hp_val = 0.0;
994  COMPLEX16 hc_val = 0.0;
995  REAL8 phasing = 0;
996  double f = freqs->data[i];
997  int j = i + offset; // shift index for frequency series if needed
998 
999  int per_thread_errcode=0;
1000 
1001  #pragma omp flush(errcode)
1002  if (errcode != XLAL_SUCCESS)
1003  goto skip;
1004 
1005  /* Generate the waveform */
1006  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
1007  if (NRTidal_version == NRTidal_V) {
1008  double window = planck_taper->data[i];
1009  REAL8 phaseTidal = phi_tidal->data[i];
1010  per_thread_errcode = PhenomPCoreOneFrequency_withTides(f, window, phaseTidal, 0.0, distance, M, phic,
1011  pAmp, pPhi, pn,
1012  &hPhenom, &phasing, &amp_prefactors, &phi_prefactors);
1013  }
1014  else if (NRTidal_version == NRTidalv2_V) {
1015  double ampTidal = amp_tidal->data[i];
1016  double window = planck_taper->data[i];
1017 
1018  /*
1019  Compute the tidal phase correction and add the 3.5PN SS and SSS contributions which
1020  are not incorporated in the TaylorF2 baseline
1021  */
1022  REAL8 phaseTidal = phi_tidal->data[i] + pn_fac*(SS_3p5PN + SSS_3p5PN)*pow(f,2./3.);
1023 
1024  per_thread_errcode = PhenomPCoreOneFrequency_withTides(f, window, phaseTidal, ampTidal, distance, M, phic,
1025  pAmp, pPhi, pn,
1026  &hPhenom, &phasing, &amp_prefactors, &phi_prefactors);
1027  }
1028  } else {
1029  per_thread_errcode = PhenomPCoreOneFrequency(f, eta, distance, M, phic,
1030  pAmp, pPhi, PCparams, pn,
1031  &hPhenom, &phasing, IMRPhenomP_version, &amp_prefactors, &phi_prefactors);
1032  }
1033 
1034  if (per_thread_errcode != XLAL_SUCCESS) {
1035  errcode = per_thread_errcode;
1036  }
1037 
1038  per_thread_errcode = PhenomPCoreTwistUp(f, hPhenom, eta, chi1_l, chi2_l, chip, M,
1039  &angcoeffs, &Y2m,
1040  alphaNNLOoffset - alpha0, epsilonNNLOoffset,
1041  &hp_val, &hc_val, IMRPhenomP_version);
1042 
1043  if (per_thread_errcode != XLAL_SUCCESS) {
1044  errcode = per_thread_errcode;
1045  #pragma omp flush(errcode)
1046  }
1047 
1048  ((*hptilde)->data->data)[j] = hp_val;
1049  ((*hctilde)->data->data)[j] = hc_val;
1050 
1051 // phis[i] = phasing;
1052 
1053  skip: /* this statement intentionally left blank */;
1054  }
1055 
1056 
1057  /* The next part computes and applies a time shift correction to make the waveform peak at t=0 */
1058 
1059  /* Set up spline for phase on fixed grid */
1060  const int n_fixed = 10;
1061  freqs_fixed = XLALCreateREAL8Sequence(n_fixed);
1062  XLAL_CHECK(freqs_fixed != NULL, XLAL_EFAULT);
1063  phase_fixed = XLALCreateREAL8Sequence(n_fixed);
1064  XLAL_CHECK(phase_fixed != NULL, XLAL_EFAULT);
1065 
1066  /* For BNS waveforms, ending frequency is f_merger; putting f_final to f_merger for IMRPhenomPv2_NRTidal and IMRPhenomPv2_NRTidalv2 */
1067  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
1068  REAL8 kappa2T = XLALSimNRTunedTidesComputeKappa2T(m1_SI, m2_SI, lambda1, lambda2);
1069  REAL8 f_merger = XLALSimNRTunedTidesMergerFrequency(M, kappa2T, q);
1070  f_final = f_merger;
1071  }
1072 
1073  /* Set up fixed frequency grid around ringdown frequency for BBH (IMRPhenomPv2) waveforms,
1074  for IMRPhenomPv2_NRTidal and IMRPhenomPv2_NRTidalv2 waveforms, the grid is set up around the merger frequency */
1075  REAL8 freqs_fixed_start = 0.8*f_final;
1076  REAL8 freqs_fixed_stop = 1.2*f_final;
1077  if (freqs_fixed_stop > fCut)
1078  freqs_fixed_stop = fCut;
1079  REAL8 delta_freqs_fixed = (freqs_fixed_stop - freqs_fixed_start) / (n_fixed - 1);
1080  for (int i=0; i < n_fixed; i++)
1081  freqs_fixed->data[i] = freqs_fixed_start + i*delta_freqs_fixed;
1082 
1083  /* Recompute tidal corrections if needed */
1084  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
1085  /* Generating the NR tidal amplitude and phase for the fixed grid */
1086  /* Get FD tidal phase correction and amplitude factor from arXiv:1706.02969 */
1087  phi_tidal_fixed = XLALCreateREAL8Sequence(n_fixed);
1088  amp_tidal_fixed = XLALCreateREAL8Sequence(n_fixed);
1089  planck_taper_fixed = XLALCreateREAL8Sequence(n_fixed);
1090  ret = XLALSimNRTunedTidesFDTidalPhaseFrequencySeries(phi_tidal_fixed, amp_tidal_fixed, planck_taper_fixed,
1091  freqs_fixed, m1_SI, m2_SI, lambda1, lambda2, NRTidal_version);
1092  XLAL_CHECK(XLAL_SUCCESS == ret, ret, "XLALSimNRTunedTidesFDTidalPhaseFrequencySeries Failed.");
1093  if (NRTidal_version == NRTidalv2_V)
1094  XLALSimInspiralGetHOSpinTerms(&SS_3p5PN_n, &SSS_3p5PN_n, X_A, X_B, chi1_l, chi2_l, quadparam1, quadparam2);
1095  }
1096 
1097  /* We need another loop to generate the phase values on the fixed grid; no need for OpenMP here */
1098  for (int i=0; i<n_fixed; i++) { // loop over frequency points in sequence
1099  COMPLEX16 hPhenom = 0.0; // IMRPhenom waveform (before precession) at a given frequency point
1100  REAL8 phasing = 0;
1101  double f = freqs_fixed->data[i];
1102  int per_thread_errcode = 0;
1103 
1104  if (errcode != XLAL_SUCCESS)
1105  goto skip_fixed;
1106 
1107  /* Generate the waveform */
1108  if (IMRPhenomP_version == IMRPhenomPv2NRTidal_V) {
1109  if (NRTidal_version == NRTidal_V) {
1110 
1111  double window = planck_taper_fixed->data[i];
1112  REAL8 phaseTidal = phi_tidal_fixed->data[i];
1113  per_thread_errcode = PhenomPCoreOneFrequency_withTides(f, window, phaseTidal, 0.0, distance, M, phic,
1114  pAmp, pPhi, pn,
1115  &hPhenom, &phasing, &amp_prefactors, &phi_prefactors);
1116  }
1117  else if (NRTidal_version == NRTidalv2_V) {
1118  double ampTidal = amp_tidal_fixed->data[i];
1119  double window = planck_taper_fixed->data[i];
1120  REAL8 phaseTidal = phi_tidal_fixed->data[i] + pn_fac*(SS_3p5PN_n + SSS_3p5PN_n)*pow(f,2./3.);
1121  per_thread_errcode = PhenomPCoreOneFrequency_withTides(f, window, phaseTidal, ampTidal, distance, M, phic,
1122  pAmp, pPhi, pn,
1123  &hPhenom, &phasing, &amp_prefactors, &phi_prefactors);
1124 
1125  }
1126  } else {
1127  per_thread_errcode = PhenomPCoreOneFrequency(f, eta, distance, M, phic,
1128  pAmp, pPhi, PCparams, pn,
1129  &hPhenom, &phasing, IMRPhenomP_version, &amp_prefactors, &phi_prefactors);
1130  }
1131 
1132  if (per_thread_errcode != XLAL_SUCCESS) {
1133  errcode = per_thread_errcode;
1134  }
1135 
1136  phase_fixed->data[i] = phasing;
1137 
1138  skip_fixed: /* this statement intentionally left blank */;
1139  }
1140 
1141  /* Correct phasing so we coalesce at t=0 (with the definition of the epoch=-1/deltaF above) */
1142  /* We apply the same time shift to hptilde and hctilde based on the overall phasing returned by PhenomPCoreOneFrequency */
1143 
1144  /* Set up spline for phase */
1145  acc_fixed = gsl_interp_accel_alloc();
1146  phiI_fixed = gsl_spline_alloc(gsl_interp_cspline, n_fixed);
1147  XLAL_CHECK(phiI_fixed, XLAL_ENOMEM, "Failed to allocate GSL spline with %d points for phase.", n_fixed);
1148  REAL8 t_corr_fixed = 0.0;
1149  gsl_spline_init(phiI_fixed, freqs_fixed->data, phase_fixed->data, n_fixed);
1150  /* Time correction is t(f_final) = 1/(2pi) dphi/df (f_final) */
1151  t_corr_fixed = gsl_spline_eval_deriv(phiI_fixed, f_final, acc_fixed) / (2*LAL_PI);
1152  // XLAL_PRINT_INFO("t_corr (fixed grid): %g\n", t_corr_fixed);
1153 
1154  /* Now correct phase */
1155  for (UINT4 i=0; i<L_fCut; i++) { // loop over frequency points in user-specified sequence
1156  double f = freqs->data[i];
1157  COMPLEX16 phase_corr = (cos(2*LAL_PI * f * t_corr_fixed) - I*sin(2*LAL_PI * f * t_corr_fixed));
1158  int j = i + offset; // shift index for frequency series if needed
1159  ((*hptilde)->data->data)[j] *= phase_corr;
1160  ((*hctilde)->data->data)[j] *= phase_corr;
1161  }
1162 
1163 
1164  cleanup:
1165  if (freqs_fixed) XLALDestroyREAL8Sequence(freqs_fixed);
1166  if (phase_fixed) XLALDestroyREAL8Sequence(phase_fixed);
1167  if (phiI_fixed) gsl_spline_free(phiI_fixed);
1168  if (acc_fixed) gsl_interp_accel_free(acc_fixed);
1169 
1170  if(PCparams) XLALFree(PCparams);
1171  if(pAmp) XLALFree(pAmp);
1172  if(pPhi) XLALFree(pPhi);
1173  if(pn) XLALFree(pn);
1174 
1175  /* If extraParams was allocated in this function and not passed in
1176  * we need to free it to prevent a leak */
1177  if(extraParams && !extraParams_in) XLALDestroyDict(extraParams);
1178 
1179  if(freqs) XLALDestroyREAL8Sequence(freqs);
1180 
1181  if (phi_tidal) XLALDestroyREAL8Sequence(phi_tidal);
1182  if (amp_tidal) XLALDestroyREAL8Sequence(amp_tidal);
1183  if (planck_taper) XLALDestroyREAL8Sequence(planck_taper);
1184  if (amp_tidal_fixed) XLALDestroyREAL8Sequence(amp_tidal_fixed);
1185  if (phi_tidal_fixed) XLALDestroyREAL8Sequence(phi_tidal_fixed);
1186  if (planck_taper_fixed) XLALDestroyREAL8Sequence(planck_taper_fixed);
1187 
1188  if( errcode != XLAL_SUCCESS ) {
1189  if(*hptilde) {
1191  *hptilde=NULL;
1192  }
1193  if(*hctilde) {
1195  *hctilde=NULL;
1196  }
1197  XLAL_ERROR(errcode);
1198  }
1199  else {
1200  return XLAL_SUCCESS;
1201  }
1202 }
1203 
1204 /* ***************************** PhenomP internal functions *********************************/
1205 
1206 /**
1207  * \f[
1208  * \newcommand{\hP}{h^\mathrm{P}}
1209  * \newcommand{\PAmp}{A^\mathrm{P}}
1210  * \newcommand{\PPhase}{\phi^\mathrm{P}}
1211  * \newcommand{\chieff}{\chi_\mathrm{eff}}
1212  * \newcommand{\chip}{\chi_\mathrm{p}}
1213  * \f]
1214  * Internal core function to calculate
1215  * plus and cross polarizations of the PhenomP model
1216  * for a single frequency.
1217  *
1218  * The general expression for the modes \f$\hP_{2m}(t)\f$
1219  * is given by Eq. 1 of arXiv:1308.3271.
1220  * We calculate the frequency domain l=2 plus and cross polarizations separately
1221  * for each m = -2, ... , 2.
1222  *
1223  * The expression of the polarizations times the \f$Y_{lm}\f$
1224  * in code notation are:
1225  * \f{equation*}{
1226  * \left(\tilde{h}_{2m}\right)_+ = e^{-2i \epsilon}
1227  * \left(e^{-i m \alpha} d^2_{-2,m} (-2Y_{2m})
1228  * + e^{+i m \alpha} d^2_{2,m} (-2Y_{2m})^*\right) \cdot \hP / 2 \,,
1229  * \f}
1230  * \f{equation*}{
1231  * \left(\tilde{h}_{2m}\right)_x = e^{-2i \epsilon}
1232  * \left(e^{-i m \alpha} d^2_{-2,m} (-2Y_{2m})
1233  * - e^{+i m \alpha} d^2_{2,m} (-2Y_{2m})^*\right) \cdot \hP / 2 \,,
1234  * \f}
1235  * where the \f$d^l_{m',m}\f$ are Wigner d-matrices evaluated at \f$-\beta\f$,
1236  * and \f$\hP\f$ is the Phenom[C,D] frequency domain model:
1237  * \f{equation*}{
1238  * \hP(f) = \PAmp(f) e^{-i \PPhase(f)} \,.
1239  * \f}
1240  *
1241  * Note that in arXiv:1308.3271, the angle \f$\beta\f$ (beta) is called iota.
1242  *
1243  * For IMRPhenomP(v1) we put all spin on the larger BH,
1244  * convention: \f$m_2 \geq m_1\f$.
1245  * Hence:
1246  * \f{eqnarray*}{
1247  * \chieff &=& \left( m_1 \cdot \chi_1 + m_2 \cdot \chi_2 \right)/M \,,\\
1248  * \chi_l &=& \chieff / m_2 \quad (\text{for } M=1) \,,\\
1249  * S_L &=& m_2^2 \chi_l = m_2 \cdot M \cdot \chieff
1250  * = \frac{q}{1+q} \cdot \chieff \quad (\text{for } M=1) \,.
1251  * \f}
1252  *
1253  * For IMRPhenomPv2 we use both aligned spins:
1254  * \f{equation*}{
1255  * S_L = \chi_1 \cdot m_1^2 + \chi_2 \cdot m_2^2 \,.
1256  * \f}
1257  *
1258  * For both IMRPhenomP(v1) and IMRPhenomPv2 we put the in-plane spin on the larger BH:
1259  * \f{equation*}{
1260  * S_\mathrm{perp} = \chip \cdot m_2^2
1261  * \f}
1262  * (perpendicular spin).
1263  */
1265  const REAL8 fHz, /**< Frequency (Hz) */
1266  const REAL8 eta, /**< Symmetric mass ratio */
1267  const REAL8 distance, /**< Distance of source (m) */
1268  const REAL8 M, /**< Total mass (Solar masses) */
1269  const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
1270  IMRPhenomDAmplitudeCoefficients *pAmp, /**< Internal IMRPhenomD amplitude coefficients */
1271  IMRPhenomDPhaseCoefficients *pPhi, /**< Internal IMRPhenomD phase coefficients */
1272  BBHPhenomCParams *PCparams, /**< Internal PhenomC parameters */
1273  PNPhasingSeries *PNparams, /**< PN inspiral phase coefficients */
1274  COMPLEX16 *hPhenom, /**< IMRPhenom waveform (before precession) */
1275  REAL8 *phasing, /**< [out] overall phasing */
1276  IMRPhenomP_version_type IMRPhenomP_version, /**< IMRPhenomP(v1) uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
1277  AmpInsPrefactors *amp_prefactors, /**< pre-calculated (cached for saving runtime) coefficients for amplitude. See LALSimIMRPhenomD_internals.c*/
1278  PhiInsPrefactors *phi_prefactors /**< pre-calculated (cached for saving runtime) coefficients for phase. See LALSimIMRPhenomD_internals.*/
1279 )
1280 {
1281  XLAL_CHECK(hPhenom != NULL, XLAL_EFAULT);
1282  XLAL_CHECK(phasing != NULL, XLAL_EFAULT);
1283 
1284  REAL8 f = fHz*LAL_MTSUN_SI*M; /* Frequency in geometric units */
1285 
1286  REAL8 aPhenom = 0.0;
1287  REAL8 phPhenom = 0.0;
1288  int errcode = XLAL_SUCCESS;
1289  UsefulPowers powers_of_f;
1290 
1291  /* Calculate Phenom amplitude and phase for a given frequency. */
1292  switch (IMRPhenomP_version) {
1293  case IMRPhenomPv1_V:
1294  XLAL_CHECK(PCparams != NULL, XLAL_EFAULT);
1295  errcode = IMRPhenomCGenerateAmpPhase( &aPhenom, &phPhenom, fHz, eta, PCparams );
1296  break;
1297  case IMRPhenomPv2_V:
1298  XLAL_CHECK(pAmp != NULL, XLAL_EFAULT);
1299  XLAL_CHECK(pPhi != NULL, XLAL_EFAULT);
1300  XLAL_CHECK(PNparams != NULL, XLAL_EFAULT);
1301  XLAL_CHECK(amp_prefactors != NULL, XLAL_EFAULT);
1302  XLAL_CHECK(phi_prefactors != NULL, XLAL_EFAULT);
1303  errcode = init_useful_powers(&powers_of_f, f);
1304  XLAL_CHECK(errcode == XLAL_SUCCESS, errcode, "init_useful_powers failed for f");
1305  aPhenom = IMRPhenDAmplitude(f, pAmp, &powers_of_f, amp_prefactors);
1306  phPhenom = IMRPhenDPhase(f, pPhi, PNparams, &powers_of_f, phi_prefactors, 1.0, 1.0);
1307  break;
1308  case IMRPhenomPv2NRTidal_V:
1309  XLAL_ERROR( XLAL_EINVAL, "Only v1 and v2 are valid IMRPhenomP versions here! The tidal version of IMRPhenomPv2 uses a separate internal function function to generate polarizations and phasing." );
1310  break;
1311  default:
1312  XLAL_ERROR( XLAL_EINVAL, "Unknown IMRPhenomP version!\nAt present only v1 and v2 are available." );
1313  break;
1314  }
1315 
1316  phPhenom -= 2.*phic; /* Note: phic is orbital phase */
1317  REAL8 amp0 = M * LAL_MRSUN_SI * M * LAL_MTSUN_SI / distance;
1318  *hPhenom = amp0 * aPhenom * (cos(phPhenom) - I*sin(phPhenom));//cexp(-I*phPhenom); /* Assemble IMRPhenom waveform. */
1319 
1320  // Return phasing for time-shift correction
1321  *phasing = -phPhenom; // ignore alpha and epsilon contributions
1322 
1323  return XLAL_SUCCESS;
1324 }
1325 
1327  const REAL8 fHz, /**< Frequency (Hz) */
1328  const REAL8 window, /**< Planck_taper */
1329  const REAL8 phaseTidal, /**< tidal phasing at a frequency sample from NRTidal infrastructure*/
1330  const REAL8 ampTidal, /**< tidal amplitude added to BBH amplitude, before Planck tapering*/
1331  const REAL8 distance, /**< Distance of source (m) */
1332  const REAL8 M, /**< Total mass (Solar masses) */
1333  const REAL8 phic, /**< Orbital phase at the peak of the underlying non precessing model (rad) */
1334  IMRPhenomDAmplitudeCoefficients *pAmp, /**< Internal IMRPhenomD amplitude coefficients */
1335  IMRPhenomDPhaseCoefficients *pPhi, /**< Internal IMRPhenomD phase coefficients */
1336  PNPhasingSeries *PNparams, /**< PN inspiral phase coefficients */
1337  COMPLEX16 *hPhenom, /**< [out] IMRPhenom waveform (before precession) */
1338  REAL8 *phasing, /**< [out] overall phasing */
1339  AmpInsPrefactors *amp_prefactors, /**< pre-calculated (cached for saving runtime) coefficients for amplitude. See LALSimIMRPhenomD_internals.c*/
1340  PhiInsPrefactors *phi_prefactors /**< pre-calculated (cached for saving runtime) coefficients for phase. See LALSimIMRPhenomD_internals.*/
1341 )
1342 {
1343 
1344  XLAL_CHECK(hPhenom != NULL, XLAL_EFAULT);
1345  XLAL_CHECK(phasing != NULL, XLAL_EFAULT);
1346 
1347  REAL8 f = fHz*LAL_MTSUN_SI*M; /* Frequency in geometric units */
1348 
1349  REAL8 aPhenom = 0.0;
1350  REAL8 phPhenom = 0.0;
1351  int errcode = XLAL_SUCCESS;
1352  UsefulPowers powers_of_f;
1353 
1354  /* Calculate Phenom amplitude and phase for a given frequency. */
1355  XLAL_CHECK(pAmp != NULL, XLAL_EFAULT);
1356  XLAL_CHECK(pPhi != NULL, XLAL_EFAULT);
1357  XLAL_CHECK(PNparams != NULL, XLAL_EFAULT);
1358  XLAL_CHECK(amp_prefactors != NULL, XLAL_EFAULT);
1359  XLAL_CHECK(phi_prefactors != NULL, XLAL_EFAULT);
1360  errcode = init_useful_powers(&powers_of_f, f);
1361  XLAL_CHECK(errcode == XLAL_SUCCESS, errcode, "init_useful_powers failed for f");
1362  aPhenom = IMRPhenDAmplitude(f, pAmp, &powers_of_f, amp_prefactors);
1363  phPhenom = IMRPhenDPhase(f, pPhi, PNparams, &powers_of_f, phi_prefactors, 1.0, 1.0);
1364 
1365  phPhenom -= 2.*phic; /* Note: phic is orbital phase */
1366  REAL8 amp0 = M * LAL_MRSUN_SI * M * LAL_MTSUN_SI / distance;
1367  *hPhenom = amp0 * (aPhenom + 2*sqrt(LAL_PI/5.) * ampTidal) * cexp(-I*(phPhenom+phaseTidal)) * window; /* Assemble IMRPhenom waveform. */
1368 
1369  // Return phasing for time-shift correction
1370  phPhenom += phaseTidal;
1371  *phasing = -phPhenom; // ignore alpha and epsilon contributions
1372 
1373  return XLAL_SUCCESS;
1374 }
1375 
1377  const REAL8 fHz, /**< Frequency (Hz) */
1378  COMPLEX16 hPhenom, /**< [in] IMRPhenom waveform (before precession) */
1379  const REAL8 eta, /**< Symmetric mass ratio */
1380  const REAL8 chi1_l, /**< Dimensionless aligned spin on companion 1 */
1381  const REAL8 chi2_l, /**< Dimensionless aligned spin on companion 2 */
1382  const REAL8 chip, /**< Dimensionless spin in the orbital plane */
1383  const REAL8 M, /**< Total mass (Solar masses) */
1384  NNLOanglecoeffs *angcoeffs, /**< Struct with PN coeffs for the NNLO angles */
1385  SpinWeightedSphericalHarmonic_l2 *Y2m, /**< Struct of l=2 spherical harmonics of spin weight -2 */
1386  const REAL8 alphaoffset, /**< f_ref dependent offset for alpha angle (azimuthal precession angle) */
1387  const REAL8 epsilonoffset, /**< f_ref dependent offset for epsilon angle */
1388  COMPLEX16 *hp, /**< [out] plus polarization \f$\tilde h_+\f$ */
1389  COMPLEX16 *hc, /**< [out] cross polarization \f$\tilde h_x\f$ */
1390  IMRPhenomP_version_type IMRPhenomP_version /**< IMRPhenomP(v1) uses IMRPhenomC, IMRPhenomPv2 uses IMRPhenomD, IMRPhenomPv2_NRTidal uses NRTidal framework with IMRPhenomPv2 */
1391 )
1392 {
1393 
1394  XLAL_CHECK(angcoeffs != NULL, XLAL_EFAULT);
1395  XLAL_CHECK(hp != NULL, XLAL_EFAULT);
1396  XLAL_CHECK(hc != NULL, XLAL_EFAULT);
1397  XLAL_CHECK(Y2m != NULL, XLAL_EFAULT);
1398 
1399  REAL8 f = fHz*LAL_MTSUN_SI*M; /* Frequency in geometric units */
1400 
1401  const REAL8 q = (1.0 + sqrt(1.0 - 4.0*eta) - 2.0*eta)/(2.0*eta);
1402  const REAL8 m1 = 1.0/(1.0+q); /* Mass of the smaller BH for unit total mass M=1. */
1403  const REAL8 m2 = q/(1.0+q); /* Mass of the larger BH for unit total mass M=1. */
1404  const REAL8 Sperp = chip*(m2*m2); /* Dimensionfull spin component in the orbital plane. S_perp = S_2_perp */
1405  REAL8 SL; /* Dimensionfull aligned spin. */
1406  const REAL8 chi_eff = (m1*chi1_l + m2*chi2_l); /* effective spin for M=1 */
1407 
1408  /* Calculate dimensionfull spins */
1409  switch (IMRPhenomP_version) {
1410  case IMRPhenomPv1_V:
1411  SL = chi_eff*m2; /* Dimensionfull aligned spin of the largest BH. SL = m2^2 chil = m2*M*chi_eff */
1412  break;
1413  case IMRPhenomPv2_V:
1414  case IMRPhenomPv2NRTidal_V:
1415  SL = chi1_l*m1*m1 + chi2_l*m2*m2; /* Dimensionfull aligned spin. */
1416  break; default:
1417  XLAL_ERROR( XLAL_EINVAL, "Unknown IMRPhenomP version!\nAt present only v1 and v2 and tidal are available." );
1418  break;
1419  }
1420 
1421  /* Compute PN NNLO angles */
1422  const REAL8 omega = LAL_PI * f;
1423  const REAL8 logomega = log(omega);
1424  const REAL8 omega_cbrt = cbrt(omega);
1425  const REAL8 omega_cbrt2 = omega_cbrt*omega_cbrt;
1426 
1427  REAL8 alpha = (angcoeffs->alphacoeff1/omega
1428  + angcoeffs->alphacoeff2/omega_cbrt2
1429  + angcoeffs->alphacoeff3/omega_cbrt
1430  + angcoeffs->alphacoeff4*logomega
1431  + angcoeffs->alphacoeff5*omega_cbrt) - alphaoffset;
1432 
1433  REAL8 epsilon = (angcoeffs->epsiloncoeff1/omega
1434  + angcoeffs->epsiloncoeff2/omega_cbrt2
1435  + angcoeffs->epsiloncoeff3/omega_cbrt
1436  + angcoeffs->epsiloncoeff4*logomega
1437  + angcoeffs->epsiloncoeff5*omega_cbrt) - epsilonoffset;
1438 
1439  /* Calculate intermediate expressions cos(beta/2), sin(beta/2) and powers thereof for Wigner d's. */
1440  REAL8 cBetah, sBetah; /* cos(beta/2), sin(beta/2) */
1441  switch (IMRPhenomP_version) {
1442  case IMRPhenomPv1_V:
1443  WignerdCoefficients_SmallAngleApproximation(&cBetah, &sBetah, omega_cbrt, SL, eta, Sperp);
1444  break;
1445  case IMRPhenomPv2_V:
1446  case IMRPhenomPv2NRTidal_V:
1447  WignerdCoefficients(&cBetah, &sBetah, omega_cbrt, SL, eta, Sperp);
1448  break;
1449  default:
1450  XLAL_ERROR( XLAL_EINVAL, " Unknown IMRPhenomP version!\nAt present only v1 and v2 and tidal are available." );
1451  break;
1452  }
1453 
1454  const REAL8 cBetah2 = cBetah*cBetah;
1455  const REAL8 cBetah3 = cBetah2*cBetah;
1456  const REAL8 cBetah4 = cBetah3*cBetah;
1457  const REAL8 sBetah2 = sBetah*sBetah;
1458  const REAL8 sBetah3 = sBetah2*sBetah;
1459  const REAL8 sBetah4 = sBetah3*sBetah;
1460 
1461  /* Compute Wigner d coefficients
1462  The expressions below agree with refX [Goldstein?] and Mathematica
1463  d2 = Table[WignerD[{2, mp, 2}, 0, -\[Beta], 0], {mp, -2, 2}]
1464  dm2 = Table[WignerD[{2, mp, -2}, 0, -\[Beta], 0], {mp, -2, 2}]
1465  */
1466  COMPLEX16 d2[5] = {sBetah4, 2*cBetah*sBetah3, sqrt_6*sBetah2*cBetah2, 2*cBetah3*sBetah, cBetah4};
1467  COMPLEX16 dm2[5] = {d2[4], -d2[3], d2[2], -d2[1], d2[0]}; /* Exploit symmetry d^2_{-2,-m} = (-1)^m d^2_{2,m} */
1468 
1469  COMPLEX16 Y2mA[5] = {Y2m->Y2m2, Y2m->Y2m1, Y2m->Y20, Y2m->Y21, Y2m->Y22};
1470  COMPLEX16 hp_sum = 0;
1471  COMPLEX16 hc_sum = 0;
1472 
1473  /* Sum up contributions to \tilde h+ and \tilde hx */
1474  /* Precompute powers of e^{i m alpha} */
1475  COMPLEX16 cexp_i_alpha = cexp(+I*alpha);
1476  COMPLEX16 cexp_2i_alpha = cexp_i_alpha*cexp_i_alpha;
1477  COMPLEX16 cexp_mi_alpha = 1.0/cexp_i_alpha;
1478  COMPLEX16 cexp_m2i_alpha = cexp_mi_alpha*cexp_mi_alpha;
1479  COMPLEX16 cexp_im_alpha[5] = {cexp_m2i_alpha, cexp_mi_alpha, 1.0, cexp_i_alpha, cexp_2i_alpha};
1480  for(int m=-2; m<=2; m++) {
1481  COMPLEX16 T2m = cexp_im_alpha[-m+2] * dm2[m+2] * Y2mA[m+2]; /* = cexp(-I*m*alpha) * dm2[m+2] * Y2mA[m+2] */
1482  COMPLEX16 Tm2m = cexp_im_alpha[m+2] * d2[m+2] * conj(Y2mA[m+2]); /* = cexp(+I*m*alpha) * d2[m+2] * conj(Y2mA[m+2]) */
1483  hp_sum += T2m + Tm2m;
1484  hc_sum += +I*(T2m - Tm2m);
1485  }
1486 
1487  COMPLEX16 eps_phase_hP = cexp(-2*I*epsilon) * hPhenom / 2.0;
1488  *hp = eps_phase_hP * hp_sum;
1489  *hc = eps_phase_hP * hc_sum;
1490 
1491  return XLAL_SUCCESS;
1492 
1493 }
1494 
1495 /**
1496  * Next-to-next-to-leading order PN coefficients
1497  * for Euler angles \f$\alpha\f$ and \f$\epsilon\f$.
1498  */
1500  NNLOanglecoeffs *angcoeffs, /**< [out] Structure to store results */
1501  const REAL8 q, /**< Mass-ratio (convention q>1) */
1502  const REAL8 chil, /**< Dimensionless aligned spin of the largest BH */
1503  const REAL8 chip) /**< Dimensionless spin component in the orbital plane */
1504 {
1505  const REAL8 m2 = q/(1. + q);
1506  const REAL8 m1 = 1./(1. + q);
1507  const REAL8 dm = m1 - m2;
1508  const REAL8 mtot = 1.;
1509  const REAL8 eta = m1*m2; /* mtot = 1 */
1510  const REAL8 eta2 = eta*eta;
1511  const REAL8 eta3 = eta2*eta;
1512  const REAL8 eta4 = eta3*eta;
1513  const REAL8 mtot2 = mtot*mtot;
1514  const REAL8 mtot4 = mtot2*mtot2;
1515  const REAL8 mtot6 = mtot4*mtot2;
1516  const REAL8 mtot8 = mtot6*mtot2;
1517  const REAL8 chil2 = chil*chil;
1518  const REAL8 chip2 = chip*chip;
1519  const REAL8 chip4 = chip2*chip2;
1520  const REAL8 dm2 = dm*dm;
1521  const REAL8 dm3 = dm2*dm;
1522  const REAL8 m2_2 = m2*m2;
1523  const REAL8 m2_3 = m2_2*m2;
1524  const REAL8 m2_4 = m2_3*m2;
1525  const REAL8 m2_5 = m2_4*m2;
1526  const REAL8 m2_6 = m2_5*m2;
1527  const REAL8 m2_7 = m2_6*m2;
1528  const REAL8 m2_8 = m2_7*m2;
1529 
1530  XLAL_CHECK_VOID(angcoeffs != NULL, XLAL_EFAULT);
1531 
1532  angcoeffs->alphacoeff1 = (-0.18229166666666666 - (5*dm)/(64.*m2));
1533 
1534  angcoeffs->alphacoeff2 = ((-15*dm*m2*chil)/(128.*mtot2*eta) - (35*m2_2*chil)/(128.*mtot2*eta));
1535 
1536  angcoeffs->alphacoeff3 = (-1.7952473958333333 - (4555*dm)/(7168.*m2) -
1537  (15*chip2*dm*m2_3)/(128.*mtot4*eta2) -
1538  (35*chip2*m2_4)/(128.*mtot4*eta2) - (515*eta)/384. - (15*dm2*eta)/(256.*m2_2) -
1539  (175*dm*eta)/(256.*m2));
1540 
1541  angcoeffs->alphacoeff4 = - (35*LAL_PI)/48. - (5*dm*LAL_PI)/(16.*m2) +
1542  (5*dm2*chil)/(16.*mtot2) + (5*dm*m2*chil)/(3.*mtot2) +
1543  (2545*m2_2*chil)/(1152.*mtot2) -
1544  (5*chip2*dm*m2_5*chil)/(128.*mtot6*eta3) -
1545  (35*chip2*m2_6*chil)/(384.*mtot6*eta3) + (2035*dm*m2*chil)/(21504.*mtot2*eta) +
1546  (2995*m2_2*chil)/(9216.*mtot2*eta);
1547 
1548  angcoeffs->alphacoeff5 = (4.318908476114694 + (27895885*dm)/(2.1676032e7*m2) -
1549  (15*chip4*dm*m2_7)/(512.*mtot8*eta4) -
1550  (35*chip4*m2_8)/(512.*mtot8*eta4) -
1551  (485*chip2*dm*m2_3)/(14336.*mtot4*eta2) +
1552  (475*chip2*m2_4)/(6144.*mtot4*eta2) +
1553  (15*chip2*dm2*m2_2)/(256.*mtot4*eta) + (145*chip2*dm*m2_3)/(512.*mtot4*eta) +
1554  (575*chip2*m2_4)/(1536.*mtot4*eta) + (39695*eta)/86016. + (1615*dm2*eta)/(28672.*m2_2) -
1555  (265*dm*eta)/(14336.*m2) + (955*eta2)/576. + (15*dm3*eta2)/(1024.*m2_3) +
1556  (35*dm2*eta2)/(256.*m2_2) + (2725*dm*eta2)/(3072.*m2) - (15*dm*m2*LAL_PI*chil)/(16.*mtot2*eta) -
1557  (35*m2_2*LAL_PI*chil)/(16.*mtot2*eta) + (15*chip2*dm*m2_7*chil2)/(128.*mtot8*eta4) +
1558  (35*chip2*m2_8*chil2)/(128.*mtot8*eta4) +
1559  (375*dm2*m2_2*chil2)/(256.*mtot4*eta) + (1815*dm*m2_3*chil2)/(256.*mtot4*eta) +
1560  (1645*m2_4*chil2)/(192.*mtot4*eta));
1561 
1562  angcoeffs->epsiloncoeff1 = (-0.18229166666666666 - (5*dm)/(64.*m2));
1563 
1564  angcoeffs->epsiloncoeff2 = ((-15*dm*m2*chil)/(128.*mtot2*eta) - (35*m2_2*chil)/(128.*mtot2*eta));
1565 
1566  angcoeffs->epsiloncoeff3 = (-1.7952473958333333 - (4555*dm)/(7168.*m2) - (515*eta)/384. -
1567  (15*dm2*eta)/(256.*m2_2) - (175*dm*eta)/(256.*m2));
1568 
1569  angcoeffs->epsiloncoeff4 = - (35*LAL_PI)/48. - (5*dm*LAL_PI)/(16.*m2) +
1570  (5*dm2*chil)/(16.*mtot2) + (5*dm*m2*chil)/(3.*mtot2) +
1571  (2545*m2_2*chil)/(1152.*mtot2) + (2035*dm*m2*chil)/(21504.*mtot2*eta) +
1572  (2995*m2_2*chil)/(9216.*mtot2*eta);
1573 
1574  angcoeffs->epsiloncoeff5 = (4.318908476114694 + (27895885*dm)/(2.1676032e7*m2) + (39695*eta)/86016. +
1575  (1615*dm2*eta)/(28672.*m2_2) - (265*dm*eta)/(14336.*m2) + (955*eta2)/576. +
1576  (15*dm3*eta2)/(1024.*m2_3) + (35*dm2*eta2)/(256.*m2_2) +
1577  (2725*dm*eta2)/(3072.*m2) - (15*dm*m2*LAL_PI*chil)/(16.*mtot2*eta) - (35*m2_2*LAL_PI*chil)/(16.*mtot2*eta) +
1578  (375*dm2*m2_2*chil2)/(256.*mtot4*eta) + (1815*dm*m2_3*chil2)/(256.*mtot4*eta) +
1579  (1645*m2_4*chil2)/(192.*mtot4*eta));
1580 }
1581 
1582 /**
1583  * Simple 2PN version of the orbital angular momentum L,
1584  * without any spin terms expressed as a function of v.
1585  * For IMRPhenomP(v2).
1586  *
1587  * Reference:
1588  * - Boh&eacute; et al, 1212.5520v2 Eq 4.7 first line
1589  */
1590 static REAL8 L2PNR(
1591  const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
1592  const REAL8 eta) /**< Symmetric mass-ratio */
1593 {
1594  const REAL8 eta2 = eta*eta;
1595  const REAL8 x = v*v;
1596  const REAL8 x2 = x*x;
1597  return (eta*(1.0 + (1.5 + eta/6.0)*x + (3.375 - (19.0*eta)/8. - eta2/24.0)*x2)) / sqrt(x);
1598 }
1599 
1600 /**
1601  * Simple 2PN version of the orbital angular momentum L,
1602  * without any spin terms expressed as a function of v.
1603  * For IMRPhenomP(v1).
1604  *
1605  * Reference:
1606  * - Kidder, Phys. Rev. D 52, 821–847 (1995), Eq. 2.9
1607  */
1609  const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
1610  const REAL8 eta) /**< Symmetric mass-ratio */
1611 {
1612  const REAL8 mu = eta; /* M=1 */
1613  const REAL8 v2 = v*v;
1614  const REAL8 v3 = v2*v;
1615  const REAL8 v4 = v3*v;
1616  const REAL8 eta2 = eta*eta;
1617  const REAL8 b = (4.75 + eta/9.)*eta*v4;
1618 
1619 
1620  return mu*sqrt((1 - ((3 - eta)*v2)/3. + b)/v2)*
1621  (1 + ((1 - 3*eta)*v2)/2. + (3*(1 - 7*eta + 13*eta2)*v4)/8. +
1622  ((14 - 41*eta + 4*eta2)*v4)/(4.*pow_2_of(1 - ((3 - eta)*v2)/3. + b)) +
1623  ((3 + eta)*v2)/(1 - ((3 - eta)*v2)/3. + b) +
1624  ((7 - 10*eta - 9*eta2)*v4)/(2.*(1 - ((3 - eta)*v2)/3. + b)));
1625 }
1626 
1627 /** Expressions used for the WignerD symbol
1628  * with full expressions for the angles.
1629  * Used for IMRPhenomP(v2):
1630  * \f{equation}{
1631  * \cos(\beta) = \hat J . \hat L
1632  * = \left( 1 + \left( S_\mathrm{p} / (L + S_L) \right)^2 \right)^{-1/2}
1633  * = \left( L + S_L \right) / \sqrt{ \left( L + S_L \right)^2 + S_p^2 }
1634  * = \mathrm{sign}\left( L + S_L \right) \cdot \left( 1 + \left( S_p / \left(L + S_L\right)\right)^2 \right)^{-1/2}
1635  * \f}
1636  */
1638  REAL8 *cos_beta_half, /**< [out] cos(beta/2) */
1639  REAL8 *sin_beta_half, /**< [out] sin(beta/2) */
1640  const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
1641  const REAL8 SL, /**< Dimensionfull aligned spin */
1642  const REAL8 eta, /**< Symmetric mass-ratio */
1643  const REAL8 Sp) /**< Dimensionfull spin component in the orbital plane */
1644 {
1645  XLAL_CHECK_VOID(cos_beta_half != NULL, XLAL_EFAULT);
1646  XLAL_CHECK_VOID(sin_beta_half != NULL, XLAL_EFAULT);
1647  /* We define the shorthand s := Sp / (L + SL) */
1648  const REAL8 L = L2PNR(v, eta);
1649  // We ignore the sign of L + SL below.
1650  REAL8 s = Sp / (L + SL); /* s := Sp / (L + SL) */
1651  REAL8 s2 = s*s;
1652  REAL8 cos_beta = 1.0 / sqrt(1.0 + s2);
1653  *cos_beta_half = + sqrt( (1.0 + cos_beta) / 2.0 ); /* cos(beta/2) */
1654  *sin_beta_half = + sqrt( (1.0 - cos_beta) / 2.0 ); /* sin(beta/2) */
1655 }
1656 
1657 /** Expressions used for the WignerD symbol
1658  * with small angle approximation.
1659  * Used for IMRPhenomP(v1):
1660  * \f{equation}{
1661  * \cos(\beta) = \hat J . \hat L
1662  * = \left(1 + \left( S_\mathrm{p} / (L + S_L)\right)^2 \right)^{-1/2}
1663  * \f}
1664  * We use the expression
1665  * \f{equation}{
1666  * \cos(\beta/2) \approx (1 + s^2 / 4 )^{-1/2} \,,
1667  * \f}
1668  * where \f$s := S_p / (L + S_L)\f$.
1669  */
1671  REAL8 *cos_beta_half, /**< Output: cos(beta/2) */
1672  REAL8 *sin_beta_half, /**< Output: sin(beta/2) */
1673  const REAL8 v, /**< Cubic root of (Pi * Frequency (geometric)) */
1674  const REAL8 SL, /**< Dimensionfull aligned spin */
1675  const REAL8 eta, /**< Symmetric mass-ratio */
1676  const REAL8 Sp) /**< Dimensionfull spin component in the orbital plane */
1677 {
1678  XLAL_CHECK_VOID(cos_beta_half != NULL, XLAL_EFAULT);
1679  XLAL_CHECK_VOID(sin_beta_half != NULL, XLAL_EFAULT);
1680  REAL8 s = Sp / (L2PNR_v1(v, eta) + SL); /* s := Sp / (L + SL) */
1681  REAL8 s2 = s*s;
1682  *cos_beta_half = 1.0/sqrt(1.0 + s2/4.0); /* cos(beta/2) */
1683  *sin_beta_half = sqrt(1.0 - 1.0/(1.0 + s2/4.0)); /* sin(beta/2) */
1684 }
1685 
1686 /**
1687  * In this helper function we check whether the maximum opening angle during the evolution
1688  * becomes larger than pi/2 or pi/4, in which case a warning is issued.
1689  */
1691  const REAL8 m1, /**< Mass of companion 1 (solar masses) */
1692  const REAL8 m2, /**< Mass of companion 2 (solar masses) */
1693  const REAL8 chi1_l, /**< Aligned spin of BH 1 */
1694  const REAL8 chi2_l, /**< Aligned spin of BH 2 */
1695  const REAL8 chip /**< Dimensionless spin in the orbital plane */
1696 ) {
1697  REAL8 M = m1+m2;
1698  REAL8 m1_normalized = m1/M;
1699  REAL8 m2_normalized = m2/M;
1700  REAL8 eta = m1_normalized*m2_normalized;
1701  REAL8 v_at_max_beta = sqrt( 2*(9. + eta - sqrt(1539 - 1008*eta - 17*eta*eta)) / 3. / (eta*eta + 57.*eta - 81.) );
1702  // The coefficients above come from finding the roots of the derivative of L2PN
1703  REAL8 SL = m1_normalized*m1_normalized*chi1_l + m2_normalized*m2_normalized*chi2_l;
1704  REAL8 Sperp = m2_normalized*m2_normalized * chip;
1705  REAL8 cBetah, sBetah;
1706  WignerdCoefficients(&cBetah, &sBetah, v_at_max_beta, SL, eta, Sperp);
1707  REAL8 L_min = L2PNR(v_at_max_beta,eta);
1708  REAL8 max_beta = 2*acos(cBetah);
1709  /** If L+SL becomes <0, WignerdCoefficients does not track the angle between J and L anymore (see tech doc, choice of + sign
1710  * so that the Wigner coefficients are OK in the aligned spin limit) and the model may become pathological as one moves away from
1711  * the aligned spin limit.
1712  * If this does not happen, then max_beta is the actual maximum opening angle as predicted by the model.
1713  */
1714  if (L_min + SL < 0. && chip > 0.)
1715  XLAL_PRINT_WARNING("The maximum opening angle exceeds Pi/2.\nThe model may be pathological in this regime.");
1716  else if (max_beta > LAL_PI_4)
1717  XLAL_PRINT_WARNING("The maximum opening angle %g is larger than Pi/4.\nThe model has not been tested against NR in this regime.", max_beta);
1718 }
1719 
1720 /**
1721  * Wrapper for final-spin formula based on:
1722  * - IMRPhenomD's FinalSpin0815() for aligned spins.
1723  *
1724  * We use their convention m1>m2
1725  * and put <b>all in-plane spin on the larger BH</b>.
1726  *
1727  * In the aligned limit return the FinalSpin0815 value.
1728  */
1730  const REAL8 m1, /**< Mass of companion 1 (solar masses) */
1731  const REAL8 m2, /**< Mass of companion 2 (solar masses) */
1732  const REAL8 chi1_l, /**< Aligned spin of BH 1 */
1733  const REAL8 chi2_l, /**< Aligned spin of BH 2 */
1734  const REAL8 chip) /**< Dimensionless spin in the orbital plane */
1735 {
1736  const REAL8 M = m1+m2;
1737  REAL8 eta = m1*m2/(M*M);
1738  if (eta > 0.25) nudge(&eta, 0.25, 1e-6);
1739 
1740  REAL8 af_parallel, q_factor;
1741  if (m1 >= m2) {
1742  q_factor = m1/M;
1743  af_parallel = FinalSpin0815(eta, chi1_l, chi2_l);
1744  }
1745  else {
1746  q_factor = m2/M;
1747  af_parallel = FinalSpin0815(eta, chi2_l, chi1_l);
1748  }
1749 
1750  REAL8 Sperp = chip * q_factor*q_factor;
1751  REAL8 af = copysign(1.0, af_parallel) * sqrt(Sperp*Sperp + af_parallel*af_parallel);
1752  return af;
1753 }
1754 
1755 /**
1756  * Wrapper for final-spin formula based on:
1757  * - Barausse \& Rezzolla, Astrophys.J.Lett.704:L40-L44, 2009,
1758  * arXiv:0904.2577
1759  *
1760  * We use their convention m1>m2
1761  * and put <b>all spin on the larger BH</b>:
1762  *
1763  * a1 = (chip, 0, chi), a2 = (0,0,0), L = (0,0,1)
1764  */
1766  const REAL8 nu, /**< Symmetric mass-ratio */
1767  const REAL8 chi, /**< Effective aligned spin of the binary: chi = (m1*chi1 + m2*chi2)/M */
1768  const REAL8 chip) /**< Dimensionless spin in the orbital plane */
1769 {
1770 
1771  const REAL8 a1_x = chip;
1772  const REAL8 a1_y = 0;
1773  const REAL8 a1_z = chi;
1774  const REAL8 a2_x = 0;
1775  const REAL8 a2_y = 0;
1776  const REAL8 a2_z = 0;
1777 
1778  const REAL8 a1 = sqrt(a1_x*a1_x + a1_y*a1_y + a1_z*a1_z);
1779  const REAL8 a2 = sqrt(a2_x*a2_x + a2_y*a2_y + a2_z*a2_z);
1780 
1781  const REAL8 cos_alpha = (a1*a2 == 0) ? 0.0 : a1_z*a2_z/(a1*a2); /* cos(alpha) = \hat a1 . \hat a2 (Eq. 7) */
1782  const REAL8 cos_beta_tilde = (a1 == 0) ? 0.0 : a1_z/a1; /* \cos(\tilde \beta) = \hat a1 . \hat L (Eq. 9) */
1783  const REAL8 cos_gamma_tilde = (a2 == 0) ? 0.0 : a2_z/a2; /* \cos(\tilde \gamma) = \hat a2 . \hat L (Eq. 9) */
1784 
1785  return FinalSpinBarausse2009(nu, a1, a2, cos_alpha, cos_beta_tilde, cos_gamma_tilde);
1786 }
1787 
1788 /**
1789  * Final-spin formula based on:
1790  * - Barausse \& Rezzolla, Astrophys.J.Lett.704:L40-L44, 2009,
1791  * arXiv:0904.2577
1792  *
1793  * We use their convention m1>m2.
1794  */
1796  const REAL8 nu, /**< Symmetric mass-ratio */
1797  const REAL8 a1, /**< |a_1| norm of dimensionless spin vector for BH 1 */
1798  const REAL8 a2, /**< |a_2| norm of dimensionless spin vector for BH 2 */
1799  const REAL8 cos_alpha, /**< \f$\cos(\alpha) = \hat a_1 . \hat a_2\f$ (Eq. 7) */
1800  const REAL8 cos_beta_tilde, /**< \f$\cos(\tilde \beta) = \hat a_1 . \hat L\f$ (Eq. 9) */
1801  const REAL8 cos_gamma_tilde) /**< \f$\cos(\tilde \gamma) = \hat a_2 . \hat L\f$ (Eq. 9)*/
1802 {
1803  REAL8 q = (2*nu)/(1 + sqrt(1 - 4*nu) - 2*nu);
1804 
1805  /* These parameters are defined in eq. 3. */
1806  const REAL8 s4 = -0.1229;
1807  const REAL8 s5 = 0.4537;
1808  const REAL8 t0 = -2.8904;
1809  const REAL8 t2 = -3.5171;
1810  const REAL8 t3 = 2.5763;
1811 
1812  /* shorthands */
1813  const REAL8 nu2 = nu*nu;
1814  const REAL8 q2 = q*q;
1815  const REAL8 q4 = q2*q2;
1816  const REAL8 q2p = 1 + q2;
1817  const REAL8 q2p2 = q2p*q2p;
1818  const REAL8 qp = 1 + q;
1819  const REAL8 qp2 = qp*qp;
1820  const REAL8 a1_2 = a1*a1;
1821  const REAL8 a2_2 = a2*a2;
1822 
1823  /* l = \tilde l/(m1*m2), where \tilde l = S_fin - (S1 + S2) = L - J_rad (Eq. 4) */
1824  const REAL8 l = 2*sqrt(3.0) + t2*nu + t3*nu2
1825  + (s4 / q2p2) * (a1_2 + a2_2*q4 + 2*a1*a2*q2*cos_alpha)
1826  + ((s5*nu + t0 + 2)/q2p) * (a1*cos_beta_tilde + a2*cos_gamma_tilde*q2); /* |l| (Eq. 10) */
1827  const REAL8 l2 = l*l;
1828 
1829  /* a_fin = S_fin/M^2 (Eq. 6) */
1830  const REAL8 a_fin = (1.0 / qp2) * sqrt(a1_2 + a2_2*q4 + 2*a1*a2*q2*cos_alpha + 2*(a1*cos_beta_tilde + a2*q2*cos_gamma_tilde)*l*q + l2*q2);
1831  return a_fin;
1832 }
1833 
1834 
1835 /**
1836  * PhenomC parameters for modified ringdown,
1837  * uses final spin formula of:
1838  * - Barausse \& Rezzolla, Astrophys.J.Lett.704:L40-L44, 2009,
1839  * arXiv:0904.2577
1840  */
1842  const REAL8 m1, /**< Mass of companion 1 (solar masses) */
1843  const REAL8 m2, /**< Mass of companion 2 (solar masses) */
1844  const REAL8 chi, /**< Reduced aligned spin of the binary chi = (m1*chi1 + m2*chi2)/M */
1845  const REAL8 chip, /**< Dimensionless spin in the orbital plane */
1846  LALDict *extraParams) /**< linked list that may contain the extra testing GR parameters and/or tidal parameters */
1847 {
1848 
1849  BBHPhenomCParams *p = NULL;
1850  p = ComputeIMRPhenomCParams(m1, m2, chi, extraParams); /* populate parameters with the original PhenomC setup */
1851  if( !p )
1853 
1854  const REAL8 M = m1 + m2;
1855  const REAL8 eta = m1 * m2 / (M * M);
1856 
1857  REAL8 finspin = FinalSpinBarausse2009_all_spin_on_larger_BH(eta, chi, chip);
1858  if( fabs(finspin) > 1.0 ) {
1859  XLAL_PRINT_WARNING("Warning: final spin magnitude %g > 1. Setting final spin magnitude = 1.", finspin);
1860  finspin = copysign(1.0, finspin);
1861  }
1862 
1863  p->afin = finspin;
1864 
1865  /* Get the Ringdown frequency */
1866  REAL8 prefac = (1./(2.*LAL_PI)) * LAL_C_SI * LAL_C_SI * LAL_C_SI / (LAL_G_SI * M * LAL_MSUN_SI);
1867  REAL8 k1 = 1.5251;
1868  REAL8 k2 = -1.1568;
1869  REAL8 k3 = 0.1292;
1870 
1871  p->fRingDown = (prefac * (k1 + k2 * pow(1. - fabs(finspin), k3)));
1872  p->MfRingDown = p->m_sec * p->fRingDown;
1873 
1874  /* Get the quality factor of ring-down, using Eq (5.6) of Main paper (arxiv.org/pdf/1005.3306v3.pdf) */
1875  p->Qual = (0.7000 + (1.4187 * pow(1.0 - fabs(finspin), -0.4990)) );
1876 
1877  /* Get the transition frequencies, at which the model switches phase and
1878  * amplitude prescriptions, as used in Eq.(5.9), (5.13) of the Main paper */
1879  p->f1 = 0.1 * p->fRingDown;
1880  p->f2 = p->fRingDown;
1881  p->f0 = 0.98 * p->fRingDown;
1882  p->d1 = 0.005;
1883  p->d2 = 0.005;
1884  p->d0 = 0.015;
1885 
1886  /* Get the coefficients beta1, beta2, defined in Eq 5.7 of the main paper */
1887  REAL8 Mfrd = p->MfRingDown;
1888 
1889  p->b2 = ((-5./3.)* p->a1 * pow(Mfrd,(-8./3.)) - p->a2/(Mfrd*Mfrd) -
1890  (p->a3/3.)*pow(Mfrd,(-4./3.)) + (2./3.)* p->a5 * pow(Mfrd,(-1./3.)) + p->a6)/eta;
1891 
1892  REAL8 psiPMrd = IMRPhenomCGeneratePhasePM( p->fRingDown, eta, p );
1893 
1894  p->b1 = psiPMrd - (p->b2 * Mfrd);
1895 
1896  /* Taking the upper cut-off frequency as 0.15M */
1897  /*p->fCut = (1.7086 * eta * eta - 0.26592 * eta + 0.28236) / p->piM;*/
1898  p->fCut = 0.15 / p->m_sec;
1899 
1900  return p;
1901 }
1902 
1903 
1904 // This function determines whether x and y are approximately equal to a relative accuracy epsilon.
1905 // Note that x and y are compared to relative accuracy, so this function is not suitable for testing whether a value is approximately zero.
1906 static bool approximately_equal(REAL8 x, REAL8 y, REAL8 epsilon) {
1907  return !gsl_fcmp(x, y, epsilon);
1908 }
1909 
1910 // If x and X are approximately equal to relative accuracy epsilon then set x = X.
1911 // If X = 0 then use an absolute comparison.
1912 static void nudge(REAL8 *x, REAL8 X, REAL8 epsilon) {
1913  if (X != 0.0) {
1914  if (approximately_equal(*x, X, epsilon)) {
1915  XLAL_PRINT_INFO("Nudging value %.15g to %.15g\n", *x, X);
1916  *x = X;
1917  }
1918  }
1919  else {
1920  if (fabs(*x - X) < epsilon)
1921  *x = X;
1922  }
1923 }
1924 
void XLALDestroyDict(LALDict *dict)
LALDict * XLALCreateDict(void)
void XLALSimInspiralGetHOSpinTerms(REAL8 *SS_3p5PN, REAL8 *SSS_3p5PN, REAL8 X_A, REAL8 X_B, REAL8 chi1, REAL8 chi2, REAL8 quadparam1, REAL8 quadparam2)
Function to add 3.5PN spin-squared and 3.5PN spin-cubed terms.
int XLALSimNRTunedTidesFDTidalPhaseFrequencySeries(const REAL8Sequence *phi_tidal, const REAL8Sequence *amp_tidal, const REAL8Sequence *planck_taper, const REAL8Sequence *fHz, REAL8 m1_SI, REAL8 m2_SI, REAL8 lambda1, REAL8 lambda2, NRTidal_version_type NRTidal_version)
Function to call the frequency domain tidal correction over an array of input frequencies.
double XLALSimNRTunedTidesMergerFrequency(const REAL8 mtot_MSUN, const REAL8 kappa2T, const REAL8 q)
compute the merger frequency of a BNS system.
double XLALSimNRTunedTidesComputeKappa2T(REAL8 m1_SI, REAL8 m2_SI, REAL8 lambda1, REAL8 lambda2)
convenient function to compute tidal coupling constant.
static size_t NextPow2(const size_t n)
static int IMRPhenomCGenerateAmpPhase(REAL8 *amplitude, REAL8 *phasing, REAL8 f, REAL8 eta, const BBHPhenomCParams *params)
static REAL8 IMRPhenomCGeneratePhasePM(REAL8 f, REAL8 eta, const BBHPhenomCParams *params)
main functions
static BBHPhenomCParams * ComputeIMRPhenomCParams(const REAL8 m1, const REAL8 m2, const REAL8 chi, LALDict *extraParams)
UsefulPowers powers_of_pi
useful powers of LAL_PI, calculated once and kept constant - to be initied with a call to init_useful...
#define f_CUT
Dimensionless frequency (Mf) at which define the end of the waveform.
Internal function for IMRPhenomD phenomenological waveform model. See LALSimIMRPhenom....
static void ComputeIMRPhenDPhaseConnectionCoefficients(IMRPhenomDPhaseCoefficients *p, PNPhasingSeries *pn, PhiInsPrefactors *prefactors, double Rholm, double Taulm)
This function aligns the three phase parts (inspiral, intermediate and merger-rindown) such that they...
static double IMRPhenDPhase(double f, IMRPhenomDPhaseCoefficients *p, PNPhasingSeries *pn, UsefulPowers *powers_of_f, PhiInsPrefactors *prefactors, double Rholm, double Taulm)
This function computes the IMR phase given phenom coefficients.
static int init_phi_ins_prefactors(PhiInsPrefactors *prefactors, IMRPhenomDPhaseCoefficients *p, PNPhasingSeries *pn)
static double FinalSpin0815(double eta, double chi1, double chi2)
Wrapper function for FinalSpin0815_s.
static void ComputeIMRPhenomDPhaseCoefficients(IMRPhenomDPhaseCoefficients *p, double eta, double chi1, double chi2, double finspin, LALDict *extraParams)
A struct containing all the parameters that need to be calculated to compute the phenomenological pha...
static void ComputeIMRPhenomDAmplitudeCoefficients(IMRPhenomDAmplitudeCoefficients *p, double eta, double chi1, double chi2, double finspin)
A struct containing all the parameters that need to be calculated to compute the phenomenological amp...
static double Subtract3PNSS(double m1, double m2, double M, double eta, double chi1, double chi2)
Subtract 3PN spin-spin term below as this is in LAL's TaylorF2 implementation (LALSimInspiralPNCoeffi...
static double IMRPhenDAmplitude(double f, IMRPhenomDAmplitudeCoefficients *p, UsefulPowers *powers_of_f, AmpInsPrefactors *prefactors)
This function computes the IMR amplitude given phenom coefficients.
static int init_amp_ins_prefactors(AmpInsPrefactors *prefactors, IMRPhenomDAmplitudeCoefficients *p)
static int init_useful_powers(UsefulPowers *p, REAL8 number)
static double pow_2_of(double number)
calc square of number without floating point 'pow'
double PhenomInternal_OrbAngMom3PN(const double f_orb_hz, const double m1, const double m2, const double s1x, const double s1y, const double s1z, const double s2x, const double s2y, const double s2z, const double f_0, const int ExpansionOrder)
Wrapper function for XLALOrbitalAngMom3PNSpinning.
static void ComputeNNLOanglecoeffs(NNLOanglecoeffs *angcoeffs, const REAL8 q, const REAL8 chil, const REAL8 chip)
Next-to-next-to-leading order PN coefficients for Euler angles and .
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_in, double deltaF, IMRPhenomP_version_type IMRPhenomP_version, NRTidal_version_type NRTidal_version, LALDict *extraParams)
Internal core function to calculate plus and cross polarizations of the PhenomP model for a set of fr...
#define ROTATEZ(angle, vx, vy, vz)
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)
Expressions used for the WignerD symbol with full expressions for the angles.
static void CheckMaxOpeningAngle(const REAL8 m1, const REAL8 m2, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip)
In this helper function we check whether the maximum opening angle during the evolution becomes large...
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)
const double sqrt_6
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)
Simple 2PN version of the orbital angular momentum L, without any spin terms expressed as a function ...
static UNUSED BBHPhenomCParams * ComputeIMRPhenomCParamsRDmod(const REAL8 m1, const REAL8 m2, const REAL8 chi, const REAL8 chip, LALDict *extraParams)
PhenomC parameters for modified ringdown, uses final spin formula of:
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)
Final-spin formula based on:
static REAL8 L2PNR(const REAL8 v, const REAL8 eta)
Simple 2PN version of the orbital angular momentum L, without any spin terms expressed as a function ...
static REAL8 FinalSpinBarausse2009_all_spin_on_larger_BH(const REAL8 nu, const REAL8 chi, const REAL8 chip)
Wrapper for final-spin formula based on:
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)
Wrapper for final-spin formula based on:
#define ROTATEY(angle, vx, vy, vz)
static void WignerdCoefficients_SmallAngleApproximation(REAL8 *cos_beta_half, REAL8 *sin_beta_half, const REAL8 v, const REAL8 SL, const REAL8 eta, const REAL8 Sp)
Expressions used for the WignerD symbol with small angle approximation.
static int PhenomPCoreOneFrequency_withTides(const REAL8 fHz, const REAL8 window, const REAL8 phaseTidal, const REAL8 ampTidal, 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 bool approximately_equal(REAL8 x, REAL8 y, REAL8 epsilon)
#define MAX_TOL_ATAN
Tolerance used below which numbers are treated as zero for the calculation of atan2.
int XLALSimInspiralSetQuadMonParamsFromLambdas(LALDict *LALparams)
if you do NOT provide a quadparam[1,2] term and you DO provide lamdba[1,2] then we calculate quad-mon...
#define c
static REAL8 UNUSED q2(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q4(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
REAL8 XLALSimInspiralWaveformParamsLookupdQuadMon1(LALDict *params)
int XLALSimInspiralWaveformParamsInsertTidalLambda1(LALDict *params, REAL8 value)
REAL8 XLALSimInspiralWaveformParamsLookupTidalLambda2(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupTidalLambda1(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPNSpinOrder(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsInsertTidalLambda2(LALDict *params, REAL8 value)
REAL8 XLALSimInspiralWaveformParamsLookupdQuadMon2(LALDict *params)
int XLALSimInspiralWaveformParamsInsertdQuadMon2(LALDict *params, REAL8 value)
int XLALSimInspiralWaveformParamsInsertdQuadMon1(LALDict *params, REAL8 value)
REAL8 tmp1
int s
Definition: bh_qnmode.c:137
int l
Definition: bh_qnmode.c:135
REAL8 M
Definition: bh_qnmode.c:133
double i
Definition: bh_ringdown.c:118
double e
Definition: bh_ringdown.c:117
const double pn
const double a2
COMPLEX16FrequencySeries * XLALCreateCOMPLEX16FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
void XLALDestroyCOMPLEX16FrequencySeries(COMPLEX16FrequencySeries *series)
#define LAL_C_SI
#define LAL_MSUN_SI
#define LAL_PI
#define LAL_MTSUN_SI
#define LAL_PI_4
#define LAL_G_SI
#define LAL_MRSUN_SI
#define LIGOTIMEGPSZERO
double complex COMPLEX16
double REAL8
uint32_t UINT4
void * XLALMalloc(size_t n)
void XLALFree(void *p)
IMRPhenomP_version_type
Definition: LALSimIMR.h:73
NRTidal_version_type
Definition: LALSimIMR.h:80
@ IMRPhenomPv3_V
version 3: based on IMRPhenomD and the precession angles from Katerina Chatziioannou PhysRevD....
Definition: LALSimIMR.h:77
@ IMRPhenomPv1_V
version 1: based on IMRPhenomC
Definition: LALSimIMR.h:74
@ IMRPhenomPv2_V
version 2: based on IMRPhenomD
Definition: LALSimIMR.h:75
@ IMRPhenomPv2NRTidal_V
version Pv2_NRTidal: based on IMRPhenomPv2; NRTides added before precession; can be used with both NR...
Definition: LALSimIMR.h:76
@ NRTidal_V
version NRTidal: based on https://arxiv.org/pdf/1706.02969.pdf
Definition: LALSimIMR.h:81
@ NRTidalv2_V
version NRTidalv2: https://arxiv.org/abs/1905.06011
Definition: LALSimIMR.h:82
static REAL8 atan2tol(REAL8 a, REAL8 b, REAL8 tol)
int XLALSimIMRPhenomP(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip, const REAL8 thetaJ, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 distance, const REAL8 alpha0, const REAL8 phic, const REAL8 deltaF, const REAL8 f_min, const REAL8 f_max, const REAL8 f_ref, IMRPhenomP_version_type IMRPhenomP_version, NRTidal_version_type NRTidal_version, LALDict *extraParams)
Driver routine to compute the precessing inspiral-merger-ringdown phenomenological waveform IMRPhenom...
int XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame(REAL8 *chi1_l, REAL8 *chi2_l, REAL8 *chip, REAL8 *thetaJN, REAL8 *alpha0, REAL8 *phi_aligned, REAL8 *zeta_polariz, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 f_ref, const REAL8 phiRef, const REAL8 incl, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 s2x, const REAL8 s2y, const REAL8 s2z, IMRPhenomP_version_type IMRPhenomP_version)
Function to map LAL parameters (masses, 6 spin components, phiRef and inclination at f_ref) (assumed ...
int XLALSimIMRPhenomPCalculateModelParametersOld(REAL8 *chi1_l, REAL8 *chi2_l, REAL8 *chip, REAL8 *thetaJ, REAL8 *alpha0, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 f_ref, const REAL8 lnhatx, const REAL8 lnhaty, const REAL8 lnhatz, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 s2x, const REAL8 s2y, const REAL8 s2z, IMRPhenomP_version_type IMRPhenomP_version)
Deprecated : used the old convention (view frame for the spins) Function to map LAL parameters (masse...
int XLALSimIMRPhenomPFrequencySequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip, const REAL8 thetaJ, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 distance, const REAL8 alpha0, const REAL8 phic, const REAL8 f_ref, IMRPhenomP_version_type IMRPhenomP_version, NRTidal_version_type NRTidal_version, LALDict *extraParams)
Driver routine to compute the precessing inspiral-merger-ringdown phenomenological waveform IMRPhenom...
@ LAL_SIM_INSPIRAL_SPIN_ORDER_35PN
int XLALSimInspiralTaylorF2AlignedPhasing(PNPhasingSeries **pfa, const REAL8 m1, const REAL8 m2, const REAL8 chi1, const REAL8 chi2, LALDict *extraPars)
Returns structure containing TaylorF2 phasing coefficients for given physical parameters.
static const INT4 m
static const INT4 q
static const INT4 a
void XLALDestroyREAL8Sequence(REAL8Sequence *sequence)
REAL8Sequence * XLALCreateREAL8Sequence(size_t length)
COMPLEX16 XLALSpinWeightedSphericalHarmonic(REAL8 theta, REAL8 phi, int s, int l, int m)
const LALUnit lalStrainUnit
const LALUnit lalSecondUnit
LALUnit * XLALUnitMultiply(LALUnit *output, const LALUnit *unit1, const LALUnit *unit2)
#define XLAL_PRINT_INFO(...)
#define XLAL_ERROR_NULL(...)
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
#define XLAL_CHECK_VOID(assertion,...)
#define XLAL_PRINT_WARNING(...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
XLAL_EINVAL
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
list x
list p
list y
list mu
double alpha
Definition: sgwb.c:106
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.
REAL8 * data
useful powers in GW waveforms: 1/6, 1/3, 2/3, 4/3, 5/3, 2, 7/3, 8/3, -1, -1/6, -7/6,...
double f_min
Definition: unicorn.c:22
double f_max
Definition: unicorn.c:23