LALSimulation  5.4.0.1-fe68b98
LALSimTEOBResumS.h File Reference

This file is part of TEOBResumS. More...

Prototypes

INT4 XLALSetup_TEOB_mode_array (LALValue *ModeArray, INT4 modeType)
 Should not be needed in LAL (all parameters passed to interface function) More...
 
INT4 XLALCheck_TEOB_mode_array (LALValue *ModeArray, UINT4 use_tidal)
 
REAL8 q_to_nu (const REAL8 q)
 
REAL8 nu_to_X1 (const REAL8 nu)
 
REAL8 Eulerlog (const REAL8 x, const INT4 m)
 
void interp_spline (REAL8 *t, REAL8 *y, INT4 n, REAL8 *ti, INT4 ni, REAL8 *yi)
 
int find_point_bisection (REAL8 x, INT4 n, REAL8 *xp, INT4 o)
 
REAL8 baryc_f (REAL8 xx, INT4 n, REAL8 *f, REAL8 *x)
 
void baryc_weights (INT4 n, REAL8 *x, REAL8 *omega)
 
REAL8 baryc_f_weights (REAL8 xx, INT4 n, REAL8 *f, REAL8 *x, REAL8 *omega)
 
REAL8 interp1d (const INT4 order, REAL8 xx, INT4 nx, REAL8 *f, REAL8 *x)
 
REAL8 find_max (const INT4 n, REAL8 dx, REAL8 x0, REAL8 *f, REAL8 *fmax)
 
int D0 (REAL8 *f, REAL8 dx, INT4 n, REAL8 *df)
 XLALWignerdMatrix, XLALWignerDMatrix in <lal/SphericalHarmonics.h> More...
 
int D2 (REAL8 *f, REAL8 dx, INT4 n, REAL8 *d2f)
 
REAL8 cumint3 (REAL8 *f, REAL8 *x, const INT4 n, REAL8 *sum)
 
void unwrap (REAL8 *p, const INT4 size)
 
void unwrap_proxy (REAL8 *p, REAL8 *r, const INT4 size, const INT4 shift0)
 
int get_uniform_size (const REAL8 tf, const REAL8 t0, const REAL8 dt)
 
void Waveform_lm_t_alloc (LALTEOBResumSWaveformModeSingleTime **wav)
 Replaced by LAL COMPLEX16TimeSeries functions. More...
 
void Waveform_lm_t_free (LALTEOBResumSWaveformModeSingleTime *wav)
 
void XLALTEOBDynamicsInit (LALTEOBResumSDynamics **dyn, INT4 size, const CHAR *name)
 
void XLALTEOBDynamicsPush (LALTEOBResumSDynamics **dyn, INT4 size)
 
void XLALFreeTEOBDynamics (LALTEOBResumSDynamics *dyn)
 
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 NQCdata_alloc (NQCdata **nqc)
 
void NQCdata_free (NQCdata *nqc)
 
REAL8 time_units_factor (REAL8 M)
 
REAL8 time_units_conversion (REAL8 M, REAL8 t)
 
REAL8 radius0 (REAL8 M, REAL8 fHz)
 
REAL8 get_mrg_timestep (REAL8 q, REAL8 chi1, REAL8 chi2)
 
void XLALSimIMRComputePolarisations (REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi)
 
void XLALSimIMRComputePolarisationsPolar (REAL8Sequence *hplus_out, REAL8Sequence *hcross_out, SphHarmPolarTimeSeries *hlm, LALValue *modeArray, REAL8 amplitude_prefactor, REAL8 theta, REAL8 phi)
 
REAL8 eob_nqc_timeshift (REAL8 nu, REAL8 chi1)
 
REAL8 logQ (REAL8 x)
 
REAL8 Yagi13_fit_barlamdel (REAL8 barlam2, INT4 ell)
 
REAL8 Yagi13_fit_barsigmalambda (REAL8 barlam2)
 
REAL8 Yagi14_fit_Coct (REAL8 C_Q)
 
REAL8 Yagi14_fit_Chex (REAL8 C_Q)
 
REAL8 JFAPG_fit_Sigma_Static (REAL8 barlam2)
 
REAL8 JFAPG_fit_Sigma_Irrotational (REAL8 barlam2)
 
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 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)
 
int eob_dyn_rhs (REAL8 t, const REAL8 y[], REAL8 dy[], void *params)
 
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)
 
int eob_dyn_rhs_s (REAL8 t, const REAL8 y[], REAL8 dy[], void *params)
 
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)
 
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_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 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)
 
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_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 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_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_fLR (REAL8 r, void *params)
 
int eob_dyn_adiabLR (LALTEOBResumSDynamics *dyn, REAL8 *rLR, INT4 tidesFlag)
 
int eob_dyn_Npostadiabatic (LALTEOBResumSDynamics *dyn, REAL8 r0)
 
void eob_dyn_ic (REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[])
 
void eob_dyn_ic_s (REAL8 r0, LALTEOBResumSDynamics *dyn, REAL8 y_init[])
 
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 eob_dyn_DHeff0 (REAL8 x, void *params)
 
void eob_metric_A5PNlog (REAL8 r, REAL8 nu, REAL8 *A, REAL8 *dA, REAL8 *d2A)
 
void eob_metric_Atidal (REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *AT, REAL8 *dAT, REAL8 *d2AT)
 
void eob_metric (REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB)
 
void eob_metric_s (REAL8 r, LALTEOBResumSDynamics *dyn, REAL8 *A, REAL8 *B, REAL8 *dA, REAL8 *d2A, REAL8 *dB)
 
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)
 
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)
 
void eob_flx_Tlm (REAL8 w, REAL8 *MTlm)
 
void eob_flx_FlmNewt (REAL8 x, REAL8 nu, REAL8 *Nlm)
 
REAL8 eob_flx_HorizonFlux (REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu)
 
REAL8 eob_flx_HorizonFlux_s (REAL8 x, REAL8 Heff, REAL8 jhat, REAL8 nu, REAL8 X1, REAL8 X2, REAL8 chi1, REAL8 chi2)
 
void eob_wav_hlm (LALTEOBResumSDynamics *dyn, LALTEOBResumSWaveformModeSingleTime *hlm_t)
 
