23#include <gsl/gsl_odeiv.h>
25#include <lal/LALStdlib.h>
26#include <lal/LALConstants.h>
27#include <lal/SimulateCoherentGW.h>
28#include <lal/GeneratePPNInspiral.h>
29#include <lal/LIGOMetadataTables.h>
30#include <lal/LALDatatypes.h>
31#include <lal/LALSimInspiral.h>
56#define LALINSPIRALH_ENULL 1
57#define LALINSPIRALH_EMEM 2
58#define LALINSPIRALH_EDIV0 3
59#define LALINSPIRALH_ESIZE 4
60#define LALINSPIRALH_ECHOICE 5
61#define LALINSPIRALH_EORDER 6
62#define LALINSPIRALH_EAPPROXIMANT 7
63#define LALINSPIRALH_EPSI0 8
64#define LALINSPIRALH_EPSI3 9
65#define LALINSPIRALH_EALPHA 10
66#define LALINSPIRALH_EFCUTOFF 11
67#define LALINSPIRALH_ENOWAVEFORM 12
68#define LALINSPIRALH_ESTOPPED 13
69#define LALINSPIRALH_EROOTINIT 14
70#define LALINSPIRALH_EFLOWER 15
71#define LALINSPIRALH_EVECTOR 16
72#define LALINSPIRALH_EFLOWERINJ 17
73#define LALINSPIRALH_EORDERMISSING 18
74#define LALINSPIRALH_EBPERR 19
75#define LALINSPIRALH_ESWITCH 20
76#define LALINSPIRALH_EMASSCHOICE 21
80#define LALINSPIRALH_MSGENULL "Arguments contained an unexpected null pointer"
81#define LALINSPIRALH_MSGEMEM "Memory allocation error"
82#define LALINSPIRALH_MSGEDIV0 "Division by zero"
83#define LALINSPIRALH_MSGESIZE "Invalid input range"
84#define LALINSPIRALH_MSGECHOICE "Invalid choice for an input parameter"
85#define LALINSPIRALH_MSGEORDER "unknown order specified"
86#define LALINSPIRALH_MSGEAPPROXIMANT "Invalid model"
87#define LALINSPIRALH_MSGEPSI0 "psi0 must be > 0"
88#define LALINSPIRALH_MSGEPSI3 "psi3 must be < 0"
89#define LALINSPIRALH_MSGEALPHA "alpha must be defined positive"
90#define LALINSPIRALH_MSGEFCUTOFF "fcutoff must be defined and > 0"
91#define LALINSPIRALH_MSGENOWAVEFORM "No Waveform generated"
92#define LALINSPIRALH_MSGESTOPPED "Waveform generation stopped"
93#define LALINSPIRALH_MSGEROOTINIT "Can't find good bracket for BisectionFindRoot"
94#define LALINSPIRALH_MSGEFLOWER "fLower too low in comparison to flso"
95#define LALINSPIRALH_MSGEVECTOR "Attempting to write beyond the end of vector"
96#define LALINSPIRALH_MSGEFLOWERINJ "flower for the injection must be greater than zero"
97#define LALINSPIRALH_MSGEORDERMISSING "The PN order requested is not implemented for this approximant"
98#define LALINSPIRALH_MSGEBPERR "Error in band passing signal"
99#define LALINSPIRALH_MSGESWITCH "Unknown case in switch"
100#define LALINSPIRALH_MSGEMASSCHOICE "Improper choice for massChoice"
103#define LAL_INSPIRAL_INTERACTION_DEFAULT LAL_INSPIRAL_INTERACTION_ALL
116typedef enum tagLALInspiralInteraction {
178typedef enum tagInputMasses {
347 struct tagInspiralTemplate *
next;
348 struct tagInspiralTemplate *
fine;
424 REAL8 FTaN,
FTa1, FTa2, FTa3, FTa4, FTa5, FTa6, FTa7, FTa8, FTl6, FTl8;
429 REAL8 fTaN,
fTa1, fTa2, fTa3, fTa4, fTa5, fTa6, fTa7, fTa8;
434 REAL8 fPaN,
fPa1, fPa2, fPa3, fPa4, fPa5, fPa6, fPa7, fPa8;
459 REAL8 pfaN,
pfa2, pfa3, pfa4, pfa5, pfa6, pfa7, pfl5, pfl6;
606tagInspiralDerivativesIn
638 const gsl_odeiv_step_type *
type;
1212REAL8 XLALInspiralPhasing2_1PN (
1249REAL8 XLALInspiralPhasing3_1PN (
1293REAL8 XLALInspiralTiming2_1PN (
1327void LALInspiralFrequency3_1PN (
1333REAL8 XLALInspiralFrequency3_1PN (
REAL8 XLALInspiralPhasing2_5PN(REAL8 v, expnCoeffs *ak)
int XLALInspiralTDWaveformFromSimInspiral(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, SimInspiralTable *thisRow, REAL8 deltaT)
Generate the plus and cross polarizations for a conditioned waveform form a row of the sim_inspiral t...
REAL8 XLALInspiralFrequency3_3PN(REAL8 td, expnCoeffs *ak)
REAL8 XLALInspiralPhasing2_7PN(REAL8 v, expnCoeffs *ak)
void LALBBHPhenWaveTimeDom(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *insp_template)
REAL8 XLALInspiralPhasing1(REAL8 v, void *params)
void LALInspiralWave3(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
void LALInspiralWaveTemplates(LALStatus *status, REAL4Vector *filter1, REAL4Vector *filter2, InspiralTemplate *params)
REAL8 XLALEtaTau04(REAL8 eta, void *in)
REAL8 XLALInspiralPhasing2_0PN(REAL8 v, expnCoeffs *ak)
REAL8 XLALInspiralPhasing3_3PN(REAL8 td, expnCoeffs *ak)
void LALBCVWaveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralParameterCalc(InspiralTemplate *params)
REAL4 LALInspiralHCrossPolarization(REAL8 phase, REAL8 v, InspiralTemplate *params)
int XLALGetAdaptiveIntFromString(const CHAR *inString)
XLAL function to determine adaptive integration flag from a string.
int XLALEOBPPWaveformForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
void XLALRungeKutta4Free(rk4GSLIntegrator *integrator)
void LALSTPNWaveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
REAL8 XLALInspiralFrequency3_4PN(REAL8 td, expnCoeffs *ak)
int XLALInspiralWave2ForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALInspiralRestrictedAmplitude(InspiralTemplate *params)
int XLALBBHPhenWaveTimeDomTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *insp_template)
int XLALInspiralGetApproximantString(CHAR *output, UINT4 length, Approximant approx, LALPNOrder order)
int XLALInspiralStationaryPhaseApprox2(REAL4Vector *signalvec, InspiralTemplate *params)
INT4 XLALPSpinGenerateQNMFreq(COMPLEX8Vector *modefreqs, InspiralTemplate *params, UINT4 l, INT4 m, UINT4 nmodes, REAL8 finalMass, REAL8 finalSpin)
INT4 XLALInspiralAttachRingdownWave(REAL4Vector *Omega, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
REAL8 XLALInspiralPhasing3_6PN(REAL8 td, expnCoeffs *ak)
INT4 XLALGenerateQNMFreqV2(COMPLEX8Vector *modefreqs, InspiralTemplate *params, UINT4 l, UINT4 m, UINT4 nmodes)
As with the above function, this generates the quasinormal mode frequencies for a black hole ringdown...
int XLALSTPNFramelessWaveform(REAL4Vector *signalvec, InspiralTemplate *params)
void LALInspiralEccentricity(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
void LALBBHPhenWaveTimeDomTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *insp_template)
REAL8 XLALInspiralPhasing2_4PN(REAL8 v, expnCoeffs *ak)
void LALInspiralSetup(LALStatus *status, expnCoeffs *ak, InspiralTemplate *params)
int XLALSimInspiralChooseWaveformFromInspiralTemplate(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, InspiralTemplate *params)
Generate the plus and cross polarizations for a waveform form a row of the InspiralTemplate structure...
void LALInspiralAmplitudeCorrectedWaveTemplates(LALStatus *status, REAL4Vector *filter1, REAL4Vector *filter2, InspiralTemplate *params)
void LALInspiralWaveLength(LALStatus *status, UINT4 *n, InspiralTemplate params)
int XLALInspiralStationaryPhaseApprox1(REAL4Vector *signalvec, InspiralTemplate *params)
void LALInspiralSpinModulatedWave(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *in)
void LALSTPNFramelessWaveformForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALInspiralWave1(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralWave1Templates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALInspiralAmplitudeCorrectedWaveForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
REAL8 XLALInspiralPhasing3_0PN(REAL8 td, expnCoeffs *ak)
int XLALBBHPhenWaveAFreqDomTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
int XLALInspiralWave2Templates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
REAL8 XLALInspiralFrequency3_0PN(REAL8 td, expnCoeffs *ak)
int XLALTaylorNWaveform(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralWave1ForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALGetInspiralOnlyFromString(const CHAR *inString)
XLAL function to determine inspiral-only flag from a string.
void LALSTPNWaveformEngine(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, REAL4Vector *a, REAL4Vector *ff, REAL8Vector *phi, REAL4Vector *shift, UINT4 *countback, InspiralTemplate *params, InspiralInit *paramsInit)
int XLALBBHPhenWaveBFreqDom(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralWave2(REAL4Vector *signalvec, InspiralTemplate *params)
void LALInspiralAmplitudeCorrectedWave(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
void LALInspiralWaveForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALNRInjectionFromSimInspiral(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, SimInspiralTable *thisRow, REAL8 deltaT)
int XLALInspiralSetup(expnCoeffs *ak, InspiralTemplate *params)
INT4 XLALGenerateHybridWaveDerivatives(REAL4Vector *rwave, REAL4Vector *dwave, REAL4Vector *ddwave, REAL8Vector *timeVec, REAL4Vector *wave, REAL8Vector *matchrange, InspiralTemplate *params)
int XLALInspiralGenerateIIRSetFourierTransform(int j, int jmax, COMPLEX16 a1, COMPLEX16 b0, int delay, COMPLEX16 *hfcos, COMPLEX16 *hfsin)
REAL8 XLALInspiralPhasing2_2PN(REAL8 v, expnCoeffs *ak)
INT4 XLALInspiralHybridRingdownWave(REAL4Vector *rdwave1, REAL4Vector *rdwave2, InspiralTemplate *params, REAL4VectorSequence *inspwave1, REAL4VectorSequence *inspwave2, COMPLEX8Vector *modefreqs, REAL8Vector *matchrange)
int XLALBandPassInspiralTemplate(REAL4Sequence *sequence, REAL4 fLow, REAL4 fHigh, REAL4 fSampling)
void LALTaylorNWaveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
int XLALBBHPhenWaveBFreqDomTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALInspiralInit(LALStatus *status, InspiralTemplate *params, InspiralInit *paramsInit)
REAL8 XLALInspiralPhasing2_3PN(REAL8 v, expnCoeffs *ak)
REAL8 XLALInspiralPhasing3_7PN(REAL8 td, expnCoeffs *ak)
int XLALInspiralInit(InspiralTemplate *params, InspiralInit *paramsInit)
REAL8 XLALInspiralPhasing2_6PN(REAL8 v, expnCoeffs *ak)
REAL8 XLALInspiralTiming2_5PN(REAL8 f, void *params)
void LALInspiralDerivatives(REAL8Vector *vec1, REAL8Vector *vec2, void *params)
int XLALSTPNFramelessWaveformForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALEOBPPWaveformTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALInspiralChooseModel(LALStatus *status, expnFunc *func, expnCoeffs *ak, InspiralTemplate *params)
REAL8 XLALInspiralFrequency3_2PN(REAL8 td, expnCoeffs *ak)
int XLALBBHPhenTimeDomEngine(REAL4Vector *signalvec1, REAL4Vector *signalvec2, REAL4Vector *h, REAL4Vector *aVec, REAL4Vector *freqVec, REAL8Vector *phiVec, UINT4 *countback, InspiralTemplate *params)
REAL8 XLALInspiralTiming2_2PN(REAL8 f, void *params)
INT4 XLALPSpinFinalMassSpin(REAL8 *finalMass, REAL8 *finalSpin, InspiralTemplate *params, REAL8 energy, REAL8 *LNhvec)
int XLALInspiralCalculateIIRSetInnerProduct(COMPLEX16Vector *a1, COMPLEX16Vector *b0, INT4Vector *delay, REAL8Vector *psd, double *ip)
INT4 XLALGenerateWaveDerivatives(REAL4Vector *dwave, REAL4Vector *ddwave, REAL4Vector *wave, InspiralTemplate *params)
REAL8 XLALInspiralFrequency3_7PN(REAL8 td, expnCoeffs *ak)
REAL8 XLALInspiralFrequency3_6PN(REAL8 td, expnCoeffs *ak)
int XLALEOBWaveform(REAL4Vector *signalvec, InspiralTemplate *params)
REAL8 XLALInspiralTofV(REAL8, void *)
REAL8 XLALChirpTimeReducedSpin(REAL8 v, REAL8 m1, REAL8 m2, REAL8 spin1, REAL8 spin2, UINT4 pnOrder)
Compute the chirp time of the "reduced-spin" templates.
int XLALEOBWaveformForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
INT4 XLALGenerateWaveDerivative(REAL8Vector *dwave, REAL8Vector *wave, REAL8 dt)
void LALBBHPhenWaveFreqDomTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALTaylorT4Waveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
int XLALEOBWaveformTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
int XLALInspiralWave3Templates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
REAL8 XLALInspiralVelocity(TofVIn *params)
REAL8 XLALInspiralTofVIntegrand(REAL8 v, void *params)
int XLALTaylorF2ReducedSpin(REAL4Vector *signalvec, InspiralTemplate *params)
Generate the "reduced-spin templates" proposed in http://arxiv.org/abs/1107.1267.
INT4 XLALGenerateQNMFreq(COMPLEX8Vector *modefreqs, InspiralTemplate *params, UINT4 l, UINT4 m, UINT4 nmodes)
void LALBBHPhenTimeDomEngine(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, REAL4Vector *h, REAL4Vector *aVec, REAL4Vector *freqVec, REAL8Vector *phiVec, UINT4 *countback, InspiralTemplate *params)
REAL8 XLALInspiralTiming2_6PN(REAL8 f, void *params)
void LALEtaTau04(LALStatus *status, REAL8 *x, REAL8 eta, void *in)
INT4 XLALPSpinInspiralAttachRingdownWave(REAL8Vector *signalvec, InspiralTemplate *params, UINT4 *attpos, UINT4 nmodes, UINT4 l, INT4 m, REAL8 finalMass, REAL8 finalSpin)
int XLALInspiralChooseModel(expnFunc *func, expnCoeffs *ak, InspiralTemplate *params)
void LALEtaTau02(LALStatus *status, REAL8 *x, REAL8 eta, void *in)
REAL4 LALInspiralHPlusPolarization(REAL8 phase, REAL8 v, InspiralTemplate *params)
int XLALBBHPhenWaveTimeDomForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALEOBPPWaveform(REAL4Vector *signalvec, InspiralTemplate *params)
REAL8 XLALInspiralTiming2_3PN(REAL8 f, void *params)
void LALBCVSpinWaveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
int XLALTaylorF2ReducedSpinTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
Generate two orthogonal "reduced-spin" templates.
int XLALBBHPhenWaveAFreqDom(REAL4Vector *signalvec, InspiralTemplate *params)
void LALRungeKutta4(LALStatus *, REAL8Vector *, rk4GSLIntegrator *, void *)
REAL8 XLALInspiralPhiofVIntegrand(REAL8, void *)
REAL8 XLALInspiralTiming2_4PN(REAL8 f, void *params)
void LALBBHPhenWaveFreqDom(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralWave3(REAL4Vector *signalvec, InspiralTemplate *params)
INT4 XLALInspiralHybridAttachRingdownWave(REAL4Vector *signalvec1, REAL4Vector *signalvec2, INT4 l, INT4 m, REAL8Vector *timeVec, REAL8Vector *matchrange, InspiralTemplate *params)
REAL8 XLALInspiralFrequency3_5PN(REAL8 td, expnCoeffs *ak)
void LALSTPNWaveformTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
int XLALTaylorEtWaveform(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALInspiralWave3ForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
REAL8 XLALInspiralPhasing3_4PN(REAL8 td, expnCoeffs *ak)
void LALTaylorT4WaveformTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
REAL8 XLALEtaTau02(REAL8 eta, void *in)
void LALBBHPhenWaveTimeDomForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALBBHPhenWaveTimeDom(REAL4Vector *signalvec, InspiralTemplate *insp_template)
void XLALSimInjectNinjaSignals(REAL4TimeSeries *chan, const char *ifo, REAL8 dynRange, SimInspiralTable *events)
int XLALInspiralGenerateIIRSet(REAL8Vector *amp, REAL8Vector *phase, double epsilon, double alpha, double beta, double padding, COMPLEX16Vector **a1, COMPLEX16Vector **b0, INT4Vector **delay)
void LALInspiralParameterCalc(LALStatus *status, InspiralTemplate *params)
INT4 XLALPSpinInspiralRingdownWave(REAL8Vector *rdwave, InspiralTemplate *params, REAL8Vector *inspwave, COMPLEX8Vector *modefreqs, UINT4 nmodes)
INT4 XLALInspiralRingdownWave(REAL4Vector *rdwave1, REAL4Vector *rdwave2, InspiralTemplate *params, REAL4VectorSequence *inspwave1, REAL4VectorSequence *inspwave2, COMPLEX8Vector *modefreqs, UINT4 nmodes)
void LALSTPNWaveformForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALInspiralIIRSetResponse(COMPLEX16Vector *a1, COMPLEX16Vector *b0, INT4Vector *delay, COMPLEX16Vector *response)
REAL8 XLALInspiralTiming2_7PN(REAL8 f, void *params)
void LALInspiralEccentricityTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALSTPNFramelessWaveformTemplates(LALStatus *status, REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
int XLALSTPNFramelessWaveformTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
int XLALSimInspiralChooseWaveformFromSimInspiral(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, SimInspiralTable *thisRow, REAL8 deltaT)
Generate the plus and cross polarizations for a waveform form a row of the sim_inspiral table.
void LALSTPNFramelessWaveform(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
REAL8 XLALInspiralPhasing3_2PN(REAL8 td, expnCoeffs *ak)
REAL8 XLALInspiralPhasing3_5PN(REAL8 td, expnCoeffs *ak)
void LALTaylorT4WaveformForInjection(LALStatus *status, CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
int XLALTaylorEtWaveformTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
void LALInspiralWave(LALStatus *status, REAL4Vector *signalvec, InspiralTemplate *params)
rk4GSLIntegrator * XLALRungeKutta4Init(INT4 n, rk4In *input)
REAL8 XLALInspiralTiming2_0PN(REAL8 f, void *params)
int XLALRungeKutta4(REAL8Vector *yout, rk4GSLIntegrator *integrator, void *params)
REAL8() InspiralPhasing2(REAL8 v, expnCoeffs *ak)
REAL8 EnergyFunction(REAL8 v, expnCoeffs *ak)
REAL8() InspiralTiming2(REAL8 f, void *params)
InputMasses
This structure is one of the members of the InspiralTemplate structure.
int XLALGetInteractionFromString(const CHAR *inString)
XLAL function to determine LALInspiralInteraction from a string.
void() TestFunction(REAL8Vector *vector1, REAL8Vector *vector2, void *params)
REAL8() InspiralPhasing3(REAL8 td, expnCoeffs *ak)
REAL8() InspiralFrequency3(REAL8 td, expnCoeffs *ak)
REAL8 FluxFunction(REAL8 v, expnCoeffs *ak)
LALInspiralInteraction
Enumeration to specify which interaction will be used in the waveform generation.
@ massesAndSpin
UNDOCUMENTED.
@ t01
unused; shouldn't be used
@ fixedPsi
The two psi values are given by the input parameter structure.
@ totalMassAndMu
total mass and reduced mass
@ fixedMasses
The two masses are given by the input parameter structure.
@ totalMassAndEta
total mass and symmetric mass ratio
@ psi0Andpsi3
BCV parameters and .
@ t03
chirptimes and , and
@ bhns
One of the mass is a Neutron star and the other a black hole (m1 [minMass-3] and m2 [3-maxMass])
@ totalMassUAndEta
total mass and eta but uniform distribution in totalMass
@ m1Andm2
component masses
@ minmaxTotalMass
UNDOCUMENTED.
@ fixedTau
The two tau values are given by the input parameter structure.
@ LAL_INSPIRAL_INTERACTION_TIDAL_5PN
Leading-order tidal interaction.
@ LAL_INSPIRAL_INTERACTION_ALL
all spin and tidal interactions
@ LAL_INSPIRAL_INTERACTION_TIDAL_6PN
Next-to-leading-order tidal interaction.
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_25PN
Next-to-leading-order spin-orbit interaction.
@ LAL_INSPIRAL_INTERACTION_QUAD_MONO_2PN
Quadrupole-monopole interaction.
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_3PN
Spin-spin interaction.
@ LAL_INSPIRAL_INTERACTION_NONE
No spin, tidal or other interactions.
@ LAL_INSPIRAL_INTERACTION_SPIN_SPIN_2PN
Spin-spin interaction.
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_15PN
Leading order spin-orbit interaction.
@ LAL_INSPIRAL_INTERACTION_ALL_SPIN
all spin interactions, no tidal interactions
@ LAL_INSPIRAL_INTERACTION_SPIN_SPIN_SELF_2PN
Spin-spin-self interaction.
int XLALPSpinInspiralRDTemplates(REAL4Vector *signalvec1, REAL4Vector *signalvec2, InspiralTemplate *params)
Function to produce waveform templates.
int XLALPSpinInspiralRD(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALPSpinInspiralRDFreqDom(REAL4Vector *signalvec, InspiralTemplate *params)
int XLALPSpinInspiralRDForInjection(CoherentGW *waveform, InspiralTemplate *params, PPNParamStruc *ppnParams)
Function Module to produce injection waveforms.
These are the input structures needed to solve for the mass ratio given the chirptimes or .
This structure is needed to solve the differential equation giving the evolution of the orbital angul...
REAL8 fourM1Plus
fourM1Plus: = (all masses expressed in seconds)
REAL8 v
parameter of 'integration': v=sqrt(M/r)
REAL8 thirtytwoBy5etc
thirtytwoBy5etc:=
REAL8 M
Total mass of the binary (in seconds)
REAL8 magS1
The constant spin magnitude of the primary.
REAL8 fourM2Plus
fourM2Plus: = (all masses expressed in seconds)
REAL8 magS2
The constant spin magnitude of the secondary.
REAL8 oneBy2Mcube
oneBy2Mcube: =
REAL8 threeBy2Mcube
threeBy2Mcube: =
Structure used as an input to compute the derivatives needed in solving the phasing formula when the ...
EOBNonQCCoeffs * nqcCoeffs
Structures used to compute the phase of the signal from the ‘beginning’, when the veolcity parameter ...
The inspiral waveform parameter structure containing information about the waveform to be generated.
REAL8 alpha4
UNDOCUMENTED.
REAL8 signalAmplitude
dimensionless amplitude of the signal (input, currently unused)
REAL8 startTime
starting time of the waveform (in sec); if different from zero, the waveform will start with an insta...
REAL8 alpha6
UNDOCUMENTED.
Approximant approximant
Post-Newtonain approximant to be used in generating the wave (input)
REAL8 t5
2.5 post-Newtonian chirp time in seconds (output)
INT4 ieta
parameter that tells whether the symmetric mass ratio should be set to zero in the PN expansions of ...
INT4Vector * segmentIdVec
Vector of segment that have been filtered against this template needed for the LDAS implementation of...
REAL8 psi3
BCV parameter .
REAL8 mu
reduced mass (in solar mass) (input/output)
INT4 nStartPad
Number of leading elements in the signal generation to be set to zero (input) If template is requeste...
REAL8 t2
first post-Newtonian chirp time in seconds (input/output)
struct tagInspiralTemplate * fine
Linked list to the next fine bank template (currently not filled by inspiral or bank codes)
REAL8 inclination
Inclination of the orbit (currently not in use)
LALPNOrder ampOrder
UNDOCUMENTED.
REAL8 eta
symmetric mass ratio (input/output)
REAL8 alpha5
UNDOCUMENTED.
REAL8 totalMass
total mass of the binary in solar mass (input/output)
REAL8 alpha
BCV amplitude correction factor .
REAL8 alpha1
UNDOCUMENTED.
REAL8 mass1
Mass of the primary in solar mass (input/output)
REAL8 distance
Distance to the binary in seconds.
REAL8 Theta
The 3PN unknown flux parameter; likely to be around unity; most waveform generation routines take the...
REAL8 startPhase
starting phase of the waveform in radians (input)
REAL8 eccentricity
initial eccentricity of the orbit (currently not in use)
REAL8 fCutoff
upper frequency cutoff in Hz to be used in generating the waveform; If the last stable orbit frequenc...
REAL8 sourcePhi
Azimuth angle in the direction to the source.
long event_id
UNDOCUMENTED.
REAL8 t4
second post-Newtonian chirp time in seconds (output)
REAL8 t3
1.5 post-Newtonian chirp time in seconds (input/output)
UINT4 inspiralOnly
UNDOCUMENTED.
REAL8 psi0
BCV parameter .
REAL8 alpha2
UNDOCUMENTED.
INT4 nEndPad
Number of trailing bins to be set to zero, the resulting waveform will have at least this many bins z...
REAL8 t7
3.5 post-Newtonian chirp time in seconds (output)
REAL8 alpha3
UNDOCUMENTED.
REAL8 orbitPhi0
Initial azimuth angle of the orbit.
REAL8 fFinal
final frequency reached, in units of Hz (output)
REAL8 t0
Newtonain chirp time in seconds (input/output)
REAL8 chirpMass
chirp mass of the binary in solar mass (output)
REAL8 tC
total chirp time seconds (output)
REAL8 mass2
Mass of the secondary in solar mass (mass1 need not be larger than mass2 (input/output)
UINT4 fixedStep
UNDOCUMENTED.
LIGOTimeGPS end_time
UNDOCUMENTED.
REAL8 orbitTheta0
Initial co-latitute of the orbit.
REAL4 minMatch
The minimal match specified by the user when the bank that contains this template was created.
LALPNOrder order
Post-Newtonain order to be used in generating the wave (input)
INT4 level
Flag used in heirarical serached to indicate if this is a coarse or a fine template.
REAL8 vFinal
final velocity parameter, in units of the speed of light (output)
InputMasses massChoice
The pair of (mass) parameters given (see structure defining this member for more details) (input)
struct tagInspiralTemplate * next
Linked list to the next coarse bank template (currently not filled by inspiral or bank codes)
REAL8 fLower
lower frequency cutoff of the detector in Hz (input)
REAL8 chi
dimensionless spin of black hole (ie mass1)
REAL8 OmegaS
The 3PN (unknown) parameter; calculated to be equal to zero by Damour, Jaranowski and Schaffer (input...
INT4 number
Unique ID number for this template needed for the LDAS implementation of the inspiral search.
REAL8 tSampling
Sampling rate in Hz (input)
REAL8 t6
third post-Newtonian chirp time in seconds (output)
REAL8 sourceTheta
Co-latitute in the direction to the source.
LALInspiralInteraction interaction
UNDOCUMENTED.
REAL8 kappa
cosine of angle between spin of mass1 and orb ang mom
LALSimInspiralFrameAxis axisChoice
UNDOCUMENTED.
This structure stores the parameters for constructing a restricted post-Newtonian waveform.
Structures used to compute the phase of the signal from the ‘beginning’, when the veolcity parameter ...
This structure contains various post-Newtonian and P-approximant expansion coefficients; the meanings...
Structure to hold the pointers to the generic functions defined above.
InspiralPhasing2 * phasing2
InspiralPhasing3 * phasing3
InspiralTiming2 * timing2
InspiralFrequency3 * frequency3
Structure containing steps and controls for the GSL Runge-Kutta solver.
gsl_odeiv_control * control
const gsl_odeiv_step_type * type
gsl_odeiv_evolve * evolve
Structure used as an input to Runge-Kutta solver.