LALSimulation  5.4.0.1-fe68b98
Module LALSimInspiralSpinTaylor.c

Detailed Description

Routines for generating spin precessing Taylor waveforms.

Deprecated:
Use XLALSimInspiralChooseTDWaveform() instead

Prototypes

int XLALSimInspiralTransformPrecessingObsoleteInitialConditions (REAL8 *incl, REAL8 *S1x, REAL8 *S1y, REAL8 *S1z, REAL8 *S2x, REAL8 *S2y, REAL8 *S2z, REAL8 thetaJN, REAL8 phiJL, REAL8 theta1, REAL8 theta2, REAL8 phi12, REAL8 chi1, REAL8 chi2, REAL8 m1, REAL8 m2, REAL8 fRef)
 Function to specify the desired orientation of a precessing binary in terms of several angles and then compute the vector components in the so-called "radiation frame" (with the z-axis along the direction of propagation) as needed to specify binary configuration for ChooseTDWaveform. More...
 
int XLALSimInspiralInitialConditionsPrecessingApproxs (REAL8 *inc, REAL8 *S1x, REAL8 *S1y, REAL8 *S1z, REAL8 *S2x, REAL8 *S2y, REAL8 *S2z, const REAL8 inclIn, const REAL8 S1xIn, const REAL8 S1yIn, const REAL8 S1zIn, const REAL8 S2xIn, const REAL8 S2yIn, const REAL8 S2zIn, const REAL8 m1, const REAL8 m2, const REAL8 fRef, const REAL8 phiRef, LALSimInspiralFrameAxis frameChoice)
 Function to specify the desired orientation of the spin components of a precessing binary. More...
 
int XLALSimInspiralSpinTaylorPNEvolveOrbit (REAL8TimeSeries **V, REAL8TimeSeries **Phi, REAL8TimeSeries **S1x, REAL8TimeSeries **S1y, REAL8TimeSeries **S1z, REAL8TimeSeries **S2x, REAL8TimeSeries **S2y, REAL8TimeSeries **S2z, REAL8TimeSeries **LNhatx, REAL8TimeSeries **LNhaty, REAL8TimeSeries **LNhatz, REAL8TimeSeries **E1x, REAL8TimeSeries **E1y, REAL8TimeSeries **E1z, const REAL8 deltaT, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 fStart, const REAL8 fEnd, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 s2x, const REAL8 s2y, const REAL8 s2z, const REAL8 lnhatx, const REAL8 lnhaty, const REAL8 lnhatz, const REAL8 e1x, const REAL8 e1y, const REAL8 e1z, const REAL8 lambda1, const REAL8 lambda2, const REAL8 quadparam1, const REAL8 quadparam2, const LALSimInspiralSpinOrder spinO, const LALSimInspiralTidalOrder tideO, const INT4 phaseO, const INT4 lscorr, const Approximant approx)
 This function evolves the orbital equations for a precessing binary using the "TaylorT1/T5/T4" approximant for solving the orbital dynamics (see arXiv:0907.0700 for a review of the various PN approximants). More...
 
int XLALSimInspiralSpinTaylorPNEvolveOrbitOnlyFinal (REAL8TimeSeries **V, REAL8TimeSeries **Phi, REAL8TimeSeries **S1x, REAL8TimeSeries **S1y, REAL8TimeSeries **S1z, REAL8TimeSeries **S2x, REAL8TimeSeries **S2y, REAL8TimeSeries **S2z, REAL8TimeSeries **LNhatx, REAL8TimeSeries **LNhaty, REAL8TimeSeries **LNhatz, REAL8TimeSeries **E1x, REAL8TimeSeries **E1y, REAL8TimeSeries **E1z, const REAL8 deltaT, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 fStart, const REAL8 fEnd, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 s2x, const REAL8 s2y, const REAL8 s2z, const REAL8 lnhatx, const REAL8 lnhaty, const REAL8 lnhatz, const REAL8 e1x, const REAL8 e1y, const REAL8 e1z, const REAL8 lambda1, const REAL8 lambda2, const REAL8 quadparam1, const REAL8 quadparam2, const LALSimInspiralSpinOrder spinO, const LALSimInspiralTidalOrder tideO, const INT4 phaseO, const INT4 lscorr, const Approximant approx)
 This function evolves the orbital equations for a precessing binary using the "TaylorT1/T5/T4" approximant for solving the orbital dynamics (see arXiv:0907.0700 for a review of the various PN approximants). More...
 
int XLALSimInspiralSpinTaylorT4OLD (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 UNUSED v0, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, REAL8 lambda1, REAL8 lambda2, REAL8 quadparam1, REAL8 quadparam2, LALDict *LALparams, int phaseO, int amplitudeO)
 Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T4" method. More...
 
int XLALSimInspiralSpinTaylorT4 (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, LALDict *LALparams)
 
