LALSimulation  5.4.0.1-89842e6
LALSimIMRPhenomX_PNR.c
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2022 Cardiff University
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 #ifdef __cplusplus
21 extern "C"
22 {
23 #endif
24 
25 #include <lal/SphericalHarmonics.h>
26 #include "LALSimIMR.h"
27 
28 #include "LALSimIMRPhenomX.h"
29 #include "LALSimIMRPhenomXPHM.h"
33 
34 #include "LALSimIMRPhenomX_PNR.h"
37 
38 #ifndef _OPENMP
39 #define omp ignore
40 #endif
41 
42 #ifndef PHENOMXHMDEBUG
43 #define DEBUG 0
44 #else
45 #define DEBUG 1
46 #endif
47 
48  /**
49  * @addtogroup LALSimIMRPhenomX_c
50  * @{
51  *
52  * @name Routines for PNR angles
53  * @{
54  *
55  * @author Eleanor Hamilton, Sebastian Khan, Jonathan E. Thompson
56  *
57  * @brief C code for routines to implement PNR angles in IMRPhenomXP/IMRPhenomXPHM.
58  *
59  * This is a C-code impementation of the PNR angle model
60  * outlined in arXiv:2107.08876.
61  *
62  * Available flags:
63  * PhenomXPNRUseTunedAngles:
64  * - 0 : Disable PNR angles. Use the default NNLO/MSA precession angles in IMRPhenomXP/HM without tuning.
65  * - 1 : Enable PNR angles.
66  */
67 
68  /**
69  * This is an external wrapper to generate the (2,2) PNR angles,
70  * following the prescription outlined in arXiv:2107.08876,
71  * given the standard inputs given to generate FD waveforms.
72  */
74  REAL8Sequence **alphaPNR, /**< [out] Alpha precession angle (rad) */
75  REAL8Sequence **betaPNR, /**< [out] Beta precession angle (rad) */
76  REAL8Sequence **gammaPNR, /**< [out] Gamma precession angle (rad) */
77  REAL8Sequence **freqs, /**< [out] Frequency array (Hz) */
78  REAL8 *alphaPNR_ref, /**< [out] reference value of alpha (rad) */
79  REAL8 *gammaPNR_ref, /**< [out] reference value of gamma (rad) */
80  REAL8 m1_SI, /**< mass of companion 1 (kg) */
81  REAL8 m2_SI, /**< mass of companion 2 (kg) */
82  REAL8 chi1x, /**< x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
83  REAL8 chi1y, /**< y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
84  REAL8 chi1z, /**< z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
85  REAL8 chi2x, /**< x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
86  REAL8 chi2y, /**< y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
87  REAL8 chi2z, /**< z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
88  REAL8 inclination, /**< Angle between orbital angular momentum and line-of-sight vector at reference frequency (rad) */
89  REAL8 deltaF, /**< Frequency spacing (Hz) */
90  REAL8 f_min, /**< Starting GW frequency (Hz) */
91  REAL8 f_max, /**< Ending GW frequency (Hz) */
92  REAL8 fRef_In, /**< Reference frequency (Hz) */
93  LALDict *lalParams /**< LAL Dictionary struct */
94  )
95  {
96  /* Simple check on masses and spins */
97  UINT4 status = 0;
98  status = XLALIMRPhenomXPCheckMassesAndSpins(&m1_SI, &m2_SI, &chi1x, &chi1y, &chi1z, &chi2x, &chi2y, &chi2z);
99  XLAL_CHECK(
100  XLAL_SUCCESS == status,
101  XLAL_EFUNC,
102  "Error: XLALIMRPhenomXPCheckMassesAndSpins failed in XLALSimIMRPhenomX_PNR_GeneratePNRAngles.\n");
103 
104  /* Ensure we have a dictionary */
105  LALDict *lalParams_aux;
106  if (lalParams == NULL)
107  {
108  lalParams_aux = XLALCreateDict();
109  }
110  else
111  {
112  lalParams_aux = XLALDictDuplicate(lalParams);
113  }
114 
115  /* Make sure we're calling the tuned angles */
117  if (!UsePNR)
118  {
119  UsePNR = 1;
121  }
122 
123  /* Map fRef to the start frequency if need be, then make sure it's within the frequency range */
124  REAL8 fRef = (fRef_In == 0.0) ? f_min : fRef_In;
125  XLAL_CHECK(
126  ((f_min <= fRef) && (fRef <= f_max)),
127  XLAL_EFUNC,
128  "Error: fRef needs to be within the specified minimum and maximum frequency values!\n");
129 
130  /* Check that the passed deltaF is non-zero */
131  XLAL_CHECK(
132  deltaF > 0,
133  XLAL_EFUNC,
134  "Error: deltaF needs to be a positive number!\n");
135 
136  /* Generate a uniformly sampled frequency grid of spacing deltaF. */
137  /* Frequencies will be set using only the lower and upper bounds that we passed */
138  size_t iStart = (size_t)(f_min / deltaF);
139  size_t iStop = (size_t)(f_max / deltaF) + 1;
140 
141  XLAL_CHECK(
142  (iStart <= iStop),
143  XLAL_EDOM,
144  "Error: the starting frequency index is greater than the stopping index! Please ensure that f_min <= f_max.\n");
145 
146  /* Allocate memory for frequency array and terminate if this fails */
147  *freqs = XLALCreateREAL8Sequence(iStop - iStart);
148  if (!(*freqs))
149  {
150  XLAL_ERROR(XLAL_EFUNC, "Frequency array allocation failed.");
151  }
152  *alphaPNR = XLALCreateREAL8Sequence(iStop - iStart);
153  if (!(*alphaPNR))
154  {
155  XLAL_ERROR(XLAL_EFUNC, "Alpha array allocation failed.");
156  }
157  *betaPNR = XLALCreateREAL8Sequence(iStop - iStart);
158  if (!(*betaPNR))
159  {
160  XLAL_ERROR(XLAL_EFUNC, "Beta array allocation failed.");
161  }
162  *gammaPNR = XLALCreateREAL8Sequence(iStop - iStart);
163  if (!(*gammaPNR))
164  {
165  XLAL_ERROR(XLAL_EFUNC, "Gamma array allocation failed.");
166  }
167 
168  /* Populate frequency array */
169  for (UINT4 i = iStart; i < iStop; i++)
170  {
171  (*freqs)->data[i - iStart] = i * deltaF;
172  }
173 
174  /* Specify arbitrary parameters for the angle generation */
175  REAL8 distance = 1.0;
176  REAL8 phiRef = 0.0;
177 
178  /* Use the true minimum and maximum frequency values */
179  REAL8 f_min_eval = (*freqs)->data[0];
180  REAL8 f_max_eval = (*freqs)->data[(*freqs)->length - 1];
181 
182  /* Initialize PhenomX Waveform struct and check that it initialized correctly */
184  pWF = XLALMalloc(sizeof(IMRPhenomXWaveformStruct));
185  status = IMRPhenomXSetWaveformVariables(pWF, m1_SI, m2_SI, chi1z, chi2z, deltaF, fRef, phiRef, f_min_eval, f_max_eval, distance, inclination, lalParams_aux, DEBUG);
186  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetWaveformVariables failed.\n");
187 
188  /* Initialize PhenomX Precession struct and check that it generated successfully */
190  pPrec = XLALMalloc(sizeof(IMRPhenomXPrecessionStruct));
192  pWF,
193  pPrec,
194  m1_SI,
195  m2_SI,
196  chi1x,
197  chi1y,
198  chi1z,
199  chi2x,
200  chi2y,
201  chi2z,
202  lalParams_aux,
203  DEBUG);
204  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetPrecessionVariables failed.\n");
205 
206  /* See if PNR angles were turned off in pPrec check */
207  UsePNR = pPrec->IMRPhenomXPNRUseTunedAngles;
208 
209  if(UsePNR){
210  /* Generate the tuned precession angles */
211  status = IMRPhenomX_PNR_GeneratePNRAngles(*alphaPNR, *betaPNR, *gammaPNR, *freqs, pWF, pPrec, lalParams_aux);
212  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomX_PNR_GeneratePNRAngles failed.\n");
213 
214  /* Pass through the reference values of alpha and gamma
215  * NOTE: gamma = - epsilon */
216  *alphaPNR_ref = pPrec->alpha_offset;
217  *gammaPNR_ref = -(pPrec->epsilon_offset);
218  }
219  else{
220  /* Generate PN angles */
222  /* MSA angles generated using built-in function */
223 
224  REAL8Sequence *cosBeta = XLALCreateREAL8Sequence(iStop - iStart);
225  if (!cosBeta)
226  {
227  XLAL_ERROR(XLAL_EFUNC, "cosBeta array allocation failed.");
228  }
229 
231  alphaPNR, gammaPNR, &cosBeta, *freqs,
232  m1_SI,m2_SI,
233  chi1x,chi1y,chi1z,
234  chi2x,chi2y,chi2z,
235  inclination,
236  fRef_In, 2, lalParams_aux
237  );
238 
239  for (UINT4 i = 0; i < (*freqs)->length; i++)
240  {
241  (*betaPNR)->data[i] = acos(cosBeta->data[i]);
242  }
243 
244  XLALDestroyREAL8Sequence(cosBeta);
245 
246  /* Reference values already applied in function */
247  *alphaPNR_ref = 0.0;
248  *gammaPNR_ref = 0.0;
249  }
250  else{
251  /* NNLO angles */
252 
253  REAL8Sequence *cosBeta = XLALCreateREAL8Sequence(iStop - iStart);
254  if (!cosBeta)
255  {
256  XLAL_ERROR(XLAL_EFUNC, "cosBeta array allocation failed.");
257  }
258 
260  alphaPNR, gammaPNR, &cosBeta, *freqs,
261  m1_SI,m2_SI,
262  chi1x,chi1y,chi1z,
263  chi2x,chi2y,chi2z,
264  inclination,
265  fRef_In, 2, lalParams_aux
266  );
267 
268  for (UINT4 i = 0; i < (*freqs)->length; i++)
269  {
270  (*betaPNR)->data[i] = acos(cosBeta->data[i]);
271  }
272 
273  XLALDestroyREAL8Sequence(cosBeta);
274 
275  /* Pass through the reference values of alpha and gamma
276  * NOTE: gamma = - epsilon */
277  *alphaPNR_ref = pPrec->alpha_offset;
278  *gammaPNR_ref = -(pPrec->epsilon_offset);
279  }
280  }
281 
282  /* Clean up memory allocation */
283  LALFree(pPrec);
284  LALFree(pWF);
285  XLALDestroyDict(lalParams_aux);
286 
287  return XLAL_SUCCESS;
288  }
289 
290  /**
291  * This is an external wrapper to generate the (l,m) PNR angles,
292  * following the prescriptions outlined in arXiv:2107.08876
293  * and arXiv:##### FIXME: add reference,
294  * given the standard inputs given to generate FD waveforms.
295  */
297  REAL8Sequence **alphaPNR, /**< [out] Alpha precession angle (rad) */
298  REAL8Sequence **betaPNR, /**< [out] Beta precession angle (rad) */
299  REAL8Sequence **gammaPNR, /**< [out] Gamma precession angle (rad) */
300  REAL8Sequence **freqs, /**< [out] Frequency array (Hz) */
301  REAL8 *alphaPNR_ref, /**< [out] Reference value of alpha (rad) */
302  REAL8 *gammaPNR_ref, /**< [out] Reference value of gamma (rad) */
303  REAL8 m1_SI, /**< mass of companion 1 (kg) */
304  REAL8 m2_SI, /**< mass of companion 2 (kg) */
305  REAL8 chi1x, /**< x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
306  REAL8 chi1y, /**< y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
307  REAL8 chi1z, /**< z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
308  REAL8 chi2x, /**< x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
309  REAL8 chi2y, /**< y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
310  REAL8 chi2z, /**< z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
311  REAL8 inclination, /**< Angle between orbital angular momentum and line-of-sight vector at reference frequency (rad) */
312  REAL8 deltaF, /**< Frequency spacing (Hz) */
313  REAL8 f_min, /**< Starting GW frequency (Hz) */
314  REAL8 f_max, /**< Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. */
315  REAL8 fRef_In, /**< Reference frequency (Hz) */
316  INT4 ell, /**< Orbital index (int) */
317  INT4 emmprime, /**< Azimuthal index (int) */
318  LALDict *lalParams /**< LAL Dictionary struct */
319  )
320  {
321  /* Simple check on masses and spins */
322  UINT4 status = 0;
323  status = XLALIMRPhenomXPCheckMassesAndSpins(&m1_SI, &m2_SI, &chi1x, &chi1y, &chi1z, &chi2x, &chi2y, &chi2z);
324  XLAL_CHECK(
325  XLAL_SUCCESS == status,
326  XLAL_EFUNC,
327  "Error: XLALIMRPhenomXPCheckMassesAndSpins failed in XLALSimIMRPhenomX_PNR_GeneratePNRAnglesHM.\n");
328 
329  /* Ensure we have a dictionary */
330  LALDict *lalParams_aux;
331  if (lalParams == NULL)
332  {
333  lalParams_aux = XLALCreateDict();
334  }
335  else
336  {
337  lalParams_aux = XLALDictDuplicate(lalParams);
338  }
339 
340  /* Make sure we're calling the tuned angles */
342  if (!UsePNR)
343  {
344  UsePNR = 1;
346  }
347 
348  /* make sure the HM multipoles are activated
349  * we need both (2,2) and (l,m) */
350  LALValue *ModeArray = XLALSimInspiralWaveformParamsLookupModeArray(lalParams_aux);
351 
352  /* If the mode array is empty, populate using a default choice of modes */
353  if (ModeArray == NULL)
354  {
355  ModeArray = XLALSimInspiralCreateModeArray();
356  }
357  XLALSimInspiralModeArrayActivateMode(ModeArray, 2, 2);
358  XLALSimInspiralModeArrayActivateMode(ModeArray, 2, -2);
359  XLALSimInspiralModeArrayActivateMode(ModeArray, ell, emmprime);
360  XLALSimInspiralModeArrayActivateMode(ModeArray, ell, -emmprime);
361 
362  XLALSimInspiralWaveformParamsInsertModeArray(lalParams_aux, ModeArray);
363  XLALDestroyValue(ModeArray);
364 
365  /* Map fRef to the start frequency if need be, then make sure it's within the frequency range */
366  REAL8 fRef = (fRef_In == 0.0) ? f_min : fRef_In;
367  XLAL_CHECK(
368  ((f_min <= fRef) && (fRef <= f_max)),
369  XLAL_EFUNC,
370  "Error: fRef needs to be within the specified minimum and maximum frequency values!\n");
371 
372  /* Check that the passed deltaF is non-zero */
373  XLAL_CHECK(
374  deltaF > 0,
375  XLAL_EFUNC,
376  "Error: deltaF needs to be a positive number!\n");
377 
378  /* Generate a uniformly sampled frequency grid of spacing deltaF. */
379  /* Frequencies will be set using only the lower and upper bounds that we passed */
380  size_t iStart = (size_t)(f_min / deltaF);
381  size_t iStop = (size_t)(f_max / deltaF) + 1;
382 
383  XLAL_CHECK(
384  (iStart <= iStop),
385  XLAL_EDOM,
386  "Error: the starting frequency index is greater than the stopping index! Please ensure that f_min <= f_max.\n");
387 
388  /* Allocate memory for frequency array and terminate if this fails */
389  *freqs = XLALCreateREAL8Sequence(iStop - iStart);
390  if (!(*freqs))
391  {
392  XLAL_ERROR(XLAL_EFUNC, "Frequency array allocation failed.");
393  }
394  *alphaPNR = XLALCreateREAL8Sequence(iStop - iStart);
395  if (!(*alphaPNR))
396  {
397  XLAL_ERROR(XLAL_EFUNC, "Alpha array allocation failed.");
398  }
399  *betaPNR = XLALCreateREAL8Sequence(iStop - iStart);
400  if (!(*betaPNR))
401  {
402  XLAL_ERROR(XLAL_EFUNC, "Beta array allocation failed.");
403  }
404  *gammaPNR = XLALCreateREAL8Sequence(iStop - iStart);
405  if (!(*gammaPNR))
406  {
407  XLAL_ERROR(XLAL_EFUNC, "Gamma array allocation failed.");
408  }
409 
410  /* Populate frequency array */
411  for (UINT4 i = iStart; i < iStop; i++)
412  {
413  (*freqs)->data[i - iStart] = i * deltaF;
414  }
415 
416  /* Specify arbitrary parameters for the angle generation */
417  REAL8 distance = 1.0;
418  REAL8 phiRef = 0.0;
419 
420  /* Use the true minimum and maximum frequency values */
421  REAL8 f_min_eval = (*freqs)->data[0];
422  REAL8 f_max_eval = (*freqs)->data[(*freqs)->length - 1];
423 
424  /* Initialize PhenomX Waveform struct and check that it initialized correctly */
426  pWF = XLALMalloc(sizeof(IMRPhenomXWaveformStruct));
427  status = IMRPhenomXSetWaveformVariables(pWF, m1_SI, m2_SI, chi1z, chi2z, deltaF, fRef, phiRef, f_min_eval, f_max_eval, distance, inclination, lalParams_aux, DEBUG);
428  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetWaveformVariables failed.\n");
429 
430  /* Initialize PhenomX Precession struct and check that it generated successfully */
432  pPrec = XLALMalloc(sizeof(IMRPhenomXPrecessionStruct));
434  pWF,
435  pPrec,
436  m1_SI,
437  m2_SI,
438  chi1x,
439  chi1y,
440  chi1z,
441  chi2x,
442  chi2y,
443  chi2z,
444  lalParams_aux,
445  DEBUG);
446  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetPrecessionVariables failed.\n");
447 
448  /* See if PNR angles were turned off in pPrec check */
449  UsePNR = pPrec->IMRPhenomXPNRUseTunedAngles;
450 
451  if(UsePNR){
452  /* First generate (2,2)-angle interpolants */
454  status = IMRPhenomX_PNR_GeneratePNRAngleInterpolants(hm_angle_spline, pWF, pPrec, lalParams_aux);
455  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomX_PNR_GeneratePNRAngleInterpolants failed.\n");
456 
457  REAL8 M = pWF->Mtot;
458 
459  /* Interpolate the angles */
460  /* If the function is called with (l,m)=(2,2) then just interpolate onto the frequency grid*/
461  if ((ell == 2) && (emmprime == 2))
462  {
463  for (size_t i = 0; i < (*freqs)->length; i++)
464  {
465  (*alphaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->alpha_spline, (*freqs)->data[i], hm_angle_spline->alpha_acc);
466  (*betaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->beta_spline, (*freqs)->data[i], hm_angle_spline->beta_acc);
467  (*gammaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->gamma_spline, (*freqs)->data[i], hm_angle_spline->gamma_acc);
468  }
469  }
470  else
471  { /* Called with l!=2 or m!=2, so we map the (2,2) angles onto a rescaled frequency grid */
472 
473  /* Collect the (2,2) and (l,m) ringdown frequencies for use in the frequency map */
474  REAL8 Mf_RD_22 = pWF->fRING;
475  REAL8 Mf_RD_lm = IMRPhenomXHM_GenerateRingdownFrequency(ell, emmprime, pWF);
476 
477  /* Generate interpolation transition frequencies */
478  REAL8 Mf_high = 0.0;
479  REAL8 Mf_low = 0.0;
480 
481  status = IMRPhenomX_PNR_LinearFrequencyMapTransitionFrequencies(&Mf_low, &Mf_high, emmprime, Mf_RD_22, Mf_RD_lm, pPrec);
482  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomX_PNR_LinearFrequencyMapTransitionFrequencies failed.\n");
483 
484  UINT4 toggleInspiralScaling = pPrec->PNRInspiralScaling;
485 
486  #if DEBUG == 1
487  // Save angles into a file
488  FILE *fileangle;
489  char fileSpecII[40];
490  sprintf(fileSpecII, "interpolation_frequencies.dat");
491 
492  fileangle = fopen(fileSpecII, "w");
493 
494  fprintf(fileangle, "Mf_RD_22: %.16e\n", Mf_RD_22);
495  fprintf(fileangle, "Mf_RD_lm: %.16e\n", Mf_RD_lm);
496 
497  fprintf(fileangle, "Mf_low: %.16e\n", Mf_low);
498  fprintf(fileangle, "Mf_high: %.16e\n", Mf_high);
499 
500  fclose(fileangle);
501  #endif
502  for (size_t i = 0; i < (*freqs)->length; i++)
503  {
504  REAL8 Mf = XLALSimIMRPhenomXUtilsHztoMf((*freqs)->data[i], M);
505  /* Calculate the mapped frequency */
506  REAL8 Mf_mapped = IMRPhenomX_PNR_LinearFrequencyMap(Mf, ell, emmprime, Mf_low, Mf_high, Mf_RD_22, Mf_RD_lm, toggleInspiralScaling);
507  REAL8 f_mapped = XLALSimIMRPhenomXUtilsMftoHz(Mf_mapped, M);
508 
509  (*alphaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->alpha_spline, f_mapped, hm_angle_spline->alpha_acc);
510  (*betaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->beta_spline, f_mapped, hm_angle_spline->beta_acc);
511  (*gammaPNR)->data[i] = gsl_spline_eval(hm_angle_spline->gamma_spline, f_mapped, hm_angle_spline->gamma_acc);
512  }
513  }
514 
515  /* Here we assign the reference values of alpha and gamma to their values in the precession struct */
516  /* NOTE: the contribution from pPrec->alpha0 is assigned in IMRPhenomX_PNR_RemapThetaJSF */
517  pPrec->alpha_offset = gsl_spline_eval(hm_angle_spline->alpha_spline, pWF->fRef, hm_angle_spline->alpha_acc);
518  /* NOTE: the sign is flipped between gamma and epsilon */
519  pPrec->epsilon_offset = -gsl_spline_eval(hm_angle_spline->gamma_spline, pWF->fRef, hm_angle_spline->gamma_acc) - pPrec->epsilon0; // note the sign difference between gamma and epsilon
520 
521  /* Remap the J-frame sky location to use beta instead of ThetaJN */
522  REAL8 betaPNR_ref = gsl_spline_eval(hm_angle_spline->beta_spline, pWF->fRef, hm_angle_spline->beta_acc);
523  status = IMRPhenomX_PNR_RemapThetaJSF(betaPNR_ref, pWF, pPrec, lalParams_aux);
524  XLAL_CHECK(
525  XLAL_SUCCESS == status,
526  XLAL_EFUNC,
527  "Error: IMRPhenomX_PNR_RemapThetaJSF failed in IMRPhenomX_PNR_GeneratePNRAnglesHM.");
528 
529  /* Pass through the reference values of alpha and gamma */
530  *alphaPNR_ref = pPrec->alpha_offset;
531  *gammaPNR_ref = -(pPrec->epsilon_offset);
532 
533  gsl_spline_free(hm_angle_spline->alpha_spline);
534  gsl_spline_free(hm_angle_spline->beta_spline);
535  gsl_spline_free(hm_angle_spline->gamma_spline);
536 
537  gsl_interp_accel_free(hm_angle_spline->alpha_acc);
538  gsl_interp_accel_free(hm_angle_spline->beta_acc);
539  gsl_interp_accel_free(hm_angle_spline->gamma_acc);
540 
541  LALFree(hm_angle_spline);
542  }
543  else{
544  /* Generate PN angles */
546  /* MSA angles generated using built-in function */
547 
548  REAL8Sequence *cosBeta = XLALCreateREAL8Sequence(iStop - iStart);
549  if (!cosBeta)
550  {
551  XLAL_ERROR(XLAL_EFUNC, "cosBeta array allocation failed.");
552  }
553 
555  alphaPNR, gammaPNR, &cosBeta, *freqs,
556  m1_SI,m2_SI,
557  chi1x,chi1y,chi1z,
558  chi2x,chi2y,chi2z,
559  inclination,
560  fRef_In, emmprime, lalParams_aux
561  );
562 
563  for (UINT4 i = 0; i < (*freqs)->length; i++)
564  {
565  (*betaPNR)->data[i] = acos(cosBeta->data[i]);
566  }
567 
568  XLALDestroyREAL8Sequence(cosBeta);
569 
570  /* Reference values already applied in function */
571  *alphaPNR_ref = 0.0;
572  *gammaPNR_ref = 0.0;
573  }
574  else{
575  /* NNLO angles */
576 
577  REAL8Sequence *cosBeta = XLALCreateREAL8Sequence(iStop - iStart);
578  if (!cosBeta)
579  {
580  XLAL_ERROR(XLAL_EFUNC, "cosBeta array allocation failed.");
581  }
582 
584  alphaPNR, gammaPNR, &cosBeta, *freqs,
585  m1_SI,m2_SI,
586  chi1x,chi1y,chi1z,
587  chi2x,chi2y,chi2z,
588  inclination,
589  fRef_In, emmprime, lalParams_aux
590  );
591 
592  for (UINT4 i = 0; i < (*freqs)->length; i++)
593  {
594  (*betaPNR)->data[i] = acos(cosBeta->data[i]);
595  }
596 
597  XLALDestroyREAL8Sequence(cosBeta);
598 
599  /* Pass through the reference values of alpha and gamma
600  * NOTE: gamma = - epsilon */
601  *alphaPNR_ref = pPrec->alpha_offset;
602  *gammaPNR_ref = -(pPrec->epsilon_offset);
603  }
604  }
605 
606  /* Clean up memory allocation */
607  LALFree(pPrec);
608  LALFree(pWF);
609  XLALDestroyDict(lalParams_aux);
610 
611  return XLAL_SUCCESS;
612  }
613 
614  /**
615  * Generate the tuned precession angles outlined in arXiv:2107.08876.
616  * The general flow is as follows:
617  *
618  * - First check if deltaF > 0:
619  *
620  * --> if it is, then we generate angles sampled on a uniform frequency grid
621  * and then compute the values of the angles at fRef using linear interpolation.
622  *
623  * --> otherwise we expect non-uniform frequencies, so we generate interpolants
624  * and then evaluate them at fRef.
625  *
626  * - The reference values of alpha and gamma are stored in the pPrec struct for use later,
627  * and then we re-map the J-frame sky position using the reference value of beta.
628  *
629  * - NOTE: alpha0 is recomputed in IMRPhenomX_PNR_RemapThetaJSF, hence it is not used here
630  * unlike epsilon0. Instead it is applied inside IMRPhenomX_PNR_RemapThetaJSF.
631  *
632  */
634  REAL8Sequence *alphaPNR, /**< alpha precession angle (rad) */
635  REAL8Sequence *betaPNR, /**< beta precession angle (rad) */
636  REAL8Sequence *gammaPNR, /**< gamma precession angle (rad) */
637  const REAL8Sequence *freqs, /**< input frequency array (Hz) */
638  IMRPhenomXWaveformStruct *pWF, /**< waveform struct */
639  IMRPhenomXPrecessionStruct *pPrec, /**< precession struct */
640  LALDict *lalParams /**< LAL dictionary struct */
641  )
642  {
643  /* Make sure the incoming pointers lead to something initialized */
644  XLAL_CHECK(alphaPNR != NULL, XLAL_EFAULT);
645  XLAL_CHECK(betaPNR != NULL, XLAL_EFAULT);
646  XLAL_CHECK(gammaPNR != NULL, XLAL_EFAULT);
647  XLAL_CHECK(freqs != NULL, XLAL_EFAULT);
648  XLAL_CHECK(pWF != NULL, XLAL_EFAULT);
649  XLAL_CHECK(pPrec != NULL, XLAL_EFAULT);
650  XLAL_CHECK(lalParams != NULL, XLAL_EFAULT);
651 
652  REAL8 f_ref = pWF->fRef;
653  REAL8 deltaF = pWF->deltaF;
654 
655  /* Make sure we're supposed to be here */
656  int UsePNR = pPrec->IMRPhenomXPNRUseTunedAngles;
657  XLAL_CHECK(
658  UsePNR,
659  XLAL_EFUNC,
660  "Error: PNR angles called without being activated!\n");
661 
662  INT4 status;
663 
664  /* Check for uniform frequency series */
665  if (deltaF > 0.0)
666  {
667 
668  /* Generate PNR structs */
669  IMRPhenomXWaveformStruct *pWF_SingleSpin = NULL;
670  IMRPhenomXPrecessionStruct *pPrec_SingleSpin = NULL;
671  IMRPhenomX_PNR_alpha_parameters *alphaParams = NULL;
672  IMRPhenomX_PNR_beta_parameters *betaParams = NULL;
673 
675  &pWF_SingleSpin,
676  &pPrec_SingleSpin,
677  &alphaParams,
678  &betaParams,
679  pWF,
680  pPrec,
681  lalParams);
682  XLAL_CHECK(
683  XLAL_SUCCESS == status,
684  XLAL_EFUNC,
685  "Error: IMRPhenomX_PNR_PopulateStructs failed!\n");
686 
687  /* generate angles */
689  alphaPNR,
690  betaPNR,
691  gammaPNR,
692  freqs,
693  pWF_SingleSpin,
694  pPrec_SingleSpin,
695  alphaParams,
696  betaParams,
697  pWF,
698  pPrec,
699  lalParams);
700  XLAL_CHECK(
701  XLAL_SUCCESS == status,
702  XLAL_EFUNC,
703  "Error: IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies failed in IMRPhenomX_PNR_GeneratePNRAngles.");
704 
705  /* Here we assign the reference values of alpha and gamma to their values in the precession struct.
706  * For the uniformly-sampled arrays, we linearly interpolate the angles between the
707  * next-to-nearest frequencies to get the reference values. */
708 
709  /* NOTE: the contribution from pPrec->alpha0 is assigned in IMRPhenomX_PNR_RemapThetaJSF */
710  pPrec->alpha_offset = IMRPhenomX_PNR_AngleAtFRef(alphaPNR, f_ref, freqs, deltaF);
711  /* NOTE: the sign is flipped between gamma and epsilon */
712  pPrec->epsilon_offset = -IMRPhenomX_PNR_AngleAtFRef(gammaPNR, f_ref, freqs, deltaF) - pPrec->epsilon0;
713 
714  /* Remap the J-frame sky location to use beta instead of ThetaJN */
715  REAL8 betaPNR_ref = IMRPhenomX_PNR_AngleAtFRef(betaPNR, f_ref, freqs, deltaF);
716  status = IMRPhenomX_PNR_RemapThetaJSF(betaPNR_ref, pWF, pPrec, lalParams);
717  XLAL_CHECK(
718  XLAL_SUCCESS == status,
719  XLAL_EFUNC,
720  "Error: IMRPhenomX_PNR_RemapThetaJSF failed in IMRPhenomX_PNR_GeneratePNRAngles.");
721 
722  /* clean this up */
724  &pWF_SingleSpin,
725  &pPrec_SingleSpin,
726  &alphaParams,
727  &betaParams);
728  }
729  else
730  { /* Non-uniform frequency array, so we generate interpolants to evaluate */
732  /* Generate the angle interpolants */
734  angle_spline,
735  pWF, pPrec,
736  lalParams);
737  XLAL_CHECK(
738  XLAL_SUCCESS == status,
739  XLAL_EFUNC,
740  "Error: IMRPhenomX_PNR_GeneratePNRAngleInterpolants failed in IMRPhenomX_PNR_GeneratePNRAngles");
741 
742  /* Fill the angle arrays with the values of the angles at the frequncy points */
743  REAL8 fval = 0;
744  for (UINT4 i = 0; i < freqs->length; i++)
745  {
746  fval = freqs->data[i];
747  alphaPNR->data[i] = gsl_spline_eval(angle_spline->alpha_spline, fval, angle_spline->alpha_acc);
748  betaPNR->data[i] = gsl_spline_eval(angle_spline->beta_spline, fval, angle_spline->beta_acc);
749  gammaPNR->data[i] = gsl_spline_eval(angle_spline->gamma_spline, fval, angle_spline->gamma_acc);
750  }
751 
752  /* Here we assign the reference values of alpha and gamma to their values in the precession struct */
753  /* NOTE: the contribution from pPrec->alpha0 is assigned in IMRPhenomX_PNR_RemapThetaJSF */
754  pPrec->alpha_offset = gsl_spline_eval(angle_spline->alpha_spline, f_ref, angle_spline->alpha_acc);
755  /* NOTE: the sign is flipped between gamma and epsilon */
756  pPrec->epsilon_offset = -gsl_spline_eval(angle_spline->gamma_spline, f_ref, angle_spline->gamma_acc) - pPrec->epsilon0;
757 
758  /* Remap the J-frame sky location to use beta instead of ThetaJN */
759  REAL8 betaPNR_ref = gsl_spline_eval(angle_spline->beta_spline, f_ref, angle_spline->beta_acc);
760  status = IMRPhenomX_PNR_RemapThetaJSF(betaPNR_ref, pWF, pPrec, lalParams);
761  XLAL_CHECK(
762  XLAL_SUCCESS == status,
763  XLAL_EFUNC,
764  "Error: IMRPhenomX_PNR_RemapThetaJSF failed in IMRPhenomX_PNR_GeneratePNRAngles.");
765 
766  /* Free up memory */
767  gsl_spline_free(angle_spline->alpha_spline);
768  gsl_spline_free(angle_spline->beta_spline);
769  gsl_spline_free(angle_spline->gamma_spline);
770 
771  gsl_interp_accel_free(angle_spline->alpha_acc);
772  gsl_interp_accel_free(angle_spline->beta_acc);
773  gsl_interp_accel_free(angle_spline->gamma_acc);
774 
775  LALFree(angle_spline);
776  }
777 
778  return XLAL_SUCCESS;
779  }
780 
781  /**
782  * Generate the tuned precession angles outlined in arXiv:2107.08876
783  * specifically on a uniform frequency grid.
784  *
785  * - First, we populate the required alpha and beta parameter structs,
786  * along with the two-spin structs if needed.
787  *
788  * - Next we store the two connection frequencies possibly needed for the
789  * HM angle frequency mapping.
790  *
791  * - Check if we should be attaching the MR contributions to beta: this is
792  * determined by calling the function IMRPhenomX_PNR_AttachMRBeta.
793  *
794  * --> If yes, loop through the frequencies and generate full alpha and beta
795  *
796  * --> If no, loop through the frequencies and generate full alpha but only
797  * the inspiral beta
798  *
799  * - Finally, if an initialized struct is passed for gamma, compute it.
800  *
801  */
803  REAL8Sequence *alphaPNR,
804  REAL8Sequence *betaPNR,
805  REAL8Sequence *gammaPNR,
806  const REAL8Sequence *freqs,
807  IMRPhenomXWaveformStruct *pWF_SingleSpin,
808  IMRPhenomXPrecessionStruct *pPrec_SingleSpin,
809  IMRPhenomX_PNR_alpha_parameters *alphaParams,
810  IMRPhenomX_PNR_beta_parameters *betaParams,
813  LALDict *lalParams)
814  {
815  /* Make sure the incoming pointers lead to something initialized */
816  /* Gamma could be a NULL pointer, so we don't check here */
817  XLAL_CHECK(alphaPNR != NULL, XLAL_EFAULT);
818  XLAL_CHECK(betaPNR != NULL, XLAL_EFAULT);
819  XLAL_CHECK(freqs != NULL, XLAL_EFAULT);
820  XLAL_CHECK(pWF != NULL, XLAL_EFAULT);
821  XLAL_CHECK(pPrec != NULL, XLAL_EFAULT);
822  XLAL_CHECK(lalParams != NULL, XLAL_EFAULT);
823 
824  REAL8 M = pWF->Mtot;
825 
826  /* Make sure we're supposed to be here */
827  int UsePNR = pPrec->IMRPhenomXPNRUseTunedAngles;
828  XLAL_CHECK(
829  UsePNR,
830  XLAL_EFUNC,
831  "Error: PNR angles called without being activated!\n");
832 
833  INT4 status;
834 
835  #if DEBUG == 1
836  IMRPhenomX_PNR_AngleParameterDebugPrint(alphaParams, betaParams);
837  #endif
838 
839  REAL8 Mf = 0.0;
840  /* generate PNR angles */
841  REAL8 q = pWF->q;
842  REAL8 chi = pPrec->chi_singleSpin;
843  /* inside calibration region */
844  if ((q <= pPrec->PNR_q_window_lower) && (chi <= pPrec->PNR_chi_window_lower))
845  {
846  /* First check to see if we attach the MR tuning to beta */
847  UINT4 attach_MR_beta = IMRPhenomX_PNR_AttachMRBeta(betaParams);
848  if (attach_MR_beta) /* yes we do! */
849  {
850  for (size_t i = 0; i < freqs->length; i++)
851  {
852  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
853  /* generate alpha and beta with MR tuning */
854  alphaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
855  betaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRBetaAtMf(Mf, betaParams, pWF, pPrec, pWF_SingleSpin, pPrec_SingleSpin);
856  }
857  }
858  else /* don't attach MR tuning to beta */
859  {
860  for (size_t i = 0; i < freqs->length; i++)
861  {
862  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
863 
864  /* generate alpha, generate beta with no MR tuning */
865  alphaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
866  betaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRBetaNoMR(Mf, pWF, pPrec);
867  }
868  }
869  }
870  /* inside transition region */
871  else if ((q <= pPrec->PNR_q_window_upper) && (chi <= pPrec->PNR_chi_window_upper))
872  {
873  /* First check to see if we attach the MR tuning to beta */
874  UINT4 attach_MR_beta = IMRPhenomX_PNR_AttachMRBeta(betaParams);
875  if (attach_MR_beta) /* yes we do! */
876  {
877  for (size_t i = 0; i < freqs->length; i++)
878  {
879  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
880  /* generate alpha and beta with MR tuning */
881  alphaPNR->data[i] = IMRPhenomX_PNR_GenerateMergedPNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
882  betaPNR->data[i] = IMRPhenomX_PNR_GenerateMergedPNRBetaAtMf(Mf, betaParams, pWF, pPrec, pWF_SingleSpin, pPrec_SingleSpin);
883  }
884  }
885  else /* don't attach MR tuning to beta */
886  {
887  for (size_t i = 0; i < freqs->length; i++)
888  {
889  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
890 
891  /* generate alpha, generate beta with no MR tuning */
892  alphaPNR->data[i] = IMRPhenomX_PNR_GenerateMergedPNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
893  betaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRBetaNoMR(Mf, pWF, pPrec);
894  }
895  }
896 
897  }
898  /* fully in outside calibration region */
899  else{
900  for (size_t i = 0; i < freqs->length; i++)
901  {
902  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
903 
904  /* generate MSA alpha, generate beta with no MR tuning */
905  alphaPNR->data[i] = IMRPhenomX_PNR_GetPNAlphaAtFreq(Mf, pWF, pPrec);
906  betaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRBetaNoMR(Mf, pWF, pPrec);
907  }
908  }
909 
910  /* If gamma is a NULL pointer, we assume that we don't want to
911  * compute it. This allows us to skip computing gamma for the
912  * interpolation code to reduce redundant steps */
913  if (gammaPNR != NULL)
914  {
915  /* generate gamma */
916  status = IMRPhenomX_PNR_GeneratePNRGamma(gammaPNR, freqs, alphaPNR, betaPNR);
917  XLAL_CHECK(
918  XLAL_SUCCESS == status,
919  XLAL_EFUNC,
920  "Error: IMRPhenomX_PNR_GeneratePNRGamma failed");
921  }
922 
923  return XLAL_SUCCESS;
924  }
925 
926  /**
927  * Internal helper function to generate PNR alpha for the
928  * antisymmetric waveform.
929  */
931  REAL8Sequence *alphaPNR,
932  const REAL8Sequence *freqs,
935  LALDict *lalParams)
936  {
937  /* Make sure the incoming pointers lead to something initialized */
938  /* Gamma could be a NULL pointer, so we don't check here */
939  XLAL_CHECK(alphaPNR != NULL, XLAL_EFAULT);
940  XLAL_CHECK(freqs != NULL, XLAL_EFAULT);
941  XLAL_CHECK(pWF != NULL, XLAL_EFAULT);
942  XLAL_CHECK(pPrec != NULL, XLAL_EFAULT);
943  XLAL_CHECK(lalParams != NULL, XLAL_EFAULT);
944 
945  REAL8 M = pWF->Mtot;
946 
947  INT4 status;
948 
949  /* generate alpha parameters */
951  status = IMRPhenomX_PNR_precompute_alpha_coefficients(alphaParams, pWF, pPrec);
952  XLAL_CHECK(
953  XLAL_SUCCESS == status,
954  XLAL_EFUNC,
955  "Error: IMRPhenomX_PNR_precompute_alpha_coefficients failed.\n");
956 
957  status = IMRPhenomX_PNR_alpha_connection_parameters(alphaParams, pWF, pPrec);
958  XLAL_CHECK(
959  XLAL_SUCCESS == status,
960  XLAL_EFUNC,
961  "Error: IMRPhenomX_PNR_alpha_connection_parameters failed.\n");
962 
963  REAL8 Mf = 0.0;
964  /* generate PNR angles */
965  REAL8 q = pWF->q;
966  REAL8 chi = pPrec->chi_singleSpin;
967  /* inside calibration region */
968  if ((q <= pPrec->PNR_q_window_lower) && (chi <= pPrec->PNR_chi_window_lower))
969  {
970  for (size_t i = 0; i < freqs->length; i++)
971  {
972  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
973 
974  /* generate alpha */
975  alphaPNR->data[i] = IMRPhenomX_PNR_GeneratePNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
976  }
977  }
978  /* inside transition region */
979  else if ((q <= pPrec->PNR_q_window_upper) && (chi <= pPrec->PNR_chi_window_upper))
980  {
981  for (size_t i = 0; i < freqs->length; i++)
982  {
983  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
984 
985  /* generate alpha */
986  alphaPNR->data[i] = IMRPhenomX_PNR_GenerateMergedPNRAlphaAtMf(Mf, alphaParams, pWF, pPrec);
987  }
988  }
989  /* fully in outside calibration region */
990  else
991  {
992  for (size_t i = 0; i < freqs->length; i++)
993  {
994  Mf = XLALSimIMRPhenomXUtilsHztoMf(freqs->data[i], M);
995 
996  /* generate MSA or NNLO alpha */
997  alphaPNR->data[i] = IMRPhenomX_PNR_GetPNAlphaAtFreq(Mf, pWF, pPrec);
998  }
999  }
1000 
1001  LALFree(alphaParams);
1002 
1003  return XLAL_SUCCESS;
1004  }
1005 
1006  /**
1007  * Generate the tuned precession angles outlined in arXiv:2107.08876
1008  * specifically populate the angle interpolant struct.
1009  *
1010  * - First, we check for an activated mode array.
1011  *
1012  * --> If it exists, then we're likely computing HM angles with these interpolants
1013  * and we need to extend the frequency range to capture the HM frequency
1014  * map. Use the activated modes to determine this.
1015  *
1016  * --> Otherwise we use the default f_min and f_max.
1017  *
1018  * - Get the required deltaF from IMRPhenomX_PNR_HMInterpolationDeltaF.
1019  *
1020  * - Compute uniform alpha and beta
1021  *
1022  * - Populate interpolants for alpha and beta, then compute gamma using these interpolants
1023  *
1024  */
1026  IMRPhenomX_PNR_angle_spline *angle_spline,
1029  LALDict *lalparams)
1030  {
1031  /* check for initialization */
1032  XLAL_CHECK(angle_spline != NULL, XLAL_EFAULT);
1033  XLAL_CHECK(pWF != NULL, XLAL_EFAULT);
1034  XLAL_CHECK(pPrec != NULL, XLAL_EFAULT);
1035  XLAL_CHECK(lalparams != NULL, XLAL_EFAULT);
1036 
1037  /* Generate PNR structs */
1038  IMRPhenomXWaveformStruct *pWF_SingleSpin = NULL;
1039  IMRPhenomXPrecessionStruct *pPrec_SingleSpin = NULL;
1040  IMRPhenomX_PNR_alpha_parameters *alphaParams = NULL;
1041  IMRPhenomX_PNR_beta_parameters *betaParams = NULL;
1042 
1043  UINT4 status = 0;
1045  &pWF_SingleSpin,
1046  &pPrec_SingleSpin,
1047  &alphaParams,
1048  &betaParams,
1049  pWF,
1050  pPrec,
1051  lalparams);
1052  XLAL_CHECK(
1053  XLAL_SUCCESS == status,
1054  XLAL_EFUNC,
1055  "Error: IMRPhenomX_PNR_PopulateStructs failed!\n");
1056 
1057  /* Store connection frequencies */
1058  pPrec->PNR_HM_Mfhigh = betaParams->Mf_beta_lower;
1059  pPrec->PNR_HM_Mflow = alphaParams->Mf_alpha_lower;
1060 
1061  /* Check for connection frequencies potentially being poorly defined for frequency interpolation */
1062  /* Check that the lower and upper connection frequencies are well behaved relative to fCut, fRING, and each other.
1063  * fCutDef is Mf = 0.3 is except for large chiEff. */
1064  if((pPrec->PNR_HM_Mfhigh > pWF->fCutDef) || (pPrec->PNR_HM_Mfhigh < 0.1 * pWF->fRING)){
1065  pPrec->PNR_HM_Mfhigh = pWF->fRING;
1066  }
1067  if((pPrec->PNR_HM_Mflow > pWF->fCutDef) || (pPrec->PNR_HM_Mfhigh < pPrec->PNR_HM_Mflow)){
1068  pPrec->PNR_HM_Mflow = pPrec->PNR_HM_Mfhigh / 2.0;
1069  }
1070 
1071  REAL8 f_min_22 = pWF->fMin;
1072  REAL8 f_max_22 = pWF->f_max_prime;
1073 
1074  if(pWF->deltaF != 0){
1075  /* need to account for finite frequency spacing in fMin */
1076  size_t iStart = (size_t) (f_min_22 / pWF->deltaF);
1077 
1078  REAL8 fMinFromDeltaF = iStart*pWF->deltaF;
1079  f_min_22 = (fMinFromDeltaF < f_min_22) ? fMinFromDeltaF : f_min_22;
1080  }
1081 
1082  /* Grab the mode array to find the minimum and maximum
1083  * frequency values that need to be interpolated. */
1084  LALValue *ModeArray = XLALSimInspiralWaveformParamsLookupModeArray(lalparams);
1085 
1086  if (ModeArray != NULL)
1087  {
1088  /* Assume we possibly have higher multipoles. */
1089  for (UINT4 ell = 2; ell <= 4; ell++)
1090  {
1091  for (UINT4 emmprime = 1; emmprime <= ell; emmprime++)
1092  {
1093  /* First check for (2,2) and multibanding */
1094  if(XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, emmprime) && ((ell == 2)&&(emmprime == 2)))
1095  {
1096  if(pPrec->MBandPrecVersion != 0)
1097  {
1098  /* For (2,2) with angle MB, find max (2,2) frequency */
1099  REAL8Sequence *coarseFreqs;
1100  XLALSimIMRPhenomXPHMMultibandingGrid(&coarseFreqs, ell, emmprime, pWF, lalparams);
1101 
1102  REAL8 MfCoarseFMax = coarseFreqs->data[coarseFreqs->length-1];
1103  REAL8 coarseFMax = XLALSimIMRPhenomXUtilsMftoHz(MfCoarseFMax, pWF->Mtot);
1104  if(coarseFMax > f_max_22){
1105  f_max_22 = coarseFMax;
1106  }
1107 
1108  XLALDestroyREAL8Sequence(coarseFreqs);
1109  }
1110  }
1111  else if(XLALSimInspiralModeArrayIsModeActive(ModeArray, ell, emmprime) && !((ell == 2)&&(emmprime == 2)))
1112  {
1113  if((pWF->q == 1) && (pWF->chi1L == pWF->chi2L) && (emmprime % 2 != 0))
1114  {
1115  continue;
1116  }
1117  REAL8 Mf_RD_22 = pWF->fRING;
1118  REAL8 Mf_RD_lm = IMRPhenomXHM_GenerateRingdownFrequency(ell, emmprime, pWF);
1119 
1120  REAL8 Mf_low = 0.0;
1121  REAL8 Mf_high = 0.0;
1122 
1123  IMRPhenomX_PNR_LinearFrequencyMapTransitionFrequencies(&Mf_low, &Mf_high, emmprime, Mf_RD_22, Mf_RD_lm, pPrec);
1124 
1125  UINT4 toggleInspiralScaling = pPrec->PNRInspiralScaling;
1126 
1127  if(pPrec->MBandPrecVersion != 0)
1128  {
1129  /* For (l,m) with angle MB, find max mapped (2,2) frequency */
1130  REAL8Sequence *coarseFreqs;
1131  XLALSimIMRPhenomXPHMMultibandingGrid(&coarseFreqs, ell, emmprime, pWF, lalparams);
1132 
1133  REAL8 MfCoarseFMin = coarseFreqs->data[0];
1134  REAL8 MfCoarseFMax = coarseFreqs->data[coarseFreqs->length-1];
1135 
1136  REAL8 MfMinMapped = IMRPhenomX_PNR_LinearFrequencyMap(MfCoarseFMin, ell, emmprime, Mf_low, Mf_high, Mf_RD_22, Mf_RD_lm, toggleInspiralScaling);
1137  REAL8 MfMaxMapped = IMRPhenomX_PNR_LinearFrequencyMap(MfCoarseFMax, ell, emmprime, Mf_low, Mf_high, Mf_RD_22, Mf_RD_lm, toggleInspiralScaling);
1138 
1139  REAL8 coarseFMin = XLALSimIMRPhenomXUtilsMftoHz(MfMinMapped, pWF->Mtot);
1140  REAL8 coarseFMax = XLALSimIMRPhenomXUtilsMftoHz(MfMaxMapped, pWF->Mtot);
1141 
1142  if(coarseFMin < f_min_22)
1143  {
1144  f_min_22 = coarseFMin;
1145  }
1146 
1147  if(coarseFMax > f_max_22)
1148  {
1149  f_max_22 = coarseFMax;
1150  }
1151 
1152  XLALDestroyREAL8Sequence(coarseFreqs);
1153  }
1154  else
1155  {
1156  REAL8 MfFMin = XLALSimIMRPhenomXUtilsHztoMf(f_min_22,pWF->Mtot);
1157  REAL8 MfFMax = XLALSimIMRPhenomXUtilsHztoMf(f_max_22,pWF->Mtot);
1158 
1159  REAL8 MfMinMapped = IMRPhenomX_PNR_LinearFrequencyMap(MfFMin, ell, emmprime, Mf_low, Mf_high, Mf_RD_22, Mf_RD_lm, toggleInspiralScaling);
1160  REAL8 MfMaxMapped = IMRPhenomX_PNR_LinearFrequencyMap(MfFMax, ell, emmprime, Mf_low, Mf_high, Mf_RD_22, Mf_RD_lm, toggleInspiralScaling);
1161 
1162  REAL8 FMin = XLALSimIMRPhenomXUtilsMftoHz(MfMinMapped, pWF->Mtot);
1163  REAL8 FMax = XLALSimIMRPhenomXUtilsMftoHz(MfMaxMapped, pWF->Mtot);
1164 
1165  if(FMin < f_min_22)
1166  {
1167  f_min_22 = FMin;
1168  }
1169 
1170  if(FMax > f_max_22)
1171  {
1172  f_max_22 = FMax;
1173  }
1174  }
1175  }
1176  }
1177  }
1178  } /* Otherwise f_min and f_max are already set to the (2,2) values */
1179 
1180  REAL8 f_min = f_min_22;
1181  REAL8 f_max = f_max_22;
1182 
1183  /* Get appropriate frequency spacing */
1184  REAL8 deltaF_int = IMRPhenomX_PNR_HMInterpolationDeltaF(f_min, pWF, pPrec);
1185 
1186  /* guarantee additional space to ensure no extrapolation and to reduce effects
1187  * of natural boundary conditions */
1188  f_min = (f_min - 2.0 * deltaF_int < 0) ? f_min / 2.0 : f_min - 2.0 * deltaF_int;
1189  f_max += 2.0 * deltaF_int;
1190 
1191  /* Frequencies will be set using the lower and upper bounds and computed deltaF */
1192  size_t iStart = (size_t)(f_min / deltaF_int);
1193  size_t iStop = (size_t)(f_max / deltaF_int) + 1;
1194 
1195 #if DEBUG == 1
1196  /* Save interpolation parameters into a file */
1197  FILE *fileangle;
1198  char fileSpecII[40];
1199  sprintf(fileSpecII, "interpolation_parameters.dat");
1200 
1201  fileangle = fopen(fileSpecII, "w");
1202 
1203  fprintf(fileangle, "f_min_22: %.16e\n", f_min_22);
1204  fprintf(fileangle, "f_max_22: %.16e\n", f_max_22);
1205  fprintf(fileangle, "f_min: %.16e\n", f_min);
1206  fprintf(fileangle, "f_max: %.16e\n", f_max);
1207  fprintf(fileangle, "fCut: %.16e\n", pWF->fCut);
1208  fprintf(fileangle, "f_min from df: %.16e\n", iStart * deltaF_int);
1209  fprintf(fileangle, "f_max from df: %.16e\n", iStop * deltaF_int);
1210  fprintf(fileangle, "deltaF: %.16e\n\n", deltaF_int);
1211 
1212  fclose(fileangle);
1213 #endif
1214 
1215  XLAL_CHECK((iStart <= iStop), XLAL_EDOM,
1216  "minimum freq index %zu and maximum freq index %zu do not fulfill 0<=ind_min<=ind_max.\n", iStart, iStop);
1217 
1218  /* Allocate memory for frequency array and angles, and terminate if this fails */
1219  REAL8Sequence *alphaPNR_22 = NULL;
1220  REAL8Sequence *betaPNR_22 = NULL;
1221  REAL8Sequence *gammaPNR_22 = NULL;
1222  REAL8Sequence *freqs_22 = NULL;
1223 
1224  freqs_22 = XLALCreateREAL8Sequence(iStop - iStart);
1225  if (!freqs_22)
1226  {
1227  XLAL_ERROR(XLAL_EFUNC, "freqs_22 array allocation failed.");
1228  }
1229  alphaPNR_22 = XLALCreateREAL8Sequence(iStop - iStart);
1230  if (!alphaPNR_22)
1231  {
1232  XLAL_ERROR(XLAL_EFUNC, "alphaPNR_22 array allocation failed.");
1233  }
1234  betaPNR_22 = XLALCreateREAL8Sequence(iStop - iStart);
1235  if (!betaPNR_22)
1236  {
1237  XLAL_ERROR(XLAL_EFUNC, "gammaPNR_22 array allocation failed.");
1238  }
1239  gammaPNR_22 = XLALCreateREAL8Sequence(iStop - iStart);
1240  if (!gammaPNR_22)
1241  {
1242  XLAL_ERROR(XLAL_EFUNC, "gammaPNR_22 array allocation failed.");
1243  }
1244 
1245  /* Populate frequency array */
1246  for (UINT4 i = iStart; i < iStop; i++)
1247  {
1248  freqs_22->data[i - iStart] = i * deltaF_int;
1249  }
1250 
1251  /* Generate the (2,2) angles on this uniform frequency grid.
1252  * Pass gamma as a NULL pointer, since we can save a step by
1253  * re-using the gsl splines constructed to compute gamma later on */
1255  alphaPNR_22,
1256  betaPNR_22,
1257  NULL,
1258  freqs_22,
1259  pWF_SingleSpin,
1260  pPrec_SingleSpin,
1261  alphaParams,
1262  betaParams,
1263  pWF,
1264  pPrec,
1265  lalparams);
1266  XLAL_CHECK(
1267  XLAL_SUCCESS == status,
1268  XLAL_EFUNC,
1269  "Error: IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies failed in IMRPhenomX_PNR_GeneratePNRAngleInterpolants.\n");
1270 
1271  /* Construct the splines for interpolation and populate the alpha and beta splines */
1272  UINT4 flen = freqs_22->length;
1273 
1274  gsl_spline *a_spline = gsl_spline_alloc(gsl_interp_cspline, flen);
1275  gsl_spline *b_spline = gsl_spline_alloc(gsl_interp_cspline, flen);
1276  gsl_spline *g_spline = gsl_spline_alloc(gsl_interp_cspline, flen);
1277 
1278  gsl_interp_accel *a_acc = gsl_interp_accel_alloc();
1279  gsl_interp_accel *b_acc = gsl_interp_accel_alloc();
1280  gsl_interp_accel *g_acc = gsl_interp_accel_alloc();
1281 
1282  gsl_spline_init(a_spline, freqs_22->data, alphaPNR_22->data, flen);
1283  gsl_spline_init(b_spline, freqs_22->data, betaPNR_22->data, flen);
1284 
1285  angle_spline->alpha_spline = a_spline;
1286  angle_spline->beta_spline = b_spline;
1287 
1288  angle_spline->alpha_acc = a_acc;
1289  angle_spline->beta_acc = b_acc;
1290 
1291  /* Now use the alpha and beta splines to compute gamma */
1293  gammaPNR_22,
1294  freqs_22,
1295  angle_spline);
1296  XLAL_CHECK(
1297  XLAL_SUCCESS == status,
1298  XLAL_EFUNC,
1299  "Error: IMRPhenomX_PNR_GeneratePNRGamma_FromInterpolants failed in IMRPhenomX_PNR_GeneratePNRAngleInterpolants.\n");
1300 
1301  /* Reset the gsl accelerators for alpha and beta
1302  * Not sure if this is helpful =\_(~~)_/= */
1303  gsl_interp_accel_reset(angle_spline->alpha_acc);
1304  gsl_interp_accel_reset(angle_spline->beta_acc);
1305 
1306  /* Populate the gamma spline */
1307  gsl_spline_init(g_spline, freqs_22->data, gammaPNR_22->data, flen);
1308 
1309  angle_spline->gamma_spline = g_spline;
1310  angle_spline->gamma_acc = g_acc;
1311 
1312  /* Clean up memory allocation
1313  * In particular, the splines copy over the data
1314  * so we can get rid of the 22 angle arrays */
1315  XLALDestroyValue(ModeArray);
1316 
1317  XLALDestroyREAL8Sequence(alphaPNR_22);
1318  XLALDestroyREAL8Sequence(betaPNR_22);
1319  XLALDestroyREAL8Sequence(gammaPNR_22);
1320  XLALDestroyREAL8Sequence(freqs_22);
1321 
1322  /* clean this up */
1324  &pWF_SingleSpin,
1325  &pPrec_SingleSpin,
1326  &alphaParams,
1327  &betaParams);
1328 
1329  return XLAL_SUCCESS;
1330  }
1331 
1332  /**
1333  * External wrapper for IMRPhenomX_PNR_GenerateEffectiveRingdownFreq
1334  *
1335  * Computes the effective ringdown frequency for a given (l,m) multipole
1336  * following the prescription given in arXiv:##### FIXME: add citation
1337  *
1338  */
1340  REAL8 m1_SI, /**< mass of companion 1 (kg) */
1341  REAL8 m2_SI, /**< mass of companion 2 (kg) */
1342  REAL8 chi1x, /**< x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1343  REAL8 chi1y, /**< y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1344  REAL8 chi1z, /**< z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1345  REAL8 chi2x, /**< x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1346  REAL8 chi2y, /**< y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1347  REAL8 chi2z, /**< z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1348  REAL8 f_min, /**< Starting GW frequency (Hz) */
1349  REAL8 f_max, /**< Ending GW frequency (Hz) */
1350  REAL8 fRef, /**< Reference frequency (Hz) */
1351  UINT4 ell, /**< Orbital index */
1352  UINT4 emmprime, /**< azimuthal index */
1353  LALDict *lalParams /**< LAL Dictionary struct */
1354  )
1355  {
1356  /* Ensure we have a dictionary */
1357  LALDict *lalParams_aux;
1358  if (lalParams == NULL)
1359  {
1360  lalParams_aux = XLALCreateDict();
1361  }
1362  else
1363  {
1364  lalParams_aux = XLALDictDuplicate(lalParams);
1365  }
1366 
1367  UINT4 status = 0;
1368  status = XLALIMRPhenomXPCheckMassesAndSpins(&m1_SI, &m2_SI, &chi1x, &chi1y, &chi1z, &chi2x, &chi2y, &chi2z);
1369  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: XLALIMRPhenomXPCheckMassesAndSpins failed.\n");
1370 
1372  pWF = XLALMalloc(sizeof(IMRPhenomXWaveformStruct));
1373  status = IMRPhenomXSetWaveformVariables(pWF, m1_SI, m2_SI, chi1z, chi2z, 0.0, fRef, 0.0, f_min, f_max, 1.0, 0.0, lalParams_aux, 0);
1374  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetWaveformVariables failed.\n");
1375 
1376  /* Initialize IMR PhenomX Precession struct and check that it generated successfully */
1378  pPrec = XLALMalloc(sizeof(IMRPhenomXPrecessionStruct));
1380  pWF,
1381  pPrec,
1382  m1_SI,
1383  m2_SI,
1384  chi1x,
1385  chi1y,
1386  chi1z,
1387  chi2x,
1388  chi2y,
1389  chi2z,
1390  lalParams_aux,
1391  DEBUG);
1392  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetPrecessionVariables failed.\n");
1393 
1394  /* get effective ringdown frequency */
1395  REAL8 effRD = IMRPhenomX_PNR_GenerateEffectiveRingdownFreq(pWF, ell, emmprime, lalParams_aux);
1396 
1397  /* clean up memory allocation */
1398  LALFree(pWF);
1399  LALFree(pPrec);
1400  XLALDestroyDict(lalParams_aux);
1401 
1402  return effRD;
1403  }
1404 
1405  /**
1406  * External wrapper for IMRPhenomX_PNR_GenerateRingdownPNRBeta.
1407  *
1408  * Generate PNR beta as described in Eq. 60 or arXiv:2107.08876 and
1409  * evaluate at the upper connection frequency f_f described in Sec. 8C.
1410  */
1412  REAL8 m1_SI, /**< mass of companion 1 (kg) */
1413  REAL8 m2_SI, /**< mass of companion 2 (kg) */
1414  REAL8 chi1x, /**< x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1415  REAL8 chi1y, /**< y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1416  REAL8 chi1z, /**< z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1417  REAL8 chi2x, /**< x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1418  REAL8 chi2y, /**< y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1419  REAL8 chi2z, /**< z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1420  REAL8 deltaF, /**< Sampling Frequency (Hz) */
1421  REAL8 f_min, /**< Starting GW frequency (Hz) */
1422  REAL8 f_max, /**< Ending GW frequency (Hz) */
1423  REAL8 fRef_In, /**< Reference frequency (Hz) */
1424  LALDict *lalParams /**< LAL Dictionary struct */
1425  )
1426  {
1427  /* Ensure we have a dictionary */
1428  LALDict *lalParams_aux;
1429  if (lalParams == NULL)
1430  {
1431  lalParams_aux = XLALCreateDict();
1432  }
1433  else
1434  {
1435  lalParams_aux = XLALDictDuplicate(lalParams);
1436  }
1437 
1439  if (!UsePNR)
1440  {
1441  UsePNR = 1;
1443  }
1444 
1445  REAL8 fRef = (fRef_In == 0.0) ? f_min : fRef_In;
1446 
1447  UINT4 status;
1448 
1449  REAL8 distance = 1.0;
1450  REAL8 inclination = 0.0;
1451  REAL8 phiRef = 0.0;
1452 
1453  /* Initialize IMR PhenomX Waveform struct and check that it initialized correctly */
1455  pWF = XLALMalloc(sizeof(IMRPhenomXWaveformStruct));
1456  status = IMRPhenomXSetWaveformVariables(pWF, m1_SI, m2_SI, chi1z, chi2z, deltaF, fRef, phiRef, f_min, f_max, distance, inclination, lalParams_aux, DEBUG);
1457  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetWaveformVariables failed.\n");
1458 
1459  /* Initialize IMR PhenomX Precession struct and check that it generated successfully */
1461  pPrec = XLALMalloc(sizeof(IMRPhenomXPrecessionStruct));
1463  pWF,
1464  pPrec,
1465  m1_SI,
1466  m2_SI,
1467  chi1x,
1468  chi1y,
1469  chi1z,
1470  chi2x,
1471  chi2y,
1472  chi2z,
1473  lalParams_aux,
1474  DEBUG);
1475  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetPrecessionVariables failed.\n");
1476 
1477  /* get the RD beta */
1478  REAL8 finalBeta = IMRPhenomX_PNR_GenerateRingdownPNRBeta(pWF, pPrec);
1479 
1480  /* free up memory allocation */
1481  LALFree(pPrec);
1482  LALFree(pWF);
1483  XLALDestroyDict(lalParams_aux);
1484 
1485  return finalBeta;
1486  }
1487 
1488  /* Function for testing deltaF spacing for interpolation */
1490  INT4 *twospin, /**< UNDOCUMENTED */
1491  INT4 *precversion, /**< UNDOCUMENTED */
1492  REAL8 *deltaF, /**< UNDOCUMENTED */
1493  REAL8 m1_SI, /**< mass of companion 1 (kg) */
1494  REAL8 m2_SI, /**< mass of companion 2 (kg) */
1495  REAL8 chi1x, /**< x-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1496  REAL8 chi1y, /**< y-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1497  REAL8 chi1z, /**< z-component of the dimensionless spin of object 1 w.r.t. Lhat = (0,0,1) */
1498  REAL8 chi2x, /**< x-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1499  REAL8 chi2y, /**< y-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1500  REAL8 chi2z, /**< z-component of the dimensionless spin of object 2 w.r.t. Lhat = (0,0,1) */
1501  REAL8 f_min, /**< Starting GW frequency (Hz) */
1502  REAL8 f_max, /**< Ending GW frequency (Hz); Defaults to Mf = 0.3 if no f_max is specified. */
1503  REAL8 fRef, /**< Reference frequency (Hz) */
1504  LALDict *lalParams /**< LAL Dictionary struct */)
1505  {
1506  /* Ensure we have a dictionary */
1507  LALDict *lalParams_aux;
1508  if (lalParams == NULL)
1509  {
1510  lalParams_aux = XLALCreateDict();
1511  }
1512  else
1513  {
1514  lalParams_aux = XLALDictDuplicate(lalParams);
1515  }
1516 
1517  UINT4 status = 0;
1518  status = XLALIMRPhenomXPCheckMassesAndSpins(&m1_SI, &m2_SI, &chi1x, &chi1y, &chi1z, &chi2x, &chi2y, &chi2z);
1519  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: XLALIMRPhenomXPCheckMassesAndSpins failed.\n");
1520 
1522  pWF = XLALMalloc(sizeof(IMRPhenomXWaveformStruct));
1523  status = IMRPhenomXSetWaveformVariables(pWF, m1_SI, m2_SI, chi1z, chi2z, 0.0, fRef, 0.0, f_min, f_max, 1.0, 0.0, lalParams_aux, 0);
1524  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetWaveformVariables failed.\n");
1525 
1526  /* Initialize IMR PhenomX Precession struct and check that it generated successfully */
1528  pPrec = XLALMalloc(sizeof(IMRPhenomXPrecessionStruct));
1529 
1530  status = IMRPhenomXGetAndSetPrecessionVariables( // needed to adjust ringdown frequency in pWF
1531  pWF,
1532  pPrec,
1533  m1_SI,
1534  m2_SI,
1535  chi1x,
1536  chi1y,
1537  chi1z,
1538  chi2x,
1539  chi2y,
1540  chi2z,
1541  lalParams_aux,
1542  DEBUG);
1543  XLAL_CHECK(XLAL_SUCCESS == status, XLAL_EFUNC, "Error: IMRPhenomXSetPrecessionVariables failed.\n");
1544 
1545  *deltaF = IMRPhenomX_PNR_HMInterpolationDeltaF(f_min, pWF, pPrec);
1546  *twospin = IMRPhenomX_PNR_CheckTwoSpin(pPrec);
1547  *precversion = pPrec->IMRPhenomXPrecVersion;
1548 
1549  LALFree(pWF);
1550  LALFree(pPrec);
1551  XLALDestroyDict(lalParams_aux);
1552 
1553  return XLAL_SUCCESS;
1554  }
1555 
1556  /** @} */
1557  /** @} */
1558 
1559 #ifdef __cplusplus
1560 }
1561 #endif
void XLALDestroyDict(LALDict *dict)
LALDict * XLALDictDuplicate(LALDict *old)
LALDict * XLALCreateDict(void)
#define LALFree(p)
#define DEBUG
int IMRPhenomX_PNR_precompute_alpha_coefficients(IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
This function evaluates the coefficients outlined in Sec 7A of arXiv:2107.08876 for alpha.
REAL8 IMRPhenomX_PNR_GetPNAlphaAtFreq(REAL8 Mf, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
Here we write a wrapper function to produce either MSA or NNLO alpha for use in IMRPhenomX_PNR_Genera...
REAL8 IMRPhenomX_PNR_GeneratePNRAlphaAtMf(REAL8 Mf, const IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
This function evaluates Eq.
int IMRPhenomX_PNR_alpha_connection_parameters(IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
This function evaluates the connection frequencies Mf1 and Mf2 detailed in Sec.
REAL8 IMRPhenomX_PNR_GenerateMergedPNRAlphaAtMf(REAL8 Mf, const IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
This function generates the blended PNR and PN expressions for alpha for the transition region of par...
UINT4 IMRPhenomX_PNR_AttachMRBeta(const IMRPhenomX_PNR_beta_parameters *betaParams)
Determine whether to attach the MR contributions to beta.
REAL8 IMRPhenomX_PNR_GeneratePNRBetaNoMR(REAL8 Mf, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
This function evaluates only the rescaled inspiral beta given in Eq.
REAL8 IMRPhenomX_PNR_GenerateRingdownPNRBeta(IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
We evaluate beta at the final Mf_beta_upper connection frequency to approximate the final value of be...
REAL8 IMRPhenomX_PNR_GenerateMergedPNRBetaAtMf(REAL8 Mf, const IMRPhenomX_PNR_beta_parameters *betaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, IMRPhenomXWaveformStruct *pWF_SingleSpin, IMRPhenomXPrecessionStruct *pPrec_SingleSpin)
This function generates beta with the tuned angles and PN expressions blended during merger-ringdown.
REAL8 IMRPhenomX_PNR_GeneratePNRBetaAtMf(REAL8 Mf, const IMRPhenomX_PNR_beta_parameters *betaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, IMRPhenomXWaveformStruct *pWF_SingleSpin, IMRPhenomXPrecessionStruct *pPrec_SingleSpin)
This function evaluates Eqs.
void IMRPhenomX_PNR_FreeStructs(IMRPhenomXWaveformStruct **pWF_SingleSpin, IMRPhenomXPrecessionStruct **pPrec_SingleSpin, IMRPhenomX_PNR_alpha_parameters **alphaParams, IMRPhenomX_PNR_beta_parameters **betaParams)
void IMRPhenomX_PNR_AngleParameterDebugPrint(IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomX_PNR_beta_parameters *betaParams)
Print various parameters in the angle structs.
INT4 IMRPhenomX_PNR_CheckTwoSpin(IMRPhenomXPrecessionStruct *pPrec)
This function quickly checks to see if we expect two-spin effects to be present in the inspiral prece...
REAL8 IMRPhenomX_PNR_GenerateEffectiveRingdownFreq(IMRPhenomXWaveformStruct *pWF, UINT4 ell, UINT4 emmprime, LALDict *lalParams)
This code recomputes the skymapped locations in the J-frame using the new value of beta computed from...
INT4 IMRPhenomX_PNR_PopulateStructs(IMRPhenomXWaveformStruct **pWF_SingleSpin, IMRPhenomXPrecessionStruct **pPrec_SingleSpin, IMRPhenomX_PNR_alpha_parameters **alphaParams, IMRPhenomX_PNR_beta_parameters **betaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
These two functions create and populate the required parameter structs for PNR.
REAL8 IMRPhenomX_PNR_HMInterpolationDeltaF(REAL8 f_min, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec)
Here we compute an appropriate deltaF to be used when generating the (2,2) angle interpolants and map...
REAL8 IMRPhenomX_PNR_LinearFrequencyMap(REAL8 Mf, REAL8 ell, REAL8 emm, REAL8 Mf_lower, REAL8 Mf_upper, REAL8 Mf_RD_22, REAL8 Mf_RD_lm, UINT4 INSPIRAL)
Computes a linear frequency map for the HM PNR angles based on the linear frequency mapping used orig...
INT4 IMRPhenomX_PNR_GeneratePNRGamma(REAL8Sequence *gamma, const REAL8Sequence *freqs, const REAL8Sequence *alpha, const REAL8Sequence *beta)
This function computes the frequency integral outlined in Eq.
INT4 IMRPhenomX_PNR_RemapThetaJSF(REAL8 beta_ref, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
This code recomputes the skymapped locations in the J-frame using the new value of beta computed from...
int IMRPhenomX_PNR_GeneratePNRGamma_FromInterpolants(REAL8Sequence *gamma, const REAL8Sequence *freqs, IMRPhenomX_PNR_angle_spline *ab_splines)
This function computes the frequency integral outlined in Eq.
INT4 IMRPhenomX_PNR_LinearFrequencyMapTransitionFrequencies(REAL8 *Mf_low, REAL8 *Mf_high, REAL8 emmprime, REAL8 Mf_RD_22, REAL8 Mf_RD_lm, IMRPhenomXPrecessionStruct *pPrec)
Compute the transition frequencies for the HM PNR angle mapping.
REAL8 IMRPhenomX_PNR_AngleAtFRef(const REAL8Sequence *angle, const REAL8 f_ref, const REAL8Sequence *freqs, const REAL8 deltaF)
Evaluates a function at two points and interpolates between them.
int IMRPhenomXSetWaveformVariables(IMRPhenomXWaveformStruct *wf, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 chi1L_In, const REAL8 chi2L_In, const REAL8 deltaF, const REAL8 fRef, const REAL8 phi0, const REAL8 f_min, const REAL8 f_max, const REAL8 distance, const REAL8 inclination, LALDict *LALParams, UNUSED const UINT4 debug)
int IMRPhenomXGetAndSetPrecessionVariables(IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, LALDict *lalParams, INT4 debug_flag)
Function to populate the IMRPhenomXPrecessionStruct:
REAL8 IMRPhenomXHM_GenerateRingdownFrequency(UINT4 ell, UINT4 emm, IMRPhenomXWaveformStruct *wf22)
INT4 XLALSimIMRPhenomXPHMMultibandingGrid(REAL8Sequence **coarseFreqs, UINT4 ell, UINT4 emmprime, IMRPhenomXWaveformStruct *pWF, LALDict *lalParams)
INT4 XLALIMRPhenomXPCheckMassesAndSpins(REAL8 *m1, REAL8 *m2, REAL8 *chi1x, REAL8 *chi1y, REAL8 *chi1z, REAL8 *chi2x, REAL8 *chi2y, REAL8 *chi2z)
Check if m1 > m2.
REAL8 XLALSimIMRPhenomXUtilsHztoMf(REAL8 fHz, REAL8 Mtot_Msun)
Convert from frequency in Hz to geometric frequency.
REAL8 XLALSimIMRPhenomXUtilsMftoHz(REAL8 Mf, REAL8 Mtot_Msun)
Convert from geometric frequency to Hz.
int XLALSimInspiralWaveformParamsInsertModeArray(LALDict *params, LALValue *value)
LALValue * XLALSimInspiralWaveformParamsLookupModeArray(LALDict *params)
INT4 XLALSimInspiralWaveformParamsLookupPhenomXPrecVersion(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXPNRUseTunedAngles(LALDict *params, INT4 value)
INT4 XLALSimInspiralWaveformParamsLookupPhenomXPNRUseTunedAngles(LALDict *params)
void XLALDestroyValue(LALValue *value)
#define fprintf
REAL8 M
Definition: bh_qnmode.c:133
double i
Definition: bh_ringdown.c:118
double REAL8
uint32_t UINT4
int32_t INT4
void * XLALMalloc(size_t n)
int XLALSimIMRPhenomX_PNR_GeneratePNRAnglesHM(REAL8Sequence **alphaPNR, REAL8Sequence **betaPNR, REAL8Sequence **gammaPNR, REAL8Sequence **freqs, REAL8 *alphaPNR_ref, REAL8 *gammaPNR_ref, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, INT4 ell, INT4 emmprime, LALDict *lalParams)
This is an external wrapper to generate the (l,m) PNR angles, following the prescriptions outlined in...
int XLALSimIMRPhenomX_PNR_GeneratePNRAngles(REAL8Sequence **alphaPNR, REAL8Sequence **betaPNR, REAL8Sequence **gammaPNR, REAL8Sequence **freqs, REAL8 *alphaPNR_ref, REAL8 *gammaPNR_ref, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, LALDict *lalParams)
This is an external wrapper to generate the (2,2) PNR angles, following the prescription outlined in ...
REAL8 XLALSimIMRPhenomX_PNR_GenerateEffectiveRingdownFreq(REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 f_min, REAL8 f_max, REAL8 fRef, UINT4 ell, UINT4 emmprime, LALDict *lalParams)
External wrapper for IMRPhenomX_PNR_GenerateEffectiveRingdownFreq.
int IMRPhenomX_PNR_GeneratePNRAlphaForAntisymmetry(REAL8Sequence *alphaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
Internal helper function to generate PNR alpha for the antisymmetric waveform.
int XLALSimIMRPhenomXPMSAAngles(REAL8Sequence **alpha_of_f, REAL8Sequence **gamma_of_f, REAL8Sequence **cosbeta_of_f, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 fRef_In, INT4 mprime, LALDict *lalParams)
int IMRPhenomX_PNR_GeneratePNRAngles_UniformFrequencies(REAL8Sequence *alphaPNR, REAL8Sequence *betaPNR, REAL8Sequence *gammaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF_SingleSpin, IMRPhenomXPrecessionStruct *pPrec_SingleSpin, IMRPhenomX_PNR_alpha_parameters *alphaParams, IMRPhenomX_PNR_beta_parameters *betaParams, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
Generate the tuned precession angles outlined in arXiv:2107.08876 specifically on a uniform frequency...
int XLALSimIMRPhenomX_PNR_InterpolationDeltaF(INT4 *twospin, INT4 *precversion, REAL8 *deltaF, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 f_min, REAL8 f_max, REAL8 fRef, LALDict *lalParams)
int IMRPhenomX_PNR_GeneratePNRAngles(REAL8Sequence *alphaPNR, REAL8Sequence *betaPNR, REAL8Sequence *gammaPNR, const REAL8Sequence *freqs, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalParams)
Generate the tuned precession angles outlined in arXiv:2107.08876.
int XLALSimIMRPhenomXPPNAngles(REAL8Sequence **alpha_of_f, REAL8Sequence **gamma_of_f, REAL8Sequence **cosbeta_of_f, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 inclination, REAL8 fRef_In, INT4 mprime, LALDict *lalParams)
REAL8 XLALSimIMRPhenomX_PNR_GenerateRingdownPNRBeta(REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 fRef_In, LALDict *lalParams)
External wrapper for IMRPhenomX_PNR_GenerateRingdownPNRBeta.
int IMRPhenomX_PNR_GeneratePNRAngleInterpolants(IMRPhenomX_PNR_angle_spline *angle_spline, IMRPhenomXWaveformStruct *pWF, IMRPhenomXPrecessionStruct *pPrec, LALDict *lalparams)
Generate the tuned precession angles outlined in arXiv:2107.08876 specifically populate the angle int...
LALValue * XLALSimInspiralCreateModeArray(void)
Create a LALValue pointer to store the mode array.
int XLALSimInspiralModeArrayIsModeActive(LALValue *modes, unsigned l, int m)
LALValue * XLALSimInspiralModeArrayActivateMode(LALValue *modes, unsigned l, int m)
static const INT4 q
void XLALDestroyREAL8Sequence(REAL8Sequence *sequence)
REAL8Sequence * XLALCreateREAL8Sequence(size_t length)
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
string status
REAL8 Mf_alpha_lower
connection frequency
gsl_spline * beta_spline
beta cubic spline
gsl_interp_accel * beta_acc
beta cubic spline accelerator
gsl_interp_accel * alpha_acc
alpha cubic spline accelerator
gsl_interp_accel * gamma_acc
gamma cubic spline accelerator
gsl_spline * gamma_spline
gamma cubic spline
gsl_spline * alpha_spline
alpha cubic spline
REAL8 PNR_HM_Mflow
Mf_alpha_lower stored from alphaParams struct, 2 A4 / 7 from arXiv:2107.08876.
INT4 MBandPrecVersion
Flag to control multibanding for Euler angles.
INT4 IMRPhenomXPrecVersion
Flag to set version of Euler angles used.
UINT4 PNRInspiralScaling
Enforce inpsiral scaling for HM angles outside of calibration window.
REAL8 chi_singleSpin
Magnitude of effective single spin used for tapering two-spin angles, Eq.
REAL8 PNR_HM_Mfhigh
Mf_beta_lower stored from betaParams struct, Eq.
REAL8 * data
double f_min
Definition: unicorn.c:22
double f_max
Definition: unicorn.c:23