void eob_wav_deltalm (REAL8 Hreal, REAL8 Omega, REAL8 nu, REAL8 *dlm)
 
void eob_wav_hhatlmTail (REAL8 Omega, REAL8 Hreal, REAL8 bphys, LALTEOBResumSWaveformModeSingleTime *tlm)
 
void eob_wav_speedyTail (REAL8 Omega, REAL8 Hreal, REAL8 bphys, LALTEOBResumSWaveformModeSingleTime *tlm)
 
void eob_wav_hlmNewt (REAL8 r, REAL8 Omega, REAL8 phi, REAL8 nu, LALTEOBResumSWaveformModeSingleTime *hNewt)
 
void eob_wav_hlmTidal (REAL8 x, LALTEOBResumSDynamics *dyn, REAL8 *hTidallm)
 
void eob_wav_flm_coeffs (REAL8 nu, REAL8 clm[KMAX][6])
 
void eob_wav_flm (REAL8 x, REAL8 nu, REAL8 clm[KMAX][6], 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)
 
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_hlmNQC_find_a1a2a3 (LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *h, SphHarmPolarTimeSeries *hnqc)
 
void eob_wav_hlmNQC_find_a1a2a3_mrg (LALTEOBResumSDynamics *dyn_mrg, SphHarmPolarTimeSeries *hlm_mrg, SphHarmPolarTimeSeries *hnqc, LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *hlm)
 
void eob_wav_hlmNQC (REAL8 nu, REAL8 r, REAL8 prstar, REAL8 Omega, REAL8 ddotr, NQCcoefs *nqc, LALTEOBResumSWaveformModeSingleTime *psilmnqc)
 
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 eob_wav_ringdown (LALTEOBResumSDynamics *dyn, SphHarmPolarTimeSeries *hlm)
 
REAL8 eob_dyn_bisecOmegaorb0 (LALTEOBResumSDynamics *dyn, REAL8 omg_orb0, REAL8 r0_kepl)
 
REAL8 eob_dyn_r0_Kepler (REAL8 f0)
 
REAL8 eob_dyn_r0_eob (REAL8 f0, LALTEOBResumSDynamics *dyn)
 
REAL8 eob_dyn_Omegaorb0 (REAL8 r, void *params)
 
void XLALTEOBDynamicsInterp (LALTEOBResumSDynamics *dyn, const INT4 size, const REAL8 t0, const REAL8 dt, const char *name)
 
void XLALTEOBDynamicsExtract (LALTEOBResumSDynamics *dyna, const REAL8 to, const REAL8 tn, LALTEOBResumSDynamics **dynb, const char *name)
 
void XLALTEOBDynamicsJoin (LALTEOBResumSDynamics *dyna, LALTEOBResumSDynamics *dynb, REAL8 to)
 
void eob_wav_hlmNQC_nospin201602 (REAL8 nu, REAL8 r, REAL8 prstar, REAL8 Omega, REAL8 ddotr, LALTEOBResumSWaveformModeSingleTime *hlmnqc)
 
void eob_nqc_setcoefs (LALTEOBResumSDynamics *dyn)
 
void eob_nqc_setcoefs_nospin201602 (REAL8 nu, NQCcoefs *nqc)
 
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)
 

Detailed Description

This file is part of TEOBResumS.

Copyright (C) 2017-2018 Alessandro Nagar, Sebastiano Bernuzzi, Sarp Ackay, Gregorio Carullo, Walter Del Pozzo, Ka Wa Tsang, Michalis Agathos LALSimulation implementation by Michalis Agathos

TEOBResumS is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.

TEOBResumS is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.

Header file of the TEOBResumS C code

This file contains all the macros, typdef, and routine prototype. Doxygen documentation should go here.

Definition in file LALSimTEOBResumS.h.

Go to the source code of this file.

Data Structures

struct  LALTEOBResumSWaveformModeSingleTime
 Multipolar waveform at given time, comes at handy. More...
 
struct  NQCcoefs
 Multipolar coefficients for NQC waveform. More...
 
struct  NQCdata
 NQC data for flux and waveform. More...
 
struct  LALTEOBResumSDynamics
 Dynamics data type. More...
 

Macros

#define POSTPOSTCIRCULAR   1 /* use post-post-circular initial conditions by default */
 
#define RC_EXCLUDESPINSPINTIDES   0 /* use tidally deformed centr. radius with self-spin and tides by default */
 
#define USEBTIDALPOTENTIAL   1 /* add B LO tidal potential */
 
#define ODESTOPAFTERNDT   (4) /* UNUSED: stop ODE integration after this many timesteps beyond the peak */
 
#define MAX_ITER   (200)
 
#define TOLERANCE   (1e-14)
 
#define POSTADIABATIC_RMIN_BNS   (14)
 
#define POSTADIABATIC_RMIN_BBH   (14)
 
#define POSTADIABATIC_DR   (0.1)
 
#define POSTADIABATIC_NSTEP_MIN   (10)
 
#define POSTADIABATIC_N   (8)
 
#define TEOB_R0_THRESHOLD   (14)
 
#define TEOB_LAMBDA_TOL   (1.0)
 
#define TEOB_GE_TIDES_DEFAULT   LAL_SIM_INSPIRAL_GETIDES_GSF23
 
#define TEOB_GM_TIDES_DEFAULT   LAL_SIM_INSPIRAL_GMTIDES_PN
 
#define TEOB_MODES_BNS_DEFAULT   TEOB_MODES_22
 
#define TEOB_MODES_BBH_DEFAULT   TEOB_MODES_22
 
#define TEOB_MODES_NSBH_DEFAULT   TEOB_MODES_22
 
#define TEOB_MODES_BHNS_DEFAULT   TEOB_MODES_22
 
#define RINGDOWN_EXTEND_ARRAY   (1500)
 
#define INTERP_UNIFORM_GRID_DEFAULT   INTERP_UNIFORM_GRID_HLM
 
#define DEBUG   0
 Macros. More...
 
#define TEOB_STRLEN   (128) /** Standard string length */
 
#define TEOBResumS_Info   "TEOBResumS code (C) 2017\n"
 
#define TEOBResumS_Usage(x)   {printf("%sUSAGE:\t%s <parfile>\n", TEOBResumS_Info, x);}
 