int XLALSimInspiralSpinTaylorT1OLD (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 UNUSED v0, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, REAL8 lambda1, REAL8 lambda2, REAL8 quadparam1, REAL8 quadparam2, LALDict *LALparams, int phaseO, int amplitudeO)
 Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T1" method. More...
 
int XLALSimInspiralSpinTaylorT1 (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, LALDict *LALparams)
 
int XLALSimInspiralSpinTaylorT5 (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, LALDict *LALparams)
 Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T5" method. More...
 
int XLALSimInspiralPrecessingPTFQWaveforms (REAL8TimeSeries **Q1, REAL8TimeSeries **Q2, REAL8TimeSeries **Q3, REAL8TimeSeries **Q4, REAL8TimeSeries **Q5, REAL8TimeSeries *V, REAL8TimeSeries *Phi, REAL8TimeSeries *S1x, REAL8TimeSeries *S1y, REAL8TimeSeries *S1z, REAL8TimeSeries *S2x, REAL8TimeSeries *S2y, REAL8TimeSeries *S2z, REAL8TimeSeries *LNhatx, REAL8TimeSeries *LNhaty, REAL8TimeSeries *LNhatz, REAL8TimeSeries *E1x, REAL8TimeSeries *E1y, REAL8TimeSeries *E1z, REAL8 m1, REAL8 m2, REAL8 r)
 Compute the physical template family "Q" vectors for a spinning, precessing binary when provided time series of all the dynamical quantities. More...
 
int XLALSimInspiralSpinTaylorT4PTFQVecs (REAL8TimeSeries **Q1, REAL8TimeSeries **Q2, REAL8TimeSeries **Q3, REAL8TimeSeries **Q4, REAL8TimeSeries **Q5, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 chi1, REAL8 kappa1, REAL8 fStart, REAL8 lambda1, REAL8 lambda2, LALSimInspiralSpinOrder spinO, LALSimInspiralTidalOrder tideO, int phaseO)
 Driver routine to compute the physical template family "Q" vectors using the "T4" method. More...
 
int XLALSimInspiralSpinTaylorT4Fourier (COMPLEX16FrequencySeries **hplus, COMPLEX16FrequencySeries **hcross, REAL8 fMin, REAL8 fMax, REAL8 deltaF, INT4 kMax, REAL8 phiRef, REAL8 v0, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, REAL8 lambda1, REAL8 lambda2, REAL8 quadparam1, REAL8 quadparam2, LALDict *LALparams, INT4 phaseO, INT4 amplitudeO, INT4 phiRefAtEnd)
 Driver routine to compute a precessing post-Newtonian inspiral waveform in the Fourier domain with phasing computed from energy balance using the so-called "T4" method. More...
 
int XLALSimInspiralSpinTaylorT5Fourier (COMPLEX16FrequencySeries **hplus, COMPLEX16FrequencySeries **hcross, REAL8 fMin, REAL8 fMax, REAL8 deltaF, INT4 kMax, REAL8 phiRef, REAL8 v0, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 fRef, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhatx, REAL8 lnhaty, REAL8 lnhatz, REAL8 e1x, REAL8 e1y, REAL8 e1z, REAL8 lambda1, REAL8 lambda2, REAL8 quadparam1, REAL8 quadparam2, LALDict *LALparams, INT4 phaseO, INT4 amplitudeO, INT4 phiRefAtEnd)
 Driver routine to compute a precessing post-Newtonian inspiral waveform in the Fourier domain with phasing computed from energy balance using the so-called "T2" method see arXiv: 0907.0700 for its definition, but in its "T5" variant, see arXiv: 1107.1267. More...
 
int XLALSimInspiralSpinTaylorF2 (COMPLEX16FrequencySeries **hplus_out, COMPLEX16FrequencySeries **hcross_out, const REAL8 phi_ref, const REAL8 deltaF, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 lnhatx, const REAL8 lnhaty, const REAL8 lnhatz, const REAL8 fStart, const REAL8 fEnd, const REAL8 f_ref, const REAL8 r, LALDict *moreParams, const INT4 phaseO, const INT4 amplitudeO)
 Computes the stationary phase approximation to the Fourier transform of a chirp waveform with phase given by Eq. More...
 
int XLALSimInspiralSpinTaylorT5duplicate (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 fStart, REAL8 r, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 incAngle, int phaseO, int amplitudeO)
 Generate time-domain generic spinning PN waveforms in the SpinTaylorT5 approximaton. More...
 

Function Documentation

◆ XLALSimInspiralTransformPrecessingObsoleteInitialConditions()

int XLALSimInspiralTransformPrecessingObsoleteInitialConditions ( REAL8 incl,
REAL8 S1x,
REAL8 S1y,
REAL8 S1z,
REAL8 S2x,
REAL8 S2y,
REAL8 S2z,
REAL8  thetaJN,
REAL8  phiJL,
REAL8  theta1,
REAL8  theta2,
REAL8  phi12,
REAL8  chi1,
REAL8  chi2,
REAL8  m1,
REAL8  m2,
REAL8  fRef 
)

