LALSimulation  5.4.0.1-fe68b98
LALSimTEOBResumS.h
Go to the documentation of this file.
1 
2 /**
3  * This file is part of TEOBResumS
4  *
5  * Copyright (C) 2017-2018 Alessandro Nagar, Sebastiano Bernuzzi,
6  * Sarp Ackay, Gregorio Carullo, Walter Del Pozzo, Ka Wa Tsang, Michalis Agathos
7  * LALSimulation implementation by Michalis Agathos
8  *
9  * TEOBResumS is free software; you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation; either version 2 of the License, or
12  * (at your option) any later version.
13  *
14  * TEOBResumS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program. If not, see http://www.gnu.org/licenses/.
21  *
22  */
23 
24 #ifndef _LALSIM_TEOB_RESUMS_H
25 #define _LALSIM_TEOB_RESUMS_H
26 
27 #if defined(__cplusplus)
28 extern "C" {
29 #elif 0
30 } /* so that editors will match preceding brace */
31 #endif
32 
33 /**
34  * @file LALSimTEOBResumS.h
35  * @brief Header file of the TEOBResumS C code
36  *
37  * This file contains all the macros, typdef, and routine prototype.
38  * Doxygen documentation should go here.
39  *
40  */
41 
42 #include <stdlib.h>
43 #include <stdio.h>
44 #include <stdbool.h>
45 #include <complex.h>
46 #include <math.h>
47 #include <string.h>
48 #include <time.h>
49 
50 #include <gsl/gsl_math.h>
51 #include <gsl/gsl_sf.h>
52 #include <gsl/gsl_complex.h>
53 #include <gsl/gsl_complex_math.h>
54 #include <gsl/gsl_roots.h>
55 #include <gsl/gsl_errno.h>
56 #include <gsl/gsl_complex.h>
57 #include <gsl/gsl_spline.h>
58 #include <gsl/gsl_matrix.h>
59 #include <gsl/gsl_odeiv2.h>
60 
61 #include <lal/Date.h>
62 #include <lal/LALString.h>
63 #include <lal/LALDatatypes.h>
64 #include <lal/Units.h>
65 #include <lal/BandPassTimeSeries.h>
66 #include <lal/SphericalHarmonics.h>
67 #include <lal/LALSimInspiral.h>
68 #include <lal/LALConstants.h>
69 #include <lal/LALStdlib.h>
70 #include <lal/Sequence.h>
71 #include <lal/TimeSeries.h>
72 #include <lal/FrequencySeries.h>
73 #include <lal/TimeFreqFFT.h>
74 #include "LALSimSphHarmSeries.h"
75 #include "LALSimSphHarmMode.h"
76 #include "LALSimIMR.h"
77 
78 
79 #ifndef POSTPOSTCIRCULAR
80 #define POSTPOSTCIRCULAR 1 /* use post-post-circular initial conditions by default */
81 #endif
82 
83 #ifndef RC_EXCLUDESPINSPINTIDES
84 #define RC_EXCLUDESPINSPINTIDES 0 /* use tidally deformed centr. radius with self-spin and tides by default */
85 #endif
86 
87 #ifndef USEBTIDALPOTENTIAL
88 #define USEBTIDALPOTENTIAL 1 /* add B LO tidal potential */
89 #endif
90 
91 #ifndef ODESTOPAFTERNDT
92 #define ODESTOPAFTERNDT (4) /* UNUSED: stop ODE integration after this many timesteps beyond the peak */
93 #endif
94 
95 #define MAX_ITER (200)
96 #define TOLERANCE (1e-14)
97 
98 /* Macros for postadiabatic stage and ID */
99 #define POSTADIABATIC_RMIN_BNS (14)
100 #define POSTADIABATIC_RMIN_BBH (14)
101 #define POSTADIABATIC_DR (0.1)
102 #define POSTADIABATIC_NSTEP_MIN (10)
103 #define POSTADIABATIC_N (8)
104 #define TEOB_R0_THRESHOLD (14)
105 
106 #define TEOB_LAMBDA_TOL (1.0)
107 
108 #define TEOB_GE_TIDES_DEFAULT LAL_SIM_INSPIRAL_GETIDES_GSF23
109 #define TEOB_GM_TIDES_DEFAULT LAL_SIM_INSPIRAL_GMTIDES_PN
110 
111 #define TEOB_MODES_BNS_DEFAULT TEOB_MODES_22
112 #define TEOB_MODES_BBH_DEFAULT TEOB_MODES_22
113 #define TEOB_MODES_NSBH_DEFAULT TEOB_MODES_22
114 #define TEOB_MODES_BHNS_DEFAULT TEOB_MODES_22
115 
116 #define RINGDOWN_EXTEND_ARRAY (1500)
117 // #define DT_MERGER_INTERP (0.5)
118 
119 #define INTERP_UNIFORM_GRID_DEFAULT INTERP_UNIFORM_GRID_HLM
120 
121 /** Macros */
122 #ifndef DEBUG
123 #define DEBUG 0
124 #endif
125 
126 #define TEOB_STRLEN (128) /** Standard string length */
127 #define TEOBResumS_Info "TEOBResumS code (C) 2017\n"
128 #define TEOBResumS_Usage(x) {printf("%sUSAGE:\t%s <parfile>\n", TEOBResumS_Info, x);}
129 #define SIGN(x,y) ((y) >= 0.0 ? fabs(x) : -fabs(x))
130 #define typeof __typeof__
131 #define MAX(a,b) \
132  ({ typeof (a) _a = (a); \
133  typeof (b) _b = (b); \
134  _a > _b ? _a : _b; })
135 #define MIN(a,b) \
136  ({ typeof (a) _a = (a); \
137  typeof (b) _b = (b); \
138  _a < _b ? _a : _b; })
139 #define MAX3(a,b,c) (((a) > (b)) ? MAX(a,c) : MAX(b,c))
140 #define MIN3(a,b,c) (((a) < (b)) ? MIN(a,c) : MIN(b,c))
141 #define SQ(a) ((a)*(a))
142 #define DEQUAL(a,b,eps) (fabs((a)-(b))<(eps)) /** double compare */
143 #define DUNEQUAL(a,b,eps) (fabs((a)-(b))>(eps))
144 #define STREQUAL(s,t) ((strcmp((s),(t))==0)) /** string compare */
145 #define SWAPTRS(a,b) \
146  ({ \
147  typeof(a) temp; \
148  temp = a; \
149  a = b; \
150  b = temp; \
151  })
152 
153 /* Useful constants */
154 // LAL_SQRT2
155 // #define Sqrt2 (1.41421356237309504880168872420969808)
156 #define Sqrt3 (1.73205080756887729352744634150587237)
157 // LAL_SQRT1_2
158 //#define ooSqrt2 (0.707106781186547524400844362104849039284836)
159 #define Log1 (0.)
160 // LAL_LN2
161 #define Log2 (0.693147180559945309417232)
162 #define Log3 (1.09861228866810969139525)
163 #define Log4 (1.38629436111989061883446)
164 #define Log5 (1.60943791243410037460076)
165 #define Log6 (1.79175946922805500081248)
166 #define Log7 (1.94591014905531330510535)
167 //#define EulerGamma_Log2 (1.27036284546147817002374) /** EulerGamma + Log2 */
168 
169 /** Index list of EOB evolved variables */
170 enum{
176 };
177 // static const char* LALTEOBResumSEvolutionVariables[] = {"r","phi","Prstar","Pphi"};
178 
179 /** Index list of EOB variables for initial data */
180 enum{
190 };
191 // static const char* LALTEOBResumSVariables[] = {"r","phi","Pphi","Prstar","Pr","j","E0","Omega"};
192 
193 /** Index list of EOB dynamical variables (to be stored in arrays) */
194 enum{
204 };
205 
206 /** Index list of TEOB mode-list options */
207 enum{
211 };
212 
213 // static const char* LALTEOBResumSDynamicalVariables[] = {"r","phi","Pphi","MOmega","ddor","Prstar","MOmega_orb","E"};
214 
215 #define KMAX (35) /** Multipolar linear index, max value */
216 #define PMTERMS_eps (1) /** Switch on Fujita-Iyer point-mass terms. This is hard-coded here */
217 
218 #define OUTPUT_LM_LEN (35)
219 
220 
221 /** List of options for centrifugal radius */
222 enum{
229  RC_NOPT
230 };
231 
232 enum{
235  SS_NOPT
236 };
237 
238 /** List of options for ODE timestepping */
239 enum{
244 };
245 static const char* const LALTEOBResumSODETimeStep[] = {"uniform","adaptive","adaptive+uniform_after_LSO","undefined"};
246 
247 /** List of options for postadiabatic inspiral */
248 enum{
253 };
254 
255 /** List of options for interp_uniform_grid */
256 enum{
261 };
262 
263 /** List of options for NQC waveform nqc_coeffs_hlm and nqc_coeffs_flx */
264 enum{
269  NQC_NOPT
270 };
271 
272 /** Error handler for root finders */
273 enum{
279 };
280 static const char* const LALTEOBResumSRootErrors[] = {"none","root is not bracketed.","root finder did not converged.","root finder failed."};
281 #define ROOTFINDER(i, x) {if ( ((i) = (x)) && ((i)>ROOT_ERRORS_NO) ) { exit(LALTEOBResumSRootErrors[(i)]); }} //TODO: CHECK THIS MACRO (LOGIC INVOLVED)
282 
283 /** Maps between linear index and the corresponding (l, m) multipole indices */
284 extern const INT4 TEOB_LINDEX[KMAX]; /* defined in LALSimIMRTEOBResumS.c */
285 extern const INT4 TEOB_MINDEX[KMAX]; /* defined in LALSimIMRTEOBResumS.c */
286 
287 
288 /** Multipolar waveform at given time, comes at handy */
289 typedef struct tagLALTEOBResumSWaveformModeSingleTime
290 {
292  REAL8 ampli[KMAX]; /* amplitude */
293  REAL8 phase[KMAX]; /* phase */
294  /* kmask is currently unused in the LAL version */
295  // INT4 kmask[KMAX]; /* mask for multipoles */
297 
298 /** Func pointer types */
299 typedef void (*EOBWavFlmSFunc)(REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8 [KMAX][6], int, REAL8 *, REAL8 *);
301 
302 /** Multipolar coefficients for NQC waveform */
303 typedef struct tagNQCcoefs
304 {
305  double a1[KMAX];
306  double a2[KMAX];
307  double a3[KMAX];
308  double b1[KMAX];
309  double b2[KMAX];
310  double b3[KMAX];
311  double n[KMAX][6];
312  int activemode[KMAX]; /* mask for modes with nonzero NQC */
313  int maxk; /* index of maximum active mode */
314  int add; /* flag */
315 } NQCcoefs;
316 
317 /** NQC data for flux and waveform */
318 typedef struct tagNQCdata
319 {
322 } NQCdata;
323 
324 /** Dynamics data type */
325 typedef struct tagLALTEOBResumSDynamics
326 {
328  /* various pointwise variables */
329  INT4 store; /* store following values? */
330  INT4 noflux; /* compute rhs without flux */
331  REAL8 t, r, phi, pphi, prstar, ddotr, Omg, Omg_orb;
332  REAL8 H, Heff, Heff_orb, E, jhat, r_omega, psi, v_phi;
333  REAL8 A, dA, d2A, B, dB;
334  REAL8 MOmg, MOmg_prev, tMOmgpeak;
335  /* stuff for ODE solver */
336  REAL8 y[TEOB_EVOLVE_NVARS]; /* rhs storage */
338  REAL8 y0[TEOB_ID_NVARS]; /* ID storage */
339  REAL8 dt, t_stop, ti;
342  bool ode_stop, ode_stop_MOmgpeak, ode_stop_radius;
343  INT4 nqc_coeffs_hlm, nqc_coeffs_flx;
345  REAL8 clm[KMAX][6];
346  /* arrays */
350  /* func pointers */
353  /* key parameters for quick access */
354  REAL8 M, nu, q, X1, X2;
355  REAL8 chi1, chi2, S1, S2, S, Sstar, a1, a2, aK2, C_Q1, C_Q2, C_Oct1, C_Oct2, C_Hex1, C_Hex2, cN3LO;
356  REAL8 rLR, rLSO;
357  REAL8 kapA2,kapA3,kapA4, kapB2,kapB3,kapB4, kapT2,kapT3,kapT4, khatA2,khatB2;
358  REAL8 bar_alph2_1, bar_alph2_2, bar_alph3_1, bar_alph3_2, bar_alph2j_1;
359  REAL8 kapA2j, kapB2j, kapT2j;
360  REAL8 rLR_tidal, pGSF_tidal;
361  REAL8 Mbhf, abhf; /* final BH */
362  INT4 use_tidal, use_spins, use_tidal_gravitomagnetic;
363  INT4 bhns; /* 0:BNS, 1:BHNS, 2:NSBH */
365 
366 /* Function protoypes grouped based on file */
367 
368 /** Should not be needed in LAL (all parameters passed to interface function) */
369 /* TEOBResumSPars.c */
370 //void par_db_init (void);
371 //void par_db_free (void);
372 //void par_db_default (void);
373 //void par_file_parse (const char *fname);
374 //void par_file_parse_merge (const char *fname);
375 //void par_db_write_file (const char *fname);
376 //void par_db_screen (void);
377 //void par_set_i(const char *key, int val);
378 //void par_set_b(const char *key, int val);
379 //void par_set_d(const char *key, double val);
380 //void par_set_s(const char *key, const char *val);
381 //void par_set_arrayi (const char *key, int *array, int n);
382 //void par_set_arrayd (const char *key, double *array, int n);
383 //int par_get_i(const char *key);
384 //int par_get_b(const char *key);
385 //double par_get_d(const char *key);
386 //const char * par_get_s(const char *key);
387 //int * par_get_arrayi(const char *key, int *n);
388 //double * par_get_arrayd(const char *key, int *n);
389 //void eob_set_params(char *s, int n);
390 //void eob_free_params(void);
391 
392 INT4 XLALSetup_TEOB_mode_array(LALValue *ModeArray, INT4 modeType);
393 INT4 XLALCheck_TEOB_mode_array(LALValue *ModeArray, UINT4 use_tidal);
394 
395 /* TEOBResumSUtil.c */
396 REAL8 q_to_nu(const REAL8 q);
397 REAL8 nu_to_X1(const REAL8 nu);
398 REAL8 Eulerlog(const REAL8 x,const INT4 m);
399 void interp_spline(REAL8 *t, REAL8 *y, INT4 n, REAL8 *ti, INT4 ni, REAL8 *yi);
400 int find_point_bisection(REAL8 x, INT4 n, REAL8 *xp, INT4 o);
401 REAL8 baryc_f(REAL8 xx, INT4 n, REAL8 *f, REAL8 *x);
402 void baryc_weights(INT4 n, REAL8 *x, REAL8 *omega);
403 REAL8 baryc_f_weights(REAL8 xx, INT4 n, REAL8 *f, REAL8 *x, REAL8 *omega);
404 REAL8 interp1d (const INT4 order, REAL8 xx, INT4 nx, REAL8 *f, REAL8 *x);
405 REAL8 find_max (const INT4 n, REAL8 dx, REAL8 x0, REAL8 *f, REAL8 *fmax);
406 //double fact(int n);
407 /** XLALWignerdMatrix, XLALWignerDMatrix in <lal/SphericalHarmonics.h> */
408 // double wigner_d_function(int l, int m, int s, double i);
409 //int spinsphericalharm(double *rY, double *iY, int s, int l, int m, double phi, double i);
410 
411 //void compute_hpc(SphHarmPolarTimeSeries *hlm, double nu, double M, double distance, double amplitude_prefactor, double psi, double iota, COMPLEX16TimeSeries *hpc);
412 
413 int D0(REAL8 *f, REAL8 dx, INT4 n, REAL8 *df);
414 int D2(REAL8 *f, REAL8 dx, INT4 n, REAL8 *d2f);
415 //int D0_x(double *f, double *x, int n, double *df);
416 
417 //double cumtrapz(double *f, double *x, const int n, double *sum);
418 REAL8 cumint3(REAL8 *f, REAL8 *x, const INT4 n, REAL8 *sum);
419 
420 void unwrap(REAL8 *p, const INT4 size);
421 void unwrap_proxy(REAL8 *p, REAL8 *r, const INT4 size, const INT4 shift0);
422 
423 int get_uniform_size(const REAL8 tf, const REAL8 t0, const REAL8 dt);
424 
425 /** Replaced by LAL COMPLEX16TimeSeries functions */
426 //void Waveform_alloc (Waveform **wav, int size, const char *name);
427 //void Waveform_push (Waveform **wav, int size);
428 //void Waveform_output (Waveform *wav);
429 //void Waveform_free (Waveform *wav);
430 //void Waveform_lm_alloc (LALTEOBResumSWaveformMode **wav, int size, const char *name);
431 //void Waveform_lm_push (LALTEOBResumSWaveformMode **wav, int size);
432 //void Waveform_lm_output (LALTEOBResumSWaveformMode *wav);
433 //void Waveform_lm_output_reim (LALTEOBResumSWaveformMode *wav);
434 //void Waveform_lm_free (LALTEOBResumSWaveformMode *wav);
437 
438 void XLALTEOBDynamicsInit (LALTEOBResumSDynamics **dyn, INT4 size, const CHAR *name);
442  LALValue **pModeArray,
443  const REAL8 m1,
444  const REAL8 m2,
445  const REAL8 S1x,
446  const REAL8 S1y,
447  const REAL8 S1z,
448  const REAL8 S2x,
449  const REAL8 S2y,
450  const REAL8 S2z,
451  const REAL8 deltaT,
452  const REAL8 lambda1,
453  const REAL8 lambda2,
454  const REAL8 lambda1oct,
455  const REAL8 lambda2oct,
456  const REAL8 lambda1hex,
457  const REAL8 lambda2hex,
458  LALDict *LALparams);
459 
460 void NQCdata_alloc (NQCdata **nqc);
461 void NQCdata_free (NQCdata *nqc);
462 
465 REAL8 radius0(REAL8 M, REAL8 fHz);
466 REAL8 get_mrg_timestep(REAL8 q, REAL8 chi1, REAL8 chi2);
467 
468 void XLALSimIMRComputePolarisations(REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi);
469 
470 void XLALSimIMRComputePolarisationsPolar(REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmPolarTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi);
471 
472 /* TEOBResumSFits.c */
473 //double eob_c3_fit_global(double nu, double chi1, double chi2, double X1, double X2, double a1, double a2);
474 //double eob_nqc_dtfit(const double chi, const double chi0);
476 REAL8 logQ(REAL8 x);
483 void HealyBBHFitRemnant(REAL8 chi1, REAL8 chi2, REAL8 q, REAL8 *mass, REAL8 *spin);
484 REAL8 JimenezFortezaRemnantSpin(REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2);
485 
486 void QNMHybridFitCab(REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 aK, REAL8 Mbh, REAL8 abh, REAL8 *a1, REAL8 *a2, REAL8 *a3, REAL8 *a4, REAL8 *b1, REAL8 *b2, REAL8 *b3, REAL8 *b4, REAL8 *sigmar, REAL8 *sigmai, INT4 usespins);
487 //double eob_approxLR(const double nu);
488 
489 /* TEOBResumSDynamics.c */
490 int eob_dyn_rhs(REAL8 t, const REAL8 y[], REAL8 dy[], void *params);
491 void eob_ham(REAL8 nu, REAL8 r, REAL8 pph, REAL8 prstar, REAL8 A, REAL8 dA,
492  REAL8 *H, REAL8 *Heff, REAL8 *dHeff_dr, REAL8 *dHeff_dprstar, REAL8 *dHeff_dpphi);
493 int eob_dyn_rhs_s(REAL8 t, const REAL8 y[], REAL8 dy[], void *params);
494 void eob_ham_s(REAL8 nu, REAL8 r, REAL8 rc, REAL8 drc_dr, REAL8 pphi, REAL8 prstar, REAL8 S, REAL8 Sstar, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 aK2, REAL8 c3, REAL8 A, REAL8 dA, REAL8 *H, REAL8 *Heff, REAL8 *Heff_orb, REAL8 *dHeff_dr, REAL8 *dHeff_dprstar, REAL8 *dHeff_dpphi, REAL8 *d2Heff_dprstar20);
495 void eob_dyn_s_GS(REAL8 r, REAL8 rc, REAL8 drc_dr, REAL8 aK2, REAL8 prstar, REAL8 pph, REAL8 nu, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 cN3LO, REAL8 *ggm);
496 //void eob_dyn_s_get_rc(REAL8 r, REAL8 nu, REAL8 at1,REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
497 
498 void eob_dyn_s_get_rc_LO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
499 void eob_dyn_s_get_rc_NLO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
500 void eob_dyn_s_get_rc_NNLO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
501 void eob_dyn_s_get_rc_NNLO_S4(REAL8 r, REAL8 nu, REAL8 at1,REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 C_Oct1, REAL8 C_Oct2, REAL8 C_Hex1, REAL8 C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
502 void eob_dyn_s_get_rc_NOSPIN(REAL8 r, REAL8 UNUSED nu, REAL8 UNUSED at1, REAL8 UNUSED at2, REAL8 UNUSED aK2, REAL8 UNUSED C_Q1, REAL8 UNUSED C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 UNUSED usetidal, REAL8 UNUSED *rc, REAL8 UNUSED *drc_dr, REAL8 UNUSED *d2rc_dr2);
503 void eob_dyn_s_get_rc_NOTIDES(REAL8 r, REAL8 nu, REAL8 at1,REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2);
504 REAL8 eob_dyn_fLR(REAL8 r, void * params);
505 int eob_dyn_adiabLR(LALTEOBResumSDynamics *dyn, REAL8 *rLR, INT4 tidesFlag);
506 //double eob_dyn_fLSO(double r, void * params);
507 //int eob_dyn_adiabLSO(LALTEOBResumSDynamics *dyn, REAL8 *rLSO);
508 
509 /* TEOBResumSPostAdiabatic.c */
511 
512 /* TEOBResumSInitialCondition.c */
513 void eob_dyn_ic(REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[]);
514 void eob_dyn_ic_s(REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[]);
515 REAL8 eob_dyn_bisecHeff0_s(REAL8 nu, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 c3, REAL8 pph, REAL8 rorb, REAL8 A, REAL8 dA, REAL8 rc, REAL8 drc_dr, REAL8 ak2, REAL8 S, REAL8 Ss);
516 REAL8 eob_dyn_DHeff0(REAL8 x, void *params);
517 
518 /* TEOBResumSMetric.c */
519 void eob_metric_A5PNlog(REAL8 r, REAL8 nu, REAL8 *A, REAL8 *dA, REAL8 *d2A);
520 void eob_metric_Atidal(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *AT, REAL8 *dAT, REAL8 *d2AT);
521 void eob_metric(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB);
522 void eob_metric_s(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB);
523 
524 /* TEOBResumSFlux.c */
525 REAL8 eob_flx_Flux(REAL8 x, REAL8 Omega, REAL8 r_omega, REAL8 E, REAL8 Heff, REAL8 jhat, REAL8 r, REAL8 pr_star, REAL8 ddotr, LALTEOBResumSDynamics *dyn);
526 REAL8 eob_flx_Flux_s(REAL8 x, REAL8 Omega, REAL8 r_omega, REAL8 E, REAL8 Heff, REAL8 jhat, REAL8 r, REAL8 pr_star, REAL8 ddotr, LALTEOBResumSDynamics *dyn);
527 void eob_flx_Tlm(REAL8 w, REAL8 *MTlm);
528 void eob_flx_FlmNewt(REAL8 x, REAL8 nu, REAL8 *Nlm);
529 REAL8 eob_flx_HorizonFlux(REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu);
530 REAL8 eob_flx_HorizonFlux_s(REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2);
531 
532 /* TEOBResumSWaveform.c */
534 
535 void eob_wav_deltalm(REAL8 Hreal,REAL8 Omega,REAL8 nu, REAL8 *dlm);
539 void eob_wav_hlmTidal(REAL8 x, LALTEOBResumSDynamics *dyn, REAL8 *hTidallm);
540 void eob_wav_flm_coeffs(REAL8 nu, REAL8 clm[KMAX][6]);
541 void eob_wav_flm(REAL8 x,REAL8 nu, REAL8 clm[KMAX][6], REAL8 *rholm, REAL8 *flm);
542 //void eob_wav_flm_old(double x,double nu, double *rholm, double *flm);
543 
544 void eob_wav_flm_s_SSLO(REAL8 x, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 a1, REAL8 a2, REAL8 C_Q1, REAL8 C_Q2, REAL8 clm[KMAX][6], int usetidal, REAL8 *rholm, REAL8 *flm);
545 void eob_wav_flm_s_SSNLO(REAL8 x, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 a1, REAL8 a2, REAL8 C_Q1, REAL8 C_Q2, REAL8 clm[KMAX][6], int usetidal, REAL8 *rholm, REAL8 *flm);
546 //void eob_wav_flm_s_old(double x, double nu, double X1, double X2, double chi1, double chi2, double a1, double a2, double C_Q1, double C_Q2, int usetidal, double *rholm, double *flm);
549 void eob_wav_hlmNQC(REAL8 nu, REAL8 r, REAL8 prstar, REAL8 Omega, REAL8 ddotr, NQCcoefs *nqc, LALTEOBResumSWaveformModeSingleTime *psilmnqc);
550 void eob_wav_ringdown_template(REAL8 x, REAL8 a1, REAL8 a2, REAL8 a3, REAL8 a4, REAL8 b1, REAL8 b2, REAL8 b3, REAL8 b4, REAL8 sigmar, REAL8 sigmai, REAL8 *psi);
552 
557 void unwrap(REAL8 *p, const INT4 size);
558 INT4 get_uniform_size(const REAL8 tN, const REAL8 t0, const REAL8 dt);
559 void XLALTEOBDynamicsInterp (LALTEOBResumSDynamics *dyn, const INT4 size, const REAL8 t0, const REAL8 dt, const char *name);
560 void XLALTEOBDynamicsExtract (LALTEOBResumSDynamics *dyna, const REAL8 to, const REAL8 tn, LALTEOBResumSDynamics **dynb, const char *name);
563  REAL8 r,
564  REAL8 prstar,
565  REAL8 Omega,
566  REAL8 ddotr,
571 
573 
574 
575 #if 0
576 { /* so that editors will match succeeding brace */
577 #elif defined(__cplusplus)
578 }
579 #endif
580 
581 #endif // of #ifndef _LALSIM_TEOB_RESUMS_H
const double b1
const double b2
const char * name
@ SS_NLO
@ SS_NOPT
@ SS_LO
void baryc_weights(INT4 n, REAL8 *x, REAL8 *omega)
void eob_wav_hlm(LALTEOBResumSDynamics *dyn, LALTEOBResumSWaveformModeSingleTime *hlm_t)
void eob_dyn_s_GS(REAL8 r, REAL8 rc, REAL8 drc_dr, REAL8 aK2, REAL8 prstar, REAL8 pph, REAL8 nu, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 cN3LO, REAL8 *ggm)
void eob_ham(REAL8 nu, REAL8 r, REAL8 pph, REAL8 prstar, REAL8 A, REAL8 dA, REAL8 *H, REAL8 *Heff, REAL8 *dHeff_dr, REAL8 *dHeff_dprstar, REAL8 *dHeff_dpphi)
void eob_dyn_s_get_rc_NOSPIN(REAL8 r, REAL8 UNUSED nu, REAL8 UNUSED at1, REAL8 UNUSED at2, REAL8 UNUSED aK2, REAL8 UNUSED C_Q1, REAL8 UNUSED C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 UNUSED usetidal, REAL8 UNUSED *rc, REAL8 UNUSED *drc_dr, REAL8 UNUSED *d2rc_dr2)
void eob_dyn_s_get_rc_NNLO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2)
void eob_nqc_setcoefs(LALTEOBResumSDynamics *dyn)
REAL8 cumint3(REAL8 *f, REAL8 *x, const INT4 n, REAL8 *sum)
void XLALTEOBDynamicsPush(LALTEOBResumSDynamics **dyn, INT4 size)
void Waveform_lm_t_free(LALTEOBResumSWaveformModeSingleTime *wav)
static const char *const LALTEOBResumSODETimeStep[]
void eob_dyn_s_get_rc_LO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2)
void unwrap(REAL8 *p, const INT4 size)
void XLALFreeTEOBDynamics(LALTEOBResumSDynamics *dyn)
void eob_dyn_s_get_rc_NLO(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2)
REAL8 eob_flx_Flux_s(REAL8 x, REAL8 Omega, REAL8 r_omega, REAL8 E, REAL8 Heff, REAL8 jhat, REAL8 r, REAL8 pr_star, REAL8 ddotr, LALTEOBResumSDynamics *dyn)
static const char *const LALTEOBResumSRootErrors[]
void eob_wav_flm_s_SSNLO(REAL8 x, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 a1, REAL8 a2, REAL8 C_Q1, REAL8 C_Q2, REAL8 clm[KMAX][6], int usetidal, REAL8 *rholm, REAL8 *flm)
void eob_wav_flm_s_SSLO(REAL8 x, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 a1, REAL8 a2, REAL8 C_Q1, REAL8 C_Q2, REAL8 clm[KMAX][6], int usetidal, REAL8 *rholm, REAL8 *flm)
REAL8 Yagi13_fit_barsigmalambda(REAL8 barlam2)
@ INTERP_UNIFORM_GRID_OFF
@ INTERP_UNIFORM_GRID_HPC
@ INTERP_UNIFORM_GRID_HLM
@ INTERP_UNIFORM_GRID_NOPT
void eob_metric(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB)
REAL8 eob_flx_HorizonFlux_s(REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2)
REAL8 time_units_factor(REAL8 M)
void XLALSimIMRComputePolarisationsPolar(REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmPolarTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi)
REAL8 eob_dyn_bisecOmegaorb0(LALTEOBResumSDynamics *dyn, REAL8 omg_orb0, REAL8 r0_kepl)
void eob_dyn_ic_s(REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[])
void eob_wav_flm(REAL8 x, REAL8 nu, REAL8 clm[KMAX][6], REAL8 *rholm, REAL8 *flm)
void XLALTEOBDynamicsInterp(LALTEOBResumSDynamics *dyn, const INT4 size, const REAL8 t0, const REAL8 dt, const char *name)
void eob_dyn_s_get_rc_NNLO_S4(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 C_Oct1, REAL8 C_Oct2, REAL8 C_Hex1, REAL8 C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2)
void XLALSimIMRComputePolarisations(REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi)
const INT4 TEOB_MINDEX[KMAX]
int eob_dyn_Npostadiabatic(LALTEOBResumSDynamics *dyn, REAL8 r0)
void HealyBBHFitRemnant(REAL8 chi1, REAL8 chi2, REAL8 q, REAL8 *mass, REAL8 *spin)
REAL8 JimenezFortezaRemnantSpin(REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2)
void eob_wav_deltalm(REAL8 Hreal, REAL8 Omega, REAL8 nu, REAL8 *dlm)
int eob_dyn_adiabLR(LALTEOBResumSDynamics *dyn, REAL8 *rLR, INT4 tidesFlag)
@ TEOB_MODES_22
@ TEOB_MODES_ALL
@ TEOB_MODES_NOPT
void eob_ham_s(REAL8 nu, REAL8 r, REAL8 rc, REAL8 drc_dr, REAL8 pphi, REAL8 prstar, REAL8 S, REAL8 Sstar, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 aK2, REAL8 c3, REAL8 A, REAL8 dA, REAL8 *H, REAL8 *Heff, REAL8 *Heff_orb, REAL8 *dHeff_dr, REAL8 *dHeff_dprstar, REAL8 *dHeff_dpphi, REAL8 *d2Heff_dprstar20)
INT4 XLALSetup_TEOB_mode_array(LALValue *ModeArray, INT4 modeType)
Should not be needed in LAL (all parameters passed to interface function)
void XLALTEOBDynamicsSetParams(LALTEOBResumSDynamics *dyn, LALValue **pModeArray, const REAL8 m1, const REAL8 m2, const REAL8 S1x, const REAL8 S1y, const REAL8 S1z, const REAL8 S2x, const REAL8 S2y, const REAL8 S2z, const REAL8 deltaT, const REAL8 lambda1, const REAL8 lambda2, const REAL8 lambda1oct, const REAL8 lambda2oct, const REAL8 lambda1hex, const REAL8 lambda2hex, LALDict *LALparams)
void eob_wav_ringdown_template(REAL8 x, REAL8 a1, REAL8 a2, REAL8 a3, REAL8 a4, REAL8 b1, REAL8 b2, REAL8 b3, REAL8 b4, REAL8 sigmar, REAL8 sigmai, REAL8 *psi)
void XLALTEOBDynamicsExtract(LALTEOBResumSDynamics *dyna, const REAL8 to, const REAL8 tn, LALTEOBResumSDynamics **dynb, const char *name)
void eob_wav_speedyTail(REAL8 Omega, REAL8 Hreal, REAL8 bphys, LALTEOBResumSWaveformModeSingleTime *tlm)
REAL8 baryc_f_weights(REAL8 xx, INT4 n, REAL8 *f, REAL8 *x, REAL8 *omega)
#define TEOB_STRLEN
@ RC_NOPT
@ RC_NLO
@ RC_LO
@ RC_NNLOS4
@ RC_NNLO
@ RC_NOSPIN
@ RC_NOTIDES
REAL8 eob_flx_Flux(REAL8 x, REAL8 Omega, REAL8 r_omega, REAL8 E, REAL8 Heff, REAL8 jhat, REAL8 r, REAL8 pr_star, REAL8 ddotr, LALTEOBResumSDynamics *dyn)
void eob_flx_FlmNewt(REAL8 x, REAL8 nu, REAL8 *Nlm)
REAL8 get_mrg_timestep(REAL8 q, REAL8 chi1, REAL8 chi2)
void eob_metric_Atidal(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *AT, REAL8 *dAT, REAL8 *d2AT)
void NQCdata_free(NQCdata *nqc)
void XLALTEOBDynamicsJoin(LALTEOBResumSDynamics *dyna, LALTEOBResumSDynamics *dynb, REAL8 to)
REAL8 q_to_nu(const REAL8 q)
#define KMAX
void eob_wav_hlmNQC_nospin201602(REAL8 nu, REAL8 r, REAL8 prstar, REAL8 Omega, REAL8 ddotr, LALTEOBResumSWaveformModeSingleTime *hlmnqc)
REAL8 baryc_f(REAL8 xx, INT4 n, REAL8 *f, REAL8 *x)
void eob_wav_flm_coeffs(REAL8 nu, REAL8 clm[KMAX][6])
REAL8 Eulerlog(const REAL8 x, const INT4 m)
@ TEOB_EVOLVE_RAD
@ TEOB_EVOLVE_NVARS
@ TEOB_EVOLVE_PHI
@ TEOB_EVOLVE_PPHI
@ TEOB_EVOLVE_PRSTAR
int D2(REAL8 *f, REAL8 dx, INT4 n, REAL8 *d2f)
void eob_wav_hlmNQC_find_a1a2a3_mrg(LALTEOBResumSDynamics *dyn_mrg, SphHarmPolarTimeSeries *hlm_mrg, SphHarmPolarTimeSeries *hnqc, LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *hlm)
REAL8 eob_dyn_bisecHeff0_s(REAL8 nu, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 c3, REAL8 pph, REAL8 rorb, REAL8 A, REAL8 dA, REAL8 rc, REAL8 drc_dr, REAL8 ak2, REAL8 S, REAL8 Ss)
REAL8 JFAPG_fit_Sigma_Irrotational(REAL8 barlam2)
REAL8 eob_c3_fit_global(REAL8 nu, REAL8 chi1, REAL8 chi2, REAL8 X1, REAL8 X2, REAL8 a1, REAL8 a2)
void XLALSphHarmPolarJoin(SphHarmPolarTimeSeries *hlma, SphHarmPolarTimeSeries *hlmb, REAL8 to)
void eob_metric_s(REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB)
REAL8 eob_dyn_r0_eob(REAL8 f0, LALTEOBResumSDynamics *dyn)
REAL8 nu_to_X1(const REAL8 nu)
int get_uniform_size(const REAL8 tf, const REAL8 t0, const REAL8 dt)
void(* EOBDynSGetRCFunc)(REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, INT4, REAL8 *, REAL8 *, REAL8 *)
@ TEOB_PA_OFF
@ TEOB_PA_ONLY
@ TEOB_PA_NOPT
@ TEOB_PA_ON
void eob_nqc_setcoefs_nospin201602(REAL8 nu, NQCcoefs *nqc)
REAL8 Yagi14_fit_Chex(REAL8 C_Q)
REAL8 Yagi13_fit_barlamdel(REAL8 barlam2, INT4 ell)
void eob_wav_hlmNewt(REAL8 r, REAL8 Omega, REAL8 phi, REAL8 nu, LALTEOBResumSWaveformModeSingleTime *hNewt)
void QNMHybridFitCab(REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2, REAL8 aK, REAL8 Mbh, REAL8 abh, REAL8 *a1, REAL8 *a2, REAL8 *a3, REAL8 *a4, REAL8 *b1, REAL8 *b2, REAL8 *b3, REAL8 *b4, REAL8 *sigmar, REAL8 *sigmai, INT4 usespins)
void eob_wav_ringdown(LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *hlm)
void eob_wav_hlmNQC(REAL8 nu, REAL8 r, REAL8 prstar, REAL8 Omega, REAL8 ddotr, NQCcoefs *nqc, LALTEOBResumSWaveformModeSingleTime *psilmnqc)
@ NQC_NOPT
@ NQC_NR_NOSPIN
@ NQC_FILE
@ NQC_OFF
@ NQC_COMPUTE
@ TEOB_ID_PRSTAR
@ TEOB_ID_PR
@ TEOB_ID_E0
@ TEOB_ID_PHI
@ TEOB_ID_RAD
@ TEOB_ID_PPHI
@ TEOB_ID_NVARS
@ TEOB_ID_OMGJ
@ TEOB_ID_J
@ ODE_TSTEP_UNIFORM
@ ODE_TSTEP_NOPT
@ ODE_TSTEP_ADAPTIVE_UNIFORM_AFTER_LSO
@ ODE_TSTEP_ADAPTIVE
@ ROOT_ERRORS_NOSUCC
@ ROOT_ERRORS_MAXITS
@ ROOT_ERRORS_NO
@ ROOT_ERRORS
@ ROOT_ERRORS_BRACKET
void unwrap_proxy(REAL8 *p, REAL8 *r, const INT4 size, const INT4 shift0)
REAL8 eob_nqc_timeshift(REAL8 nu, REAL8 chi1)
void eob_wav_hhatlmTail(REAL8 Omega, REAL8 Hreal, REAL8 bphys, LALTEOBResumSWaveformModeSingleTime *tlm)
REAL8 eob_flx_HorizonFlux(REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu)
void eob_metric_A5PNlog(REAL8 r, REAL8 nu, REAL8 *A, REAL8 *dA, REAL8 *d2A)
int eob_dyn_rhs(REAL8 t, const REAL8 y[], REAL8 dy[], void *params)
REAL8 logQ(REAL8 x)
REAL8 interp1d(const INT4 order, REAL8 xx, INT4 nx, REAL8 *f, REAL8 *x)
void Waveform_lm_t_alloc(LALTEOBResumSWaveformModeSingleTime **wav)
Replaced by LAL COMPLEX16TimeSeries functions.
REAL8 time_units_conversion(REAL8 M, REAL8 t)
int D0(REAL8 *f, REAL8 dx, INT4 n, REAL8 *df)
XLALWignerdMatrix, XLALWignerDMatrix in <lal/SphericalHarmonics.h>
REAL8 radius0(REAL8 M, REAL8 fHz)
INT4 XLALCheck_TEOB_mode_array(LALValue *ModeArray, UINT4 use_tidal)
void NQCdata_alloc(NQCdata **nqc)
@ TEOB_PRSTAR
@ TEOB_OMGORB
@ TEOB_DYNAMICS_NVARS
@ TEOB_DDOTR
@ TEOB_E0
@ TEOB_RAD
@ TEOB_MOMG
@ TEOB_PHI
@ TEOB_PPHI
int eob_dyn_rhs_s(REAL8 t, const REAL8 y[], REAL8 dy[], void *params)
void eob_wav_hlmNQC_find_a1a2a3(LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *h, SphHarmPolarTimeSeries *hnqc)
const INT4 TEOB_LINDEX[KMAX]
Maps between linear index and the corresponding (l, m) multipole indices.
REAL8 eob_dyn_Omegaorb0(REAL8 r, void *params)
void XLALTEOBDynamicsInit(LALTEOBResumSDynamics **dyn, INT4 size, const CHAR *name)
REAL8 Yagi14_fit_Coct(REAL8 C_Q)
void eob_dyn_s_get_rc_NOTIDES(REAL8 r, REAL8 nu, REAL8 at1, REAL8 at2, REAL8 aK2, REAL8 C_Q1, REAL8 C_Q2, REAL8 UNUSED C_Oct1, REAL8 UNUSED C_Oct2, REAL8 UNUSED C_Hex1, REAL8 UNUSED C_Hex2, INT4 usetidal, REAL8 *rc, REAL8 *drc_dr, REAL8 *d2rc_dr2)
REAL8 eob_dyn_r0_Kepler(REAL8 f0)
void(* EOBWavFlmSFunc)(REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8[KMAX][6], int, REAL8 *, REAL8 *)
Func pointer types.
REAL8 eob_dyn_fLR(REAL8 r, void *params)
void eob_flx_Tlm(REAL8 w, REAL8 *MTlm)
void interp_spline(REAL8 *t, REAL8 *y, INT4 n, REAL8 *ti, INT4 ni, REAL8 *yi)
void eob_wav_hlmTidal(REAL8 x, LALTEOBResumSDynamics *dyn, REAL8 *hTidallm)
REAL8 JFAPG_fit_Sigma_Static(REAL8 barlam2)
REAL8 eob_dyn_DHeff0(REAL8 x, void *params)
void eob_dyn_ic(REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[])
int find_point_bisection(REAL8 x, INT4 n, REAL8 *xp, INT4 o)
REAL8 find_max(const INT4 n, REAL8 dx, REAL8 x0, REAL8 *f, REAL8 *fmax)
static REAL8 sum(REAL8 *arr1, REAL8 *arr2, int k)
REAL8 Hreal
REAL8 M
Definition: bh_qnmode.c:133
double dt
Definition: bh_ringdown.c:113
double theta
Definition: bh_sphwf.c:118
const double H
const double a4
const double a2
const double w
sigmaKerr data[0]
const double B
const double nx
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
static const INT4 r
static const INT4 m
static const INT4 q
list y
Dynamics data type.
EOBDynSGetRCFunc eob_dyn_s_get_rc
EOBWavFlmSFunc eob_wav_flm_s
Multipolar waveform at given time, comes at handy.
Multipolar coefficients for NQC waveform.
NQC data for flux and waveform.
NQCcoefs * hlm
NQCcoefs * flx
Structure to carry a collection of spherical harmonic modes in COMPLEX16 time series.
Definition: burst.c:245
double deltaT
Definition: unicorn.c:24