#define SIGN(x, y)   ((y) >= 0.0 ? fabs(x) : -fabs(x))
 
#define typeof   __typeof__
 
#define MAX(a, b)
 
#define MIN(a, b)
 
#define MAX3(a, b, c)   (((a) > (b)) ? MAX(a,c) : MAX(b,c))
 
#define MIN3(a, b, c)   (((a) < (b)) ? MIN(a,c) : MIN(b,c))
 
#define SQ(a)   ((a)*(a))
 
#define DEQUAL(a, b, eps)   (fabs((a)-(b))<(eps)) /** double compare */
 
#define DUNEQUAL(a, b, eps)   (fabs((a)-(b))>(eps))
 
#define STREQUAL(s, t)   ((strcmp((s),(t))==0)) /** string compare */
 
#define SWAPTRS(a, b)
 
#define Sqrt3   (1.73205080756887729352744634150587237)
 
#define Log1   (0.)
 
#define Log2   (0.693147180559945309417232)
 
#define Log3   (1.09861228866810969139525)
 
#define Log4   (1.38629436111989061883446)
 
#define Log5   (1.60943791243410037460076)
 
#define Log6   (1.79175946922805500081248)
 
#define Log7   (1.94591014905531330510535)
 
#define KMAX   (35) /** Multipolar linear index, max value */
 
#define PMTERMS_eps   (1) /** Switch on Fujita-Iyer point-mass terms. This is hard-coded here */
 
#define OUTPUT_LM_LEN   (35)
 
#define ROOTFINDER(i, x)   {if ( ((i) = (x)) && ((i)>ROOT_ERRORS_NO) ) { exit(LALTEOBResumSRootErrors[(i)]); }}
 

Typedefs

typedef void(* EOBWavFlmSFunc) (REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8[KMAX][6], int, REAL8 *, REAL8 *)
 Func pointer types. More...
 
typedef void(* EOBDynSGetRCFunc) (REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, INT4, REAL8 *, REAL8 *, REAL8 *)
 

Enumerations

enum  {
  TEOB_EVOLVE_RAD , TEOB_EVOLVE_PHI , TEOB_EVOLVE_PRSTAR , TEOB_EVOLVE_PPHI ,
  TEOB_EVOLVE_NVARS
}
 Index list of EOB evolved variables. More...
 
enum  {
  TEOB_ID_RAD , TEOB_ID_PHI , TEOB_ID_PPHI , TEOB_ID_PRSTAR ,
  TEOB_ID_PR , TEOB_ID_J , TEOB_ID_E0 , TEOB_ID_OMGJ ,
  TEOB_ID_NVARS
}
 Index list of EOB variables for initial data. More...
 
enum  {
  TEOB_RAD , TEOB_PHI , TEOB_PPHI , TEOB_MOMG ,
  TEOB_DDOTR , TEOB_PRSTAR , TEOB_OMGORB , TEOB_E0 ,
  TEOB_DYNAMICS_NVARS
}
 Index list of EOB dynamical variables (to be stored in arrays) More...
 
enum  { TEOB_MODES_ALL , TEOB_MODES_22 , TEOB_MODES_NOPT }
 Index list of TEOB mode-list options. More...
 
enum  {
  RC_LO , RC_NLO , RC_NNLO , RC_NNLOS4 ,
  RC_NOSPIN , RC_NOTIDES , RC_NOPT
}
 List of options for centrifugal radius. More...
 
enum  { SS_LO , SS_NLO , SS_NOPT }
 
enum  { ODE_TSTEP_UNIFORM , ODE_TSTEP_ADAPTIVE , ODE_TSTEP_ADAPTIVE_UNIFORM_AFTER_LSO , ODE_TSTEP_NOPT }
 List of options for ODE timestepping. More...
 
enum  { TEOB_PA_OFF , TEOB_PA_ON , TEOB_PA_ONLY , TEOB_PA_NOPT }
 List of options for postadiabatic inspiral. More...
 
enum  { INTERP_UNIFORM_GRID_OFF , INTERP_UNIFORM_GRID_HPC , INTERP_UNIFORM_GRID_HLM , INTERP_UNIFORM_GRID_NOPT }
 List of options for interp_uniform_grid. More...
 
enum  {
  NQC_OFF , NQC_COMPUTE , NQC_NR_NOSPIN , NQC_FILE ,
  NQC_NOPT
}
 List of options for NQC waveform nqc_coeffs_hlm and nqc_coeffs_flx. More...
 
enum  {
  ROOT_ERRORS_NO , ROOT_ERRORS_BRACKET , ROOT_ERRORS_MAXITS , ROOT_ERRORS_NOSUCC ,
  ROOT_ERRORS
}
 Error handler for root finders. More...
 

Variables

static const char *const LALTEOBResumSODETimeStep [] = {"uniform","adaptive","adaptive+uniform_after_LSO","undefined"}
 
static const char *const LALTEOBResumSRootErrors [] = {"none","root is not bracketed.","root finder did not converged.","root finder failed."}
 
const INT4 TEOB_LINDEX [KMAX]
 Maps between linear index and the corresponding (l, m) multipole indices. More...
 
const INT4 TEOB_MINDEX [KMAX]
 

Macro Definition Documentation

◆ POSTPOSTCIRCULAR

#define POSTPOSTCIRCULAR   1 /* use post-post-circular initial conditions by default */

Definition at line 80 of file LALSimTEOBResumS.h.

◆ RC_EXCLUDESPINSPINTIDES

#define RC_EXCLUDESPINSPINTIDES   0 /* use tidally deformed centr. radius with self-spin and tides by default */

Definition at line 84 of file LALSimTEOBResumS.h.

◆ USEBTIDALPOTENTIAL

#define USEBTIDALPOTENTIAL   1 /* add B LO tidal potential */

Definition at line 88 of file LALSimTEOBResumS.h.

◆ ODESTOPAFTERNDT

#define ODESTOPAFTERNDT   (4) /* UNUSED: stop ODE integration after this many timesteps beyond the peak */

Definition at line 92 of file LALSimTEOBResumS.h.

◆ MAX_ITER

#define MAX_ITER   (200)