Function to specify the desired orientation of a precessing binary in terms of several angles and then compute the vector components in the so-called "radiation frame" (with the z-axis along the direction of propagation) as needed to specify binary configuration for ChooseTDWaveform.

Input: thetaJN is the inclination between total angular momentum (J) and the direction of propagation (N) theta1 and theta2 are the inclinations of S1 and S2 measured from the Newtonian orbital angular momentum (L_N) phi12 is the difference in azimuthal angles of S1 and S2. chi1, chi2 are the dimensionless spin magnitudes ( \(0 \le chi1,2 \le 1\)) phiJL is the azimuthal angle of L_N on its cone about J. m1, m2, f_ref are the component masses and reference GW frequency, they are needed to compute the magnitude of L_N, and thus J.

Output: incl - inclination angle of L_N relative to N x, y, z components S1 and S2 (unit spin vectors times their dimensionless spin magnitudes - i.e. they have unit magnitude for extremal BHs and smaller magnitude for slower spins).

NOTE: Here the "total" angular momentum is computed as J = L_N + S1 + S2 where L_N is the Newtonian orbital angular momentum. In fact, there are PN corrections to L which contribute to J that are NOT ACCOUNTED FOR in this function. This is done so the function does not need to know about the PN order of the system and to avoid subtleties with spin-orbit contributions to L. Also, it is believed that the difference in Jhat with or without these PN corrections to L is quite small.

NOTE: fRef = 0 is not a valid choice. If you will pass fRef=0 into ChooseWaveform, then here pass in f_min, the starting GW frequency

The various rotations in this transformation are described in more detail in a Mathematica notebook available here: https://www.lsc-group.phys.uwm.edu/ligovirgo/cbcnote/Waveforms/TransformPrecessingInitialConditions

!!!UNREVIEWED!!!

Parameters
inclInclination angle of L_N (returned)
S1xS1 x component (returned)
S1yS1 y component (returned)
S1zS1 z component (returned)
S2xS2 x component (returned)
S2yS2 y component (returned)
S2zS2 z component (returned)
thetaJNzenith angle between J and N (rad)
phiJLazimuthal angle of L_N on its cone about J (rad)
theta1zenith angle between S1 and LNhat (rad)
theta2zenith angle between S2 and LNhat (rad)
phi12difference in azimuthal angle btwn S1, S2 (rad)
chi1dimensionless spin of body 1
chi2dimensionless spin of body 2
m1mass of body 1 (kg)
m2mass of body 2 (kg)
fRefreference GW frequency (Hz)

Definition at line 3829 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralInitialConditionsPrecessingApproxs()

int XLALSimInspiralInitialConditionsPrecessingApproxs ( REAL8 inc,
REAL8 S1x,
REAL8 S1y,
REAL8 S1z,
REAL8 S2x,
REAL8 S2y,
REAL8 S2z,
const REAL8  inclIn,
const REAL8  S1xIn,
const REAL8  S1yIn,
const REAL8  S1zIn,
const REAL8  S2xIn,
const REAL8  S2yIn,
const REAL8  S2zIn,
const REAL8  m1,
const REAL8  m2,
const REAL8  fRef,
const REAL8  phiRef,
LALSimInspiralFrameAxis  frameChoice 
)

Function to specify the desired orientation of the spin components of a precessing binary.

Input: Depending on the values of the variable frameChoice input parameters may represent different quantities. incl is the angle between in the X-Y-Z frame specified by

  • frameChoice = LAL_SIM_INSPIRAL_FRAME_AXIS_TOTAL_J J // Z, N = (0, sin(inc), cos(inc))
  • frameChoice = LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L (default) L // Z, N = (0, sin(inc), cos(inc))
  • frameChoice = LAL_SIM_INSPIRAL_FRAME_AXIS_VIEW
  • N // Z, L = (sin(inc), 0, cos(inc))

Cartesian components of S1 and S2 are given wrt axes x-y-Z (ORBITAL_L and VIEW cases), being x the vector pointing from body 2 to body 1. x-y axes are rotated wrt X-Y by phiRef counterclockwise, i.e. if S1=(a,b,0) in fame x-y-Z, in frame X-Y-Z it will have components S1'=(a cos(phiRef) - b sin(phiRef), a sin(phiRef) + b cos(phiRef), 0). In the TOTAL_J case spin components are given wrt J, with the x-Z plane spanned by the line connecting the two bodies and J.

m1, m2, f_ref are the component masses and reference GW frequency, they are needed to compute the magnitude of L_N and. phi_ref is needed since the spin components in the plane of the orbit are defined wrt assuming the x axis is parallel to the line joining the binary consistuents, from body 2 to body 1.

