Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALInspiral 5.0.3.1-b246709
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALInspiral.h
Go to the documentation of this file.
1/*
2* Copyright (C) 2007 Stas Babak, David Churches, Drew Keppel, Duncan Brown, Jolien Creighton, David McKechan, Patrick Brady, Peter Shawhan, Reinhard Prix, B.S. Sathyaprakash, Anand Sengupta, Craig Robinson , Sean Seader, Thomas Cokelaer, Riccardo Sturani, Laszlo Vereb
3*
4* This program is free software; you can redistribute it and/or modify
5* it under the terms of the GNU General Public License as published by
6* the Free Software Foundation; either version 2 of the License, or
7* (at your option) any later version.
8*
9* This program is distributed in the hope that it will be useful,
10* but WITHOUT ANY WARRANTY; without even the implied warranty of
11* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12* GNU General Public License for more details.
13*
14* You should have received a copy of the GNU General Public License
15* along with with program; see the file COPYING. If not, write to the
16* Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17* MA 02110-1301 USA
18*/
19
20#ifndef _LALINSPIRAL_H
21#define _LALINSPIRAL_H
22
23#include <gsl/gsl_odeiv.h>
24
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>
32
33#ifdef __cplusplus
34extern "C" {
35#endif
36
37/**
38 * \addtogroup LALInspiral_h
39 * \author Churches, D. K , B. S. Sathyaprakash, T. Cokelaer.
40 *
41 * \brief %Header file for the template generation codes.
42 *
43 * ### Synopsis ###
44 *
45 * \code
46 * #include <lal/LALInspiral.h>
47 * \endcode
48 *
49 * This header file covers routines that are used in template generation.
50 *
51 */
52/** @{ */
53
54/**\name Error Codes */
55/** @{ */
56#define LALINSPIRALH_ENULL 1 /**< Arguments contained an unexpected null pointer */
57#define LALINSPIRALH_EMEM 2 /**< Memory allocation error */
58#define LALINSPIRALH_EDIV0 3 /**< Division by zero */
59#define LALINSPIRALH_ESIZE 4 /**< Invalid input range */
60#define LALINSPIRALH_ECHOICE 5 /**< Invalid choice for an input parameter */
61#define LALINSPIRALH_EORDER 6 /**< unknown order specified */
62#define LALINSPIRALH_EAPPROXIMANT 7 /**< Invalid model */
63#define LALINSPIRALH_EPSI0 8 /**< psi0 must be > 0 */
64#define LALINSPIRALH_EPSI3 9 /**< psi3 must be < 0 */
65#define LALINSPIRALH_EALPHA 10 /**< alpha must be defined positive */
66#define LALINSPIRALH_EFCUTOFF 11 /**< fcutoff must be defined and > 0 */
67#define LALINSPIRALH_ENOWAVEFORM 12 /**< No Waveform generated */
68#define LALINSPIRALH_ESTOPPED 13 /**< Waveform generation stopped */
69#define LALINSPIRALH_EROOTINIT 14 /**< Can't find good bracket for BisectionFindRoot */
70#define LALINSPIRALH_EFLOWER 15 /**< fLower too low in comparison to flso */
71#define LALINSPIRALH_EVECTOR 16 /**< Attempting to write beyond the end of vector */
72#define LALINSPIRALH_EFLOWERINJ 17 /**< flower for the injection must be greater than zero */
73#define LALINSPIRALH_EORDERMISSING 18 /**< The PN order requested is not implemented for this approximant */
74#define LALINSPIRALH_EBPERR 19 /**< Error in band passing signal */
75#define LALINSPIRALH_ESWITCH 20 /**< Unknown case in switch */
76#define LALINSPIRALH_EMASSCHOICE 21 /**< Improper choice for massChoice */
77/** @} */
78
79/** \cond DONT_DOXYGEN */
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"
101/** \endcond */
102
103#define LAL_INSPIRAL_INTERACTION_DEFAULT LAL_INSPIRAL_INTERACTION_ALL
104
105/**
106 * XLAL function to determine LALInspiralInteraction from a string.
107 * DEPRECATED: USE LALSimInspiralSpinOrder, LALSimInspiralTidalOrder INSTEAD
108 */
109int XLALGetInteractionFromString(const CHAR *inString);
110
111/**
112 * Enumeration to specify which interaction will be used in the waveform
113 * generation. Their combination also can be used by the bitwise or.
114 * DEPRECATED: USE LALSimInspiralSpinOrder, LALSimInspiralTidalOrder INSTEAD
115 */
116typedef enum tagLALInspiralInteraction {
117 LAL_INSPIRAL_INTERACTION_NONE = 0, /**< No spin, tidal or other interactions */
118 LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_15PN = 1, /**< Leading order spin-orbit interaction */
119 LAL_INSPIRAL_INTERACTION_SPIN_SPIN_2PN = 1 << 1, /**< Spin-spin interaction */
120 LAL_INSPIRAL_INTERACTION_SPIN_SPIN_SELF_2PN = 1 << 2, /**< Spin-spin-self interaction */
121 LAL_INSPIRAL_INTERACTION_QUAD_MONO_2PN = 1 << 3, /**< Quadrupole-monopole interaction */
122 LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_25PN = 1 << 4, /**< Next-to-leading-order spin-orbit interaction */
123 LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_3PN = 1 << 5, /**< Spin-spin interaction */
124 LAL_INSPIRAL_INTERACTION_TIDAL_5PN = 1 << 6, /**< Leading-order tidal interaction */
125 LAL_INSPIRAL_INTERACTION_TIDAL_6PN = 1 << 7, /**< Next-to-leading-order tidal interaction */
126 LAL_INSPIRAL_INTERACTION_ALL_SPIN = (1 << 6) - 1, /**< all spin interactions, no tidal interactions */
127 LAL_INSPIRAL_INTERACTION_ALL = (1 << 8) - 1 /**< all spin and tidal interactions */
129
130
131/**
132 * These are the input structures needed to solve for the mass
133 * ratio \f$\eta\f$ given the chirptimes \f$(\tau_0,\, \tau_2)\f$ or
134 * \f$(\tau_0, \, \tau_4).\f$
135 *
136 * Here, \c t2\f$ = \tau_2,\f$ \c A2 \f$ = A_2 ({\tau_0}/{A_0})^{3/5},\f$ and \c B2 \f$=B_2\f$,
137 * where \f$A_0 = 5/[256 (\pi f_{s} )^{8/3}],\f$ \f$A_2 = 3715 / [64512 (\pi f_s)^2],\f$
138 * \f$B_2 = 4620/3715.\f$
139 *
140 * Similarly, \c t4 \f$ = \tau_4,\f$ \c A4 \f$ = A_4 ({\tau_0}/{A_0})^{1/5},\f$
141 * \c B4 \f$=B_4\f$ and \c C4 \f$=C_4,\f$ where
142 * where \f$A_0 = 5/[256 (\pi f_{s} )^{8/3}],\f$
143 * \f$A_4 = 5 \times 3058673/ [128 \times 1016064 (\pi f_s)^{4/3}],\f$
144 * \f$B_4 = 5429 \times 1016064 /(1008 \times 3058673),\f$ and \f$C_4 = 617 \times
145 * 1016064/(144 \times 3058673).\f$
146 */
147/** @{ */
148typedef struct
149tagEtaTau02In
150{
154} EtaTau02In;
155
156typedef struct
157tagEtaTau04In
158{
163} EtaTau04In;
164/** @} */
165
166
167/**
168 * This structure is one of the members of the \c InspiralTemplate structure.
169 * A user can specify the parameters of a binary using any of the following combination of \e masses:
170 * m1Andm2, totalMassAndEta, totalMassUAndEta, totalMassAndMu, t01, t02, t03, t04, psi0Andpsi3
171 *
172 * The LALRandomInspiralSignal uses that structure as an input. Since the injected
173 * waveform are not necessarely wanted to be random, we also provide the following
174 * options
175 * bhns, fixedMasses, fixedPsi, fixedTau
176 *
177 */
178typedef enum tagInputMasses {
179 m1Andm2, /**< component masses */
180 totalMassAndEta, /**< total mass and symmetric mass ratio */
181 totalMassUAndEta, /**< total mass and eta but uniform distribution in totalMass */
182 totalMassAndMu, /**< total mass and reduced mass */
183 t01, /**< unused; shouldn't be used */
184 t02, /**< chirptimes \f$\tau_0\f$ and \f$\tau_2\f$ */
185 t03, /**< chirptimes \f$\tau_0\f$ and \f$\tau_3\f$, and */
186 t04, /**< chirptimes \f$\tau_0\f$ and \f$\tau_4\f$ */
187 psi0Andpsi3, /**< BCV parameters \f$\psi_0\f$ and \f$\psi_3\f$ */
188
189 bhns, /**< One of the mass is a Neutron star and the other a black hole (m1 \f$\in\f$ [minMass-3] and m2 \f$\in\f$ [3-maxMass]) */
190 fixedMasses, /**< The two masses are given by the input parameter structure */
191 fixedPsi, /**< The two psi values are given by the input parameter structure */
192 fixedTau, /**< The two tau values are given by the input parameter structure */
193 massesAndSpin, /**< UNDOCUMENTED */
194 minmaxTotalMass, /**< UNDOCUMENTED */
195 spinOnly /**< UNDOCUMENTED */
197
198
199
200/**
201 * The inspiral waveform parameter structure containing information about the waveform to be generated.
202 */
203typedef struct
204tagInspiralTemplate
205{
206 /** \name Parameters needed to generate Taylor/Pade waveforms */
207 /** @{ */
208 Approximant approximant; /**< Post-Newtonain approximant to be used in generating the wave (input) */
209 LALPNOrder order; /**< Post-Newtonain order to be used in generating the wave (input) */
210 LALPNOrder ampOrder; /**< UNDOCUMENTED */
211 REAL8 mass1; /**< Mass of the primary in solar mass (input/output) */
212 REAL8 mass2; /**< Mass of the secondary in solar mass (\c mass1 need not be larger than \c mass2 (input/output) */
213 REAL8 fCutoff; /**< upper frequency cutoff in Hz to be used in generating the waveform;
214 * If the last stable orbit frequency is smaller than the upper cutoff it will be used
215 * in terminating the waveform instead of fCutoff (input)
216 */
217 REAL8 fLower; /**< lower frequency cutoff of the detector in Hz (input) */
218 REAL8 tSampling; /**< Sampling rate in Hz (input) */
219 REAL8 distance; /**< Distance to the binary in seconds */
220 REAL8 signalAmplitude; /**< dimensionless amplitude of the signal (input, currently unused) */
221 REAL8 startPhase; /**< starting phase of the waveform in radians (input) */
222 REAL8 startTime; /**< starting time of the waveform (in sec); if different from zero, the
223 * waveform will start with an instantaneous frequency different from fLower and reach
224 * fLower at time (approximately) zero (input, not used in Stationary phase approximation)
225 */
226 INT4 ieta; /**< parameter that tells whether the symmetric mass ratio \f$\eta\f$
227 * should be set to zero in the PN expansions of GW flux and binding energy;
228 * If <tt>ieta=0</tt> \f$\eta\f$ will be set to zero, otherwise the appropriate
229 * value of \f$\eta\f$ from the given parameters will be used
230 */
231 /** @} */
232
233 /** \name Additional parameters for EOB waveforms */
234 /** @{ */
235 REAL8 Theta; /**< The 3PN unknown flux parameter; likely to be around unity;
236 * most waveform generation routines take theta to be zero; Robustness of the EOB waveform
237 * has been demonstrated for \f$-2 < \f$ \c Theta \f$< 2\f$ (input)
238 */
239 REAL8 Zeta2; /**< UNDOCUMENTED */
240 /** @} */
241
242 /** \name Parameters for BCV1 template */
243 /** @{ */
244 REAL8 alpha; /**< BCV amplitude correction factor \f$\alpha f_\mathrm{cut}^{2/3}\f$ */
245 REAL8 psi0; /**< BCV parameter \f$\psi_0\f$ */
246 REAL8 psi3; /**< BCV parameter \f$\psi_3\f$ */
247 /** @} */
248
249 /** \name Additional parameters for BCV2 template */
250 /** @{ */
251 REAL8 beta; /**< UNDOCUMENTED */
252 REAL8 alpha1; /**< UNDOCUMENTED */
253 REAL8 alpha2; /**< UNDOCUMENTED */
254 REAL8 alpha3; /**< UNDOCUMENTED */
255 REAL8 alpha4; /**< UNDOCUMENTED */
256 REAL8 alpha5; /**< UNDOCUMENTED */
257 REAL8 alpha6; /**< UNDOCUMENTED */
258 /** @} */
259
260 /** \name Parameters for spinning BH waveform */
261 /** @{ */
262 REAL8 inclination; /**< Inclination of the orbit (currently not in use) */
263 REAL8 orbitTheta0; /**< Initial co-latitute of the orbit */
264 REAL8 orbitPhi0; /**< Initial azimuth angle of the orbit */
265 REAL8 spin1[3]; /**< Spin vector of the primary */
266 REAL8 spin2[3]; /**< Spin vector of the secondary */
267 REAL8 sourceTheta; /**< Co-latitute in the direction to the source */
268 REAL8 sourcePhi; /**< Azimuth angle in the direction to the source */
270 /** @} */
271
272 /** \name Spin parameters for the PTF template */\
273 /** @{ */
274 REAL8 chi; /**< dimensionless spin of black hole (ie mass1) */
275 REAL8 kappa; /**< cosine of angle between spin of mass1 and orb ang mom */
276 /** @} */
277
278 /** \name Parameters which are currently might be used */
279 /** @{ */
280 REAL8 eccentricity; /**< initial eccentricity of the orbit (currently not in use) */
281 /** @} */
282
284
285/**
286 * \name Paramters which are computed using LALInspiralParameterCalc
287 * Note that tc and fFinal are computed during waveform generation!!!
288 */
289 /** @{ */
290 REAL8 chirpMass; /**< chirp mass of the binary \f$=\eta^{3/5} m\f$ in solar mass (output) */
291 REAL8 eta; /**< symmetric mass ratio \f$\eta=m_1m_2/m^2\f$ (input/output) */
292 REAL8 totalMass; /**< total mass of the binary \f$m=m_1+m_2\f$ in solar mass (input/output) */
293 REAL8 fFinal; /**< final frequency reached, in units of Hz (output) */
294 REAL8 t0; /**< Newtonain chirp time in seconds (input/output) */
295 REAL8 t2; /**< first post-Newtonian chirp time in seconds (input/output) */
296 REAL8 t3; /**< 1.5 post-Newtonian chirp time in seconds (input/output) */
297 REAL8 t4; /**< second post-Newtonian chirp time in seconds (output) */
298 REAL8 t5; /**< 2.5 post-Newtonian chirp time in seconds (output) */
299 REAL8 t6; /**< third post-Newtonian chirp time in seconds (output) */
300 REAL8 t7; /**< 3.5 post-Newtonian chirp time in seconds (output) */
301 REAL8 tC; /**< total chirp time seconds (output) */
302
303 REAL4 minMatch; /**< The minimal match specified by the user when the bank that contains this template was created */
304 REAL8 mu; /**< reduced mass (in solar mass) (input/output) */
305 INT4 level; /**< Flag used in heirarical serached to indicate if this is a coarse or a fine template */
306 INT4 number; /**< Unique ID number for this template needed for the LDAS implementation of the inspiral search */
307 INT4 nStartPad; /**< Number of leading elements in the signal generation to be set to zero (input)
308 * If template is requested, that value must be set to zero; In the injection routines related to inject
309 * package, that nStartPad is set to zero; However, for injection performed using the inspiral package,
310 * that value can be set to non zero
311 */
312 INT4 nEndPad; /**< Number of trailing bins to be set to zero, the resulting waveform will have at least this many bins zero at the end, probably
313 * more since we always deal with an integer power of 2 array (input)
314 */
315 REAL8 OmegaS; /**< The 3PN (unknown) parameter; calculated to be equal to zero by Damour, Jaranowski and Schaffer (input) */
316 REAL8 vFinal; /**< final velocity parameter, in units of the speed of light (output) */
317
318 /* REAL8 vInitial; // initial velocity parameter, in units of the speed of light (output)
319 REAL8 rFinal; // final 'separation' between the bodies, in units of total mass (output)
320 REAL8 rInitial; // initial radial separation of the two, in units of total mass bodies (used only in EOB waveforms) (output)
321 REAL8 rLightRing;*/ // radial coordinate at the light ring, in units of total mass (output)
322
323 InputMasses massChoice; /**< The pair of (mass) parameters given (see structure defining this member for more details) (input) */
324 INT4Vector *segmentIdVec; /**< %Vector of segment that have been filtered against this template needed for the LDAS implementation of the inspiral search */
325 LIGOTimeGPS end_time; /**< UNDOCUMENTED */
326 long event_id; /**< UNDOCUMENTED */
327 CHAR ifo[LIGOMETA_IFO_MAX]; /**< UNDOCUMENTED */
328
329 REAL4 Gamma[10]; /**< Gamma[] is a vector that stores the upper triangular part of the metric in
330 * the space of parameters; For time domain searches, Gamma[0,...,5] stores
331 * the following information :<br>
332 * Gamma[0] -> (tc,tc) metric component<br>
333 * Gamma[1] -> (tc,t0) metric component<br>
334 * Gamma[2] -> (tc,t3) metric component<br>
335 * Gamma[3] -> (t0,t0) metric component<br>
336 * Gamma[4] -> (t0,t3) metric component<br>
337 * Gamma[5] -> (t3,t3) metric component<br>
338 * For spinBCV searches, (in 4 dimensions) Gamma[0,...,9] would be required
339 */
340 REAL4 qmParameter[2]; /**< UNDOCUMENTED */
342
343 UINT4 fixedStep; /**< UNDOCUMENTED */
344 UINT4 inspiralOnly; /**< UNDOCUMENTED */
345 /** @} */
346
347 struct tagInspiralTemplate *next; /**< Linked list to the next coarse bank template (currently not filled by inspiral or bank codes) */
348 struct tagInspiralTemplate *fine; /**< Linked list to the next fine bank template (currently not filled by inspiral or bank codes) */
350
351
352/**
353 * This is a structure needed by the inner workings of the inspiral wave generation code
354 */
355typedef struct
356tagInspiralToffInput
357{
370
371
372/**
373 * This structure is needed to solve the differential equation
374 * giving the evolution of the orbital angular momentum and the
375 * spin angular momenta in the case of spinning black hole binaries.
376 */
377typedef struct
378tagInspiralACSTParams
379{
380 REAL8 v; /**< parameter of 'integration': v=sqrt(M/r) */
381 REAL8 magS1; /**< The constant spin magnitude of the primary */
382 REAL8 magS2; /**< The constant spin magnitude of the secondary */
383 REAL8 NCap[3]; /**< Source direction (unit vector) in detector coordinate system */
384 REAL8 spin1[3]; /**< Spin of the larger body */
385 REAL8 M; /**< Total mass of the binary (in seconds) */
386 REAL8 fourM1Plus; /**< fourM1Plus: = \f$(4 m_1+3 m_2)/(2 m_1 M^3)\f$ (all masses expressed in seconds) */
387 REAL8 fourM2Plus; /**< fourM2Plus: = \f$(4 m_2+3 m_1)/(2 m_2 M^3)\f$ (all masses expressed in seconds) */
388 REAL8 oneBy2Mcube; /**< oneBy2Mcube: = \f$1/(2 M^3)\f$ */
389 REAL8 threeBy2Mcube; /**< threeBy2Mcube: = \f$3/(2 M^3)\f$ */
390 REAL8 thirtytwoBy5etc; /**< thirtytwoBy5etc:= \f$(32/5) \eta^2 M\f$ */
392
393/**
394 * This structure contains various post-Newtonian and P-approximant expansion
395 * coefficients; the meanings of the coefficients is indicated as comments
396 * before each list.
397 */
398typedef struct
399tagexpnCoeffs {
400 int ieta;
401
402 /** \name Coefficients in the Taylor expansion of new energy function */
403 /** @{ */
404 REAL8 eTaN, eTa1, eTa2, eTa3;
405 /** @} */
406
407 /** \name Coefficients in the Pade expression of new energy function */
408 /** @{ */
409 REAL8 ePaN, ePa1, ePa2, ePa3;
410 /** @} */
411
412 /** \name Coefficients in the Taylor expansion of usual energy function */
413 /** @{ */
414 REAL8 ETaN, ETa1, ETa2, ETa3;
415 /** @} */
416
417 /** \name Coefficients in the Taylor expansion of the derivative of the usual energy function */
418 /** @{ */
419 REAL8 dETaN, dETa1, dETa2, dETa3;
420 /** @} */
421
422 /** \name Taylor expansion coefficients of energy flux */
423 /** @{ */
424 REAL8 FTaN, FTa1, FTa2, FTa3, FTa4, FTa5, FTa6, FTa7, FTa8, FTl6, FTl8;
425 /** @} */
426
427 /** \name Taylor expansion coefficients of factored flux */
428 /** @{ */
429 REAL8 fTaN, fTa1, fTa2, fTa3, fTa4, fTa5, fTa6, fTa7, fTa8;
430 /** @} */
431
432 /* \name Coefficients of the corresponding P-approximant */
433 /** @{ */
434 REAL8 fPaN, fPa1, fPa2, fPa3, fPa4, fPa5, fPa6, fPa7, fPa8;
435 /** @} */
436
437 /** \name Taylor expansion coefficents in t(v) */
438 /** @{ */
439 REAL8 tvaN, tva2, tva3, tva4, tva5, tva6, tva7, tvl6;
440 /** @} */
441
442 /** \name Taylor expansion coefficents in phi(v) */
443 /** @{ */
444 REAL8 pvaN, pva2, pva3, pva4, pva5, pva6, pva7, pvl6;
445 /** @} */
446
447 /** \name Taylor expansion coefficents in phi(t) */
448 /** @{ */
449 REAL8 ptaN, pta2, pta3, pta4, pta5, pta6, pta7, ptl6;
450 /** @} */
451
452 /** \name Taylor expansion coefficents in f(t) */
453 /** @{ */
454 REAL8 ftaN, fta2, fta3, fta4, fta5, fta6, fta7, ftl6;
455 /** @} */
456
457 /** \name Taylor expansion coefficents in psi(f) in the Fourier phase */
458 /** @{ */
459 REAL8 pfaN, pfa2, pfa3, pfa4, pfa5, pfa6, pfa7, pfl5, pfl6;
460 /** @} */
461
462 /** \name Taylor expansion for the spinning case */
463 /** @{ */
464 REAL8 ST[9], thetahat ;
465 /** @} */
466
467 /** \name Sampling rate and interval */
468 /** @{ */
469 REAL8 samplingrate, samplinginterval;
470 /** @} */
471
472 /** \name Symmetric mass ratio, total mass, component masses */
473 /** @{ */
474 REAL8 eta, totalmass, m1, m2;
475 /** @} */
476
477 /** \name Unknown 3PN parameters, euler constant */
478 /** @{ */
479 REAL8 lambda, theta, EulerC, omegaS, zeta2;
480 /** @} */
481
482/**
483 * \name Initial and final values of frequency, time, velocity; lso
484 * values of velocity and frequency; final phase.
485 */
486 /** @{ */
487 REAL8 f0, fn, t0, tn, v0, vn, vf, vlso, flso, phiC;
488 /** @} */
489
490 /** \name Last stable orbit and pole defined by various Taylor and P-approximants */
491 /** @{ */
492 REAL8 vlsoT0, vlsoT2, vlsoT4, vlsoT6;
493 REAL8 vlsoP0, vlsoP2, vlsoP4, vlsoP6;
495 REAL8 vpoleP4, vpoleP6;
497 /** @} */
498
499} expnCoeffs;
500
501
502/**
503 * \name Energy, flux, phase, time and frequency functions.
504 * The following functions are generic function definitions that will be used in
505 * template generation. The function <tt>LALInspiralChooseModel,</tt>
506 * which is called by wave generation interface code, points these
507 * functions to the appropriate specific functions depending on the
508 * choices made by the user.
509 */
510/** @{ */
512 REAL8 v,
513 expnCoeffs *ak);
514
516 REAL8 v,
517 expnCoeffs *ak);
518
519typedef void (TestFunction)(
520 REAL8Vector *vector1,
521 REAL8Vector *vector2,
522 void *params);
523
525 REAL8 v,
526 expnCoeffs *ak);
527
529 REAL8 td,
530 expnCoeffs *ak);
531
533 REAL8 td,
534 expnCoeffs *ak);
535
537 REAL8 f,
538 void *params);
539/** @} */
540
541/**
542 * Structure to hold the pointers to the generic functions defined above
543 */
544typedef struct
545tagexpnFunc
546{
553} expnFunc;
554
555
556/**
557 * Structure needed to compute the time elapsed
558 * from/to the starting epoch of the waveform when the velocity
559 * parameter was \f$v_0,\f$ to/from the current epoch when velocity
560 * parameter is \f$v\f$
561 */
562typedef struct
563tagTofVIn
564{
565 REAL8 t;
566 REAL8 v0;
567 REAL8 t0;
568 REAL8 vlso;
569 REAL8 totalmass;
570 EnergyFunction *dEnergy;
571 FluxFunction *flux;
572 expnCoeffs *coeffs;
573} TofVIn;
574
575/**
576 * Structure needed to compute the time elapsed
577 * from/to the starting epoch of the waveform when the velocity
578 * parameter was \f$v_0,\f$ to/from the current epoch when velocity
579 * parameter is \f$v\f$
580 */
581typedef struct
582tagTofVIntegrandIn
583{
584 EnergyFunction *dEnergy;
585 FluxFunction *flux;
586 expnCoeffs *coeffs;
588
589/** UNDOCUMENTED */
590typedef struct
591tagEOBNonQCCoeffs
592{
593 REAL8 a1;
594 REAL8 a2;
595 REAL8 a3;
596 REAL8 a4;
597 REAL8 b1;
598 REAL8 b2;
600
601/**
602 * Structure used as an input to compute the derivatives needed in solving the phasing formula when the
603 * \c approximant is <tt>TaylorT1, TaylorP1</tt> or <tt>EOB</tt>
604 */
605typedef struct
606tagInspiralDerivativesIn
607{
614
615/**
616 * Structure used as an input to Runge-Kutta solver.
617 */
618typedef struct
619tagrk4In
620{
630} rk4In;
631
632/**
633 * Structure containing steps and controls for the GSL Runge-Kutta solver.
634 */
635typedef struct
636tagrk4GSLIntegrator
637{
638 const gsl_odeiv_step_type *type;
639 gsl_odeiv_step *step;
640 gsl_odeiv_control *control;
641 gsl_odeiv_evolve *evolve;
645
646
647/**
648 * Structures used to compute the phase of the signal from the `beginning', when the
649 * veolcity parameter is \f$v_0,\f$ to a time when the velocity parameter
650 * has evolved to a user input value \f$v\f$.
651 */
652typedef struct
653tagInspiralPhaseIn
654{
661
662/**
663 * Structures used to compute the phase of the signal from the `beginning', when the
664 * veolcity parameter is \f$v_0,\f$ to a time when the velocity parameter
665 * has evolved to a user input value \f$v\f$.
666 */
667typedef struct
668tagPhiofVIntegrandIn
669{
674
675typedef struct
676tagInspiralInit
677{
681
683
684/** @} */ /* end:LALInspiral_h */
685
686
687/* ---------- Function prototypes ---------- */
688
689/* --- HERE ARE SOME USEFUL PROTOTYPE FOR LENGTH, PARAMETER CALCULATION... --- */
690
694
697
700
703 UINT4 *n,
705
708 expnFunc *func,
709 expnCoeffs *ak,
711
713 expnFunc *func,
714 expnCoeffs *ak,
716
717void LALInspiralSetup (
719 expnCoeffs *ak,
721
723 expnCoeffs *ak,
725
726
727/*
728 * FIXME: The current SWIG wrappings remove LAL and XLAL prefixes from
729 * functions such that the InspiralInit struct collides with
730 * XLALInspiralInit. Since many places refer to the struct, let's
731 * make {XLAL,LAL}InspiralInit invisible to SWIG.
732 */
733#ifndef SWIG
734void
738 InspiralInit *paramsInit);
739
740int
743 InspiralInit *paramsInit);
744#endif /* SWIG */
745
746/**
747 * Generate the plus and cross polarizations for a waveform
748 * form a row of the sim_inspiral table.
749 *
750 * Parses a row from the sim_inspiral table and passes the appropriate members
751 * to XLALSimInspiralChooseWaveform().
752 *
753 * FIXME: this should eventually be moved to lalsimulation
754 * along with the appropriate string parsing functions
755 */
757 REAL8TimeSeries **hplus, /**< +-polarization waveform */
758 REAL8TimeSeries **hcross, /**< x-polarization waveform */
759 SimInspiralTable *thisRow, /**< row from the sim_inspiral table containing waveform parameters */
760 REAL8 deltaT /**< sampling interval */
761 );
762
763/**
764 * Generate the plus and cross polarizations for a conditioned waveform
765 * form a row of the sim_inspiral table.
766 *
767 * Parses a row from the sim_inspiral table and passes the appropriate members
768 * to XLALSimInspiralTD().
769 *
770 */
772 REAL8TimeSeries **hplus, /**< +-polarization waveform */
773 REAL8TimeSeries **hcross, /**< x-polarization waveform */
774 SimInspiralTable *thisRow, /**< row from the sim_inspiral table containing waveform parameters */
775 REAL8 deltaT /**< time step (s) */
776 );
777
778/**
779 * Generate the plus and cross polarizations for a waveform
780 * form a row of the InspiralTemplate structure.
781 *
782 * Parses the InspiralTemplate stucture and passes the appropriate members
783 * to XLALSimInspiralChooseWaveform().
784 */
786 REAL8TimeSeries **hplus, /**< +-polarization waveform */
787 REAL8TimeSeries **hcross, /**< x-polarization waveform */
788 InspiralTemplate *params /**< stucture containing waveform parameters */
789 );
790
791/* --- HERE ARE THE WAVEFORMS/MODELS PROTOTYPES --- */
794 REAL4Vector *signalvec,
796
799 REAL4Vector *filter1,
800 REAL4Vector *filter2,
802
803void
806 CoherentGW *waveform,
808 PPNParamStruc *ppnParams);
809
810
811void LALInspiralWave(
813 REAL4Vector *signalvec,
815
818 REAL4Vector *filter1,
819 REAL4Vector *filter2,
821
824 CoherentGW *waveform,
826 PPNParamStruc *ppnParams);
827
829 REAL4Vector *signalvec,
831
833 REAL4Vector *signalvec1,
834 REAL4Vector *signalvec2,
836
838 CoherentGW *waveform,
840 PPNParamStruc *ppnParams);
841
844 REAL4Vector *signalvec,
846
849 REAL4Vector *signalvec1,
850 REAL4Vector *signalvec2,
852
854 REAL4Vector *signalvec,
856
858 REAL4Vector *signalvec1,
859 REAL4Vector *signalvec2,
861
863 CoherentGW *waveform,
865 PPNParamStruc *ppnParams);
866
867void LALInspiralWave3 (
869 REAL4Vector *signalvec,
871
873 REAL4Vector *signalvec,
875
877 REAL4Vector *signalvec1,
878 REAL4Vector *signalvec2,
880
882 CoherentGW *waveform,
884 PPNParamStruc *ppnParams);
885
886int
888 REAL4Vector *signalvec,
890
891int
893 REAL4Vector *signalvec,
895
897 REAL4Vector *signalvec,
899
901 REAL4Vector *signalvec1,
902 REAL4Vector *signalvec2,
904
906 CoherentGW *waveform,
908 PPNParamStruc *ppnParams);
909
911 REAL4Vector *signalvec,
913
915 REAL4Vector *signalvec1,
916 REAL4Vector *signalvec2,
918
920 CoherentGW *waveform,
922 PPNParamStruc *ppnParams);
923
924void LALBCVWaveform(
926 REAL4Vector *signalvec,
928
930 REAL4Vector *signalvec,
932
934 REAL4Vector *signalvec1,
935 REAL4Vector *signalvec2,
937
940 REAL4Vector *signalvec,
942
945 REAL4Vector *signalvec1,
946 REAL4Vector *signalvec2,
948 );
949
952 CoherentGW *waveform,
954 PPNParamStruc *ppnParams);
955
958 REAL4Vector *signalvec,
960
963 REAL4Vector *signalvec,
965
967 REAL4Vector *signalvec,
969
972 REAL4Vector *signalvec,
973 InspiralTemplate *in);
974
975
976void
979 CoherentGW *waveform,
981 PPNParamStruc *ppnParams);
982
983void
986 REAL4Vector *signalvec1,
987 REAL4Vector *signalvec2,
989 ) ;
990
991void LALSTPNWaveform(
993 REAL4Vector *signalvec,
995
996void
999 REAL4Vector *signalvec1,
1000 REAL4Vector *signalvec2,
1001 REAL4Vector *a,
1002 REAL4Vector *ff,
1003 REAL8Vector *phi,
1004 REAL4Vector *shift,
1005 UINT4 *countback,
1007 InspiralInit *paramsInit
1008 );
1009
1010void
1013 REAL4Vector *signalvec,
1015
1016int
1018 REAL4Vector *signalvec,
1020
1021void
1024 REAL4Vector *signalvec1,
1025 REAL4Vector *signalvec2,
1027
1028int
1030 REAL4Vector *signalvec1,
1031 REAL4Vector *signalvec2,
1033
1034void
1037 CoherentGW *waveform,
1039 PPNParamStruc *ppnParams
1040 );
1041
1042int
1044 CoherentGW *waveform,
1046 PPNParamStruc *ppnParams
1047 );
1048
1049/* Phen-Spin waveform functions*/
1050
1051
1053 REAL4Vector *signalvec,
1055
1057 REAL4Vector *signalvec1,
1058 REAL4Vector *signalvec2,
1060 );
1061
1063 REAL4Vector * signalvec,
1065
1067 CoherentGW *waveform,
1069 PPNParamStruc *ppnParams);
1070
1071/* Phenomenological waveform generation functions */
1072
1074 REAL4Vector *signalvec,
1076
1078 REAL4Vector *signalvec,
1080
1082 REAL4Vector *signalvec1,
1083 REAL4Vector *signalvec2,
1085
1087 REAL4Vector *signalvec1,
1088 REAL4Vector *signalvec2,
1090
1092 REAL4Vector *signalvec,
1093 InspiralTemplate *insp_template);
1094
1096 REAL4Vector *signalvec1,
1097 REAL4Vector *signalvec2,
1098 InspiralTemplate *insp_template);
1099
1101 REAL4Vector *signalvec1,
1102 REAL4Vector *signalvec2,
1103 REAL4Vector *h,
1104 REAL4Vector *aVec,
1105 REAL4Vector *freqVec,
1106 REAL8Vector *phiVec,
1107 UINT4 *countback,
1109
1111 CoherentGW *waveform,
1113 PPNParamStruc *ppnParams);
1114
1115
1116/* DEPRECATED: Compatibility layer for phenomenological waveforms */
1117
1120 REAL4Vector *signalvec,
1122
1125 REAL4Vector *signalvec1,
1126 REAL4Vector *signalvec2,
1128
1131 REAL4Vector *signalvec,
1132 InspiralTemplate *insp_template);
1133
1136 REAL4Vector *signalvec1,
1137 REAL4Vector *signalvec2,
1138 InspiralTemplate *insp_template);
1139
1142 REAL4Vector *signalvec1,
1143 REAL4Vector *signalvec2,
1144 REAL4Vector *h,
1145 REAL4Vector *aVec,
1146 REAL4Vector *freqVec,
1147 REAL8Vector *phiVec,
1148 UINT4 *countback,
1150
1153 CoherentGW *waveform,
1155 PPNParamStruc *ppnParams);
1156
1157/* end DEPRECATED */
1158
1159/* Reduced-spin PN templats */
1162
1164 REAL4Vector *signalvec2,
1166
1168 REAL8 spin2, UINT4 pnOrder);
1169
1170/* --- OTHER PROTOTYPES --- */
1171void LALEtaTau02(
1173 REAL8 *x,
1174 REAL8 eta,
1175 void *in);
1176
1178 REAL8 eta,
1179 void *in);
1180
1181void LALEtaTau04(
1183 REAL8 *x,
1184 REAL8 eta,
1185 void *in);
1186
1188 REAL8 eta,
1189 void *in);
1190
1192 REAL8Vector *vec1,
1193 REAL8Vector *vec2,
1194 void *params);
1195
1197 TofVIn *params);
1198
1200 REAL8 v,
1201 void *params);
1202
1204 REAL8,
1205 void *);
1206
1208 REAL8 v,
1209 expnCoeffs *ak);
1210
1211#if 0 /* DO NOT EXIST */
1212REAL8 XLALInspiralPhasing2_1PN (
1213 REAL8 v,
1214 expnCoeffs *ak);
1215#endif
1216
1218 REAL8 v,
1219 expnCoeffs *ak);
1220
1222 REAL8 v,
1223 expnCoeffs *ak);
1224
1226 REAL8 v,
1227 expnCoeffs *ak);
1228
1230 REAL8 v,
1231 expnCoeffs *ak);
1232
1234 REAL8 v,
1235 expnCoeffs *ak);
1236
1238 REAL8 v,
1239 expnCoeffs *ak);
1240
1241
1242
1243
1245 REAL8 td,
1246 expnCoeffs *ak);
1247
1248#if 0 /* DO NOT EXIST */
1249REAL8 XLALInspiralPhasing3_1PN (
1250 REAL8 td,
1251 expnCoeffs *ak);
1252#endif
1253
1255 REAL8 td,
1256 expnCoeffs *ak);
1257
1259 REAL8 td,
1260 expnCoeffs *ak);
1261
1263 REAL8 td,
1264 expnCoeffs *ak);
1265
1267 REAL8 td,
1268 expnCoeffs *ak);
1269
1271 REAL8 td,
1272 expnCoeffs *ak);
1273
1275 REAL8 td,
1276 expnCoeffs *ak);
1277
1278
1280 REAL8,
1281 void *);
1282
1284 REAL8 v,
1285 void *params
1286 );
1287
1289 REAL8 f,
1290 void *params);
1291
1292#if 0 /* DO NOT EXIST */
1293REAL8 XLALInspiralTiming2_1PN (
1294 REAL8 f,
1295 void *params);
1296#endif
1297
1299 REAL8 f,
1300 void *params);
1301
1303 REAL8 f,
1304 void *params);
1305
1307 REAL8 f,
1308 void *params);
1309
1311 REAL8 f,
1312 void *params);
1313
1315 REAL8 f,
1316 void *params);
1317
1319 REAL8 f,
1320 void *params);
1321
1323 REAL8 td,
1324 expnCoeffs *ak);
1325
1326#if 0 /* DO NOT EXIST */
1327void LALInspiralFrequency3_1PN (
1329 REAL8 *frequency,
1330 REAL8 td,
1331 expnCoeffs *ak);
1332
1333REAL8 XLALInspiralFrequency3_1PN (
1334 REAL8 td,
1335 expnCoeffs *ak);
1336#endif
1337
1339 REAL8 td,
1340 expnCoeffs *ak);
1341
1343 REAL8 td,
1344 expnCoeffs *ak);
1345
1347 REAL8 td,
1348 expnCoeffs *ak);
1349
1351 REAL8 td,
1352 expnCoeffs *ak);
1353
1355 REAL8 td,
1356 expnCoeffs *ak);
1357
1359 REAL8 td,
1360 expnCoeffs *ak);
1361
1363 REAL8 phase,
1364 REAL8 v,
1366
1368 REAL8 phase,
1369 REAL8 v,
1371
1373 INT4 n,
1374 rk4In *input);
1375
1376void LALRungeKutta4(
1377 LALStatus *,
1378 REAL8Vector *,
1380 void *);
1381
1382int
1384 REAL8Vector *yout,
1385 rk4GSLIntegrator *integrator,
1386 void *params
1387 );
1388
1390 rk4GSLIntegrator *integrator);
1391
1392/* --- TEST PROTOTYPES --- */
1393
1395 REAL4Vector *rdwave1,
1396 REAL4Vector *rdwave2,
1398 REAL4VectorSequence *inspwave1,
1399 REAL4VectorSequence *inspwave2,
1400 COMPLEX8Vector *modefreqs,
1401 REAL8Vector *matchrange
1402 );
1403
1405 REAL4Vector *rdwave1,
1406 REAL4Vector *rdwave2,
1408 REAL4VectorSequence *inspwave1,
1409 REAL4VectorSequence *inspwave2,
1410 COMPLEX8Vector *modefreqs,
1411 UINT4 nmodes
1412 );
1414 REAL4Vector *rwave,
1415 REAL4Vector *dwave,
1416 REAL4Vector *ddwave,
1417 REAL8Vector *timeVec,
1418 REAL4Vector *wave,
1419 REAL8Vector *matchrange,
1421 );
1422
1424 REAL4Vector *dwave,
1425 REAL4Vector *ddwave,
1426 REAL4Vector *wave,
1428 );
1429
1431 COMPLEX8Vector *modefreqs,
1433 UINT4 l,
1434 UINT4 m,
1435 UINT4 nmodes
1436 );
1437
1439 COMPLEX8Vector *modefreqs,
1441 UINT4 l,
1442 UINT4 m,
1443 UINT4 nmodes
1444 );
1445
1446/* made static in LALInspiralRingdownWave.c */
1447/* INT4 XLALFinalMassSpin(
1448 REAL8 *finalMass,
1449 REAL8 *finalSpin,
1450 InspiralTemplate *params
1451 ); */
1452
1454 REAL4Vector *signalvec1,
1455 REAL4Vector *signalvec2,
1456 INT4 l,
1457 INT4 m,
1458 REAL8Vector *timeVec,
1459 REAL8Vector *matchrange,
1461
1463 REAL4Vector *Omega,
1464 REAL4Vector *signalvec1,
1465 REAL4Vector *signalvec2,
1467
1468/**
1469 * XLAL function to determine adaptive integration flag from a string. Returns
1470 * 1 if string contains 'fixedStep', otherwise returns 0 to signal
1471 * adaptive integration should be used.
1472 */
1473int XLALGetAdaptiveIntFromString(const CHAR *inString);
1474
1475/**
1476 * XLAL function to determine inspiral-only flag from a string. Returns
1477 * 1 if string contains 'inspiralOnly', otherwise returns 0 to signal
1478 * full inspiral-merger-ringdown waveform should be generated.
1479 */
1480int XLALGetInspiralOnlyFromString(const CHAR *inString);
1481
1483 REAL8Vector *rdwave,
1485 REAL8Vector *inspwave,
1486 COMPLEX8Vector *modefreqs,
1487 UINT4 nmodes
1488 );
1489
1491 REAL8Vector *dwave,
1492 REAL8Vector *wave,
1493 REAL8 dt
1494 );
1495
1497 COMPLEX8Vector *modefreqs,
1499 UINT4 l,
1500 INT4 m,
1501 UINT4 nmodes,
1502 REAL8 finalMass,
1503 REAL8 finalSpin
1504 );
1505
1507 REAL8 *finalMass,
1508 REAL8 *finalSpin,
1510 REAL8 energy,
1511 REAL8 *LNhvec
1512 );
1513
1515 REAL8Vector *signalvec,
1517 UINT4 *attpos,
1518 UINT4 nmodes,
1519 UINT4 l,
1520 INT4 m,
1521 REAL8 finalMass,
1522 REAL8 finalSpin
1523);
1524
1526 UINT4 length,
1527 Approximant approx,
1528 LALPNOrder order
1529 );
1530
1532 REAL4Sequence *sequence,
1533 REAL4 fLow,
1534 REAL4 fHigh,
1535 REAL4 fSampling
1536 );
1537
1539 REAL8Vector *amp,
1540 REAL8Vector *phase,
1541 double epsilon,
1542 double alpha,
1543 double beta,
1544 double padding,
1545 COMPLEX16Vector **a1,
1546 COMPLEX16Vector **b0,
1547 INT4Vector **delay
1548 );
1549
1551 COMPLEX16Vector *a1,
1552 COMPLEX16Vector *b0,
1553 INT4Vector *delay,
1554 COMPLEX16Vector *response
1555 );
1556
1558 int j,
1559 int jmax,
1560 COMPLEX16 a1,
1561 COMPLEX16 b0,
1562 int delay,
1563 COMPLEX16 *hfcos,
1564 COMPLEX16 *hfsin
1565 );
1566
1568 COMPLEX16Vector *a1,
1569 COMPLEX16Vector *b0,
1570 INT4Vector *delay,
1571 REAL8Vector *psd,
1572 double *ip
1573 );
1574
1575int
1577 REAL8TimeSeries **hplus,
1578 REAL8TimeSeries **hcross,
1579 SimInspiralTable *thisRow,
1580 REAL8 deltaT
1581);
1582
1584 REAL4TimeSeries* chan,
1585 const char *ifo,
1586 REAL8 dynRange,
1587 SimInspiralTable* events
1588);
1589
1590/*---------------------------------------------------------------- */
1591
1592#ifdef __cplusplus
1593}
1594#endif
1595
1596#endif /* _LALINSPIRAL_H */
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)
Definition: LALEtaTau04.c:89
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)
Definition: LALEtaTau04.c:67
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)
Definition: LALEtaTau02.c:66
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)
Definition: LALEtaTau02.c:88
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)
const double b1
const double b2
#define LIGOMETA_IFO_MAX
double theta
const double a4
const double a2
double complex COMPLEX16
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
float REAL4
REAL8() InspiralPhasing2(REAL8 v, expnCoeffs *ak)
Definition: LALInspiral.h:524
REAL8 EnergyFunction(REAL8 v, expnCoeffs *ak)
Definition: LALInspiral.h:511
REAL8() InspiralTiming2(REAL8 f, void *params)
Definition: LALInspiral.h:536
InputMasses
This structure is one of the members of the InspiralTemplate structure.
Definition: LALInspiral.h:178
int XLALGetInteractionFromString(const CHAR *inString)
XLAL function to determine LALInspiralInteraction from a string.
void() TestFunction(REAL8Vector *vector1, REAL8Vector *vector2, void *params)
Definition: LALInspiral.h:519
REAL8() InspiralPhasing3(REAL8 td, expnCoeffs *ak)
Definition: LALInspiral.h:528
REAL8() InspiralFrequency3(REAL8 td, expnCoeffs *ak)
Definition: LALInspiral.h:532
REAL8 FluxFunction(REAL8 v, expnCoeffs *ak)
Definition: LALInspiral.h:515
LALInspiralInteraction
Enumeration to specify which interaction will be used in the waveform generation.
Definition: LALInspiral.h:116
@ massesAndSpin
UNDOCUMENTED.
Definition: LALInspiral.h:193
@ t01
unused; shouldn't be used
Definition: LALInspiral.h:183
@ fixedPsi
The two psi values are given by the input parameter structure.
Definition: LALInspiral.h:191
@ totalMassAndMu
total mass and reduced mass
Definition: LALInspiral.h:182
@ fixedMasses
The two masses are given by the input parameter structure.
Definition: LALInspiral.h:190
@ totalMassAndEta
total mass and symmetric mass ratio
Definition: LALInspiral.h:180
@ psi0Andpsi3
BCV parameters and .
Definition: LALInspiral.h:187
@ t03
chirptimes and , and
Definition: LALInspiral.h:185
@ bhns
One of the mass is a Neutron star and the other a black hole (m1 [minMass-3] and m2 [3-maxMass])
Definition: LALInspiral.h:189
@ t02
chirptimes and
Definition: LALInspiral.h:184
@ totalMassUAndEta
total mass and eta but uniform distribution in totalMass
Definition: LALInspiral.h:181
@ m1Andm2
component masses
Definition: LALInspiral.h:179
@ minmaxTotalMass
UNDOCUMENTED.
Definition: LALInspiral.h:194
@ fixedTau
The two tau values are given by the input parameter structure.
Definition: LALInspiral.h:192
@ spinOnly
UNDOCUMENTED.
Definition: LALInspiral.h:195
@ t04
chirptimes and
Definition: LALInspiral.h:186
@ LAL_INSPIRAL_INTERACTION_TIDAL_5PN
Leading-order tidal interaction.
Definition: LALInspiral.h:124
@ LAL_INSPIRAL_INTERACTION_ALL
all spin and tidal interactions
Definition: LALInspiral.h:127
@ LAL_INSPIRAL_INTERACTION_TIDAL_6PN
Next-to-leading-order tidal interaction.
Definition: LALInspiral.h:125
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_25PN
Next-to-leading-order spin-orbit interaction.
Definition: LALInspiral.h:122
@ LAL_INSPIRAL_INTERACTION_QUAD_MONO_2PN
Quadrupole-monopole interaction.
Definition: LALInspiral.h:121
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_3PN
Spin-spin interaction.
Definition: LALInspiral.h:123
@ LAL_INSPIRAL_INTERACTION_NONE
No spin, tidal or other interactions.
Definition: LALInspiral.h:117
@ LAL_INSPIRAL_INTERACTION_SPIN_SPIN_2PN
Spin-spin interaction.
Definition: LALInspiral.h:119
@ LAL_INSPIRAL_INTERACTION_SPIN_ORBIT_15PN
Leading order spin-orbit interaction.
Definition: LALInspiral.h:118
@ LAL_INSPIRAL_INTERACTION_ALL_SPIN
all spin interactions, no tidal interactions
Definition: LALInspiral.h:126
@ LAL_INSPIRAL_INTERACTION_SPIN_SPIN_SELF_2PN
Spin-spin-self interaction.
Definition: LALInspiral.h:120
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.
LALSimInspiralFrameAxis
Approximant
LALPNOrder
static const INT4 m
static const INT4 a
double t0
These are the input structures needed to solve for the mass ratio given the chirptimes or .
Definition: LALInspiral.h:150
This structure is needed to solve the differential equation giving the evolution of the orbital angul...
Definition: LALInspiral.h:379
REAL8 fourM1Plus
fourM1Plus: = (all masses expressed in seconds)
Definition: LALInspiral.h:386
REAL8 v
parameter of 'integration': v=sqrt(M/r)
Definition: LALInspiral.h:380
REAL8 thirtytwoBy5etc
thirtytwoBy5etc:=
Definition: LALInspiral.h:390
REAL8 M
Total mass of the binary (in seconds)
Definition: LALInspiral.h:385
REAL8 magS1
The constant spin magnitude of the primary.
Definition: LALInspiral.h:381
REAL8 fourM2Plus
fourM2Plus: = (all masses expressed in seconds)
Definition: LALInspiral.h:387
REAL8 magS2
The constant spin magnitude of the secondary.
Definition: LALInspiral.h:382
REAL8 oneBy2Mcube
oneBy2Mcube: =
Definition: LALInspiral.h:388
REAL8 threeBy2Mcube
threeBy2Mcube: =
Definition: LALInspiral.h:389
Structure used as an input to compute the derivatives needed in solving the phasing formula when the ...
Definition: LALInspiral.h:607
EnergyFunction * dEnergy
Definition: LALInspiral.h:609
FluxFunction * flux
Definition: LALInspiral.h:610
EOBNonQCCoeffs * nqcCoeffs
Definition: LALInspiral.h:612
expnCoeffs * coeffs
Definition: LALInspiral.h:611
expnFunc func
Definition: LALInspiral.h:680
expnCoeffs ak
Definition: LALInspiral.h:679
Structures used to compute the phase of the signal from the ‘beginning’, when the veolcity parameter ...
Definition: LALInspiral.h:654
EnergyFunction * dEnergy
Definition: LALInspiral.h:657
FluxFunction * flux
Definition: LALInspiral.h:658
expnCoeffs * coeffs
Definition: LALInspiral.h:659
The inspiral waveform parameter structure containing information about the waveform to be generated.
Definition: LALInspiral.h:205
REAL8 alpha4
UNDOCUMENTED.
Definition: LALInspiral.h:255
REAL8 signalAmplitude
dimensionless amplitude of the signal (input, currently unused)
Definition: LALInspiral.h:220
REAL8 startTime
starting time of the waveform (in sec); if different from zero, the waveform will start with an insta...
Definition: LALInspiral.h:222
REAL8 alpha6
UNDOCUMENTED.
Definition: LALInspiral.h:257
Approximant approximant
Post-Newtonain approximant to be used in generating the wave (input)
Definition: LALInspiral.h:208
REAL8 t5
2.5 post-Newtonian chirp time in seconds (output)
Definition: LALInspiral.h:298
INT4 ieta
parameter that tells whether the symmetric mass ratio should be set to zero in the PN expansions of ...
Definition: LALInspiral.h:226
INT4Vector * segmentIdVec
Vector of segment that have been filtered against this template needed for the LDAS implementation of...
Definition: LALInspiral.h:324
REAL8 psi3
BCV parameter .
Definition: LALInspiral.h:246
REAL8 mu
reduced mass (in solar mass) (input/output)
Definition: LALInspiral.h:304
INT4 nStartPad
Number of leading elements in the signal generation to be set to zero (input) If template is requeste...
Definition: LALInspiral.h:307
REAL8 t2
first post-Newtonian chirp time in seconds (input/output)
Definition: LALInspiral.h:295
struct tagInspiralTemplate * fine
Linked list to the next fine bank template (currently not filled by inspiral or bank codes)
Definition: LALInspiral.h:348
REAL8 inclination
Inclination of the orbit (currently not in use)
Definition: LALInspiral.h:262
LALPNOrder ampOrder
UNDOCUMENTED.
Definition: LALInspiral.h:210
REAL8 eta
symmetric mass ratio (input/output)
Definition: LALInspiral.h:291
REAL8 alpha5
UNDOCUMENTED.
Definition: LALInspiral.h:256
REAL8 totalMass
total mass of the binary in solar mass (input/output)
Definition: LALInspiral.h:292
REAL8 alpha
BCV amplitude correction factor .
Definition: LALInspiral.h:244
REAL8 alpha1
UNDOCUMENTED.
Definition: LALInspiral.h:252
REAL8 mass1
Mass of the primary in solar mass (input/output)
Definition: LALInspiral.h:211
REAL8 distance
Distance to the binary in seconds.
Definition: LALInspiral.h:219
REAL8 Theta
The 3PN unknown flux parameter; likely to be around unity; most waveform generation routines take the...
Definition: LALInspiral.h:235
REAL8 startPhase
starting phase of the waveform in radians (input)
Definition: LALInspiral.h:221
REAL8 eccentricity
initial eccentricity of the orbit (currently not in use)
Definition: LALInspiral.h:280
REAL8 fCutoff
upper frequency cutoff in Hz to be used in generating the waveform; If the last stable orbit frequenc...
Definition: LALInspiral.h:213
REAL8 sourcePhi
Azimuth angle in the direction to the source.
Definition: LALInspiral.h:268
REAL8 Zeta2
UNDOCUMENTED.
Definition: LALInspiral.h:239
long event_id
UNDOCUMENTED.
Definition: LALInspiral.h:326
REAL8 t4
second post-Newtonian chirp time in seconds (output)
Definition: LALInspiral.h:297
REAL8 t3
1.5 post-Newtonian chirp time in seconds (input/output)
Definition: LALInspiral.h:296
UINT4 inspiralOnly
UNDOCUMENTED.
Definition: LALInspiral.h:344
REAL8 psi0
BCV parameter .
Definition: LALInspiral.h:245
REAL8 alpha2
UNDOCUMENTED.
Definition: LALInspiral.h:253
INT4 nEndPad
Number of trailing bins to be set to zero, the resulting waveform will have at least this many bins z...
Definition: LALInspiral.h:312
REAL8 t7
3.5 post-Newtonian chirp time in seconds (output)
Definition: LALInspiral.h:300
REAL8 alpha3
UNDOCUMENTED.
Definition: LALInspiral.h:254
REAL8 polarisationAngle
Definition: LALInspiral.h:269
REAL8 orbitPhi0
Initial azimuth angle of the orbit.
Definition: LALInspiral.h:264
REAL8 fFinal
final frequency reached, in units of Hz (output)
Definition: LALInspiral.h:293
REAL8 t0
Newtonain chirp time in seconds (input/output)
Definition: LALInspiral.h:294
REAL8 chirpMass
chirp mass of the binary in solar mass (output)
Definition: LALInspiral.h:290
REAL8 tC
total chirp time seconds (output)
Definition: LALInspiral.h:301
REAL8 mass2
Mass of the secondary in solar mass (mass1 need not be larger than mass2 (input/output)
Definition: LALInspiral.h:212
UINT4 fixedStep
UNDOCUMENTED.
Definition: LALInspiral.h:343
LIGOTimeGPS end_time
UNDOCUMENTED.
Definition: LALInspiral.h:325
REAL8 orbitTheta0
Initial co-latitute of the orbit.
Definition: LALInspiral.h:263
REAL4 minMatch
The minimal match specified by the user when the bank that contains this template was created.
Definition: LALInspiral.h:303
LALPNOrder order
Post-Newtonain order to be used in generating the wave (input)
Definition: LALInspiral.h:209
INT4 level
Flag used in heirarical serached to indicate if this is a coarse or a fine template.
Definition: LALInspiral.h:305
REAL8 vFinal
final velocity parameter, in units of the speed of light (output)
Definition: LALInspiral.h:316
InputMasses massChoice
The pair of (mass) parameters given (see structure defining this member for more details) (input)
Definition: LALInspiral.h:323
struct tagInspiralTemplate * next
Linked list to the next coarse bank template (currently not filled by inspiral or bank codes)
Definition: LALInspiral.h:347
REAL8 fLower
lower frequency cutoff of the detector in Hz (input)
Definition: LALInspiral.h:217
REAL8 chi
dimensionless spin of black hole (ie mass1)
Definition: LALInspiral.h:274
REAL8 OmegaS
The 3PN (unknown) parameter; calculated to be equal to zero by Damour, Jaranowski and Schaffer (input...
Definition: LALInspiral.h:315
INT4 number
Unique ID number for this template needed for the LDAS implementation of the inspiral search.
Definition: LALInspiral.h:306
REAL8 tSampling
Sampling rate in Hz (input)
Definition: LALInspiral.h:218
REAL8 beta
UNDOCUMENTED.
Definition: LALInspiral.h:251
REAL8 t6
third post-Newtonian chirp time in seconds (output)
Definition: LALInspiral.h:299
REAL8 sourceTheta
Co-latitute in the direction to the source.
Definition: LALInspiral.h:267
LALInspiralInteraction interaction
UNDOCUMENTED.
Definition: LALInspiral.h:341
REAL8 kappa
cosine of angle between spin of mass1 and orb ang mom
Definition: LALInspiral.h:275
LALSimInspiralFrameAxis axisChoice
UNDOCUMENTED.
Definition: LALInspiral.h:283
This is a structure needed by the inner workings of the inspiral wave generation code.
Definition: LALInspiral.h:357
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 ...
Definition: LALInspiral.h:669
expnCoeffs * coeffs
Definition: LALInspiral.h:672
EnergyFunction * dEnergy
Definition: LALInspiral.h:670
FluxFunction * flux
Definition: LALInspiral.h:671
This structure contains various post-Newtonian and P-approximant expansion coefficients; the meanings...
Definition: LALInspiral.h:399
REAL8 tva2
Definition: LALInspiral.h:439
REAL8 pva2
Definition: LALInspiral.h:444
REAL8 ePa1
Definition: LALInspiral.h:409
REAL8 FTa1
Definition: LALInspiral.h:424
REAL8 eTa1
Definition: LALInspiral.h:404
REAL8 EulerC
Definition: LALInspiral.h:479
REAL8 vlsoT0
Definition: LALInspiral.h:492
REAL8 vpoleP4
Definition: LALInspiral.h:495
REAL8 fTa1
Definition: LALInspiral.h:429
REAL8 ETa1
Definition: LALInspiral.h:414
REAL8 fPa1
Definition: LALInspiral.h:434
REAL8 pta2
Definition: LALInspiral.h:449
REAL8 samplinginterval
Definition: LALInspiral.h:469
REAL8 vpolePP
Definition: LALInspiral.h:496
REAL8 fta2
Definition: LALInspiral.h:454
REAL8 vlsoP0
Definition: LALInspiral.h:493
REAL8 pfa2
Definition: LALInspiral.h:459
REAL8 vlsoPP
Definition: LALInspiral.h:494
REAL8 dETa1
Definition: LALInspiral.h:419
Structure to hold the pointers to the generic functions defined above.
Definition: LALInspiral.h:546
InspiralPhasing2 * phasing2
Definition: LALInspiral.h:550
EnergyFunction * dEnergy
Definition: LALInspiral.h:547
InspiralPhasing3 * phasing3
Definition: LALInspiral.h:551
InspiralTiming2 * timing2
Definition: LALInspiral.h:549
FluxFunction * flux
Definition: LALInspiral.h:548
InspiralFrequency3 * frequency3
Definition: LALInspiral.h:552
Structure containing steps and controls for the GSL Runge-Kutta solver.
Definition: LALInspiral.h:637
gsl_odeiv_control * control
Definition: LALInspiral.h:640
gsl_odeiv_step * step
Definition: LALInspiral.h:639
const gsl_odeiv_step_type * type
Definition: LALInspiral.h:638
gsl_odeiv_evolve * evolve
Definition: LALInspiral.h:641
Structure used as an input to Runge-Kutta solver.
Definition: LALInspiral.h:620
REAL8Vector * dym
Definition: LALInspiral.h:626
REAL8 x
Definition: LALInspiral.h:622
REAL8Vector * yt
Definition: LALInspiral.h:625
REAL8Vector * dydx
Definition: LALInspiral.h:624
REAL8 h
Definition: LALInspiral.h:628
REAL8Vector * dyt
Definition: LALInspiral.h:627
INT4 n
Definition: LALInspiral.h:629
REAL8Vector * y
Definition: LALInspiral.h:623
TestFunction * function
Definition: LALInspiral.h:621