Definition at line 95 of file LALSimTEOBResumS.h.

◆ TOLERANCE

#define TOLERANCE   (1e-14)

Definition at line 96 of file LALSimTEOBResumS.h.

◆ POSTADIABATIC_RMIN_BNS

#define POSTADIABATIC_RMIN_BNS   (14)

Definition at line 99 of file LALSimTEOBResumS.h.

◆ POSTADIABATIC_RMIN_BBH

#define POSTADIABATIC_RMIN_BBH   (14)

Definition at line 100 of file LALSimTEOBResumS.h.

◆ POSTADIABATIC_DR

#define POSTADIABATIC_DR   (0.1)

Definition at line 101 of file LALSimTEOBResumS.h.

◆ POSTADIABATIC_NSTEP_MIN

#define POSTADIABATIC_NSTEP_MIN   (10)

Definition at line 102 of file LALSimTEOBResumS.h.

◆ POSTADIABATIC_N

#define POSTADIABATIC_N   (8)

Definition at line 103 of file LALSimTEOBResumS.h.

◆ TEOB_R0_THRESHOLD

#define TEOB_R0_THRESHOLD   (14)

Definition at line 104 of file LALSimTEOBResumS.h.

◆ TEOB_LAMBDA_TOL

#define TEOB_LAMBDA_TOL   (1.0)

Definition at line 106 of file LALSimTEOBResumS.h.

◆ TEOB_GE_TIDES_DEFAULT

#define TEOB_GE_TIDES_DEFAULT   LAL_SIM_INSPIRAL_GETIDES_GSF23

Definition at line 108 of file LALSimTEOBResumS.h.

◆ TEOB_GM_TIDES_DEFAULT

#define TEOB_GM_TIDES_DEFAULT   LAL_SIM_INSPIRAL_GMTIDES_PN

Definition at line 109 of file LALSimTEOBResumS.h.

◆ TEOB_MODES_BNS_DEFAULT

#define TEOB_MODES_BNS_DEFAULT   TEOB_MODES_22

Definition at line 111 of file LALSimTEOBResumS.h.

◆ TEOB_MODES_BBH_DEFAULT

#define TEOB_MODES_BBH_DEFAULT   TEOB_MODES_22

Definition at line 112 of file LALSimTEOBResumS.h.

◆ TEOB_MODES_NSBH_DEFAULT

#define TEOB_MODES_NSBH_DEFAULT   TEOB_MODES_22

Definition at line 113 of file LALSimTEOBResumS.h.

◆ TEOB_MODES_BHNS_DEFAULT

#define TEOB_MODES_BHNS_DEFAULT   TEOB_MODES_22

Definition at line 114 of file LALSimTEOBResumS.h.

◆ RINGDOWN_EXTEND_ARRAY

#define RINGDOWN_EXTEND_ARRAY   (1500)

Definition at line 116 of file LALSimTEOBResumS.h.

◆ INTERP_UNIFORM_GRID_DEFAULT

#define INTERP_UNIFORM_GRID_DEFAULT   INTERP_UNIFORM_GRID_HLM

Definition at line 119 of file LALSimTEOBResumS.h.

◆ DEBUG

#define DEBUG   0

Macros.

Definition at line 123 of file LALSimTEOBResumS.h.

◆ TEOB_STRLEN

#define TEOB_STRLEN   (128) /** Standard string length */

Definition at line 126 of file LALSimTEOBResumS.h.

◆ TEOBResumS_Info

#define TEOBResumS_Info   "TEOBResumS code (C) 2017\n"

Definition at line 127 of file LALSimTEOBResumS.h.

◆ TEOBResumS_Usage

#define TEOBResumS_Usage (   x)    {printf("%sUSAGE:\t%s <parfile>\n", TEOBResumS_Info, x);}

Definition at line 128 of file LALSimTEOBResumS.h.

◆ SIGN

#define SIGN (   x,
  y 
)    ((y) >= 0.0 ? fabs(x) : -fabs(x))

Definition at line 129 of file LALSimTEOBResumS.h.

◆ typeof

#define typeof   __typeof__

Definition at line 130 of file LALSimTEOBResumS.h.

◆ MAX

#define MAX (   a,
 
)
Value:
({ typeof (a) _a = (a); \
typeof (b) _b = (b); \
_a > _b ? _a : _b; })
#define typeof
static const INT4 a

Definition at line 131 of file LALSimTEOBResumS.h.

◆ MIN

#define MIN (   a,
 
)
Value:
({ typeof (a) _a = (a); \
typeof (b) _b = (b); \
_a < _b ? _a : _b; })

Definition at line 135 of file LALSimTEOBResumS.h.

◆ MAX3

#define MAX3 (   a,
  b,
  c 
)    (((a) > (b)) ? MAX(a,c) : MAX(b,c))

Definition at line 139 of file LALSimTEOBResumS.h.

◆ MIN3

#define MIN3 (   a,
  b,
  c 
)    (((a) < (b)) ? MIN(a,c) : MIN(b,c))

Definition at line 140 of file LALSimTEOBResumS.h.

◆ SQ

#define SQ (   a)    ((a)*(a))

Definition at line 141 of file LALSimTEOBResumS.h.

◆ DEQUAL

#define DEQUAL (   a,
  b,
  eps 
)    (fabs((a)-(b))<(eps)) /** double compare */

Definition at line 142 of file LALSimTEOBResumS.h.

◆ DUNEQUAL

#define DUNEQUAL (   a,
  b,
  eps 
)    (fabs((a)-(b))>(eps))

Definition at line 143 of file LALSimTEOBResumS.h.

◆ STREQUAL

#define STREQUAL (   s,
 
)    ((strcmp((s),(t))==0)) /** string compare */

Definition at line 144 of file LALSimTEOBResumS.h.

◆ SWAPTRS

#define SWAPTRS (   a,
 
)
Value:
({ \
typeof(a) temp; \
temp = a; \
a = b; \
b = temp; \
})

Definition at line 145 of file LALSimTEOBResumS.h.

◆ Sqrt3

#define Sqrt3   (1.73205080756887729352744634150587237)

Definition at line 156 of file LALSimTEOBResumS.h.

◆ Log1

#define Log1   (0.)

Definition at line 159 of file LALSimTEOBResumS.h.