Output: x, y, z components S1 and S2 wrt N inc angle between L and N, i.e. N=(0,0,1), LN=(sin(inc),0,cos(inc)) as required by precessing waveforms routines in SpinTaylor and IRMPhenomP codes.

NOTE: Here the "total" angular momentum is computed as J = L_N(1+l_1PN) + S1 + S2 where L_N(1+l_1PN) is the 1PN-corrected orbital angular momentum, which is parallel to Newtonian angular momentum. PN corrections to L from 1.5PN order on (wrt to Newtonian value) are NOT ACCOUNTED FOR in this function.

Parameters
incinclination angle (returned)
S1xS1x component (returned)
S1yS1y component (returned)
S1zS1z component (returned)
S2xS2x components (returned)
S2yS2y components (returned)
S2zS2z components (returned)
inclInInclination angle in input
S1xInS1 x component
S1yInS1 y component
S1zInS1 z component
S2xInS2 x component
S2yInS2 y component
S2zInS2 z component
m1mass of body 1 (kg)
m2mass of body 2 (kg)
fRefreference GW frequency (Hz)
phiRefreference GW phase
frameChoiceFlag to identify axis wrt which spin components are given. Pass in NULL (or None in python) for default (OrbitalL)

Definition at line 4010 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorPNEvolveOrbit()

int XLALSimInspiralSpinTaylorPNEvolveOrbit ( REAL8TimeSeries **  V,
REAL8TimeSeries **  Phi,
REAL8TimeSeries **  S1x,
REAL8TimeSeries **  S1y,
REAL8TimeSeries **  S1z,
REAL8TimeSeries **  S2x,
REAL8TimeSeries **  S2y,
REAL8TimeSeries **  S2z,
REAL8TimeSeries **  LNhatx,
REAL8TimeSeries **  LNhaty,
REAL8TimeSeries **  LNhatz,
REAL8TimeSeries **  E1x,
REAL8TimeSeries **  E1y,
REAL8TimeSeries **  E1z,
const REAL8  deltaT,
const REAL8  m1_SI,
const REAL8  m2_SI,
const REAL8  fStart,
const REAL8  fEnd,
const REAL8  s1x,
const REAL8  s1y,
const REAL8  s1z,
const REAL8  s2x,
const REAL8  s2y,
const REAL8  s2z,
const REAL8  lnhatx,
const REAL8  lnhaty,
const REAL8  lnhatz,
const REAL8  e1x,
const REAL8  e1y,
const REAL8  e1z,
const REAL8  lambda1,
const REAL8  lambda2,
const REAL8  quadparam1,
const REAL8  quadparam2,
const LALSimInspiralSpinOrder  spinO,
const LALSimInspiralTidalOrder  tideO,
const INT4  phaseO,
const INT4  lscorr,
const Approximant  approx 
)

This function evolves the orbital equations for a precessing binary using the "TaylorT1/T5/T4" approximant for solving the orbital dynamics (see arXiv:0907.0700 for a review of the various PN approximants).

It returns time series of the "orbital velocity", orbital phase, and components for both individual spin vectors, the "Newtonian" orbital angular momentum (which defines the instantaneous plane) and "E1", a basis vector in the instantaneous orbital plane. Note that LNhat and E1 completely specify the instantaneous orbital plane. It also returns the time and phase of the final time step

For input, the function takes the two masses, the initial orbital phase, Values of S1, S2, LNhat, E1 vectors at starting time, the desired time step size, the starting GW frequency, and PN order at which to evolve the phase,

NOTE: All vectors are given in the frame where the z-axis is set by the angular momentum at reference frequency, the x-axis is chosen orthogonal to it, and the y-axis is given by the RH rule. Initial values must be passed in this frame, and the time series of the vector components will also be returned in this frame.

Review completed on git hash ...

Parameters
Vpost-Newtonian parameter [returned]
Phiorbital phase [returned]
S1xSpin1 vector x component [returned]
S1y" " " y component [returned] @param S1z " " " z component [returned]
S2xSpin2 vector x component [returned]
S2y" " " y component [returned] @param S2z " " " z component [returned]
LNhatxunit orbital ang. mom. x [returned]
LNhaty" " " y component [returned] @param LNhatz " " " z component [returned]
E1xorb. plane basis vector x[returned]
E1y" " " y component [returned] @param E1z " " " z component [returned]
deltaTsampling interval (s)
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
fStartstarting GW frequency
fEndending GW frequency, fEnd=0 means integrate as far forward as possible
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
spinOtwice PN order of spin effects
tideOtwice PN order of tidal effects
phaseOtwice post-Newtonian order
lscorrflag to control L_S terms
approxPN approximant (SpinTaylorT1/T5/T4)

Definition at line 4270 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorPNEvolveOrbitOnlyFinal()

