LALSimulation  5.4.0.1-fe68b98
LALSimInspiralEccentricityFD.c
Go to the documentation of this file.
1 /* Copyright (C) 2014 E. A. Huerta, Prayush Kumar
2  *
3  * This program is free software; you can redistribute it and/or modify
4  * it under the terms of the GNU General Public License as published by
5  * the Free Software Foundation; either version 2 of the License, or
6  * (at your option) any later version.
7  *
8  * This program is distributed in the hope that it will be useful,
9  * but WITHOUT ANY WARRANTY; without even the implied warranty of
10  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11  * GNU General Public License for more details.
12  *
13  * You should have received a copy of the GNU General Public License
14  * along with with program; see the file COPYING. If not, write to the
15  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
16  * MA 02110-1301 USA
17  *
18  * This code describes the gravitational radiation emitted by compact binaries with low to moderate values of eccentricity [0 < e < 0.5]. The waveform phase includes post-Newtonian corrections up to 3.5 order and eccentricity corrections up to order e^8 at each post-Newtonian order. The waveform amplitude is modeled within the restricted post-Newtonian approach.
19  */
20 
21 
22 #include <stdio.h>
23 #include <stdlib.h>
24 #include <math.h>
25 #include <complex.h>
26 #include <lal/Date.h>
27 #include <lal/FrequencySeries.h>
28 #include <lal/LALConstants.h>
29 #include <lal/LALDatatypes.h>
30 #include <lal/LALSimInspiral.h>
31 #include <lal/Units.h>
32 #include <lal/XLALError.h>
33 
34 
35 #include "gsl/gsl_matrix.h"
36 #include "gsl/gsl_vector.h"
37 #include "gsl/gsl_linalg.h"
38 #include "gsl/gsl_eigen.h"
39 #include "gsl/gsl_permutation.h"
40 #include <gsl/gsl_errno.h>
41 #include <gsl/gsl_math.h>
42 #include <gsl/gsl_roots.h>
43 #include <gsl/gsl_sf_bessel.h>
44 #include <gsl/gsl_sf_gamma.h>
45 #include <gsl/gsl_complex.h>
46 #include <gsl/gsl_complex_math.h>
47 #include <gsl/gsl_errno.h>
48 #include <gsl/gsl_spline.h>
49 #include <gsl/gsl_errno.h>
50 #include <gsl/gsl_multiroots.h>
51 #include <gsl/gsl_rng.h>
52 
54 
55 #define gamma (0.577215664901532860606512090)
56 
57 #ifdef __GNUC__
58 #define UNUSED __attribute__ ((unused))
59 #else
60 #define UNUSED
61 #endif
62 
63 
64 typedef struct
65 tagexpnCoeffsEPC {
66 
67  /*Input params*/
68  REAL8 f0, e0;
69 
70  /* Coefficients to reconstruct the phase*/
71  REAL8 a0, a1,a2,a3,a4,a5,a6,a7,a8,a9,a10,a11, a12,a13,a14, a15, a16, a17, a18, a19, a20, a21, a22, a23, a24, a25, a26, a27, a28, a29, a30, a31, a32, a33, a34, a35, a36, a37, a38, a39, a40;
72 
73  /*Coefficients to reconstruct z_real*/
74  /*plus component*/
75  REAL8 b1p, b2p, b3p, b4p, b5p, b6p, b7p, b8p, b9p, b10p, b11p, b12p, b13p, b14p, b15p, b16p, b17p, b18p, b19p, b20p, b21p, b22p, b23p, b24p, b25p, b26p, b27p, b28p, b29p;
76  /*cross component*/
77  REAL8 b1c, b2c, b3c, b4c, b5c, b6c, b7c, b8c, b9c, b10c, b11c, b12c, b13c, b14c, b15c, b16c, b17c, b18c, b19c, b20c, b21c, b22c, b23c, b24c, b25c, b26c, b27c, b28c, b29c;
78 
79  /*Coefficients to reconstruct z_imaginary*/
80  /*plus component*/
81  REAL8 d1p, d2p, d3p, d4p, d5p, d6p, d7p, d8p, d9p, d10p, d11p, d12p, d13p, d14p, d15p, d16p, d17p, d18p, d19p, d20p, d21p, d22p, d23p, d24p, d25p, d26p, d27p, d28p, d29p;
82  /*cross component*/
83  REAL8 d1c, d2c, d3c, d4c, d5c, d6c, d7c, d8c, d9c, d10c, d11c, d12c, d13c, d14c, d15c, d16c, d17c, d18c, d19c, d20c, d21c, d22c, d23c, d24c, d25c, d26c, d27c, d28c, d29c;
84  /* symmetric mass ratio, total mass, component masses*/
85  //REAL8 nu,m,m1,m2,mu,Mtotal;
86 
88 
89 static REAL8 PhaseAccTaylorF2(int k, REAL8 f, expnCoeffsEPC *ak){
90 
91  REAL8 k4, k5, k19, k19ov3, k19ov9, k19ov18, f1, f2, f3, f4, f8, fact, lo1, lo2, lo3, lo4,lfact;
92  k4=f*f*f*f;
93  k5=f*f*f*f*f;
94  k19=k5*k5*k5*k4;
95  k19ov3=cbrt(k19);
96  k19ov9=cbrt(k19ov3);
97  k19ov18=sqrt(k19ov9);
98  f8=1./(k19ov18);
99  f1=f8*f8;
100  f2=f1*f1*f1;
101  f3=f1*f1;
102  f4=f3*f3;
103  fact=ak->a1*f;
104  lo1=fact;
105  lo4=cbrt(lo1);
106  lo3=lo4*lo4;
107  lo2=lo3*lo3*lo4;
108  lfact=6.*sqrt(6.);
109 
110  switch( k )
111  {
112  case 1: return( ak->a0*(1./lo2)*( 1.+ ak->a2*f1 + ak->a3*f2 + ak->a4*f3 + ak->a5*f4) );
113  break;
114 
115  case 2: return( ak->a0*(1./lo1)*(ak->a6 + ak->a7*f1 + ak->a8*f2 + ak->a9*f3 + ak->a10*f4) );
116  break;
117 
118  case 3: return( ak->a0*(1./lo3)*(ak->a11 + ak->a12*f1 + ak->a13*f2 + ak->a14*f3 + ak->a15*f4 ) );
119  break;
120 
121  case 4: return(ak->a0*(1./lo4)*(ak->a16 + ak->a17*f1 + ak->a18*f2 + ak->a19*f3 + ak->a20*f4) );
122  break;
123 
124  case 5: return( ak->a0*(ak->a21 + ak->a21*log(lfact*lo1) +ak->a22*f1+ak->a23*f2 + ak->a24*f3+ak->a25*f4 ) );
125  break;
126 
127  case 6: return(ak->a0*(lo4)*(ak->a26 + log(4.*lo4)*(ak->a27 + ak->a28*f1 + ak->a29*f2 + ak->a30*f3 + ak->a31*f4 ) + ak->a32*f1 + ak->a33*f2 + ak->a34*f3 + ak->a35*f4 ) );
128  break;
129 
130  case 7: return(ak->a0*(lo3)*( ak->a36 + ak->a37*f1+ ak->a38*f2+ ak->a39*f3+ak->a40*f4) );
131  break;
132 
133  default: return 0;
134  break;
135 
136  }
137 }
138 
139 
141  REAL8 k4, k5, k19, k19ov3, k19ov9, k19ov18, f1,f2,f3,f4,f5,f6,f7,f8;
142  k4=f*f*f*f;
143  k5=f*f*f*f*f;
144  k19=k5*k5*k5*k4;
145  k19ov3=cbrt(k19);
146  k19ov9=cbrt(k19ov3);
147  k19ov18=sqrt(k19ov9);
148  f8=1./(k19ov18);
149  f1=f8*f8;
150  f2=f1*f1*f1;
151  f3=f1*f1;
152  f4=f3*f3;
153  f5=f2*f8;
154  f6=f3*f8;
155  f7=f1*f8;
156 
157  switch( k )
158  {
159  case 1: return(ak->b1p*f5 + ak->b2p*f6 + ak->b3p*f7 + ak->b4p*f8);
160  break;
161 
162  case 2: return(ak->b5p + ak->b6p*f1 + ak->b7p*f2 + ak->b8p*f3 + ak->b9p*f4);
163  break;
164 
165  case 3: return(ak->b10p*f5 + ak->b11p*f6 + ak->b12p*f7 + ak->b13p*f8);
166  break;
167 
168  case 4: return(ak->b14p*f1 + ak->b15p*f2 + ak->b16p*f3 + ak->b17p*f4);
169  break;
170 
171  case 5: return(ak->b18p*f5 + ak->b19p*f6 + ak->b20p*f7);
172  break;
173 
174  case 6: return(ak->b21p*f2 + ak->b22p*f3 + ak->b23p*f4);
175  break;
176 
177  case 7: return(ak->b24p*f5 + ak->b25p*f6);
178  break;
179 
180  case 8: return(ak->b26p*f2 + ak->b27p*f4);
181  break;
182 
183  case 9: return(ak->b28p*f5);
184  break;
185 
186  case 10: return(ak->b29p*f4);
187  break;
188 
189  default: return 0;
190  break;
191 
192  }
193 }
194 
196  REAL8 k4, k5, k19, k19ov3, k19ov9, k19ov18, f1,f2,f3,f4,f5,f6,f7,f8;
197  k4=f*f*f*f;
198  k5=f*f*f*f*f;
199  k19=k5*k5*k5*k4;
200  k19ov3=cbrt(k19);
201  k19ov9=cbrt(k19ov3);
202  k19ov18=sqrt(k19ov9);
203  f8=1./(k19ov18);
204  f1=f8*f8;
205  f2=f1*f1*f1;
206  f3=f1*f1;
207  f4=f3*f3;
208  f5=f2*f8;
209  f6=f3*f8;
210  f7=f1*f8;
211 
212  switch( k )
213  {
214  case 1: return(ak->b1c*f5 + ak->b2c*f6 + ak->b3c*f7 + ak->b4c*f8);
215  break;
216 
217  case 2: return(ak->b5c + ak->b6c*f1 + ak->b7c*f2 + ak->b8c*f3 + ak->b9c*f4);
218  break;
219 
220  case 3: return(ak->b10c*f5 + ak->b11c*f6 + ak->b12c*f7 + ak->b13c*f8);
221  break;
222 
223  case 4: return(ak->b14c*f1 + ak->b15c*f2 + ak->b16c*f3 + ak->b17c*f4);
224  break;
225 
226  case 5: return(ak->b18c*f5 + ak->b19c*f6 + ak->b20c*f7);
227  break;
228 
229  case 6: return(ak->b21c*f2 + ak->b22c*f3 + ak->b23c*f4);
230  break;
231 
232  case 7: return(ak->b24c*f5 + ak->b25c*f6);
233  break;
234 
235  case 8: return(ak->b26c*f2 + ak->b27c*f4);
236  break;
237 
238  case 9: return(ak->b28c*f5);
239  break;
240 
241  case 10: return(ak->b29c*f4);
242  break;
243 
244  default: return 0;
245  break;
246 
247  }
248 }
249 
250 
252  REAL8 k4, k5, k19, k19ov3, k19ov9, k19ov18, f1,f2,f3,f4,f5,f6,f7,f8;
253  k4=f*f*f*f;
254  k5=f*f*f*f*f;
255  k19=k5*k5*k5*k4;
256  k19ov3=cbrt(k19);
257  k19ov9=cbrt(k19ov3);
258  k19ov18=sqrt(k19ov9);
259  f8=1./(k19ov18);
260  f1=f8*f8;
261  f2=f1*f1*f1;
262  f3=f1*f1;
263  f4=f3*f3;
264  f5=f2*f8;
265  f6=f3*f8;
266  f7=f1*f8;
267 
268  switch( k )
269  {
270  case 1: return(ak->d1p*f5 + ak->d2p*f6 + ak->d3p*f7 + ak->d4p*f8);
271  break;
272 
273  case 2: return(ak->d5p + ak->d6p*f1 + ak->d7p*f2 + ak->d8p*f3 + ak->d9p*f4);
274  break;
275 
276  case 3: return(ak->d10p*f5 + ak->d11p*f6 + ak->d12p*f7 + ak->d13p*f8);
277  break;
278 
279  case 4: return(ak->d14p*f1 + ak->d15p*f2 + ak->d16p*f3 + ak->d17p*f4);
280  break;
281 
282  case 5: return(ak->d18p*f5 + ak->d19p*f6 + ak->d20p*f7);
283  break;
284 
285  case 6: return(ak->d21p*f2 + ak->d22p*f3 + ak->d23p*f4);
286  break;
287 
288  case 7: return(ak->d24p*f5 + ak->d25p*f6);
289  break;
290 
291  case 8: return(ak->d26p*f2 + ak->d27p*f4);
292  break;
293 
294  case 9: return(ak->d28p*f5);
295  break;
296 
297  case 10: return(ak->d29p*f4);
298  break;
299 
300  default: return 0;
301  break;
302 
303  }
304 }
305 
307  REAL8 k4, k5, k19, k19ov3, k19ov9, k19ov18, f1,f2,f3,f4,f5,f6,f7,f8;
308  k4=f*f*f*f;
309  k5=f*f*f*f*f;
310  k19=k5*k5*k5*k4;
311  k19ov3=cbrt(k19);
312  k19ov9=cbrt(k19ov3);
313  k19ov18=sqrt(k19ov9);
314  f8=1./(k19ov18);
315  f1=f8*f8;
316  f2=f1*f1*f1;
317  f3=f1*f1;
318  f4=f3*f3;
319  f5=f2*f8;
320  f6=f3*f8;
321  f7=f1*f8;
322 
323  switch( k )
324  {
325  case 1: return(ak->d1c*f5 + ak->d2c*f6 + ak->d3c*f7 + ak->d4c*f8);
326  break;
327 
328  case 2: return(ak->d5c + ak->d6c*f1 + ak->d7c*f2 + ak->d8c*f3 + ak->d9c*f4);
329  break;
330 
331  case 3: return(ak->d10c*f5 + ak->d11c*f6 + ak->d12c*f7 + ak->d13c*f8);
332  break;
333 
334  case 4: return(ak->d14c*f1 + ak->d15c*f2 + ak->d16c*f3 + ak->d17c*f4);
335  break;
336 
337  case 5: return(ak->d18c*f5 + ak->d19c*f6 + ak->d20c*f7);
338  break;
339 
340  case 6: return(ak->d21c*f2 + ak->d22c*f3 + ak->d23c*f4);
341  break;
342 
343  case 7: return(ak->d24c*f5 + ak->d25c*f6);
344  break;
345 
346  case 8: return(ak->d26c*f2 + ak->d27c*f4);
347  break;
348 
349  case 9: return(ak->d28c*f5);
350  break;
351 
352  case 10: return(ak->d29c*f4);
353  break;
354 
355  default: return 0;
356  break;
357 
358  }
359 }
360 
361 
362 static int
364  expnCoeffsEPC *ak, /**< coefficients for EPC evolution [modified] */
365  const REAL8 m1, /**< Mass of companion 1 (kg) */
366  const REAL8 m2, /**< Mass of companion 2 (kg) */
367  const REAL8 fStart, /**< Start GW frequency (Hz) */
368  const REAL8 i, /**< Polar inclination of source (rad) */
369  const REAL8 inclination_azimuth, /**< Azimuthal component of inclination angles */
370  const REAL8 e_min /**< Initial eccentricity at f_min: range [0, 0.5] */
371  )
372 {
373 
374  const REAL8 m = m1 + m2;
375  const REAL8 Mtotal = m * LAL_MTSUN_SI; /* total mass in seconds */
376  const REAL8 eta = m1 * m2 / (m * m);
377  const REAL8 e0=e_min;
378  const REAL8 inc=i;
379  const REAL8 bet=inclination_azimuth;
380  const REAL8 f0=fStart;
381 
382  /* Coefficients to reconstruct the phase*/
383 
384  ak->a0 = C0(eta); ak->b1p = z1(e0, f0, inc, bet, 1., 0.);
385  ak->a1 = C1(Mtotal); ak->b2p = z2(e0, f0, inc, bet, 1., 0.);
386  ak->a2 = C2(e0, f0); ak->b3p = z3(e0, f0, inc, bet, 1., 0.);
387  ak->a3 = C3(e0, f0); ak->b4p = z4(e0, f0, inc, bet, 1., 0.);
388  ak->a4 = C4(e0, f0); ak->b5p = z5(f0, inc, bet, 1., 0.);
389  ak->a5 = C5(e0, f0); ak->b6p = z6(e0, f0, inc, bet, 1., 0.);
390  ak->a6 = C6(eta); ak->b7p = z7(e0, f0, inc, bet, 1., 0.);
391  ak->a7 = C7(eta, e0, f0); ak->b8p = z8(e0, f0, inc, bet, 1., 0.);
392  ak->a8 = C8(eta, e0, f0); ak->b9p = z9(e0, f0, inc, bet, 1., 0.);
393  ak->a9 = C9(eta, e0, f0); ak->b10p = z10(e0, f0, inc, bet, 1., 0.);
394  ak->a10 = C10(eta, e0, f0); ak->b11p = z11(e0, f0, inc, bet, 1., 0.);
395  ak->a11 = C11(eta); ak->b12p = z12(e0, f0, inc, bet, 1., 0.);
396  ak->a12 = C12(eta, e0, f0); ak->b13p = z13(e0, f0, inc, bet, 1., 0.);
397  ak->a13 = C13(eta, e0, f0); ak->b14p = z14(e0, f0, inc, bet, 1., 0.);
398  ak->a14 = C14(eta, e0, f0); ak->b15p = z15(e0, f0, inc, bet, 1., 0.);
399  ak->a15 = C15(eta, e0, f0); ak->b16p = z16(e0, f0, inc, bet, 1., 0.);
400  ak->a16 = C16(eta); ak->b17p = z17(e0, f0, inc, bet, 1., 0.);
401  ak->a17 = C17(eta, e0, f0); ak->b18p = z18(e0, f0, inc, bet, 1., 0.);
402  ak->a18 = C18(eta, e0, f0); ak->b19p = z19(e0, f0, inc, bet, 1., 0.);
403  ak->a19 = C19(eta, e0, f0); ak->b20p = z20(e0, f0, inc, bet, 1., 0.);
404  ak->a20 = C20(eta, e0, f0); ak->b21p = z21(e0, f0, inc, bet, 1., 0.);
405  ak->a21 = C21(eta); ak->b22p = z22(e0, f0, inc, bet, 1., 0.);
406  ak->a22 = C22(eta, e0, f0); ak->b23p = z23(e0, f0, inc, bet, 1., 0.);
407  ak->a23 = C23(eta, e0, f0); ak->b24p = z24(e0, f0, inc, bet, 1., 0.);
408  ak->a24 = C24(eta, e0, f0); ak->b25p = z25(e0, f0, inc, bet, 1., 0.);
409  ak->a25 = C25(eta, e0, f0); ak->b26p = z26(e0, f0, inc, bet, 1., 0.);
410  ak->a26 = C26(eta); ak->b27p = z27(e0, f0, inc, bet, 1., 0.);
411  ak->a27 = C27(eta); ak->b28p = z28(e0, f0, inc, bet, 1., 0.);
412  ak->a28 = C28(eta, e0, f0); ak->b29p = z29(e0, f0, inc, bet, 1., 0.);
413  ak->a29 = C29(eta, e0, f0); ak->d1p = q1(e0, f0, inc, bet, 1., 0.);
414  ak->a30 = C30(eta, e0, f0); ak->d2p = q2(e0, f0, inc, bet,1., 0.);
415  ak->a31 = C31(eta, e0, f0); ak->d3p = q3(e0, f0, inc, bet, 1., 0.);
416  ak->a32 = C32(eta, e0, f0); ak->d4p = q4(e0, f0, inc, bet, 1., 0.);
417  ak->a33 = C33(eta, e0, f0); ak->d5p = q5(f0, inc, bet, 1., 0.);
418  ak->a34 = C34(eta, e0, f0); ak->d6p = q6(e0, f0, inc, bet, 1., 0.);
419  ak->a35 = C35(eta, e0, f0); ak->d7p = q7(e0, f0, inc, bet, 1., 0.);
420  ak->a36 = C36(eta); ak->d8p = q8(e0, f0, inc, bet, 1., 0.);
421  ak->a37 = C37(eta, e0, f0); ak->d9p = q9(e0, f0, inc, bet, 1., 0.);
422  ak->a38 = C38(eta, e0, f0); ak->d10p = q10(e0, f0, inc, bet, 1., 0.);
423  ak->a39 = C39(eta, e0, f0); ak->d11p = q11(e0, f0, inc, bet, 1., 0.);
424  ak->a40 = C40(eta, e0, f0); ak->d12p = q12(e0, f0, inc, bet, 1., 0.);
425  ak->d13p = q13(e0, f0, inc, bet, 1., 0.); ak->d14p = q14(e0, f0, inc, bet, 1., 0.);
426  ak->d15p = q15(e0, f0, inc, bet, 1., 0.); ak->d16p = q16(e0, f0, inc, bet, 1., 0.);
427  ak->d17p = q17(e0, f0, inc, bet, 1., 0.); ak->d18p = q18(e0, f0, inc, bet, 1., 0.);
428  ak->d19p = q19(e0, f0, inc, bet, 1., 0.); ak->d20p = q20(e0, f0, inc, bet, 1., 0.);
429  ak->d21p = q21(e0, f0, inc, bet, 1., 0.); ak->d22p = q22(e0, f0, inc, bet, 1., 0.);
430  ak->d23p = q23(e0, f0, inc, bet, 1., 0.); ak->d24p = q24(e0, f0, inc, bet, 1., 0.);
431  ak->d25p = q25(e0, f0, inc, bet, 1., 0.); ak->d26p = q26(e0, f0, inc, bet, 1., 0.);
432  ak->d27p = q27(e0, f0, inc, bet, 1., 0.); ak->d28p = q28(e0, f0, inc, bet, 1., 0.);
433  ak->d29p = q29(e0, f0, inc, bet, 1., 0.);
434 
435  ak->b1c = z1(e0, f0, inc, bet, 0., 1.); ak->b2c = z2(e0, f0, inc, bet, 0., 1.);
436  ak->b3c = z3(e0, f0, inc, bet, 0., 1.); ak->b4c = z4(e0, f0, inc, bet, 0., 1.);
437  ak->b5c = z5(f0, inc, bet, 0., 1.); ak->b6c = z6(e0, f0, inc, bet, 0., 1.);
438  ak->b7c = z7(e0, f0, inc, bet, 0., 1.); ak->b8c = z8(e0, f0, inc, bet, 0., 1.);
439  ak->b9c = z9(e0, f0, inc, bet, 0., 1.); ak->b10c = z10(e0, f0, inc, bet, 0., 1.);
440  ak->b11c = z11(e0, f0, inc, bet, 0., 1.); ak->b12c = z12(e0, f0, inc, bet, 0., 1.);
441  ak->b13c = z13(e0, f0, inc, bet, 0., 1.); ak->b14c = z14(e0, f0, inc, bet, 0., 1.);
442  ak->b15c = z15(e0, f0, inc, bet, 0., 1.); ak->b16c = z16(e0, f0, inc, bet, 0., 1.);
443  ak->b17c = z17(e0, f0, inc, bet, 0., 1.); ak->b18c = z18(e0, f0, inc, bet, 0., 1.);
444  ak->b19c = z19(e0, f0, inc, bet, 0., 1.); ak->b20c = z20(e0, f0, inc, bet, 0., 1.);
445  ak->b21c = z21(e0, f0, inc, bet, 0., 1.); ak->b22c = z22(e0, f0, inc, bet, 0., 1.);
446  ak->b23c = z23(e0, f0, inc, bet, 0., 1.); ak->b24c = z24(e0, f0, inc, bet, 0., 1.);
447  ak->b25c = z25(e0, f0, inc, bet, 0., 1.); ak->b26c = z26(e0, f0, inc, bet, 0., 1.);
448  ak->b27c = z27(e0, f0, inc, bet, 0., 1.); ak->b28c = z28(e0, f0, inc, bet, 0., 1.);
449  ak->b29c = z29(e0, f0, inc, bet, 0., 1.);
450 
451  ak->d1c = q1(e0, f0, inc, bet, 0., 1.); ak->d2c = q2(e0, f0, inc, bet, 0., 1.);
452  ak->d3c = q3(e0, f0, inc, bet, 0., 1.); ak->d4c = q4(e0, f0, inc, bet, 0., 1.);
453  ak->d5c = q5(f0, inc, bet, 0., 1.); ak->d6c = q6(e0, f0, inc, bet, 0., 1.);
454  ak->d7c = q7(e0, f0, inc, bet, 0., 1.); ak->d8c = q8(e0, f0, inc, bet, 0., 1.);
455  ak->d9c = q9(e0, f0, inc, bet, 0., 1.); ak->d10c = q10(e0, f0, inc, bet, 0., 1.);
456  ak->d11c =q11(e0, f0, inc, bet, 0., 1.); ak->d12c = q12(e0, f0, inc, bet, 0., 1.);
457  ak->d13c = q13(e0, f0, inc, bet, 0., 1.); ak->d14c = q14(e0, f0, inc, bet, 0., 1.);
458  ak->d15c = q15(e0, f0, inc, bet, 0., 1.); ak->d16c = q16(e0, f0, inc, bet, 0., 1.);
459  ak->d17c = q17(e0, f0, inc, bet, 0., 1.); ak->d18c = q18(e0, f0, inc, bet, 0., 1.);
460  ak->d19c = q19(e0, f0, inc, bet, 0., 1.); ak->d20c = q20(e0, f0, inc, bet, 0., 1.);
461  ak->d21c = q21(e0, f0, inc, bet, 0., 1.); ak->d22c = q22(e0, f0, inc, bet, 0., 1.);
462  ak->d23c = q23(e0, f0, inc, bet, 0., 1.); ak->d24c = q24(e0, f0, inc, bet, 0., 1.);
463  ak->d25c = q25(e0, f0, inc, bet, 0., 1.); ak->d26c = q26(e0, f0, inc, bet, 0., 1.);
464  ak->d27c = q27(e0, f0, inc, bet, 0., 1.); ak->d28c = q28(e0, f0, inc, bet, 0., 1.);
465  ak->d29c = q29(e0, f0, inc, bet, 0., 1.);
466 
467  return 0;
468 
469 }
470 
471 
472 static REAL8 Heaviside( REAL8 x){
473 if( x >= 0 )
474  return 1;
475 
476 return 0;
477 }
478 
479 
480 
481 /**
482  * @addtogroup LALSimInspiralEccentricityFD_c
483  * @brief Routines to generate frequency-domain eccentric inspiral waveforms.
484  * @{
485  */
486 
487 
489  COMPLEX16FrequencySeries **hptilde, /**< FD plus polarization */
490  COMPLEX16FrequencySeries **hctilde, /**< FD cross polarization */
491  const REAL8 phiRef, /**< Orbital coalescence phase (rad) */
492  const REAL8 deltaF, /**< Frequency resolution */
493  const REAL8 m1_SI, /**< Mass of companion 1 (kg) */
494  const REAL8 m2_SI, /**< Mass of companion 2 (kg) */
495  const REAL8 fStart, /**< Start GW frequency (Hz) */
496  const REAL8 fEnd, /**< Highest GW frequency (Hz): end at Schwarzschild ISCO */
497  const REAL8 i, /**< Polar inclination of source (rad) */
498  const REAL8 r, /**< Distance of source (m) */
499  const REAL8 inclination_azimuth, /**< Azimuthal component of inclination angles [0, 2 M_PI]*/
500  const REAL8 e_min, /**< Initial eccentricity at frequency f_min: range [0, 0.4] */
501  const INT4 phaseO /**< Twice PN phase order */
502  )
503 
504  {
505 
506  gsl_complex cphase, exphase, czeta_FPlus,czeta_FCross, czeta_FPlus_times_exphase,czeta_FCross_times_exphase;
507 
508 
509  const REAL8 m1 = m1_SI / LAL_MSUN_SI;
510  const REAL8 m2 = m2_SI / LAL_MSUN_SI;
511  const REAL8 m = m1 + m2;
512  const REAL8 Mtotal = m * LAL_MTSUN_SI; /* total mass in seconds */
513  const REAL8 eta = m1 * m2 / (m * m);
514  const REAL8 piM = LAL_PI * Mtotal;
515  const REAL8 vISCO = 1. / sqrt(6.);
516  const REAL8 fISCO = vISCO * vISCO * vISCO / piM;
517  const REAL8 fupper = 2.*fISCO;
518  const REAL8 mchirp= pow(eta, 3./5.)*Mtotal;
519  REAL8 shft, f_max;
520  REAL8 f;
521  const REAL8 zr=0.;
522  REAL8 Amplitude, zim, Phaseorder, pc_re_hplus, pc_im_hplus, pc_re_hcross, pc_im_hcross, re_hplus, im_hplus, re_hcross, im_hcross;
523  size_t j, n, jStart;
524 
525  expnCoeffsEPC ak;
526 
527  expnCoeffsEPC* ak_ptr = &ak;
528  EPCSetup( ak_ptr, m1, m2, fStart, i, inclination_azimuth, e_min);
529 
530 
531  COMPLEX16 *data_p = NULL;
532  COMPLEX16 *data_c = NULL;
533  LIGOTimeGPS tC = {0, 0};
534 
535  COMPLEX16FrequencySeries *htilde_p;
536  COMPLEX16FrequencySeries *htilde_c;
537 
538  /* Perform some initial checks */
539  if (!hptilde) XLAL_ERROR(XLAL_EFAULT);
540  if (*hptilde) XLAL_ERROR(XLAL_EFAULT);
541  if (!hctilde) XLAL_ERROR(XLAL_EFAULT);
542  if (*hctilde) XLAL_ERROR(XLAL_EFAULT);
543  if (m1_SI <= 0) XLAL_ERROR(XLAL_EDOM);
544  if (m2_SI <= 0) XLAL_ERROR(XLAL_EDOM);
545  if (fStart <= 0) XLAL_ERROR(XLAL_EDOM);
546  if (r <= 0) XLAL_ERROR(XLAL_EDOM);
547 
548 
549  /* allocate htilde_p and htilde_c*/
550  if ( fEnd == 0. ) // End at ISCO
551  f_max = fISCO;
552  else // End at user-specified freq.
553  f_max = fEnd;
554  n = (size_t) (f_max / deltaF + 1);
555  XLALGPSAdd(&tC, -1 / deltaF); /* coalesce at t=0 */
556 
557 
558  htilde_p = XLALCreateCOMPLEX16FrequencySeries("htilde_p: FD waveform", &tC, 0.0, deltaF, &lalStrainUnit, n);
559  if (!htilde_p) XLAL_ERROR(XLAL_EFUNC);
560  memset(htilde_p->data->data, 0, n * sizeof(COMPLEX16));
561  XLALUnitDivide(&htilde_p->sampleUnits, &htilde_p->sampleUnits, &lalSecondUnit);
562 
563  htilde_c = XLALCreateCOMPLEX16FrequencySeries("htilde_c: FD waveform", &tC, 0.0, deltaF, &lalStrainUnit, n);
564  if (!htilde_c) XLAL_ERROR(XLAL_EFUNC);
565  memset(htilde_c->data->data, 0, n * sizeof(COMPLEX16));
566  XLALUnitDivide(&htilde_c->sampleUnits, &htilde_c->sampleUnits, &lalSecondUnit);
567 
568 
569  /* extrinsic parameters*/
570  Amplitude = -sqrt(5./384.)*pow(M_PI, -2./3.)*(pow(mchirp,5./6.)/r)*LAL_MRSUN_SI/LAL_MTSUN_SI;
571  shft = LAL_TWOPI * (tC.gpsSeconds + 1e-9 * tC.gpsNanoSeconds);
572 
573  jStart = (size_t) ceil(fStart / deltaF);
574  f=jStart*deltaF;
575  data_p = htilde_p->data->data;
576  data_c = htilde_c->data->data;
577 
578  /* In order to decompose the waveform in the form h = F_+ h_+ + F_x h_x we decompose the amplitude function using a two basis decomposition. Note that the functions zeta_real and zeta_im depend on several paramenters, including F_+ and F_x. Since zeta_real and zeta_im are linear function in the antenna pattern functions, czeta_FPlus and czeta_FCross are used to extract F_+ and F_x from the waveform amplitude*/
579 
580 
581  for ( j=jStart;j<n;j++) {
582  switch (phaseO) {
583  case -1:
584  case 7:
585  pc_re_hplus=0.;
586  pc_im_hplus=0.;
587  pc_re_hcross=0.;
588  pc_im_hcross=0.;
589  Phaseorder=0.;
590 
591  for(int k=1; k<8; k++){
592  Phaseorder+=-PhaseAccTaylorF2(k, f, &ak);
593  }
594 
595  for(int lm=1;lm<11;lm++){
596 
597  zim=M_PI/4. + pow(((REAL8)lm)/2.,8./3.)*Phaseorder - shft*f + ((REAL8)lm)*phiRef;
598 
599  cphase = gsl_complex_rect (zr,zim);
600 
601  czeta_FPlus = gsl_complex_rect (((double)zeta_generic_re_plus(lm, f, &ak)),((double)zeta_generic_im_plus(lm, f, &ak)));
602 
603  czeta_FCross=gsl_complex_rect (((double)zeta_generic_re_cross(lm, f, &ak)),((double)zeta_generic_im_cross(lm, f, &ak)));
604 
605  exphase=gsl_complex_exp (cphase);
606 
607  czeta_FPlus_times_exphase = gsl_complex_mul(czeta_FPlus, exphase);
608  czeta_FCross_times_exphase = gsl_complex_mul(czeta_FCross, exphase);
609 
610 
611  re_hplus=Heaviside(((double)lm)*fupper-2.*f)*pow(((double)lm)/2., 2./3.)*GSL_REAL (czeta_FPlus_times_exphase);
612  im_hplus=Heaviside(((double)lm)*fupper-2.*f)*pow(((double)lm)/2., 2./3.)*GSL_IMAG (czeta_FPlus_times_exphase);
613 
614  re_hcross=Heaviside(((double)lm)*fupper-2.*f)*pow(((double)lm)/2., 2./3.)*GSL_REAL (czeta_FCross_times_exphase);
615  im_hcross=Heaviside(((double)lm)*fupper-2.*f)*pow(((double)lm)/2., 2./3.)*GSL_IMAG (czeta_FCross_times_exphase);
616 
617 
618  pc_re_hplus+=re_hplus;
619  pc_im_hplus+=im_hplus;
620 
621  pc_re_hcross+=re_hcross;
622  pc_im_hcross+=im_hcross;
623  }
624  break;
625  default:
626  XLAL_ERROR(XLAL_ETYPE, "Invalid phase PN order %d", phaseO);
627  }
628 
629 
630 
631  /*Note that h(f)= FPlus*data_p + FCross*data_c, where the polarizations are given by:*/
632 
633  data_p[j] = Amplitude*pow(f,-7./6.)*(pc_re_hplus + pc_im_hplus*1.0j);
634  data_c[j] = Amplitude*pow(f,-7./6.)*(pc_re_hcross + pc_im_hcross*1.0j);
635 
636  f+=deltaF;
637 }
638 
639  *hptilde = htilde_p;
640  *hctilde = htilde_c;
641  return XLAL_SUCCESS;
642 
643 }
644 
645 /** @} */
static REAL8 Heaviside(REAL8 x)
static int EPCSetup(expnCoeffsEPC *ak, const REAL8 m1, const REAL8 m2, const REAL8 fStart, const REAL8 i, const REAL8 inclination_azimuth, const REAL8 e_min)
static REAL8 zeta_generic_re_plus(int k, REAL8 f, expnCoeffsEPC *ak)
static REAL8 zeta_generic_re_cross(int k, REAL8 f, expnCoeffsEPC *ak)
static REAL8 PhaseAccTaylorF2(int k, REAL8 f, expnCoeffsEPC *ak)
static REAL8 zeta_generic_im_cross(int k, REAL8 f, expnCoeffsEPC *ak)
static REAL8 zeta_generic_im_plus(int k, REAL8 f, expnCoeffsEPC *ak)
static REAL8 UNUSED z8(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C25(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C32(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q24(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q16(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C22(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C24(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z19(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q29(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z10(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C11(REAL8 UNUSED eta)
static REAL8 UNUSED z20(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q8(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z23(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q23(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C13(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C28(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C27(REAL8 UNUSED eta)
static REAL8 UNUSED q5(REAL8 UNUSED f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q3(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z16(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C39(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z18(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C15(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C19(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C36(REAL8 eta)
static REAL8 UNUSED C31(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C37(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C8(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q21(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C35(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C23(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C17(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q6(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C7(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C2(REAL8 e0, REAL8 f0)
static REAL8 UNUSED z4(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q12(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z15(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q20(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q13(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C1(REAL8 Mtotal)
static REAL8 UNUSED C4(REAL8 e0, REAL8 f0)
static REAL8 UNUSED C10(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q14(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q22(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q2(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q10(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z24(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q9(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z21(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C38(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z6(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q1(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C9(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C16(REAL8 eta)
static REAL8 UNUSED C5(REAL8 e0, REAL8 f0)
static REAL8 UNUSED C20(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z26(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C26(REAL8 eta)
static REAL8 UNUSED q11(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q4(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z12(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q18(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z29(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q17(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z14(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z17(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q15(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C3(REAL8 e0, REAL8 f0)
static REAL8 UNUSED q25(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z7(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C34(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED C14(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z13(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z25(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z11(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C21(REAL8 eta)
static REAL8 UNUSED q19(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C0(REAL8 eta)
static REAL8 UNUSED q7(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z9(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z1(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED q27(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C30(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z5(REAL8 UNUSED f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z22(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C40(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z3(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C12(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED z27(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z28(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C18(REAL8 eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q28(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C29(REAL8 UNUSED eta, REAL8 e0, REAL8 f0)
static REAL8 UNUSED q26(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED z2(REAL8 e0, REAL8 f0, REAL8 inc, REAL8 bet, REAL8 FPlus, REAL8 FCross)
static REAL8 UNUSED C6(REAL8 eta)
static REAL8 UNUSED C33(REAL8 eta, REAL8 e0, REAL8 f0)
double i
Definition: bh_ringdown.c:118
double e
Definition: bh_ringdown.c:117
const double a4
const double a2
COMPLEX16FrequencySeries * XLALCreateCOMPLEX16FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
#define LAL_MSUN_SI
#define LAL_PI
#define LAL_TWOPI
#define LAL_MTSUN_SI
#define LAL_MRSUN_SI
double complex COMPLEX16
double REAL8
int32_t INT4
int XLALSimInspiralEFD(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8 phiRef, const REAL8 deltaF, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 fStart, const REAL8 fEnd, const REAL8 i, const REAL8 r, const REAL8 inclination_azimuth, const REAL8 e_min, const INT4 phaseO)
static const INT4 r
static const INT4 m
const LALUnit lalStrainUnit
const LALUnit lalSecondUnit
LALUnit * XLALUnitDivide(LALUnit *output, const LALUnit *unit1, const LALUnit *unit2)
#define XLAL_ERROR(...)
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
XLAL_ETYPE
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
list x
COMPLEX16Sequence * data
COMPLEX16 * data
INT4 gpsNanoSeconds
double f_max
Definition: unicorn.c:23