◆ Log2

#define Log2   (0.693147180559945309417232)

Definition at line 161 of file LALSimTEOBResumS.h.

◆ Log3

#define Log3   (1.09861228866810969139525)

Definition at line 162 of file LALSimTEOBResumS.h.

◆ Log4

#define Log4   (1.38629436111989061883446)

Definition at line 163 of file LALSimTEOBResumS.h.

◆ Log5

#define Log5   (1.60943791243410037460076)

Definition at line 164 of file LALSimTEOBResumS.h.

◆ Log6

#define Log6   (1.79175946922805500081248)

Definition at line 165 of file LALSimTEOBResumS.h.

◆ Log7

#define Log7   (1.94591014905531330510535)

Definition at line 166 of file LALSimTEOBResumS.h.

◆ KMAX

#define KMAX   (35) /** Multipolar linear index, max value */

Definition at line 215 of file LALSimTEOBResumS.h.

◆ PMTERMS_eps

#define PMTERMS_eps   (1) /** Switch on Fujita-Iyer point-mass terms. This is hard-coded here */

Definition at line 216 of file LALSimTEOBResumS.h.

◆ OUTPUT_LM_LEN

#define OUTPUT_LM_LEN   (35)

Definition at line 218 of file LALSimTEOBResumS.h.

◆ ROOTFINDER

#define ROOTFINDER (   i,
  x 
)    {if ( ((i) = (x)) && ((i)>ROOT_ERRORS_NO) ) { exit(LALTEOBResumSRootErrors[(i)]); }}

Definition at line 281 of file LALSimTEOBResumS.h.

Typedef Documentation

◆ EOBWavFlmSFunc

typedef void(* EOBWavFlmSFunc) (REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8[KMAX][6], int, REAL8 *, REAL8 *)

Func pointer types.

Definition at line 299 of file LALSimTEOBResumS.h.

◆ EOBDynSGetRCFunc

typedef void(* EOBDynSGetRCFunc) (REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, REAL8, INT4, REAL8 *, REAL8 *, REAL8 *)

Definition at line 300 of file LALSimTEOBResumS.h.

Enumeration Type Documentation

◆ anonymous enum

anonymous enum

Index list of EOB evolved variables.

Enumerator
TEOB_EVOLVE_RAD 
TEOB_EVOLVE_PHI 
TEOB_EVOLVE_PRSTAR 
TEOB_EVOLVE_PPHI 
TEOB_EVOLVE_NVARS 

Definition at line 170 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

Index list of EOB variables for initial data.

Enumerator
TEOB_ID_RAD 
TEOB_ID_PHI 
TEOB_ID_PPHI 
TEOB_ID_PRSTAR 
TEOB_ID_PR 
TEOB_ID_J 
TEOB_ID_E0 
TEOB_ID_OMGJ 
TEOB_ID_NVARS 

Definition at line 180 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

Index list of EOB dynamical variables (to be stored in arrays)

Enumerator
TEOB_RAD 
TEOB_PHI 
TEOB_PPHI 
TEOB_MOMG 
TEOB_DDOTR 
TEOB_PRSTAR 
TEOB_OMGORB 
TEOB_E0 
TEOB_DYNAMICS_NVARS 

Definition at line 194 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

Index list of TEOB mode-list options.

Enumerator
TEOB_MODES_ALL 
TEOB_MODES_22 
TEOB_MODES_NOPT 

Definition at line 207 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

List of options for centrifugal radius.

Enumerator
RC_LO 
RC_NLO 
RC_NNLO 
RC_NNLOS4 
RC_NOSPIN 
RC_NOTIDES 
RC_NOPT 

Definition at line 222 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum
Enumerator
SS_LO 
SS_NLO 
SS_NOPT 

Definition at line 232 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

List of options for ODE timestepping.

Enumerator
ODE_TSTEP_UNIFORM 
ODE_TSTEP_ADAPTIVE 
ODE_TSTEP_ADAPTIVE_UNIFORM_AFTER_LSO 
ODE_TSTEP_NOPT 

Definition at line 239 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

List of options for postadiabatic inspiral.

Enumerator
TEOB_PA_OFF 
TEOB_PA_ON 
TEOB_PA_ONLY 
TEOB_PA_NOPT 

Definition at line 248 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

List of options for interp_uniform_grid.

Enumerator
INTERP_UNIFORM_GRID_OFF 
INTERP_UNIFORM_GRID_HPC 
INTERP_UNIFORM_GRID_HLM 
INTERP_UNIFORM_GRID_NOPT 

Definition at line 256 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

List of options for NQC waveform nqc_coeffs_hlm and nqc_coeffs_flx.

Enumerator
NQC_OFF 
NQC_COMPUTE 
NQC_NR_NOSPIN 
NQC_FILE 
NQC_NOPT 

Definition at line 264 of file LALSimTEOBResumS.h.

◆ anonymous enum

anonymous enum

Error handler for root finders.

Enumerator
ROOT_ERRORS_NO 
ROOT_ERRORS_BRACKET 
ROOT_ERRORS_MAXITS 
ROOT_ERRORS_NOSUCC 
ROOT_ERRORS 

Definition at line 273 of file LALSimTEOBResumS.h.

Function Documentation

◆ XLALSetup_TEOB_mode_array()

INT4 XLALSetup_TEOB_mode_array ( LALValue *  ModeArray,
INT4  modeType 
)