int XLALSimInspiralSpinTaylorPNEvolveOrbitOnlyFinal ( REAL8TimeSeries **  V,
REAL8TimeSeries **  Phi,
REAL8TimeSeries **  S1x,
REAL8TimeSeries **  S1y,
REAL8TimeSeries **  S1z,
REAL8TimeSeries **  S2x,
REAL8TimeSeries **  S2y,
REAL8TimeSeries **  S2z,
REAL8TimeSeries **  LNhatx,
REAL8TimeSeries **  LNhaty,
REAL8TimeSeries **  LNhatz,
REAL8TimeSeries **  E1x,
REAL8TimeSeries **  E1y,
REAL8TimeSeries **  E1z,
const REAL8  deltaT,
const REAL8  m1_SI,
const REAL8  m2_SI,
const REAL8  fStart,
const REAL8  fEnd,
const REAL8  s1x,
const REAL8  s1y,
const REAL8  s1z,
const REAL8  s2x,
const REAL8  s2y,
const REAL8  s2z,
const REAL8  lnhatx,
const REAL8  lnhaty,
const REAL8  lnhatz,
const REAL8  e1x,
const REAL8  e1y,
const REAL8  e1z,
const REAL8  lambda1,
const REAL8  lambda2,
const REAL8  quadparam1,
const REAL8  quadparam2,
const LALSimInspiralSpinOrder  spinO,
const LALSimInspiralTidalOrder  tideO,
const INT4  phaseO,
const INT4  lscorr,
const Approximant  approx 
)

This function evolves the orbital equations for a precessing binary using the "TaylorT1/T5/T4" approximant for solving the orbital dynamics (see arXiv:0907.0700 for a review of the various PN approximants).

It returns the final values of the "orbital velocity", orbital phase, and components for both individual spin vectors, the "Newtonian" orbital angular momentum (which defines the instantaneous plane) and "E1", a basis vector in the instantaneous orbital plane. Note that LNhat and E1 completely specify the instantaneous orbital plane. It also returns the time and phase of the final time step

For input, the function takes the two masses, the initial orbital phase, Values of S1, S2, LNhat, E1 vectors at starting time, the desired time step size, the starting GW frequency, and PN order at which to evolve the phase,

NOTE: All vectors are given in the frame where the z-axis is set by the angular momentum at reference frequency, the x-axis is chosen orthogonal to it, and the y-axis is given by the RH rule. Initial values must be passed in this frame, and the time series of the vector components will also be returned in this frame.

Review completed on git hash ...

Parameters
Vpost-Newtonian parameter [returned]
Phiorbital phase [returned]
S1xSpin1 vector x component [returned]
S1y" " " y component [returned] @param S1z " " " z component [returned]
S2xSpin2 vector x component [returned]
S2y" " " y component [returned] @param S2z " " " z component [returned]
LNhatxunit orbital ang. mom. x [returned]
LNhaty" " " y component [returned] @param LNhatz " " " z component [returned]
E1xorb. plane basis vector x[returned]
E1y" " " y component [returned] @param E1z " " " z component [returned]
deltaTsampling interval (s)
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
fStartstarting GW frequency
fEndending GW frequency, fEnd=0 means integrate as far forward as possible
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
spinOtwice PN order of spin effects
tideOtwice PN order of tidal effects
phaseOtwice post-Newtonian order
lscorrflag to control L_S terms
approxPN approximant (SpinTaylorT1/T5/T4)

Definition at line 4617 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT4OLD()

int XLALSimInspiralSpinTaylorT4OLD ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8 UNUSED  v0,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
REAL8  lambda1,
REAL8  lambda2,
REAL8  quadparam1,
REAL8  quadparam2,
LALDict *  LALparams,
int  phaseO,
int  amplitudeO 
)

Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T4" method.

This routine allows the user to specify different pN orders for the phasing and amplitude of the waveform.

The reference frequency fRef is used as follows: 1) if fRef = 0: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. The orbital phase of the last sample is set to phiRef (i.e. phiRef is the "coalescence phase", roughly speaking). THIS IS THE DEFAULT BEHAVIOR CONSISTENT WITH OTHER APPROXIMANTS

2) If fRef = fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. phiRef is used to set the orbital phase of the first sample at fStart.

3) If fRef > fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fRef. phiRef is used to set the orbital phase at fRef. The code will integrate forwards and backwards from fRef and stitch the two together to create a complete waveform. This allows one to specify the orientation of the binary in-band (or at any arbitrary point). Otherwise, the user can only directly control the initial orientation.

4) fRef < 0 or fRef >= Schwarz. ISCO are forbidden and the code will abort.

REVIEW completed on git hash 6640e79e60791d5230731acc63351676ce7ce413

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
v0tail gauge term (default = 1)
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
LALparamsLAL dictionary containing accessory parameters
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order

Definition at line 4904 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT4()

int XLALSimInspiralSpinTaylorT4 ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
LALDict *  LALparams 
)
Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
LALparamsLAL dictionary containing accessory parameters

Definition at line 4951 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT1OLD()

