23 #define UNUSED __attribute__ ((unused))
38 #include <gsl/gsl_errno.h>
39 #include <gsl/gsl_bspline.h>
40 #include <gsl/gsl_blas.h>
41 #include <gsl/gsl_min.h>
42 #include <gsl/gsl_spline.h>
43 #include <lal/Units.h>
44 #include <lal/SeqFactories.h>
45 #include <lal/LALConstants.h>
46 #include <lal/XLALError.h>
47 #include <lal/FrequencySeries.h>
49 #include <lal/StringInput.h>
50 #include <lal/Sequence.h>
51 #include <lal/LALStdio.h>
52 #include <lal/FileIO.h>
54 #include <lal/LALSimInspiral.h>
55 #include <lal/LALSimIMR.h>
58 #include <lal/LALConfig.h>
59 #ifdef LAL_PTHREAD_LOCK
68 static const double gA[] = {0.0001, 0.00011, 0.000121, 0.0001331, 0.00014641, 0.000161051, \
69 0.000177156, 0.000194872, 0.000214359, 0.000235795, 0.000259374, \
70 0.000285312, 0.000313843, 0.000345227, 0.00037975, 0.000417725, \
71 0.000459497, 0.000505447, 0.000555992, 0.000611591, 0.00067275, \
72 0.000740025, 0.000814027, 0.00089543, 0.000984973, 0.00108347, \
73 0.00119182, 0.001311, 0.0014421, 0.00158631, 0.00174494, 0.00191943, \
74 0.00211138, 0.00232252, 0.00255477, 0.00281024, 0.00309127, \
75 0.00340039, 0.00374043, 0.00411448, 0.00452593, 0.00497852, \
76 0.00547637, 0.00602401, 0.00662641, 0.00728905, 0.00801795, \
77 0.00881975, 0.00970172, 0.0106719, 0.0117391, 0.012913, 0.0142043, \
78 0.0156247, 0.0171872, 0.0189059, 0.0207965, 0.0228762, 0.0251638, \
79 0.0276801, 0.0304482, 0.033493, 0.0368423, 0.0405265, 0.0445792, \
80 0.0490371, 0.0539408, 0.0593349, 0.0652683, 0.0717952, 0.0789747, \
81 0.0868722, 0.0955594, 0.105115, 0.115627, 0.12719, 0.139908, 0.14};
83 static const double gPhi[] = {0.0001, 0.000101411, 0.000102849, 0.000104314, 0.000105806, \
84 0.000107328, 0.000108878, 0.000110459, 0.00011207, 0.000113712, \
85 0.000115387, 0.000117095, 0.000118836, 0.000120613, 0.000122424, \
86 0.000124272, 0.000126157, 0.000128081, 0.000130044, 0.000132047, \
87 0.000134091, 0.000136177, 0.000138307, 0.000140481, 0.000142701, \
88 0.000144968, 0.000147283, 0.000149648, 0.000152063, 0.000154531, \
89 0.000157052, 0.000159627, 0.00016226, 0.00016495, 0.0001677, \
90 0.000170512, 0.000173386, 0.000176325, 0.000179331, 0.000182405, \
91 0.00018555, 0.000188768, 0.000192059, 0.000195428, 0.000198876, \
92 0.000202405, 0.000206017, 0.000209716, 0.000213504, 0.000217383, \
93 0.000221357, 0.000225428, 0.000229598, 0.000233872, 0.000238253, \
94 0.000242743, 0.000247346, 0.000252066, 0.000256907, 0.000261871, \
95 0.000266965, 0.00027219, 0.000277553, 0.000283057, 0.000288707, \
96 0.000294507, 0.000300464, 0.000306582, 0.000312866, 0.000319323, \
97 0.000325958, 0.000332778, 0.000339788, 0.000346996, 0.000354409, \
98 0.000362034, 0.000369878, 0.000377949, 0.000386256, 0.000394808, \
99 0.000403612, 0.00041268, 0.00042202, 0.000431643, 0.00044156, \
100 0.000451782, 0.000462321, 0.000473188, 0.000484398, 0.000495963, \
101 0.000507897, 0.000520216, 0.000532935, 0.00054607, 0.000559639, \
102 0.000573659, 0.000588149, 0.00060313, 0.000618621, 0.000634645, \
103 0.000651225, 0.000668385, 0.00068615, 0.000704548, 0.000723607, \
104 0.000743356, 0.000763826, 0.000785052, 0.000807068, 0.000829911, \
105 0.000853621, 0.000878237, 0.000903805, 0.000930369, 0.00095798, \
106 0.000986689, 0.00101655, 0.00104762, 0.00107997, 0.00111365, \
107 0.00114874, 0.00118532, 0.00122345, 0.00126323, 0.00130475, \
108 0.00134809, 0.00139336, 0.00144067, 0.00149013, 0.00154188, \
109 0.00159603, 0.00165273, 0.00171213, 0.0017744, 0.0018397, 0.00190823, \
110 0.00198018, 0.00205578, 0.00213524, 0.00221883, 0.00230681, \
111 0.00239947, 0.00249712, 0.00260011, 0.0027088, 0.00282359, \
112 0.00294492, 0.00307324, 0.00320907, 0.00335297, 0.00350553, \
113 0.00366741, 0.00383935, 0.00402211, 0.00421656, 0.00442364, \
114 0.0046444, 0.00487997, 0.0051316, 0.00540067, 0.00568873, 0.00599744, \
115 0.0063287, 0.00668457, 0.00706737, 0.00747967, 0.00792436, \
116 0.00840463, 0.00892411, 0.00948683, 0.0100974, 0.0107608, 0.011483, \
117 0.0122706, 0.013131, 0.0140727, 0.0151056, 0.0162408, 0.0174911, \
118 0.0188713, 0.0203987, 0.0220931, 0.0239776, 0.0260796, 0.0284307, \
119 0.0310686, 0.0340377, 0.0373912, 0.0411922, 0.0455169, 0.0504574, \
120 0.0561255, 0.0626582, 0.0702238, 0.0790313, 0.0893416, 0.101483, \
121 0.115873, 0.133047, 0.14};
123 #ifdef LAL_PTHREAD_LOCK
124 static pthread_once_t SEOBNRv1ROMEffectiveSpin_is_initialized = PTHREAD_ONCE_INIT;
129 typedef struct tagSEOBNRROMdata_coeff
148 typedef struct tagSplineData
150 gsl_bspline_workspace *bwx;
151 gsl_bspline_workspace *bwy;
190 static size_t NextPow2(
const size_t n);
194 static int read_vector(
const char dir[],
const char fname[], gsl_vector *v);
195 static int read_matrix(
const char dir[],
const char fname[], gsl_matrix *
m);
197 static int load_data(
const char dir[], gsl_vector *cvec_amp, gsl_vector *cvec_phi, gsl_matrix *Bamp, gsl_matrix *Bphi, gsl_vector *cvec_amp_pre);
202 gsl_vector *cvec_amp,
203 gsl_vector *cvec_phi,
204 gsl_vector *cvec_amp_pre,
240 static int load_data(
const char dir[], gsl_vector *cvec_amp, gsl_vector *cvec_phi, gsl_matrix *Bamp, gsl_matrix *Bphi, gsl_vector *cvec_amp_pre) {
243 ret |=
read_vector(dir,
"SEOBNRv1ROM_SS_Amp_ciall.dat", cvec_amp);
244 ret |=
read_vector(dir,
"SEOBNRv1ROM_SS_Phase_ciall.dat", cvec_phi);
245 ret |=
read_matrix(dir,
"SEOBNRv1ROM_SS_Bamp_bin.dat", Bamp);
246 ret |=
read_matrix(dir,
"SEOBNRv1ROM_SS_Bphase_bin.dat", Bphi);
247 ret |=
read_vector(dir,
"SEOBNRv1ROM_SS_AmpPrefac_ci.dat", cvec_amp_pre);
253 if(!splinedata) exit(1);
260 (*splinedata)->ncx = ncx;
261 (*splinedata)->ncy = ncy;
264 double qvec[] = {1., 1.125, 1.25, 1.375, 1.5, 1.75, 2., 2.25, 2.5, 2.75, 3., 3.25, \
265 3.5, 3.75, 4., 4.25, 4.5, 4.75, 5., 5.5, 6., 6.5, 7., 7.5, 8., 8.5, \
266 9., 9.5, 10., 10.5, 11., 11.5, 12., 12.5, 13., 13.5, 14., 14.5, 15., \
267 15.5, 16., 16.5, 17., 17.5, 18., 18.5, 19., 19.5, 20., 20.5, 21., \
268 21.5, 22., 22.5, 23., 23.5, 24., 24.5, 25., 25.5, 26., 26.5, 27., \
269 27.5, 28., 28.5, 29., 29.5, 30., 30.5, 31., 31.5, 32., 32.5, 33., \
270 33.5, 34., 34.5, 35., 35.5, 36., 36.5, 37., 37.5, 38., 38.5, 39., \
271 39.5, 40., 41., 42., 43., 44., 45., 46., 47., 48., 49., 50., 51., \
272 52., 53., 54., 55., 56., 57., 58., 59., 60., 61., 62., 63., 64., 65., \
273 66., 67., 68., 69., 70., 71., 72., 73., 74., 75., 76., 77., 78., 79., \
274 80., 81., 82., 83., 84., 85., 86., 87., 88., 89., 90., 91., 92., 93., \
275 94., 95., 95.5, 96., 96.5, 97., 97.5, 98., 98.5, 98.75, 99., 99.25, \
278 double chivec[] = {-1., -0.975, -0.95, -0.925, -0.9, -0.875, -0.85, -0.825, -0.8, \
279 -0.775, -0.75, -0.725, -0.7, -0.675, -0.65, -0.625, -0.6, -0.55, \
280 -0.5, -0.45, -0.4, -0.35, -0.3, -0.25, -0.2, -0.15, -0.1, -0.05, 0., \
281 0.05, 0.1, 0.15, 0.2, 0.25, 0.3, 0.325, 0.35, 0.375, 0.4, 0.425, \
282 0.45, 0.475, 0.5, 0.525, 0.55, 0.575, 0.6};
284 const size_t nbreak_x = ncx-2;
285 const size_t nbreak_y = ncy-2;
288 gsl_bspline_workspace *bwx = gsl_bspline_alloc(4, nbreak_x);
289 gsl_bspline_workspace *bwy = gsl_bspline_alloc(4, nbreak_y);
292 gsl_vector *breakpts_x = gsl_vector_alloc(nbreak_x);
293 gsl_vector *breakpts_y = gsl_vector_alloc(nbreak_y);
295 gsl_vector_set(breakpts_x,
i, qvec[
i]);
296 for (
UINT4 j=0; j<nbreak_y; j++)
297 gsl_vector_set(breakpts_y, j, chivec[j]);
299 gsl_bspline_knots(breakpts_x, bwx);
300 gsl_bspline_knots(breakpts_y, bwy);
302 gsl_vector_free(breakpts_x);
303 gsl_vector_free(breakpts_y);
305 (*splinedata)->bwx=bwx;
306 (*splinedata)->bwy=bwy;
311 if(!splinedata)
return;
312 if(splinedata->
bwx) gsl_bspline_free(splinedata->
bwx);
313 if(splinedata->
bwy) gsl_bspline_free(splinedata->
bwy);
322 gsl_vector *cvec_amp,
323 gsl_vector *cvec_phi,
324 gsl_vector *cvec_amp_pre,
332 gsl_bspline_workspace *bwx=splinedata->
bwx;
333 gsl_bspline_workspace *bwy=splinedata->
bwy;
335 int ncx = splinedata->
ncx;
336 int ncy = splinedata->
ncy;
340 for (
int k=0; k<
nk_amp; k++) {
341 gsl_vector v = gsl_vector_subvector(cvec_amp, k*
N,
N).vector;
343 gsl_vector_set(c_amp, k, csum);
347 for (
int k=0; k<
nk_phi; k++) {
348 gsl_vector v = gsl_vector_subvector(cvec_phi, k*
N,
N).vector;
350 gsl_vector_set(c_phi, k, csum);
373 XLALPrintError(
"WARNING: You tried to setup the SEOBNRv1ROMEffectiveSpin model that was already initialised. Ignoring\n");
377 (romdata)->cvec_amp = gsl_vector_alloc(
N*
nk_amp);
378 (romdata)->cvec_phi = gsl_vector_alloc(
N*
nk_phi);
381 (romdata)->cvec_amp_pre = gsl_vector_alloc(
N);
382 ret=
load_data(dir, (romdata)->cvec_amp, (romdata)->cvec_phi, (romdata)->Bamp, (romdata)->Bphi, (romdata)->cvec_amp_pre);
393 if(romdata->cvec_amp) gsl_vector_free(romdata->cvec_amp);
394 if(romdata->cvec_phi) gsl_vector_free(romdata->cvec_phi);
395 if(romdata->Bamp) gsl_matrix_free(romdata->Bamp);
396 if(romdata->Bphi) gsl_matrix_free(romdata->Bphi);
397 if(romdata->cvec_amp_pre) gsl_vector_free(romdata->cvec_amp_pre);
404 if(!romdatacoeff) exit(1);
411 (*romdatacoeff)->c_amp = gsl_vector_alloc(
nk_amp);
412 (*romdatacoeff)->c_phi = gsl_vector_alloc(
nk_phi);
417 if(romdatacoeff->
c_amp) gsl_vector_free(romdatacoeff->
c_amp);
418 if(romdatacoeff->
c_phi) gsl_vector_free(romdatacoeff->
c_phi);
424 return 1 << (size_t) ceil(log2(n));
451 if(!hptilde || !hctilde)
454 if(*hptilde || *hctilde) {
455 XLALPrintError(
"(*hptilde) and (*hctilde) are supposed to be NULL, but got %p and %p",(*hptilde),(*hctilde));
462 if (
q > 100.0)
nudge(&
q, 100.0, 1
e-6);
463 if (chi < -1.0)
nudge(&chi, -1.0, 1
e-6);
464 if (chi > 0.6)
nudge(&chi, 0.6, 1
e-6);
467 if ( chi < -1.0 || chi > 0.6 ) {
468 XLALPrintError(
"XLAL Error - %s: chi smaller than -1 or larger than 0.6!\nSEOBNRv1ROMEffectiveSpin is only available for spins in the range -1 <= a/M <= 0.6.\n", __func__);
473 XLALPrintError(
"XLAL Error - %s: q=%lf larger than 100!\nSEOBNRv1ROMEffectiveSpin is only available for q in the range 1 <= q <= 100.\n", __func__,
q);
477 if (
q >= 20 &&
q <= 40 && chi < -0.75 && chi > -0.9) {
478 XLALPrintWarning(
"XLAL Warning - %s: q in [20,40] and chi in [-0.8]. The SEOBNRv1 model is not trustworthy in this region!\nSee Fig 15 in CQG 31 195010, 2014 for details.", __func__);
484 double fLow = freqs_in->
data[0];
485 double fHigh = freqs_in->
data[freqs_in->
length - 1];
493 double fLow_geom = fLow * Mtot_sec;
494 double fHigh_geom = fHigh * Mtot_sec;
495 double fRef_geom = fRef * Mtot_sec;
496 double deltaF_geom = deltaF * Mtot_sec;
504 XLALPrintWarning(
"Maximal frequency Mf_high=%g is greater than highest ROM frequency Mf_ROM_Max=%g. Using Mf_high=Mf_ROM_Max.", fHigh_geom,
Mf_ROM_max);
508 XLAL_ERROR(
XLAL_EDOM,
"End frequency %g is smaller than starting frequency %g!\n", fHigh_geom, fLow_geom);
510 XLALPrintWarning(
"Reference frequency Mf_ref=%g is greater than maximal frequency in ROM Mf=%g. Starting at maximal frequency in ROM.\n", fRef_geom,
Mf_ROM_max);
514 XLALPrintWarning(
"Reference frequency Mf_ref=%g is smaller than lowest frequency in ROM Mf=%g. Starting at lowest frequency in ROM.\n", fLow_geom,
Mf_ROM_min);
529 romdata->cvec_amp_pre,
530 romdata_coeff->
c_amp,
531 romdata_coeff->
c_phi,
537 XLAL_ERROR(retcode,
"Parameter-space interpolation failed.");
543 gsl_vector* amp_f = gsl_vector_alloc(
nk_amp);
544 gsl_vector* phi_f = gsl_vector_alloc(
nk_phi);
545 gsl_blas_dgemv(CblasTrans, 1.0, romdata->Bamp, romdata_coeff->
c_amp, 0.0, amp_f);
546 gsl_blas_dgemv(CblasTrans, 1.0, romdata->Bphi, romdata_coeff->
c_phi, 0.0, phi_f);
549 gsl_interp_accel *acc_amp = gsl_interp_accel_alloc();
550 gsl_spline *spline_amp = gsl_spline_alloc(gsl_interp_cspline,
nk_amp);
551 gsl_spline_init(spline_amp,
gA, gsl_vector_const_ptr(amp_f,0),
nk_amp);
553 gsl_interp_accel *acc_phi = gsl_interp_accel_alloc();
554 gsl_spline *spline_phi = gsl_spline_alloc(gsl_interp_cspline,
nk_phi);
555 gsl_spline_init(spline_phi,
gPhi, gsl_vector_const_ptr(phi_f,0),
nk_phi);
564 npts =
NextPow2(fHigh_geom / deltaF_geom) + 1;
565 if (fHigh_geom < fHigh * Mtot_sec)
566 npts =
NextPow2(fHigh * Mtot_sec / deltaF_geom) + 1;
574 double fHigh_temp = fHigh_geom / Mtot_sec;
576 UINT4 iStop = (
UINT4) ceil(fHigh_temp / deltaF);
582 freqs->
data[
i-iStart] =
i*deltaF_geom;
596 freqs->
data[
i] = freqs_in->
data[
i] * Mtot_sec;
600 if (!(*hptilde) || !(*hctilde))
603 gsl_spline_free(spline_amp);
604 gsl_spline_free(spline_phi);
605 gsl_interp_accel_free(acc_amp);
606 gsl_interp_accel_free(acc_phi);
607 gsl_vector_free(amp_f);
608 gsl_vector_free(phi_f);
612 memset((*hptilde)->data->data, 0, npts *
sizeof(
COMPLEX16));
613 memset((*hctilde)->data->data, 0, npts *
sizeof(
COMPLEX16));
621 REAL8 cosi = cos(inclination);
622 REAL8 pcoef = 0.5*(1.0 + cosi*cosi);
627 double amp0 = Mtot * amp_pre * Mtot_sec *
LAL_MRSUN_SI / (distance);
630 double phase_change = gsl_spline_eval(spline_phi, fRef_geom, acc_phi) - 2*phiRef;
634 double f = freqs->
data[
i];
637 double A = gsl_spline_eval(spline_amp, f, acc_amp);
638 double phase = gsl_spline_eval(spline_phi, f, acc_phi) - phase_change;
640 pdata[j] = pcoef * htilde;
641 cdata[j] = -I * ccoef * htilde;
651 if (Mf_final > freqs->
data[L-1])
652 Mf_final = freqs->
data[L-1];
653 if (Mf_final < freqs->
data[0])
656 gsl_spline_free(spline_amp);
657 gsl_spline_free(spline_phi);
658 gsl_interp_accel_free(acc_amp);
659 gsl_interp_accel_free(acc_phi);
660 gsl_vector_free(amp_f);
661 gsl_vector_free(phi_f);
668 REAL8 t_corr = gsl_spline_eval_deriv(spline_phi, Mf_final, acc_phi) / (2*
LAL_PI);
672 double f = freqs->
data[
i] - fRef_geom;
674 pdata[j] *= cexp(-2*
LAL_PI * I * f * t_corr);
675 cdata[j] *= cexp(-2*
LAL_PI * I * f * t_corr);
680 gsl_spline_free(spline_amp);
681 gsl_spline_free(spline_phi);
682 gsl_interp_accel_free(acc_amp);
683 gsl_interp_accel_free(acc_phi);
684 gsl_vector_free(amp_f);
685 gsl_vector_free(phi_f);
751 struct tagCOMPLEX16FrequencySeries **hptilde,
752 struct tagCOMPLEX16FrequencySeries **hctilde,
766 double Mtot = mass1+mass2;
767 double q = mass2 / mass1;
773 #ifdef LAL_PTHREAD_LOCK
784 phiRef, fRef, distance, inclination, Mtot_sec,
q, chi, freqs, 0);
798 struct tagCOMPLEX16FrequencySeries **hptilde,
799 struct tagCOMPLEX16FrequencySeries **hctilde,
815 double Mtot = mass1+mass2;
816 double q = mass2 / mass1;
825 #ifdef LAL_PTHREAD_LOCK
837 freqs->
data[0] = fLow;
838 freqs->
data[1] = fHigh;
841 phiRef, fRef, distance, inclination, Mtot_sec,
q, chi, freqs, deltaF);
859 char datafile[] =
"SEOBNRv1ROM_SS_Phase_ciall.dat";
864 char *dir = dirname(
path);
static const REAL8 Mf_ROM_max
static const REAL8 Mf_ROM_min
static void nudge(REAL8 *x, REAL8 X, REAL8 epsilon)
Auxiliary functions for SEOBNRv1/v2 reduced order modeling codes in LALSimIMRSEOBNRv2ROMDoubleSpin....
static UNUSED double SEOBNRROM_Ringdown_Mf_From_Mtot_q(const double Mtot_sec, const double q, const double chi1, const double chi2, Approximant apx)
static UNUSED REAL8 Interpolate_Coefficent_Matrix(gsl_vector *v, REAL8 eta, REAL8 chi, int ncx, int ncy, gsl_bspline_workspace *bwx, gsl_bspline_workspace *bwy)
static void SEOBNRROMdata_coeff_Cleanup(SEOBNRROMdata_coeff *romdatacoeff)
static void SplineData_Destroy(SplineData *splinedata)
static const double gPhi[]
static int SEOBNRv1ROMEffectiveSpin_Init(const char dir[])
Setup SEOBNRv1ROMEffectiveSpin model using data files installed in dir.
static void SEOBNRROMdata_coeff_Init(SEOBNRROMdata_coeff **romdatacoeff)
static void SEOBNRv1ROMEffectiveSpin_Init_LALDATA(void)
Setup SEOBNRv1ROMEffectiveSpin model using data files installed in $LAL_DATA_PATH.
static int read_vector(const char dir[], const char fname[], gsl_vector *v)
static SEOBNRROMdata __lalsim_SEOBNRv1ROMSS_data
static void SEOBNRROMdata_Cleanup(SEOBNRROMdata *romdata)
static void SplineData_Init(SplineData **splinedata)
static int TP_Spline_interpolation_2d(REAL8 q, REAL8 chi, gsl_vector *cvec_amp, gsl_vector *cvec_phi, gsl_vector *cvec_amp_pre, gsl_vector *c_amp, gsl_vector *c_phi, REAL8 *amp_pre)
static int SEOBNRv1ROMEffectiveSpinCore(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, double phiRef, double fRef, double distance, double inclination, double Mtot_sec, double q, double chi, const REAL8Sequence *freqs, double deltaF)
Core function for computing the ROM waveform.
static size_t NextPow2(const size_t n)
static bool SEOBNRv1ROMEffectiveSpin_IsSetup(void)
Helper function to check if the SEOBNRv1ROMEffectiveSpin model has been initialised.
static int load_data(const char dir[], gsl_vector *cvec_amp, gsl_vector *cvec_phi, gsl_matrix *Bamp, gsl_matrix *Bphi, gsl_vector *cvec_amp_pre)
static int SEOBNRROMdata_Init(SEOBNRROMdata *romdata, const char dir[])
static int read_matrix(const char dir[], const char fname[], gsl_matrix *m)
#define XLAL_FILE_RESOLVE_PATH(fname)
COMPLEX16FrequencySeries * XLALCreateCOMPLEX16FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
void * XLALCalloc(size_t m, size_t n)
int XLALSimIMRSEOBNRv1ROMEffectiveSpin(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, REAL8 phiRef, REAL8 deltaF, REAL8 fLow, REAL8 fHigh, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi)
Compute waveform in LAL format for the SEOBNRv1_ROM_EffectiveSpin model.
int XLALSimIMRSEOBNRv1ROMEffectiveSpinFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi)
Compute waveform in LAL format at specified frequencies for the SEOBNRv1_ROM_EffectiveSpin model.
@ SEOBNRv1
Spin-aligned EOBNR model.
void XLALDestroyREAL8Sequence(REAL8Sequence *sequence)
REAL8Sequence * XLALCreateREAL8Sequence(size_t length)
const LALUnit lalStrainUnit
const LALUnit lalSecondUnit
LALUnit * XLALUnitMultiply(LALUnit *output, const LALUnit *unit1, const LALUnit *unit2)
#define XLAL_ERROR_VOID(...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
int int XLALPrintWarning(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
int pthread_once(pthread_once_t *once_control, void(*init_routine)(void))
gsl_bspline_workspace * bwx
gsl_bspline_workspace * bwy
gsl_vector * cvec_amp_pre