Should not be needed in LAL (all parameters passed to interface function)

  • Factorial *‍/* Wigner d-function *‍/* Spin-weighted spherical harmonic

Definition at line 296 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALCheck_TEOB_mode_array()

INT4 XLALCheck_TEOB_mode_array ( LALValue *  ModeArray,
UINT4  use_tidal 
)

Definition at line 325 of file LALSimIMRTEOBResumSUtils.c.

◆ q_to_nu()

REAL8 q_to_nu ( const REAL8  q)

Definition at line 48 of file LALSimIMRTEOBResumSUtils.c.

◆ nu_to_X1()

REAL8 nu_to_X1 ( const REAL8  nu)

Definition at line 58 of file LALSimIMRTEOBResumSUtils.c.

◆ Eulerlog()

REAL8 Eulerlog ( const REAL8  x,
const INT4  m 
)

Definition at line 65 of file LALSimIMRTEOBResumSUtils.c.

◆ interp_spline()

void interp_spline ( REAL8 t,
REAL8 y,
INT4  n,
REAL8 ti,
INT4  ni,
REAL8 yi 
)

◆ find_point_bisection()

int find_point_bisection ( REAL8  x,
INT4  n,
REAL8 xp,
INT4  o 
)

Definition at line 93 of file LALSimIMRTEOBResumSUtils.c.

◆ baryc_f()

REAL8 baryc_f ( REAL8  xx,
INT4  n,
REAL8 f,
REAL8 x 
)

Definition at line 114 of file LALSimIMRTEOBResumSUtils.c.

◆ baryc_weights()

void baryc_weights ( INT4  n,
REAL8 x,
REAL8 omega 
)

Definition at line 145 of file LALSimIMRTEOBResumSUtils.c.

◆ baryc_f_weights()

REAL8 baryc_f_weights ( REAL8  xx,
INT4  n,
REAL8 f,
REAL8 x,
REAL8 omega 
)

Definition at line 165 of file LALSimIMRTEOBResumSUtils.c.

◆ interp1d()

REAL8 interp1d ( const INT4  order,
REAL8  xx,
INT4  nx,
REAL8 f,
REAL8 x 
)

Definition at line 182 of file LALSimIMRTEOBResumSUtils.c.

◆ find_max()

REAL8 find_max ( const INT4  n,
REAL8  dx,
REAL8  x0,
REAL8 f,
REAL8 fmax 
)

Definition at line 193 of file LALSimIMRTEOBResumSUtils.c.

◆ D0()

int D0 ( REAL8 f,
REAL8  dx,
INT4  n,
REAL8 df 
)

XLALWignerdMatrix, XLALWignerDMatrix in <lal/SphericalHarmonics.h>

Definition at line 473 of file LALSimIMRTEOBResumSUtils.c.

◆ D2()

int D2 ( REAL8 f,
REAL8  dx,
INT4  n,
REAL8 d2f 
)

Definition at line 493 of file LALSimIMRTEOBResumSUtils.c.

◆ cumint3()

REAL8 cumint3 ( REAL8 f,
REAL8 x,
const INT4  n,
REAL8 sum 
)

Definition at line 513 of file LALSimIMRTEOBResumSUtils.c.

◆ unwrap()

void unwrap ( REAL8 p,
const INT4  size 
)

Definition at line 559 of file LALSimIMRTEOBResumSUtils.c.

◆ unwrap_proxy()

void unwrap_proxy ( REAL8 p,
REAL8 r,
const INT4  size,
const INT4  shift0 
)

Definition at line 587 of file LALSimIMRTEOBResumSUtils.c.

◆ get_uniform_size()

int get_uniform_size ( const REAL8  tf,
const REAL8  t0,
const REAL8  dt 
)

Definition at line 624 of file LALSimIMRTEOBResumSUtils.c.

◆ Waveform_lm_t_alloc()

void Waveform_lm_t_alloc ( LALTEOBResumSWaveformModeSingleTime **  wav)

Replaced by LAL COMPLEX16TimeSeries functions.

Definition at line 639 of file LALSimIMRTEOBResumSUtils.c.

◆ Waveform_lm_t_free()

void Waveform_lm_t_free ( LALTEOBResumSWaveformModeSingleTime wav)

Definition at line 648 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALTEOBDynamicsInit()

void XLALTEOBDynamicsInit ( LALTEOBResumSDynamics **  dyn,
INT4  size,
const CHAR name 
)

Definition at line 674 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALTEOBDynamicsPush()

void XLALTEOBDynamicsPush ( LALTEOBResumSDynamics **  dyn,
INT4  size 
)

Definition at line 691 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALFreeTEOBDynamics()

void XLALFreeTEOBDynamics ( LALTEOBResumSDynamics dyn)

Definition at line 836 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALTEOBDynamicsSetParams()

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 
)

◆ NQCdata_alloc()

void NQCdata_alloc ( NQCdata **  nqc)

Definition at line 655 of file LALSimIMRTEOBResumSUtils.c.

◆ NQCdata_free()

void NQCdata_free ( NQCdata nqc)

Definition at line 665 of file LALSimIMRTEOBResumSUtils.c.

◆ time_units_factor()

REAL8 time_units_factor ( REAL8  M)

Definition at line 1204 of file LALSimIMRTEOBResumSUtils.c.

◆ time_units_conversion()

REAL8 time_units_conversion ( REAL8  M,
REAL8  t 
)

Definition at line 1209 of file LALSimIMRTEOBResumSUtils.c.

◆ radius0()

REAL8 radius0 ( REAL8  M,
REAL8  fHz 
)

Definition at line 1215 of file LALSimIMRTEOBResumSUtils.c.

◆ get_mrg_timestep()

REAL8 get_mrg_timestep ( REAL8  q,
REAL8  chi1,
REAL8  chi2 
)

Definition at line 630 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALSimIMRComputePolarisations()

void XLALSimIMRComputePolarisations ( REAL8Sequence hplus_out,
REAL8Sequence hcross_out,
SphHarmTimeSeries hlm,
LALValue *  modeArray,
REAL8  amplitude_prefactor,
REAL8  theta,
REAL8  phi 
)

Definition at line 365 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALSimIMRComputePolarisationsPolar()

void XLALSimIMRComputePolarisationsPolar ( REAL8Sequence hplus_out,
REAL8Sequence hcross_out,
SphHarmPolarTimeSeries hlm,
LALValue *  modeArray,
REAL8  amplitude_prefactor,
REAL8  theta,
REAL8  phi 
)

Definition at line 418 of file LALSimIMRTEOBResumSUtils.c.

◆ eob_nqc_timeshift()

REAL8 eob_nqc_timeshift ( REAL8  nu,
REAL8  chi1 
)

Definition at line 4433 of file LALSimIMRTEOBResumS_Internals.c.

◆ logQ()

REAL8 logQ ( REAL8  x)

Definition at line 1404 of file LALSimIMRTEOBResumSUtils.c.

◆ Yagi13_fit_barlamdel()

REAL8 Yagi13_fit_barlamdel ( REAL8  barlam2,
INT4  ell 
)

◆ Yagi13_fit_barsigmalambda()

REAL8 Yagi13_fit_barsigmalambda ( REAL8  barlam2)

Definition at line 1451 of file LALSimIMRTEOBResumSUtils.c.

◆ Yagi14_fit_Coct()

REAL8 Yagi14_fit_Coct ( REAL8  C_Q)

Definition at line 1475 of file LALSimIMRTEOBResumSUtils.c.

◆ Yagi14_fit_Chex()

REAL8 Yagi14_fit_Chex ( REAL8  C_Q)

Definition at line 1488 of file LALSimIMRTEOBResumSUtils.c.

◆ JFAPG_fit_Sigma_Static()

REAL8 JFAPG_fit_Sigma_Static ( REAL8  barlam2)

Definition at line 1518 of file LALSimIMRTEOBResumSUtils.c.

◆ JFAPG_fit_Sigma_Irrotational()

REAL8 JFAPG_fit_Sigma_Irrotational ( REAL8  barlam2)

Definition at line 1500 of file LALSimIMRTEOBResumSUtils.c.

◆ HealyBBHFitRemnant()

void HealyBBHFitRemnant ( REAL8  chi1,
REAL8  chi2,
REAL8  q,
REAL8 mass,
REAL8 spin 
)

Definition at line 1234 of file LALSimIMRTEOBResumSUtils.c.

◆ JimenezFortezaRemnantSpin()

REAL8 JimenezFortezaRemnantSpin ( REAL8  nu,
REAL8  X1,
REAL8  X2,
REAL8  chi1,
REAL8  chi2 
)

Definition at line 1336 of file LALSimIMRTEOBResumSUtils.c.

◆ QNMHybridFitCab()

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 
)

◆ eob_dyn_rhs()

int eob_dyn_rhs ( REAL8  t,
const REAL8  y[],
REAL8  dy[],
void *  params 
)

Definition at line 662 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_ham()

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 
)

Definition at line 63 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_rhs_s()

int eob_dyn_rhs_s ( REAL8  t,
const REAL8  y[],
REAL8  dy[],
void *  params 
)

Definition at line 758 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_ham_s()

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 
)

Definition at line 90 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_s_GS()

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 
)

◆ eob_dyn_s_get_rc_LO()

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 
)

◆ eob_dyn_s_get_rc_NLO()

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 
)

◆ eob_dyn_s_get_rc_NNLO()

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 
)

◆ eob_dyn_s_get_rc_NNLO_S4()

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 
)

◆ eob_dyn_s_get_rc_NOSPIN()

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 
)

◆ eob_dyn_s_get_rc_NOTIDES()

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 
)

◆ eob_dyn_fLR()

REAL8 eob_dyn_fLR ( REAL8  r,
void *  params 
)

Definition at line 4992 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_adiabLR()

int eob_dyn_adiabLR ( LALTEOBResumSDynamics dyn,
REAL8 rLR,
INT4  tidesFlag 
)

Definition at line 5005 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_Npostadiabatic()

int eob_dyn_Npostadiabatic ( LALTEOBResumSDynamics dyn,
const REAL8  r0 
)
  • Root finder for adiabatic LSO *‍/

Definition at line 5159 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_ic()

void eob_dyn_ic ( REAL8  r0,
LALTEOBResumSDynamics dyn,
REAL8  y_init[] 
)

Definition at line 411 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_ic_s()

void eob_dyn_ic_s ( REAL8  r0,
LALTEOBResumSDynamics dyn,
REAL8  y_init[] 
)

Definition at line 484 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_bisecHeff0_s()

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 
)

Definition at line 197 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_DHeff0()

REAL8 eob_dyn_DHeff0 ( REAL8  x,
void *  params 
)

Definition at line 158 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_metric_A5PNlog()

void eob_metric_A5PNlog ( REAL8  r,
REAL8  nu,
REAL8 A,
REAL8 dA,
REAL8 d2A 
)

Definition at line 1380 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_metric_Atidal()

void eob_metric_Atidal ( REAL8  r,
LALTEOBResumSDynamics dyn,
REAL8 AT,
REAL8 dAT,
REAL8 d2AT 
)

Definition at line 1471 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_metric()

void eob_metric ( REAL8  r,
LALTEOBResumSDynamics dyn,
REAL8 A,
REAL8 B,
REAL8 dA,
REAL8 d2A,
REAL8 dB 
)

Definition at line 1752 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_metric_s()

void eob_metric_s ( REAL8  r,
LALTEOBResumSDynamics dyn,
REAL8 A,
REAL8 B,
REAL8 dA,
REAL8 d2A,
REAL8 dB 
)

Definition at line 1811 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_Flux()

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 
)

Definition at line 2191 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_Flux_s()

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 
)

Definition at line 2049 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_Tlm()

void eob_flx_Tlm ( REAL8  w,
REAL8 MTlm 
)

Definition at line 1929 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_FlmNewt()

void eob_flx_FlmNewt ( REAL8  x,
REAL8  nu,
REAL8 Nlm 
)

Definition at line 1878 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_HorizonFlux()

REAL8 eob_flx_HorizonFlux ( REAL8  x,
REAL8  Heff,
REAL8  jhat,
REAL8  nu 
)

Definition at line 1952 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_flx_HorizonFlux_s()

REAL8 eob_flx_HorizonFlux_s ( REAL8  x,
REAL8  Heff,
REAL8  jhat,
REAL8  nu,
REAL8  X1,
REAL8  X2,
REAL8  chi1,
REAL8  chi2 
)

◆ eob_wav_hlm()

void eob_wav_hlm ( LALTEOBResumSDynamics dyn,
LALTEOBResumSWaveformModeSingleTime hlm_t 
)

Definition at line 4248 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_deltalm()

void eob_wav_deltalm ( REAL8  Hreal,
REAL8  Omega,
REAL8  nu,
REAL8 dlm 
)

Definition at line 2908 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_hhatlmTail()

void eob_wav_hhatlmTail ( REAL8  Omega,
REAL8  Hreal,
REAL8  bphys,
LALTEOBResumSWaveformModeSingleTime tlm 
)

Definition at line 2263 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_speedyTail()

void eob_wav_speedyTail ( REAL8  Omega,
REAL8  Hreal,
REAL8  bphys,
LALTEOBResumSWaveformModeSingleTime tlm 
)

Definition at line 2207 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_hlmNewt()

void eob_wav_hlmNewt ( REAL8  r,
REAL8  Omega,
REAL8  phi,
REAL8  nu,
LALTEOBResumSWaveformModeSingleTime hNewt 
)

Definition at line 2971 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_hlmTidal()

void eob_wav_hlmTidal ( REAL8  x,
LALTEOBResumSDynamics dyn,
REAL8 hTidallm 
)

Definition at line 3047 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_flm_coeffs()

void eob_wav_flm_coeffs ( REAL8  nu,
REAL8  clm[KMAX][6] 
)

Definition at line 2584 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_flm()

void eob_wav_flm ( REAL8  x,
REAL8  nu,
REAL8  clm[KMAX][6],
REAL8 rholm,
REAL8 flm 
)

Definition at line 2784 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_flm_s_SSLO()

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 
)

◆ eob_wav_flm_s_SSNLO()

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 
)

◆ eob_wav_hlmNQC_find_a1a2a3()

void eob_wav_hlmNQC_find_a1a2a3 ( LALTEOBResumSDynamics dyn,
SphHarmPolarTimeSeries h,
SphHarmPolarTimeSeries hnqc 
)

Definition at line 3286 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_hlmNQC_find_a1a2a3_mrg()

void eob_wav_hlmNQC_find_a1a2a3_mrg ( LALTEOBResumSDynamics dyn_mrg,
SphHarmPolarTimeSeries hlm_mrg,
SphHarmPolarTimeSeries hnqc,
LALTEOBResumSDynamics dyn,
SphHarmPolarTimeSeries hlm 
)

Definition at line 3736 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_hlmNQC()

void eob_wav_hlmNQC ( REAL8  nu,
REAL8  r,
REAL8  prstar,
REAL8  Omega,
REAL8  ddotr,
NQCcoefs nqc,
LALTEOBResumSWaveformModeSingleTime psilmnqc 
)

◆ eob_wav_ringdown_template()

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 
)

Definition at line 4784 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_wav_ringdown()

void eob_wav_ringdown ( LALTEOBResumSDynamics dyn,
SphHarmPolarTimeSeries hlm 
)

Definition at line 4805 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_bisecOmegaorb0()

REAL8 eob_dyn_bisecOmegaorb0 ( LALTEOBResumSDynamics dyn,
REAL8  omg_orb0,
REAL8  r0_kepl 
)

Definition at line 367 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_r0_Kepler()

REAL8 eob_dyn_r0_Kepler ( REAL8  f0)

Definition at line 246 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_r0_eob()

REAL8 eob_dyn_r0_eob ( REAL8  f0,
LALTEOBResumSDynamics dyn 
)

Definition at line 253 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_dyn_Omegaorb0()

REAL8 eob_dyn_Omegaorb0 ( REAL8  r,
void *  params 
)

Definition at line 266 of file LALSimIMRTEOBResumS_Internals.c.

◆ XLALTEOBDynamicsInterp()

void XLALTEOBDynamicsInterp ( LALTEOBResumSDynamics dyn,
const INT4  size,
const REAL8  t0,
const REAL8  dt,
const char name 
)

Definition at line 704 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALTEOBDynamicsExtract()

void XLALTEOBDynamicsExtract ( LALTEOBResumSDynamics dyna,
const REAL8  to,
const REAL8  tn,
LALTEOBResumSDynamics **  dynb,
const char name 
)

Definition at line 743 of file LALSimIMRTEOBResumSUtils.c.

◆ XLALTEOBDynamicsJoin()

void XLALTEOBDynamicsJoin ( LALTEOBResumSDynamics dyna,
LALTEOBResumSDynamics dynb,
REAL8  to 
)

Definition at line 795 of file LALSimIMRTEOBResumSUtils.c.

◆ eob_wav_hlmNQC_nospin201602()

void eob_wav_hlmNQC_nospin201602 ( REAL8  nu,
REAL8  r,
REAL8  prstar,
REAL8  Omega,
REAL8  ddotr,
LALTEOBResumSWaveformModeSingleTime hlmnqc 
)

Definition at line 3123 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_nqc_setcoefs()

void eob_nqc_setcoefs ( LALTEOBResumSDynamics dyn)

Definition at line 4461 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_nqc_setcoefs_nospin201602()

void eob_nqc_setcoefs_nospin201602 ( REAL8  nu,
NQCcoefs nqc 
)

Definition at line 4525 of file LALSimIMRTEOBResumS_Internals.c.

◆ eob_c3_fit_global()

REAL8 eob_c3_fit_global ( REAL8  nu,
REAL8  chi1,
REAL8  chi2,
REAL8  X1,
REAL8  X2,
REAL8  a1,
REAL8  a2 
)

◆ XLALSphHarmPolarJoin()

void XLALSphHarmPolarJoin ( SphHarmPolarTimeSeries hlma,
SphHarmPolarTimeSeries hlmb,
REAL8  to 
)

Definition at line 1536 of file LALSimIMRTEOBResumSUtils.c.

Variable Documentation

◆ LALTEOBResumSODETimeStep

const char* const LALTEOBResumSODETimeStep[] = {"uniform","adaptive","adaptive+uniform_after_LSO","undefined"}
static

Definition at line 245 of file LALSimTEOBResumS.h.

◆ LALTEOBResumSRootErrors

const char* const LALTEOBResumSRootErrors[] = {"none","root is not bracketed.","root finder did not converged.","root finder failed."}
static

Definition at line 280 of file LALSimTEOBResumS.h.

◆ TEOB_LINDEX

const INT4 TEOB_LINDEX[KMAX]
extern

Maps between linear index and the corresponding (l, m) multipole indices.

Maps between linear index and the corresponding (l, m) multipole indices.

Global vars, defined as external in header

Definition at line 59 of file LALSimIMRTEOBResumS.c.

◆ TEOB_MINDEX

const INT4 TEOB_MINDEX[KMAX]
extern

Definition at line 68 of file LALSimIMRTEOBResumS.c.