int XLALSimInspiralSpinTaylorT1OLD ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8 UNUSED  v0,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
REAL8  lambda1,
REAL8  lambda2,
REAL8  quadparam1,
REAL8  quadparam2,
LALDict *  LALparams,
int  phaseO,
int  amplitudeO 
)

Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T1" method.

This routine allows the user to specify different pN orders for the phasing and amplitude of the waveform.

The reference frequency fRef is used as follows: 1) if fRef = 0: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. The orbital phase of the last sample is set to phiRef (i.e. phiRef is the "coalescence phase", roughly speaking). THIS IS THE DEFAULT BEHAVIOR CONSISTENT WITH OTHER APPROXIMANTS

2) If fRef = fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. phiRef is used to set the orbital phase of the first sample at fStart.

3) If fRef > fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fRef. phiRef is used to set the orbital phase at fRef. The code will integrate forwards and backwards from fRef and stitch the two together to create a complete waveform. This allows one to specify the orientation of the binary in-band (or at any arbitrary point). Otherwise, the user can only directly control the initial orientation.

4) fRef < 0 or fRef >= Schwarz. ISCO are forbidden and the code will abort.

REVIEW completed on git hash 6640e79e60791d5230731acc63351676ce7ce413

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
v0tail gauge term (default = 1)
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
LALparamsLAL dictionary containing accessory parameters
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order

Definition at line 5014 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT1()

int XLALSimInspiralSpinTaylorT1 ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
LALDict *  LALparams 
)
Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
LALparamsLAL dictionary containing accessory parameters

Definition at line 5061 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT5()

int XLALSimInspiralSpinTaylorT5 ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
LALDict *  LALparams 
)

Driver routine to compute a precessing post-Newtonian inspiral waveform with phasing computed from energy balance using the so-called "T5" method.

This routine allows the user to specify different pN orders for the phasing, amplitude, spin and tidal effects of the waveform.

The reference frequency fRef is used as follows: 1) if fRef = 0: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. The orbital phase of the last sample is set to phiRef (i.e. phiRef is the "coalescence phase", roughly speaking). THIS IS THE DEFAULT BEHAVIOR CONSISTENT WITH OTHER APPROXIMANTS

2) If fRef = fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. phiRef is used to set the orbital phase of the first sample at fStart.

3) If fRef > fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fRef. phiRef is used to set the orbital phase at fRef. The code will integrate forwards and backwards from fRef and stitch the two together to create a complete waveform. This allows one to specify the orientation of the binary in-band (or at any arbitrary point). Otherwise, the user can only directly control the initial orientation.

4) fRef < 0 or fRef >= Schwarz. ISCO are forbidden and the code will abort.

REVIEW completed on git hash ...

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
LALparamsLAL dictionary containing accessory parameters

Definition at line 5124 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralPrecessingPTFQWaveforms()

int XLALSimInspiralPrecessingPTFQWaveforms ( REAL8TimeSeries **  Q1,
REAL8TimeSeries **  Q2,
REAL8TimeSeries **  Q3,
REAL8TimeSeries **  Q4,
REAL8TimeSeries **  Q5,
REAL8TimeSeries V,
REAL8TimeSeries Phi,
REAL8TimeSeries S1x,
REAL8TimeSeries S1y,
REAL8TimeSeries S1z,
REAL8TimeSeries S2x,
REAL8TimeSeries S2y,
REAL8TimeSeries S2z,
REAL8TimeSeries LNhatx,
REAL8TimeSeries LNhaty,
REAL8TimeSeries LNhatz,
REAL8TimeSeries E1x,
REAL8TimeSeries E1y,
REAL8TimeSeries E1z,
REAL8  m1,
REAL8  m2,
REAL8  r 
)

Compute the physical template family "Q" vectors for a spinning, precessing binary when provided time series of all the dynamical quantities.

These vectors always supplied to dominant order.

Based on Pan, Buonanno, Chan and Vallisneri PRD69 104017, (see also theses of Diego Fazi and Ian Harry)

NOTE: The vectors MUST be given in the so-called radiation frame where Z is the direction of propagation, X is the principal '+' axis and Y = Z x X

!!!UNREVIEWED!!!

Parameters
Q1PTF-Q1 waveform [returned]
Q2PTF-Q2 waveform [returned]
Q3PTF-Q2 waveform [returned]
Q4PTF-Q2 waveform [returned]
Q5PTF-Q2 waveform [returned]
Vpost-Newtonian parameter
Phiorbital phase
S1xSpin1 vector x component
S1ySpin1 vector y component
S1zSpin1 vector z component
S2xSpin2 vector x component
S2ySpin2 vector y component
S2zSpin2 vector z component
LNhatxunit orbital ang. mom. x comp.
LNhatyunit orbital ang. mom. y comp.
LNhatzunit orbital ang. mom. z comp.
E1xorbital plane basis vector x comp.
E1yorbital plane basis vector y comp.
E1zorbital plane basis vector z comp.
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
rdistance of source (m)

