Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-8a6b96f
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimIMREOBFactorizedWaveform.c
Go to the documentation of this file.
1/*
2* Copyright (C) 2010 Craig Robinson, Yi Pan
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
21/**
22 * \author Craig Robinson, Yi Pan
23 *
24 * \brief The functions contained within this file pre-compute the various
25 * coefficients which are required for calculating the factorized waveform
26 * in EOBNRv2. Note that for some of the higher modes, the coefficients
27 * are changed for the generation of the waveform compared to the generation
28 * of the flux. Thus we have a function which adds these additional
29 * contributions to the already computed coefficients.
30 */
31
32#include <math.h>
33#include <complex.h>
34#include "LALSimIMREOBNRv2.h"
35
36/* Include static functions */
40
41#ifndef _LALSIMIMRFACTORIZEDWAVEFORM_C
42#define _LALSIMIMRFACTORIZEDWAVEFORM_C
43
44#ifdef __GNUC__
45#define UNUSED __attribute__ ((unused))
46#else
47#define UNUSED
48#endif
49
50/**
51 * Constant which comes up in some of the EOB models. Its value is
52 * (94/3 -41/32*pi*pi)
53 */
54#define ninty4by3etc 18.687902694437592603
55
56
57static inline REAL8 XLALCalculateA5( REAL8 eta );
58
59static inline REAL8 XLALCalculateA6( REAL8 eta );
60
61static
63 REAL8 eta) UNUSED;
64
65
66/**
67 * Calculates the a5 parameter in the A potential function in EOBNRv2
68 */
69static inline
70REAL8 XLALCalculateA5( const REAL8 eta /**<< Symmetric mass ratio */
71 )
72{
73 return - 5.82827 - 143.486 * eta + 447.045 * eta * eta;
74}
75
76/**
77 * Calculates the a6 parameter in the A potential function in EOBNRv2
78 */
79static inline
80REAL8 XLALCalculateA6( const REAL8 UNUSED eta /**<< Symmetric mass ratio */
81 )
82{
83 return 184.0;
84}
85
86
87/**
88 * Function to pre-compute the coefficients in the EOB A potential function
89 */
90UNUSED static
92 EOBACoefficients * const coeffs, /**<< A coefficients (populated in function) */
93 const REAL8 eta /**<< Symmetric mass ratio */
94 )
95{
96 REAL8 eta2, eta3;
97 REAL8 a4, a5, a6;
98
99 eta2 = eta*eta;
100 eta3 = eta2 * eta;
101
102 /* Note that the definitions of a5 and a6 DO NOT correspond to those in the paper */
103 /* Therefore we have to multiply the results of our a5 and a6 finctions by eta. */
104
105 a4 = ninty4by3etc * eta;
106 a5 = XLALCalculateA5( eta ) * eta;
107 a6 = XLALCalculateA6( eta ) * eta;
108
109 coeffs->n4 = -64. + 12.*a4 + 4.*a5 + a6 + 64.*eta - 4.*eta2;
110 coeffs->n5 = 32. -4.*a4 - a5 - 24.*eta;
111 coeffs->d0 = 4.*a4*a4 + 4.*a4*a5 + a5*a5 - a4*a6 + 16.*a6
112 + (32.*a4 + 16.*a5 - 8.*a6) * eta + 4.*a4*eta2 + 32.*eta3;
113 coeffs->d1 = 4.*a4*a4 + a4*a5 + 16.*a5 + 8.*a6 + (32.*a4 - 2.*a6)*eta + 32.*eta2 + 8.*eta3;
114 coeffs->d2 = 16.*a4 + 8.*a5 + 4.*a6 + (8.*a4 + 2.*a5)*eta + 32.*eta2;
115 coeffs->d3 = 8.*a4 + 4.*a5 + 2.*a6 + 32.*eta - 8.*eta2;
116 coeffs->d4 = 4.*a4 + 2.*a5 + a6 + 16.*eta - 4.*eta2;
117 coeffs->d5 = 32. - 4.*a4 - a5 - 24. * eta;
118
119 return XLAL_SUCCESS;
120}
121
122/**
123 * This function calculates the EOB A function which using the pre-computed
124 * coefficients which should already have been calculated.
125 */
126static
127REAL8 XLALCalculateEOBA( const REAL8 r, /**<< Orbital separation (in units of total mass M) */
128 EOBACoefficients * restrict coeffs /**<< Pre-computed coefficients for the A function */
129 )
130{
131
132 REAL8 r2, r3, r4, r5;
133 REAL8 NA, DA;
134
135 /* Note that this function uses pre-computed coefficients,
136 * and assumes they have been calculated. Since this is a static function,
137 * so only used here, I assume it is okay to neglect error checking
138 */
139
140 r2 = r*r;
141 r3 = r2 * r;
142 r4 = r2*r2;
143 r5 = r4*r;
144
145
146 NA = r4 * coeffs->n4
147 + r5 * coeffs->n5;
148
149 DA = coeffs->d0
150 + r * coeffs->d1
151 + r2 * coeffs->d2
152 + r3 * coeffs->d3
153 + r4 * coeffs->d4
154 + r5 * coeffs->d5;
155
156 return NA/DA;
157}
158
159/**
160 * Calculated the derivative of the EOB A function with respect to
161 * r, using the pre-computed A coefficients
162 */
163static
164REAL8 XLALCalculateEOBdAdr( const REAL8 r, /**<< Orbital separation (in units of total mass M) */
165 EOBACoefficients * restrict coeffs /**<< Pre-computed coefficients for the A function */
166 )
167{
168 REAL8 r2, r3, r4, r5;
169
170 REAL8 NA, DA, dNA, dDA, dA;
171
172 r2 = r*r;
173 r3 = r2 * r;
174 r4 = r2*r2;
175 r5 = r4*r;
176
177 NA = r4 * coeffs->n4
178 + r5 * coeffs->n5;
179
180 DA = coeffs->d0
181 + r * coeffs->d1
182 + r2 * coeffs->d2
183 + r3 * coeffs->d3
184 + r4 * coeffs->d4
185 + r5 * coeffs->d5;
186
187 dNA = 4. * coeffs->n4 * r3
188 + 5. * coeffs->n5 * r4;
189
190 dDA = coeffs->d1
191 + 2. * coeffs->d2 * r
192 + 3. * coeffs->d3 * r2
193 + 4. * coeffs->d4 * r3
194 + 5. * coeffs->d5 * r4;
195
196 dA = dNA * DA - dDA * NA;
197
198 return dA / (DA*DA);
199}
200
201/**
202 * Calculate the EOB D function.
203 */
204static REAL8 XLALCalculateEOBD( REAL8 r, /**<< Orbital separation (in units of total mass M) */
205 REAL8 eta /**<< Symmetric mass ratio */
206 )
207{
208 REAL8 u, u2, u3;
209
210 u = 1./r;
211 u2 = u*u;
212 u3 = u2*u;
213
214 return 1./(1.+6.*eta*u2+2.*eta*(26.-3.*eta)*u3);
215}
216
217
218/**
219 * Function to calculate the EOB effective Hamiltonian for the
220 * given values of the dynamical variables. The coefficients in the
221 * A potential function should already have been computed.
222 * Note that the pr used here is the tortoise co-ordinate.
223 */
224static
225REAL8 XLALEffectiveHamiltonian( const REAL8 eta, /**<< Symmetric mass ratio */
226 const REAL8 r, /**<< Orbital separation */
227 const REAL8 pr, /**<< Tortoise co-ordinate */
228 const REAL8 pp, /**<< Momentum pphi */
229 EOBACoefficients *aCoeffs /**<< Pre-computed coefficients in A function */
230 )
231{
232
233 /* The pr used in here is the tortoise co-ordinate */
234 REAL8 r2, pr2, pp2, z3, eoba;
235
236 r2 = r * r;
237 pr2 = pr * pr;
238 pp2 = pp * pp;
239
240 eoba = XLALCalculateEOBA( r, aCoeffs );
241 z3 = 2. * ( 4. - 3. * eta ) * eta;
242 return sqrt( pr2 + eoba * ( 1. + pp2/r2 + z3*pr2*pr2/r2 ) );
243}
244
245
246/**
247 * Function which calculates the various coefficients used in the generation
248 * of the factorized waveform. These coefficients depend only on the symmetric
249 * mass ratio eta. It should be noted that this function calculates the
250 * coefficients used in calculating the flux. For generating the waveforms
251 * themselves, the coefficients have additional terms added which are calculated
252 * using XLALModifyFacWaveformCoefficients(). THe non-spinning parts of these
253 * coefficients can be found in Pan et al, arXiv:1106.1021v1 [gr-qc].
254 */
256 FacWaveformCoeffs * const coeffs, /**<< Structure containing coefficients (populated in function) */
257 const REAL8 eta /**<< Symmetric mass ratio */
258 )
259{
260
261 REAL8 eta2 = eta*eta;
262 REAL8 eta3 = eta2 * eta;
263
264 REAL8 dM, dM2; //dM3;
265
266 REAL8 a = 0;
267 REAL8 a2 = 0;
268 REAL8 a3 = 0;
269 REAL8 chiS = 0;
270 REAL8 chiA = 0;
271
272 /* Combination which appears a lot */
273 REAL8 m1Plus3eta, m1Plus3eta2, m1Plus3eta3;
274
275 dM2 = 1. - 4.*eta;
276
277 /* Check that deltaM has a reasonable value */
278 if ( dM2 < 0 )
279 {
280 XLALPrintError( "eta seems to be < 0.25 - this isn't allowed!\n" );
282 }
283
284 dM = sqrt( dM2 );
285 //dM3 = dM2 * dM;
286
287 m1Plus3eta = - 1. + 3.*eta;
288 m1Plus3eta2 = m1Plus3eta * m1Plus3eta;
289 m1Plus3eta3 = m1Plus3eta * m1Plus3eta2;
290
291 /* Initialize all coefficients to zero */
292 /* This is important, as we will not set some if dM is zero */
293 memset( coeffs, 0, sizeof( FacWaveformCoeffs ) );
294
295
296 /* l = 2 */
297
298 coeffs->delta22vh3 = 7./3.;
299 coeffs->delta22vh6 = (-4.*a)/3. + (428.*LAL_PI)/105.;
300 coeffs->delta22v8 = (20.*a)/63.;
301 coeffs->delta22vh9 = -2203./81. + (1712.*LAL_PI*LAL_PI)/315.;
302 coeffs->delta22v5 = - 24.*eta;
303
304 coeffs->rho22v2 = -43./42. + (55.*eta)/84.;
305 coeffs->rho22v3 = (-2.*(chiS + chiA*dM - chiS*eta))/3.;
306 coeffs->rho22v4 = -20555./10584. + (chiS*chiS + 2.*chiA*chiS*dM + chiA*chiA*dM2)/2.
307 - (33025.*eta)/21168. + (19583.*eta2)/42336.;
308 coeffs->rho22v5 = (-34.*a)/21.;
309 coeffs->rho22v6 = 1556919113./122245200. + (89.*a2)/252. - (48993925.*eta)/9779616.
310 - (6292061.*eta2)/3259872. + (10620745.*eta3)/39118464.
311 + (41.*eta*LAL_PI*LAL_PI)/192.;
312 coeffs->rho22v6l = - 428./105.;
313 coeffs->rho22v7 = (18733.*a)/15876. + a*a2/3.;
314 coeffs->rho22v8 = -387216563023./160190110080. + (18353.*a2)/21168. - a2*a2/8.;
315 coeffs->rho22v8l = 9202./2205.;
316 coeffs->rho22v10 = -16094530514677./533967033600.;
317 coeffs->rho22v10l = 439877./55566.;
318
319
320 if ( dM2 )
321 {
322 coeffs->delta21vh3 = 2./3.;
323 coeffs->delta21vh6 = (-17.*a)/35. + (107.*LAL_PI)/105.;
324 coeffs->delta21vh7 = (3.*a2)/140.;
325 coeffs->delta21vh9 = -272./81. + (214.*LAL_PI*LAL_PI)/315.;
326 coeffs->delta21v5 = - 493. * eta /42.;
327
328 coeffs->rho21v1 = (-3.*(chiS+chiA/dM))/(4.);
329 //coeffs->rho21v2 = -59./56 - (9.*chiAPlusChiSdM*chiAPlusChiSdM)/(32.*dM2) + (23.*eta)/84.;
330 /*coeffs->rho21v3 = (-567.*chiA*chiA*chiA - 1701.*chiA*chiA*chiS*dM
331 + chiA*(-4708. + 1701.*chiS*chiS - 2648.*eta)*(-1. + 4.*eta)
332 + chiS* dM3 *(4708. - 567.*chiS*chiS
333 + 1816.*eta))/(2688.*dM3);*/
334 coeffs->rho21v2 = -59./56. + (23.*eta)/84. - 9./32.*a2;
335 coeffs->rho21v3 = 1177./672.*a - 27./128.*a3;
336 coeffs->rho21v4 = -47009./56448.- (865.*a2)/1792. - (405.*a2*a2)/2048. - (10993.*eta)/14112.
337 + (617.*eta2)/4704.;
338 coeffs->rho21v5 = (-98635.*a)/75264. + (2031.*a*a2)/7168. - (1701.*a2*a3)/8192.;
339 coeffs->rho21v6 = 7613184941./2607897600.+ (9032393.*a2)/1806336. + (3897.*a2*a2)/16384.
340 - (15309.*a3*a3)/65536.;
341 coeffs->rho21v6l = - 107./105.;
342 coeffs->rho21v7 = (-3859374457.*a)/1159065600. - (55169.*a3)/16384.
343 + (18603.*a2*a3)/65536. - (72171.*a2*a2*a3)/262144.;
344 coeffs->rho21v7l = 107.*a/140.;
345 coeffs->rho21v8 = -1168617463883./911303737344.;
346 coeffs->rho21v8l = 6313./5880.;
347 coeffs->rho21v10 = -63735873771463./16569158860800.;
348 coeffs->rho21v10l = 5029963./5927040.;
349 }
350
351 /* l = 3 */
352 if ( dM2 )
353 {
354 coeffs->delta33vh3 = 13./10.;
355 coeffs->delta33vh6 = (-81.*a)/20. + (39.*LAL_PI)/7.;
356 coeffs->delta33vh9 = -227827./3000. + (78.*LAL_PI*LAL_PI)/7.;
357 coeffs->delta33v5 = - 80897.*eta / 2430.;
358
359 coeffs->rho33v2 = -7./6. + (2.*eta)/3.;
360 coeffs->rho33v3 = (chiS*dM*(-4. + 5.*eta) + chiA*(-4. + 19.*eta))/(6.*dM);
361 coeffs->rho33v4 = -6719./3960. + a2/2. - (1861.*eta)/990. + (149.*eta2)/330.;
362 coeffs->rho33v5 = (-4.*a)/3.;
363 coeffs->rho33v6 = 3203101567./227026800. + (5.*a2)/36.;
364 coeffs->rho33v6l = - 26./7.;
365 coeffs->rho33v7 = (5297.*a)/2970. + a*a2/3.;
366 coeffs->rho33v8 = -57566572157./8562153600.;
367 coeffs->rho33v8l = 13./3.;
368 }
369
370 coeffs->delta32vh3 = (10. + 33.*eta)/(-15.*m1Plus3eta);
371 coeffs->delta32vh4 = 4.*a;
372 coeffs->delta32vh6 = (-136.*a)/45. + (52.*LAL_PI)/21.;
373 coeffs->delta32vh9 = -9112./405. + (208.*LAL_PI*LAL_PI)/63.;
374
375 coeffs->rho32v = (4.*chiS*eta)/(-3.*m1Plus3eta);
376 coeffs->rho32v2 = (-4.*a2*eta2)/(9.*m1Plus3eta2) + (328. - 1115.*eta
377 + 320.*eta2)/(270.*m1Plus3eta);
378 coeffs->rho32v3 = (2.*(45.*a*m1Plus3eta3
379 - a*eta*(328. - 2099.*eta + 5.*(733. + 20.*a2)*eta2
380 - 960.*eta3)))/(405.*m1Plus3eta3);
381 coeffs->rho32v4 = a2/3. + (-1444528.
382 + 8050045.*eta - 4725605.*eta2 - 20338960.*eta3
383 + 3085640.*eta2*eta2)/(1603800.*m1Plus3eta2);
384 coeffs->rho32v5 = (-2788.*a)/1215.;
385 coeffs->rho32v6 = 5849948554./940355325. + (488.*a2)/405.;
386 coeffs->rho32v6l = - 104./63.;
387 coeffs->rho32v8 = -10607269449358./3072140846775.;
388 coeffs->rho32v8l = 17056./8505.;
389
390 if ( dM2 )
391 {
392 coeffs->delta31vh3 = 13./30.;
393 coeffs->delta31vh6 = (61.*a)/20. + (13.*LAL_PI)/21.;
394 coeffs->delta31vh7 = (-24.*a2)/5.;
395 coeffs->delta31vh9 = -227827./81000. + (26.*LAL_PI*LAL_PI)/63.;
396 coeffs->delta31v5 = - 17.*eta/10.;
397
398 coeffs->rho31v2 = -13./18. - (2.*eta)/9.;
399 coeffs->rho31v3 = (chiA*(-4. + 11.*eta) + chiS*dM*(-4. + 13.*eta))/(6.*dM);
400 coeffs->rho31v4 = 101./7128.
401 - (5.*a2)/6. - (1685.*eta)/1782. - (829.*eta2)/1782.;
402 coeffs->rho31v5 = (4.*a)/9.;
403 coeffs->rho31v6 = 11706720301./6129723600. - (49.*a2)/108.;
404 coeffs->rho31v6l = - 26./63.;
405 coeffs->rho31v7 = (-2579.*a)/5346. + a*a2/9.;
406 coeffs->rho31v8 = 2606097992581./4854741091200.;
407 coeffs->rho31v8l = 169./567.;
408 }
409
410 /* l = 4 */
411
412 coeffs->delta44vh3 = (112. + 219.*eta)/(-120.*m1Plus3eta);
413 coeffs->delta44vh6 = (-464.*a)/75. + (25136.*LAL_PI)/3465.;
414
415 coeffs->rho44v2 = (1614. - 5870.*eta + 2625.*eta2)/(1320.*m1Plus3eta);
416 coeffs->rho44v3 = (chiA*(10. - 39.*eta)*dM + chiS*(10. - 41.*eta
417 + 42.*eta2))/(15.*m1Plus3eta);
418 coeffs->rho44v4 = a2/2. + (-511573572.
419 + 2338945704.*eta - 313857376.*eta2 - 6733146000.*eta3
420 + 1252563795.*eta2*eta2)/(317116800.*m1Plus3eta2);
421 coeffs->rho44v5 = (-69.*a)/55.;
422 coeffs->rho44v6 = 16600939332793./1098809712000. + (217.*a2)/3960.;
423 coeffs->rho44v6l = - 12568./3465.;
424
425 if ( dM2 )
426 {
427 coeffs->delta43vh3 = (486. + 4961.*eta)/(810.*(1. - 2.*eta));
428 coeffs->delta43vh4 = (11.*a)/4.;
429 coeffs->delta43vh6 = 1571.*LAL_PI/385.;
430
431 coeffs->rho43v = (5.*(chiA - chiS*dM)*eta)/(8.*dM*(-1. + 2.*eta));
432 coeffs->rho43v2 = (222. - 547.*eta + 160.*eta2)/(176.*(-1. + 2.*eta));
433 coeffs->rho43v4 = -6894273./7047040. + (3.*a2)/8.;
434 coeffs->rho43v5 = (-12113.*a)/6160.;
435 coeffs->rho43v6 = 1664224207351./195343948800.;
436 coeffs->rho43v6l = - 1571./770.;
437 }
438
439 coeffs->delta42vh3 = (7.*(1. + 6.*eta))/(-15.*m1Plus3eta);
440 coeffs->delta42vh6 = (212.*a)/75. + (6284.*LAL_PI)/3465.;
441
442 coeffs->rho42v2 = (1146. - 3530.*eta + 285.*eta2)/(1320.*m1Plus3eta);
443 coeffs->rho42v3 = (chiA*(10. - 21.*eta)*dM + chiS*(10. - 59.*eta
444 + 78.*eta2))/(15.*m1Plus3eta);
445 coeffs->rho42v4 = a2/2. + (-114859044. + 295834536.*eta + 1204388696.*eta2 - 3047981160.*eta3
446 - 379526805.*eta2*eta2)/(317116800.*m1Plus3eta2);
447 coeffs->rho42v5 = (-7.*a)/110.;
448 coeffs->rho42v6 = 848238724511./219761942400. + (2323.*a2)/3960.;
449 coeffs->rho42v6l = - 3142./3465.;
450
451 if ( dM2 )
452 {
453 coeffs->delta41vh3 = (2. + 507.*eta)/(10.*(1. - 2.*eta));
454 coeffs->delta41vh4 = (11.*a)/12.;
455 coeffs->delta41vh6 = 1571.*LAL_PI/3465.;
456
457 coeffs->rho41v = (5.*(chiA - chiS*dM)*eta)/(8.*dM*(-1. + 2.*eta));
458 coeffs->rho41v2 = (602. - 1385.*eta + 288.*eta2)/(528.*(-1. + 2.*eta));
459 coeffs->rho41v4 = -7775491./21141120. + (3.*a2)/8.;
460 coeffs->rho41v5 = (-20033.*a)/55440. - (5*a*a2)/6.;
461 coeffs->rho41v6 = 1227423222031./1758095539200.;
462 coeffs->rho41v6l = - 1571./6930.;
463 }
464
465 /* l = 5 */
466 if ( dM2 )
467 {
468 coeffs->delta55vh3 = (96875. + 857528.*eta)/(131250.*(1 - 2*eta));
469
470 coeffs->rho55v2 = (487. - 1298.*eta + 512.*eta2)/(390.*(-1. + 2.*eta));
471 coeffs->rho55v3 = (-2.*a)/3.;
472 coeffs->rho55v4 = -3353747./2129400. + a2/2.;
473 coeffs->rho55v5 = - 241. * a / 195.;
474 }
475
476 coeffs->delta54vh3 = 8./15.;
477 coeffs->delta54vh4 = 12.*a/5.;
478
479 coeffs->rho54v2 = (-17448. + 96019.*eta - 127610.*eta2
480 + 33320.*eta3)/(13650.*(1. - 5.*eta + 5.*eta2));
481 coeffs->rho54v3 = (-2.*a)/15.;
482 coeffs->rho54v4 = -16213384./15526875. + (2.*a2)/5.;
483
484 if ( dM2 )
485 {
486 coeffs->delta53vh3 = 31./70.;
487
488 coeffs->rho53v2 = (375. - 850.*eta + 176.*eta2)/(390.*(-1. + 2.*eta));
489 coeffs->rho53v3 = (-2.*a)/3.;
490 coeffs->rho53v4 = -410833./709800. + a2/2.;
491 coeffs->rho53v5 = - 103.*a/325.;
492 }
493
494 coeffs->delta52vh3 = 4./15.;
495 coeffs->delta52vh4 = 6.*a/5.;
496
497 coeffs->rho52v2 = (-15828. + 84679.*eta - 104930.*eta2
498 + 21980.*eta3)/(13650.*(1. - 5.*eta + 5.*eta2));
499 coeffs->rho52v3 = (-2.*a)/15.;
500 coeffs->rho52v4 = -7187914./15526875. + (2.*a2)/5.;
501
502 if ( dM2 )
503 {
504 coeffs->delta51vh3 = 31./210.;
505
506 coeffs->rho51v2 = (319. - 626.*eta + 8.*eta2)/(390.*(-1. + 2.*eta));
507 coeffs->rho51v3 = (-2.*a)/3.;
508 coeffs->rho51v4 = -31877./304200. + a2/2.;
509 coeffs->rho51v5 = 139.*a/975.;
510 }
511
512 /* l = 6 */
513
514 coeffs->delta66vh3 = 43./70.;
515
516 coeffs->rho66v2 = (-106. + 602.*eta - 861.*eta2
517 + 273.*eta3)/(84.*(1. - 5.*eta + 5.*eta2));
518 coeffs->rho66v3 = (-2.*a)/3.;
519 coeffs->rho66v4 = -1025435./659736. + a2/2.;
520
521 if ( dM2 )
522 {
523 coeffs->delta65vh3 = 10./21.;
524
525 coeffs->rho65v2 = (-185. + 838.*eta - 910.*eta2
526 + 220.*eta3)/(144.*(dM2 + 3.*eta2));
527 coeffs->rho65v3 = - 2.*a/9.;
528 }
529
530 coeffs->delta64vh3 = 43./105.;
531
532 coeffs->rho64v2 = (-86. + 462.*eta - 581.*eta2
533 + 133.*eta3)/(84.*(1. - 5.*eta + 5.*eta2));
534 coeffs->rho64v3 = (-2.*a)/3.;
535 coeffs->rho64v4 = -476887./659736. + a2/2.;
536
537 if ( dM2 )
538 {
539 coeffs->delta63vh3 = 2./7.;
540
541 coeffs->rho63v2 = (-169. + 742.*eta - 750.*eta2
542 + 156.*eta3)/(144.*(dM2 + 3.*eta2));
543 coeffs->rho63v3 = - 2.*a/9.;
544 }
545
546 coeffs->delta62vh3 = 43./210.;
547
548 coeffs->rho62v2 = (-74. + 378.*eta - 413.*eta2
549 + 49.*eta3)/(84.*(1. - 5.*eta + 5.*eta2));
550 coeffs->rho62v3 = (-2.*a)/3.;
551 coeffs->rho62v4 = -817991./3298680. + a2/2.;
552
553 if ( dM2 )
554 {
555 coeffs->delta61vh3 = 2./21.;
556
557 coeffs->rho61v2 = (-161. + 694.*eta - 670.*eta2
558 + 124.*eta3)/(144.*(dM2 + 3.*eta2));
559 coeffs->rho61v3 = - 2. * a / 9.;
560 }
561
562 /* l = 7 */
563 if ( dM2 )
564 {
565 coeffs->delta77vh3 = 19./36.;
566
567 coeffs->rho77v2 = (-906. + 4246.*eta - 4963.*eta2
568 + 1380.*eta3)/(714.*(dM2 + 3.*eta2));
569 coeffs->rho77v3 = - 2.*a/3.;
570 }
571
572 coeffs->rho76v2 = (2144. - 16185.*eta + 37828.*eta2 - 29351.*eta3
573 + 6104.*eta2*eta2) / (1666.*(-1 + 7*eta - 14*eta2
574 + 7*eta3));
575
576 if ( dM2 )
577 {
578 coeffs->delta75vh3 = 95./252.;
579
580 coeffs->rho75v2 = (-762. + 3382.*eta - 3523.*eta2
581 + 804.*eta3)/(714.*(dM2 + 3.*eta2));
582 coeffs->rho75v3 = - 2.*a/3.;
583 }
584
585 coeffs->rho74v2 = (17756. - 131805.*eta + 298872.*eta2 - 217959.*eta3
586 + 41076.*eta2*eta2) / (14994.*(-1. + 7.*eta - 14.*eta2
587 + 7.*eta3));
588
589 if ( dM2 )
590 {
591 coeffs->delta73vh3 = 19./84.;
592
593 coeffs->rho73v2 = (-666. + 2806.*eta - 2563.*eta2
594 + 420.*eta3)/(714.*(dM2 + 3.*eta2));
595 coeffs->rho73v3 = - 2.*a/3.;
596 }
597
598 coeffs->rho72v2 = (16832. - 123489.*eta + 273924.*eta2 - 190239.*eta3
599 + 32760.*eta2*eta2) /(14994.*(-1. + 7.*eta - 14.*eta2
600 + 7.*eta3));
601
602 if ( dM2 )
603 {
604 coeffs->delta71vh3 = 19./252.;
605
606 coeffs->rho71v2 = (-618. + 2518.*eta - 2083.*eta2
607 + 228.*eta3)/(714.*(dM2 + 3.*eta2));
608 coeffs->rho71v3 = - 2.*a/3.;
609 }
610
611 /* l = 8 */
612
613 coeffs->rho88v2 = (3482. - 26778.*eta + 64659.*eta2 - 53445.*eta3
614 + 12243.*eta2*eta2) / (2736.*(-1. + 7.*eta - 14.*eta2
615 + 7.*eta3));
616
617 if ( dM2 )
618 {
619 coeffs->rho87v2 = (23478. - 154099.*eta + 309498.*eta2 - 207550.*eta3
620 + 38920*eta2*eta2) / (18240.*(-1 + 6*eta - 10*eta2
621 + 4*eta3));
622 }
623
624 coeffs->rho86v2 = (1002. - 7498.*eta + 17269.*eta2 - 13055.*eta3
625 + 2653.*eta2*eta2) / (912.*(-1. + 7.*eta - 14.*eta2
626 + 7.*eta3));
627
628 if ( dM2 )
629 {
630 coeffs->rho85v2 = (4350. - 28055.*eta + 54642.*eta2 - 34598.*eta3
631 + 6056.*eta2*eta2) / (3648.*(-1. + 6.*eta - 10.*eta2
632 + 4.*eta3));
633 }
634
635 coeffs->rho84v2 = (2666. - 19434.*eta + 42627.*eta2 - 28965.*eta3
636 + 4899.*eta2*eta2) / (2736.*(-1. + 7.*eta - 14.*eta2
637 + 7.*eta3));
638
639 if ( dM2 )
640 {
641 coeffs->rho83v2 = (20598. - 131059.*eta + 249018.*eta2 - 149950.*eta3
642 + 24520.*eta2*eta2) / (18240.*(-1. + 6.*eta - 10.*eta2
643 + 4.*eta3));
644 }
645
646 coeffs->rho82v2 = (2462. - 17598.*eta + 37119.*eta2 - 22845.*eta3
647 + 3063.*eta2*eta2) / (2736.*(-1. + 7.*eta - 14.*eta2
648 + 7.*eta3));
649
650 if ( dM2 )
651 {
652 coeffs->rho81v2 = (20022. - 126451.*eta + 236922.*eta2 - 138430.*eta3
653 + 21640.*eta2*eta2) / (18240.*(-1. + 6.*eta - 10.*eta2
654 + 4.*eta3));
655 }
656
657 /* All relevant coefficients should be set, so we return */
658
659 return XLAL_SUCCESS;
660}
661
662
663/**
664 * Function which adds the additional terms required for waveform generation
665 * to the factorized waveform coefficients. Note that this function only calculates
666 * additional terms not present in the flux, so the factorized waveform coefficients
667 * SHOULD ALREADY HAVE BEEN CALCULATED using XLALCalcFacWaveformCoefficients() prior
668 * to calling this function.
669 */
671 FacWaveformCoeffs * const coeffs, /**<< Structure containing coefficients */
672 const REAL8 eta /**<< Symmetric mass ratio */
673 )
674{
675
676 if ( !coeffs )
677 {
679 }
680
681 /* Tweak the relevant coefficients for the generation of the waveform */
682 coeffs->rho21v6 += -5. * eta;
683 coeffs->rho33v6 += -20. * eta;
684 coeffs->rho44v6 += -15. * eta;
685 coeffs->rho55v6 += 4. * eta;
686
687 coeffs->delta21v7 += 30. * eta;
688 coeffs->delta33v7 += -10. * eta;
689 coeffs->delta44v5 += -70. * eta;
690 coeffs->delta55v5 += 40. * eta;
691
692 return XLAL_SUCCESS;
693}
694
695/**
696 * Computes the non-Keplerian correction to the velocity as determined from the
697 * frequency obtained assuming a circular orbit. In the early stages of the evolution,
698 * this should be a number close to 1.
699 */
700static REAL8
702 REAL8Vector * restrict values, /**<< Dynamics r, phi, pr, pphi */
703 const REAL8 eta, /**<< Symmetric mass ratio */
704 EOBACoefficients *coeffs /**<< Pre-computed A coefficients */
705 )
706{
707
708 REAL8 r = values->data[0];
709 REAL8 pphi = values->data[3];
710
711 REAL8 A = XLALCalculateEOBA( r, coeffs );
712 REAL8 dA = XLALCalculateEOBdAdr( r, coeffs );
713
714 return 2. * (1. + 2. * eta * ( -1. + sqrt( (1. + pphi*pphi/(r*r)) * A ) ) )
715 / ( r*r * dA );
716}
717
718/**
719 * Computes the factorized waveform according to the prescription
720 * given in Pan et al, arXiv:1106.1021v1 [gr-qc], for a given
721 * mode l,m, for the given values of the dynamics at that point.
722 * The function returns XLAL_SUCCESS if everything works out properly,
723 * otherwise XLAL_FAILURE will be returned.
724 */
726 COMPLEX16 * restrict hlm, /**<< The value of hlm (populated by the function) */
727 REAL8Vector * restrict values, /**<< Vector containing dynamics r, phi, pr, pphi for a given point */
728 const REAL8 v, /**<< Velocity (in geometric units) */
729 const INT4 l, /**<< Mode l */
730 const INT4 m, /**<< Mode m */
731 EOBParams * restrict params /**<< Structure containing pre-computed coefficients, etc. */
732 )
733{
734
735 /* Status of function calls */
736 INT4 status;
737 INT4 i;
738
739 REAL8 eta;
740 REAL8 r, pr, pp, Omega, v2, vh, vh3, k, hathatk, eulerlogxabs;
741 REAL8 Hreal, Heff, Slm, deltalm, rholm, rholmPwrl;
742 COMPLEX16 Tlm;
743 COMPLEX16 hNewton;
744 gsl_sf_result lnr1, arg1, z2;
745
746 /* Non-Keplerian velocity */
747 REAL8 vPhi;
748
749 /* Pre-computed coefficients */
750 FacWaveformCoeffs *hCoeffs = params->hCoeffs;
751
752 if ( abs(m) > (INT4) l )
753 {
755 }
756
757
758 eta = params->eta;
759
760 /* Check our eta was sensible */
761 if ( eta > 0.25 )
762 {
763 XLALPrintError("Eta seems to be > 0.25 - this isn't allowed!\n" );
765 }
766 else if ( eta == 0.25 && m % 2 )
767 {
768 /* If m is odd and dM = 0, hLM will be zero */
769 memset( hlm, 0, sizeof( COMPLEX16 ) );
770 return XLAL_SUCCESS;
771 }
772
773 r = values->data[0];
774 pr = values->data[2];
775 pp = values->data[3];
776
777 Heff = XLALEffectiveHamiltonian( eta, r, pr, pp, params->aCoeffs );
778 Hreal = sqrt( 1.0 + 2.0 * eta * ( Heff - 1.0) );
779 v2 = v * v;
780 Omega = v2 * v;
781 vh3 = Hreal * Omega;
782 vh = cbrt(vh3);
783 eulerlogxabs = LAL_GAMMA + log( 2.0 * (REAL8)m * v );
784
785
786 /* Calculate the non-Keplerian velocity */
787 /* given by Eq. (18) of Pan et al, PRD84, 124052(2011) */
788 /* psi given by Eq. (19) of Pan et al, PRD84, 124052(2011) */
789 /* Assign temporarily to vPhi */
790 vPhi = nonKeplerianCoefficient( values, eta, params->aCoeffs );
791 /* Assign rOmega value temporarily to vPhi */
792 vPhi = r * cbrt(vPhi);
793 /* Assign rOmega * Omega to vPhi */
794 vPhi *= Omega;
795
796 /* Calculate the newtonian multipole */
797 status = XLALSimIMREOBCalculateNewtonianMultipole( &hNewton, vPhi * vPhi, vPhi/Omega,
798 values->data[1], (UINT4)l, m, params );
799 if ( status == XLAL_FAILURE )
800 {
802 }
803
804 /* Calculate the source term */
805 if ( ( (l+m)%2 ) == 0)
806 {
807 Slm = Heff;
808 }
809 else
810 {
811 Slm = v * pp;
812 }
813
814 /* Calculate the Tail term */
815 k = m * Omega;
816 hathatk = Hreal * k;
817 XLAL_CALLGSL( status = gsl_sf_lngamma_complex_e( l+1.0, -2.0*hathatk, &lnr1, &arg1 ) );
818 if (status != GSL_SUCCESS)
819 {
820 XLALPrintError("Error in GSL function\n" );
822 }
823 XLAL_CALLGSL( status = gsl_sf_fact_e( l, &z2 ) );
824 if ( status != GSL_SUCCESS)
825 {
826 XLALPrintError("Error in GSL function\n" );
828 }
829 Tlm = cexp( ( lnr1.val + LAL_PI * hathatk ) + I * (
830 arg1.val + 2.0 * hathatk * log(4.0*k/sqrt(LAL_E)) ) );
831 Tlm /= z2.val;
832
833 /* Calculate the residue phase and amplitude terms */
834 switch( l )
835 {
836 case 2:
837 switch( abs(m) )
838 {
839 case 2:
840 deltalm = vh3*(hCoeffs->delta22vh3 + vh3*(hCoeffs->delta22vh6
841 + vh*vh*(hCoeffs->delta22vh9*vh)))
842 + hCoeffs->delta22v5 *v*v2*v2 + hCoeffs->delta22v8 *v2*v2*v2*v2;
843 rholm = 1. + v2*(hCoeffs->rho22v2 + v*(hCoeffs->rho22v3
844 + v*(hCoeffs->rho22v4
845 + v*(hCoeffs->rho22v5 + v*(hCoeffs->rho22v6
846 + hCoeffs->rho22v6l*eulerlogxabs + v*(hCoeffs->rho22v7
847 + v*(hCoeffs->rho22v8 + hCoeffs->rho22v8l*eulerlogxabs
848 + (hCoeffs->rho22v10 + hCoeffs->rho22v10l * eulerlogxabs)*v2)))))));
849 break;
850 case 1:
851 deltalm = vh3*(hCoeffs->delta21vh3 + vh3*(hCoeffs->delta21vh6
852 + vh*(hCoeffs->delta21vh7 + (hCoeffs->delta21vh9)*vh*vh)))
853 + hCoeffs->delta21v5*v*v2*v2 + hCoeffs->delta21v7*v2*v2*v2*v;
854 rholm = 1. + v*(hCoeffs->rho21v1
855 + v*( hCoeffs->rho21v2 + v*(hCoeffs->rho21v3 + v*(hCoeffs->rho21v4
856 + v*(hCoeffs->rho21v5 + v*(hCoeffs->rho21v6 + hCoeffs->rho21v6l*eulerlogxabs
857 + v*(hCoeffs->rho21v7 + hCoeffs->rho21v7l * eulerlogxabs
858 + v*(hCoeffs->rho21v8 + hCoeffs->rho21v8l * eulerlogxabs
859 + (hCoeffs->rho21v10 + hCoeffs->rho21v10l * eulerlogxabs)*v2))))))));
860 break;
861 default:
863 break;
864 }
865 break;
866 case 3:
867 switch (m)
868 {
869 case 3:
870 deltalm = vh3*(hCoeffs->delta33vh3 + vh3*(hCoeffs->delta33vh6 + hCoeffs->delta33vh9*vh3))
871 + hCoeffs->delta33v5*v*v2*v2 + hCoeffs->delta33v7*v2*v2*v2*v;
872 rholm = 1. + v2*(hCoeffs->rho33v2 + v*(hCoeffs->rho33v3 + v*(hCoeffs->rho33v4
873 + v*(hCoeffs->rho33v5 + v*(hCoeffs->rho33v6 + hCoeffs->rho33v6l*eulerlogxabs
874 + v*(hCoeffs->rho33v7 + (hCoeffs->rho33v8 + hCoeffs->rho33v8l*eulerlogxabs)*v))))));
875 break;
876 case 2:
877 deltalm = vh3*(hCoeffs->delta32vh3 + vh*(hCoeffs->delta32vh4 + vh*vh*(hCoeffs->delta32vh6
878 + hCoeffs->delta32vh9*vh3)));
879 rholm = 1. + v*(hCoeffs->rho32v
880 + v*(hCoeffs->rho32v2 + v*(hCoeffs->rho32v3 + v*(hCoeffs->rho32v4 + v*(hCoeffs->rho32v5
881 + v*(hCoeffs->rho32v6 + hCoeffs->rho32v6l*eulerlogxabs
882 + (hCoeffs->rho32v8 + hCoeffs->rho32v8l*eulerlogxabs)*v2))))));
883 break;
884 case 1:
885 deltalm = vh3*(hCoeffs->delta31vh3 + vh3*(hCoeffs->delta31vh6
886 + vh*(hCoeffs->delta31vh7 + hCoeffs->delta31vh9*vh*vh)))
887 + hCoeffs->delta31v5*v*v2*v2;
888 rholm = 1. + v2*(hCoeffs->rho31v2 + v*(hCoeffs->rho31v3 + v*(hCoeffs->rho31v4
889 + v*(hCoeffs->rho31v5 + v*(hCoeffs->rho31v6 + hCoeffs->rho31v6l*eulerlogxabs
890 + v*(hCoeffs->rho31v7 + (hCoeffs->rho31v8 + hCoeffs->rho31v8l*eulerlogxabs)*v))))));
891 break;
892 default:
894 break;
895 }
896 break;
897 case 4:
898 switch (m)
899 {
900 case 4:
901 deltalm = vh3*(hCoeffs->delta44vh3 + hCoeffs->delta44vh6 *vh3)
902 + hCoeffs->delta44v5*v2*v2*v;
903 rholm = 1. + v2*(hCoeffs->rho44v2
904 + v*( hCoeffs->rho44v3 + v*(hCoeffs->rho44v4
905 + v*(hCoeffs->rho44v5 + (hCoeffs->rho44v6
906 + hCoeffs->rho44v6l*eulerlogxabs)*v))));
907 break;
908 case 3:
909 deltalm = vh3*(hCoeffs->delta43vh3 + vh*(hCoeffs->delta43vh4
910 + hCoeffs->delta43vh6*vh*vh));
911 rholm = 1. + v*(hCoeffs->rho43v
912 + v*(hCoeffs->rho43v2
913 + v2*(hCoeffs->rho43v4 + v*(hCoeffs->rho43v5
914 + (hCoeffs->rho43v6 + hCoeffs->rho43v6l*eulerlogxabs)*v))));
915 break;
916 case 2:
917 deltalm = vh3*(hCoeffs->delta42vh3 + hCoeffs->delta42vh6*vh3);
918 rholm = 1. + v2*(hCoeffs->rho42v2
919 + v*(hCoeffs->rho42v3 + v*(hCoeffs->rho42v4 + v*(hCoeffs->rho42v5
920 + (hCoeffs->rho42v6 + hCoeffs->rho42v6l*eulerlogxabs)*v))));
921 break;
922 case 1:
923 deltalm = vh3*(hCoeffs->delta41vh3 + vh*(hCoeffs->delta41vh4
924 + hCoeffs->delta41vh6*vh*vh));
925 rholm = 1. + v*(hCoeffs->rho41v
926 + v*(hCoeffs->rho41v2
927 + v2*(hCoeffs->rho41v4 + v*(hCoeffs->rho41v5
928 + (hCoeffs->rho41v6 + hCoeffs->rho41v6l*eulerlogxabs)*v))));
929 break;
930 default:
932 break;
933 }
934 break;
935 case 5:
936 switch (m)
937 {
938 case 5:
939 deltalm = hCoeffs->delta55vh3*vh3 + hCoeffs->delta55v5*v2*v2*v;
940 rholm = 1. + v2*( hCoeffs->rho55v2
941 + v*(hCoeffs->rho55v3 + v*(hCoeffs->rho55v4
942 + v*(hCoeffs->rho55v5 + hCoeffs->rho55v6*v))));
943 break;
944 case 4:
945 deltalm = vh3*(hCoeffs->delta54vh3 + hCoeffs->delta54vh4*vh);
946 rholm = 1. + v2*(hCoeffs->rho54v2 + v*(hCoeffs->rho54v3
947 + hCoeffs->rho54v4*v));
948 break;
949 case 3:
950 deltalm = hCoeffs->delta53vh3 * vh3;
951 rholm = 1. + v2*(hCoeffs->rho53v2
952 + v*(hCoeffs->rho53v3 + v*(hCoeffs->rho53v4 + hCoeffs->rho53v5*v)));
953 break;
954 case 2:
955 deltalm = vh3*(hCoeffs->delta52vh3 + hCoeffs->delta52vh4*vh);
956 rholm = 1. + v2*(hCoeffs->rho52v2 + v*(hCoeffs->rho52v3
957 + hCoeffs->rho52v4*v));
958 break;
959 case 1:
960 deltalm = hCoeffs->delta51vh3 * vh3;
961 rholm = 1. + v2*(hCoeffs->rho51v2
962 + v*(hCoeffs->rho51v3 + v*(hCoeffs->rho51v4 + hCoeffs->rho51v5*v)));
963 break;
964 default:
966 break;
967 }
968 break;
969 case 6:
970 switch (m)
971 {
972 case 6:
973 deltalm = hCoeffs->delta66vh3*vh3;
974 rholm = 1. + v2*(hCoeffs->rho66v2 + v*(hCoeffs->rho66v3
975 + hCoeffs->rho66v4*v));
976 break;
977 case 5:
978 deltalm = hCoeffs->delta65vh3*vh3;
979 rholm = 1. + v2*(hCoeffs->rho65v2 + hCoeffs->rho65v3*v);
980 break;
981 case 4:
982 deltalm = hCoeffs->delta64vh3 * vh3;
983 rholm = 1. + v2*(hCoeffs->rho64v2 + v*(hCoeffs->rho64v3
984 + hCoeffs->rho64v4*v));
985 break;
986 case 3:
987 deltalm = hCoeffs->delta63vh3 * vh3;
988 rholm = 1. + v2*(hCoeffs->rho63v2 + hCoeffs->rho63v3*v);
989 break;
990 case 2:
991 deltalm = hCoeffs->delta62vh3 * vh3;
992 rholm = 1. + v2*(hCoeffs->rho62v2 + v*(hCoeffs->rho62v3
993 + hCoeffs->rho62v4 * v));
994 break;
995 case 1:
996 deltalm = hCoeffs->delta61vh3 * vh3;
997 rholm = 1. + v2*(hCoeffs->rho61v2 + hCoeffs->rho61v3*v);
998 break;
999 default:
1001 break;
1002 }
1003 break;
1004 case 7:
1005 switch (m)
1006 {
1007 case 7:
1008 deltalm = hCoeffs->delta77vh3 * vh3;
1009 rholm = 1. + v2*(hCoeffs->rho77v2 + hCoeffs->rho77v3 * v);
1010 break;
1011 case 6:
1012 deltalm = 0.0;
1013 rholm = 1. + hCoeffs->rho76v2 * v2;
1014 break;
1015 case 5:
1016 deltalm = hCoeffs->delta75vh3 * vh3;
1017 rholm = 1. + v2*(hCoeffs->rho75v2 + hCoeffs->rho75v3*v);
1018 break;
1019 case 4:
1020 deltalm = 0.0;
1021 rholm = 1. + hCoeffs->rho74v2 * v2;
1022 break;
1023 case 3:
1024 deltalm = hCoeffs->delta73vh3 *vh3;
1025 rholm = 1. + v2*(hCoeffs->rho73v2 + hCoeffs->rho73v3 * v);
1026 break;
1027 case 2:
1028 deltalm = 0.0;
1029 rholm = 1. + hCoeffs->rho72v2 * v2;
1030 break;
1031 case 1:
1032 deltalm = hCoeffs->delta71vh3 * vh3;
1033 rholm = 1. + v2*(hCoeffs->rho71v2 +hCoeffs->rho71v3 * v);
1034 break;
1035 default:
1037 break;
1038 }
1039 break;
1040 case 8:
1041 switch (m)
1042 {
1043 case 8:
1044 deltalm = 0.0;
1045 rholm = 1. + hCoeffs->rho88v2 * v2;
1046 break;
1047 case 7:
1048 deltalm = 0.0;
1049 rholm = 1. + hCoeffs->rho87v2 * v2;
1050 break;
1051 case 6:
1052 deltalm = 0.0;
1053 rholm = 1. + hCoeffs->rho86v2 * v2;
1054 break;
1055 case 5:
1056 deltalm = 0.0;
1057 rholm = 1. + hCoeffs->rho85v2 * v2;
1058 break;
1059 case 4:
1060 deltalm = 0.0;
1061 rholm = 1. + hCoeffs->rho84v2 * v2;
1062 break;
1063 case 3:
1064 deltalm = 0.0;
1065 rholm = 1. + hCoeffs->rho83v2 * v2;
1066 break;
1067 case 2:
1068 deltalm = 0.0;
1069 rholm = 1. + hCoeffs->rho82v2 * v2;
1070 break;
1071 case 1:
1072 deltalm = 0.0;
1073 rholm = 1. + hCoeffs->rho81v2 * v2;
1074 break;
1075 default:
1077 break;
1078 }
1079 break;
1080 default:
1082 break;
1083 }
1084
1085 /* Raise rholm to the lth power */
1086 rholmPwrl = 1.0;
1087 i = l;
1088 while ( i-- )
1089 {
1090 rholmPwrl *= rholm;
1091 }
1092
1093 *hlm = Tlm * cexp(I * deltalm) * Slm * rholmPwrl;
1094 *hlm *= hNewton;
1095
1096 return XLAL_SUCCESS;
1097}
1098
1099#endif /*_LALSIMIMRFACTORIZEDWAVEFORM_C*/
static UNUSED int XLALCalculateEOBACoefficients(EOBACoefficients *const coeffs, const REAL8 eta)
Function to pre-compute the coefficients in the EOB A potential function.
static REAL8 XLALCalculateEOBdAdr(const REAL8 r, EOBACoefficients *restrict coeffs)
Calculated the derivative of the EOB A function with respect to r, using the pre-computed A coefficie...
static REAL8 XLALEffectiveHamiltonian(const REAL8 eta, const REAL8 r, const REAL8 pr, const REAL8 pp, EOBACoefficients *aCoeffs)
Function to calculate the EOB effective Hamiltonian for the given values of the dynamical variables.
static UNUSED int XLALSimIMREOBCalcFacWaveformCoefficients(FacWaveformCoeffs *const coeffs, const REAL8 eta)
Function which calculates the various coefficients used in the generation of the factorized waveform.
static REAL8 XLALCalculateA6(REAL8 eta)
static REAL8 XLALCalculateA5(REAL8 eta)
Calculates the a5 parameter in the A potential function in EOBNRv2.
static UNUSED int XLALSimIMREOBModifyFacWaveformCoefficients(FacWaveformCoeffs *const coeffs, const REAL8 eta)
Function which adds the additional terms required for waveform generation to the factorized waveform ...
static REAL8 XLALCalculateEOBD(REAL8 r, REAL8 eta) UNUSED
Calculate the EOB D function.
static REAL8 nonKeplerianCoefficient(REAL8Vector *restrict values, const REAL8 eta, EOBACoefficients *coeffs)
Computes the non-Keplerian correction to the velocity as determined from the frequency obtained assum...
static REAL8 XLALCalculateEOBA(const REAL8 r, EOBACoefficients *restrict coeffs)
This function calculates the EOB A function which using the pre-computed coefficients which should al...
static UNUSED int XLALSimIMREOBGetFactorizedWaveform(COMPLEX16 *restrict hlm, REAL8Vector *restrict values, const REAL8 v, const INT4 l, const INT4 m, EOBParams *restrict params)
Computes the factorized waveform according to the prescription given in Pan et al,...
#define ninty4by3etc
Constant which comes up in some of the EOB models.
static UNUSED int XLALSimIMREOBCalculateNewtonianMultipole(COMPLEX16 *multipole, REAL8 x, UNUSED REAL8 r, REAL8 phi, UINT4 l, INT4 m, EOBParams *params)
This function calculates the Newtonian multipole part of the factorized waveform for the EOBNRv2 mode...
static REAL8 UNUSED z3(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)
Module containing the energy and flux functions for waveform generation.
REAL8 Hreal
#define XLAL_CALLGSL(statement)
int l
Definition: bh_qnmode.c:135
double i
Definition: bh_ringdown.c:118
const double pp
const double pp2
const double pr
const double u
const double a4
const double u3
const double r2
const double a2
const double u2
#define LAL_E
#define LAL_PI
#define LAL_GAMMA
double complex COMPLEX16
double REAL8
uint32_t UINT4
int32_t INT4
static const INT4 r
static const INT4 m
static const INT4 a
#define XLAL_ERROR(...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
XLAL_SUCCESS
XLAL_EFUNC
XLAL_EINVAL
XLAL_FAILURE
string status
Structure containing the coefficients for EOBNRv2 A potential function.
Structure containing all the parameters needed for the EOB waveform.
Structure containing the coefficients for calculating the factorized waveform.
Definition: burst.c:245