Definition at line 5171 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT4PTFQVecs()

int XLALSimInspiralSpinTaylorT4PTFQVecs ( REAL8TimeSeries **  Q1,
REAL8TimeSeries **  Q2,
REAL8TimeSeries **  Q3,
REAL8TimeSeries **  Q4,
REAL8TimeSeries **  Q5,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  chi1,
REAL8  kappa1,
REAL8  fStart,
REAL8  lambda1,
REAL8  lambda2,
LALSimInspiralSpinOrder  spinO,
LALSimInspiralTidalOrder  tideO,
int  phaseO 
)

Driver routine to compute the physical template family "Q" vectors using the "T4" method.

Note that PTF describes single spin systems

This routine requires leading-order amplitude dependence but allows the user to specify the phase PN order

!!!UNREVIEWED!!!

Parameters
Q1Q1 output vector
Q2Q2 output vector
Q3Q3 output vector
Q4Q4 output vector
Q5Q5 output vector
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
chi1spin magnitude (|S1|)
kappa1L . S1 (1 if they are aligned)
fStartstart GW frequency (Hz)
lambda1(tidal deformability of mass 1) / (mass of mody 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
spinOtwice PN order of spin effects
tideOtwice PN order of tidal effects
phaseOtwice PN phase order

Definition at line 5328 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT4Fourier()

int XLALSimInspiralSpinTaylorT4Fourier ( COMPLEX16FrequencySeries **  hplus,
COMPLEX16FrequencySeries **  hcross,
REAL8  fMin,
REAL8  fMax,
REAL8  deltaF,
INT4  kMax,
REAL8  phiRef,
REAL8  v0,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
REAL8  lambda1,
REAL8  lambda2,
REAL8  quadparam1,
REAL8  quadparam2,
LALDict *  LALparams,
INT4  phaseO,
INT4  amplitudeO,
INT4  phiRefAtEnd 
)

Driver routine to compute a precessing post-Newtonian inspiral waveform in the Fourier domain with phasing computed from energy balance using the so-called "T4" method.

This routine allows the user to specify different pN orders for the phasing and amplitude of the waveform.

The reference frequency fRef is used as follows: 1) if fRef = 0: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. The orbital phase of the last sample is set to phiRef (i.e. phiRef is the "coalescence phase", roughly speaking). THIS IS THE DEFAULT BEHAVIOR CONSISTENT WITH OTHER APPROXIMANTS

2) If fRef = fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. phiRef is used to set the orbital phase of the first sample at fStart.

3) If fRef > fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fRef. phiRef is used to set the orbital phase at fRef. The code will integrate forwards and backwards from fRef and stitch the two together to create a complete waveform. This allows one to specify the orientation of the binary in-band (or at any arbitrary point). Otherwise, the user can only directly control the initial orientation.

4) fRef < 0 or fRef >= Schwarz. ISCO are forbidden and the code will abort.

It is recommended, but not necessary to set fStart slightly smaller than fMin, e.g. fStart = 9.5 for fMin = 10.

The returned Fourier series are set so that the Schwarzschild ISCO frequency corresponds to t = 0 as closely as possible.

!!!UNREVIEWED!!!

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
fMinminimum frequency of the returned series
fMaxmaximum frequency of the returned series
deltaFfrequency interval of the returned series
kMaxk_max as described in defined arXiv: 1408.5158 (min 0, max 10).
phiReforbital phase at reference pt.
v0tail gauge term (default = 1)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
LALparamsLAL dictionary containing accessory parameters
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order
phiRefAtEndwhether phiRef corresponds to the end of the inspiral

Definition at line 5434 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorT5Fourier()

int XLALSimInspiralSpinTaylorT5Fourier ( COMPLEX16FrequencySeries **  hplus,
COMPLEX16FrequencySeries **  hcross,
REAL8  fMin,
REAL8  fMax,
REAL8  deltaF,
INT4  kMax,
REAL8  phiRef,
REAL8  v0,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  fRef,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhatx,
REAL8  lnhaty,
REAL8  lnhatz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
REAL8  lambda1,
REAL8  lambda2,
REAL8  quadparam1,
REAL8  quadparam2,
LALDict *  LALparams,
INT4  phaseO,
INT4  amplitudeO,
INT4  phiRefAtEnd 
)

Driver routine to compute a precessing post-Newtonian inspiral waveform in the Fourier domain with phasing computed from energy balance using the so-called "T2" method see arXiv: 0907.0700 for its definition, but in its "T5" variant, see arXiv: 1107.1267.

This routine allows the user to specify different pN orders for the phasing and amplitude of the waveform.

The reference frequency fRef is used as follows: 1) if fRef = 0: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. The orbital phase of the last sample is set to phiRef (i.e. phiRef is the "coalescence phase", roughly speaking). THIS IS THE DEFAULT BEHAVIOR CONSISTENT WITH OTHER APPROXIMANTS

2) If fRef = fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fStart. phiRef is used to set the orbital phase of the first sample at fStart.

3) If fRef > fStart: The initial values of s1, s2, lnhat and e1 will be the values at frequency fRef. phiRef is used to set the orbital phase at fRef. The code will integrate forwards and backwards from fRef and stitch the two together to create a complete waveform. This allows one to specify the orientation of the binary in-band (or at any arbitrary point). Otherwise, the user can only directly control the initial orientation.

4) fRef < 0 or fRef >= Schwarz. ISCO are forbidden and the code will abort.

It is recommended, but not necessary to set fStart slightly smaller than fMin, e.g. fStart = 9.5 for fMin = 10.

The returned Fourier series are set so that the Schwarzschild ISCO frequency corresponds to t = 0 as closely as possible.

!!!UNREVIEWED!!!

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
fMinminimum frequency of the returned series
fMaxmaximum frequency of the returned series
deltaFfrequency interval of the returned series
kMaxk_max as described in arXiv: 1408.5158 (min 0, max 10).
phiReforbital phase at reference pt.
v0tail gauge term (default = 1)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
fRefreference GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
e1xinitial value of E1x
e1yinitial value of E1y
e1zinitial value of E1z
lambda1(tidal deformability of mass 1) / (mass of body 1)^5 (dimensionless)
lambda2(tidal deformability of mass 2) / (mass of body 2)^5 (dimensionless)
quadparam1phenom. parameter describing induced quad. moment of body 1 (=1 for BHs, ~2-12 for NSs)
quadparam2phenom. parameter describing induced quad. moment of body 2 (=1 for BHs, ~2-12 for NSs)
LALparamsLAL dictionary containing accessory parameters
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order
phiRefAtEndwhether phiRef corresponds to the end of the inspiral

Definition at line 5516 of file LALSimInspiralSpinTaylor.c.

◆ XLALSimInspiralSpinTaylorF2()

int XLALSimInspiralSpinTaylorF2 ( COMPLEX16FrequencySeries **  hplus_out,
COMPLEX16FrequencySeries **  hcross_out,
const REAL8  phi_ref,
const REAL8  deltaF,
const REAL8  m1_SI,
const REAL8  m2_SI,
const REAL8  s1x,
const REAL8  s1y,
const REAL8  s1z,
const REAL8  lnhatx,
const REAL8  lnhaty,
const REAL8  lnhatz,
const REAL8  fStart,
const REAL8  fEnd,
const REAL8  f_ref,
const REAL8  r,
LALDict *  moreParams,
const INT4  phaseO,
const INT4  amplitudeO 
)

Computes the stationary phase approximation to the Fourier transform of a chirp waveform with phase given by Eq.

 \eqref{eq_InspiralFourierPhase_f2} and amplitude given by expanding \(1/\sqrt{\dot{F}}\). If the PN order is set to -1, then the highest implemented order is used.

See arXiv:0810.5336 and arXiv:astro-ph/0504538 for spin corrections to the phasing. See arXiv:1303.7412 for spin-orbit phasing corrections at 3 and 3.5PN order

Parameters
hplus_outFD hplus waveform
hcross_outFD hcross waveform
phi_refreference orbital phase (rad)
deltaFfrequency resolution
m1_SImass of companion 1 (kg)
m2_SImass of companion 2 (kg)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
lnhatxinitial value of LNhatx
lnhatyinitial value of LNhaty
lnhatzinitial value of LNhatz
fStartstart GW frequency (Hz)
fEndhighest GW frequency (Hz) of waveform generation - if 0, end at Schwarzschild ISCO
f_refReference GW frequency (Hz) - if 0 reference point is coalescence
rdistance of source (m)
moreParamsLinked list of extra. Pass in NULL (or None in python) for standard waveform. Set "sideband",m to get a single sideband (m=-2..2)
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order

Definition at line 289 of file LALSimInspiralSpinTaylorF2.c.

◆ XLALSimInspiralSpinTaylorT5duplicate()

int XLALSimInspiralSpinTaylorT5duplicate ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  phiRef,
REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  fStart,
REAL8  r,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  incAngle,
int  phaseO,
int  amplitudeO 
)

Generate time-domain generic spinning PN waveforms in the SpinTaylorT5 approximaton.

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
phiReforbital phase at reference pt.
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
fStartstart GW frequency (Hz)
rdistance of source (m)
s1xinitial value of S1x
s1yinitial value of S1y
s1zinitial value of S1z
s2xinitial value of S2x
s2yinitial value of S2y
s2zinitial value of S2z
incAngleinclination angle with J_ini
phaseOtwice PN phase order
amplitudeOtwice PN amplitude order

Definition at line 166 of file LALSimInspiralSpinTaylorT5duplicate.c.