LALSimulation  5.4.0.1-fe68b98
LALSimIMRSEOBNRv1ROMEffectiveSpin.c
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2014 Michael Puerrer, John Veitch
3  * Reduced Order Model for SEOBNR
4  *
5  * This program is free software; you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation; either version 2 of the License, or
8  * (at your option) any later version.
9  *
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with with program; see the file COPYING. If not, write to the
17  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
18  * MA 02110-1301 USA
19  */
20 
21 
22 #ifdef __GNUC__
23 #define UNUSED __attribute__ ((unused))
24 #else
25 #define UNUSED
26 #endif
27 
28 #include <stdio.h>
29 #include <stdlib.h>
30 #include <math.h>
31 #include <complex.h>
32 #include <time.h>
33 #include <stdbool.h>
34 #include <alloca.h>
35 #include <string.h>
36 #include <libgen.h>
37 
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>
48 #include <lal/Date.h>
49 #include <lal/StringInput.h>
50 #include <lal/Sequence.h>
51 #include <lal/LALStdio.h>
52 #include <lal/FileIO.h>
53 
54 #include <lal/LALSimInspiral.h>
55 #include <lal/LALSimIMR.h>
57 
58 #include <lal/LALConfig.h>
59 #ifdef LAL_PTHREAD_LOCK
60 #include <pthread.h>
61 #endif
62 
63 
64 /********* Input data for spline basis points **************/
65 #define nk_amp 78 // number of SVD-modes == number of basis functions for amplitude
66 #define nk_phi 200 // number of SVD-modes == number of basis functions for phase
67 
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};
82 
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};
122 
123 #ifdef LAL_PTHREAD_LOCK
124 static pthread_once_t SEOBNRv1ROMEffectiveSpin_is_initialized = PTHREAD_ONCE_INIT;
125 #endif
126 
127 /*************** type definitions ******************/
128 
129 typedef struct tagSEOBNRROMdata_coeff
130 {
131  gsl_vector* c_amp;
132  gsl_vector* c_phi;
134 
136 {
138  gsl_vector* cvec_amp;
139  gsl_vector* cvec_phi;
140  gsl_matrix *Bamp;
141  gsl_matrix *Bphi;
142  gsl_vector* cvec_amp_pre;
143 };
144 typedef struct tagSEOBNRROMdata SEOBNRROMdata;
145 
146 static SEOBNRROMdata __lalsim_SEOBNRv1ROMSS_data;
147 
148 typedef struct tagSplineData
149 {
150  gsl_bspline_workspace *bwx;
151  gsl_bspline_workspace *bwy;
152  int ncx, ncy;
153 } SplineData;
154 
155 /**************** Internal functions **********************/
156 
157 static void SEOBNRv1ROMEffectiveSpin_Init_LALDATA(void);
158 static int SEOBNRv1ROMEffectiveSpin_Init(const char dir[]);
159 static bool SEOBNRv1ROMEffectiveSpin_IsSetup(void);
160 
161 static int SEOBNRROMdata_Init(SEOBNRROMdata *romdata, const char dir[]);
162 static void SEOBNRROMdata_Cleanup(SEOBNRROMdata *romdata);
163 
164 /**
165  * Core function for computing the ROM waveform.
166  * Interpolate projection coefficient data and evaluate coefficients at desired (q, chi).
167  * Construct 1D splines for amplitude and phase.
168  * Compute strain waveform from amplitude and phase.
169 */
171  COMPLEX16FrequencySeries **hptilde,
172  COMPLEX16FrequencySeries **hctilde,
173  double phiRef,
174  double fRef,
175  double distance,
176  double inclination,
177  double Mtot_sec,
178  double q,
179  double chi,
180  const REAL8Sequence *freqs, /* Frequency points at which to evaluate the waveform (Hz) */
181  double deltaF
182  /* If deltaF >0, the frequency points given in freqs are uniformly spaced with
183  * spacing deltaF. Otherwise, the frequency points are spaced non-uniformly.
184  * Then we will use deltaF = 0 to create the frequency series we return. */
185 );
186 
187 static void SEOBNRROMdata_coeff_Init(SEOBNRROMdata_coeff **romdatacoeff);
188 static void SEOBNRROMdata_coeff_Cleanup(SEOBNRROMdata_coeff *romdatacoeff);
189 
190 static size_t NextPow2(const size_t n);
191 static void SplineData_Destroy(SplineData *splinedata);
192 static void SplineData_Init(SplineData **splinedata);
193 
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);
196 
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);
198 
199 static int TP_Spline_interpolation_2d(
200  REAL8 q, // Input: q-value for which projection coefficients should be evaluated
201  REAL8 chi, // Input: chi-value for which projection coefficients should be evaluated
202  gsl_vector *cvec_amp, // Input: data for spline coefficients for amplitude
203  gsl_vector *cvec_phi, // Input: data for spline coefficients for phase
204  gsl_vector *cvec_amp_pre, // Input: data for spline coefficients for amplitude prefactor
205  gsl_vector *c_amp, // Output: interpolated projection coefficients for amplitude
206  gsl_vector *c_phi, // Output: interpolated projection coefficients for phase
207  REAL8 *amp_pre // Output: interpolated amplitude prefactor
208 );
209 
210 
211 /********************* Definitions begin here ********************/
212 
213 /** Setup SEOBNRv1ROMEffectiveSpin model using data files installed in dir
214  */
215 static int SEOBNRv1ROMEffectiveSpin_Init(const char dir[]) {
216  if(__lalsim_SEOBNRv1ROMSS_data.setup) {
217  XLALPrintError("Error: SEOBNRROMdata was already set up!");
219  }
220 
222 
223  if(__lalsim_SEOBNRv1ROMSS_data.setup) {
224  return(XLAL_SUCCESS);
225  }
226  else {
227  return(XLAL_EFAILED);
228  }
229 }
230 
231 /** Helper function to check if the SEOBNRv1ROMEffectiveSpin model has been initialised */
234  return true;
235  else
236  return false;
237 }
238 
239 // Read binary ROM data for basis functions and coefficients
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) {
241  // Load binary data for amplitude and phase spline coefficients as computed in Mathematica
242  int ret = XLAL_SUCCESS;
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);
248  return(ret);
249 }
250 
251 static void SplineData_Init( SplineData **splinedata )
252 {
253  if(!splinedata) exit(1);
254  if(*splinedata) SplineData_Destroy(*splinedata);
255 
256  (*splinedata)=XLALCalloc(1,sizeof(SplineData));
257 
258  const int ncx = 159; // points in q
259  const int ncy = 49; // points in chi
260  (*splinedata)->ncx = ncx;
261  (*splinedata)->ncy = ncy;
262 
263  // Set up B-spline basis for desired knots
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, \
276  99.5, 99.75, 100.};
277 
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};
283 
284  const size_t nbreak_x = ncx-2; // must have nbreak = n -2 for cubic splines
285  const size_t nbreak_y = ncy-2; // must have nbreak = n -2 for cubic splines
286 
287  // allocate a cubic bspline workspace (k = 4)
288  gsl_bspline_workspace *bwx = gsl_bspline_alloc(4, nbreak_x);
289  gsl_bspline_workspace *bwy = gsl_bspline_alloc(4, nbreak_y);
290 
291  // set breakpoints (and thus knots by hand)
292  gsl_vector *breakpts_x = gsl_vector_alloc(nbreak_x);
293  gsl_vector *breakpts_y = gsl_vector_alloc(nbreak_y);
294  for (UINT4 i=0; i<nbreak_x; i++)
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]);
298 
299  gsl_bspline_knots(breakpts_x, bwx);
300  gsl_bspline_knots(breakpts_y, bwy);
301 
302  gsl_vector_free(breakpts_x);
303  gsl_vector_free(breakpts_y);
304 
305  (*splinedata)->bwx=bwx;
306  (*splinedata)->bwy=bwy;
307 }
308 
309 static void SplineData_Destroy(SplineData *splinedata)
310 {
311  if(!splinedata) return;
312  if(splinedata->bwx) gsl_bspline_free(splinedata->bwx);
313  if(splinedata->bwy) gsl_bspline_free(splinedata->bwy);
314  XLALFree(splinedata);
315 }
316 
317 // Interpolate projection coefficients for amplitude and phase over the parameter space (q, chi).
318 // The multi-dimensional interpolation is carried out via a tensor product decomposition.
320  REAL8 q, // Input: q-value for which projection coefficients should be evaluated
321  REAL8 chi, // Input: chi-value for which projection coefficients should be evaluated
322  gsl_vector *cvec_amp, // Input: data for spline coefficients for amplitude
323  gsl_vector *cvec_phi, // Input: data for spline coefficients for phase
324  gsl_vector *cvec_amp_pre, // Input: data for spline coefficients for amplitude prefactor
325  gsl_vector *c_amp, // Output: interpolated projection coefficients for amplitude
326  gsl_vector *c_phi, // Output: interpolated projection coefficients for phase
327  REAL8 *amp_pre // Output: interpolated amplitude prefactor
328 ) {
329 
330  SplineData *splinedata=NULL;
331  SplineData_Init(&splinedata);
332  gsl_bspline_workspace *bwx=splinedata->bwx;
333  gsl_bspline_workspace *bwy=splinedata->bwy;
334 
335  int ncx = splinedata->ncx; // points in q
336  int ncy = splinedata->ncy; // points in chi
337  int N = ncx*ncy; // size of the data matrix for one SVD-mode
338 
339  // Evaluate the TP spline for all SVD modes - amplitude
340  for (int k=0; k<nk_amp; k++) { // For each SVD mode
341  gsl_vector v = gsl_vector_subvector(cvec_amp, k*N, N).vector; // Pick out the coefficient matrix corresponding to the k-th SVD mode.
342  REAL8 csum = Interpolate_Coefficent_Matrix(&v, q, chi, ncx, ncy, bwx, bwy);
343  gsl_vector_set(c_amp, k, csum);
344  }
345 
346  // Evaluate the TP spline for all SVD modes - phase
347  for (int k=0; k<nk_phi; k++) { // For each SVD mode
348  gsl_vector v = gsl_vector_subvector(cvec_phi, k*N, N).vector; // Pick out the coefficient matrix corresponding to the k-th SVD mode.
349  REAL8 csum = Interpolate_Coefficent_Matrix(&v, q, chi, ncx, ncy, bwx, bwy);
350  gsl_vector_set(c_phi, k, csum);
351  }
352 
353  // Evaluate the TP spline for the amplitude prefactor
354  *amp_pre = Interpolate_Coefficent_Matrix(cvec_amp_pre, q, chi, ncx, ncy, bwx, bwy);
355 
356  SplineData_Destroy(splinedata);
357 
358  return(0);
359 }
360 
361 /* Set up a new ROM model, using data contained in dir */
362 static int SEOBNRROMdata_Init(SEOBNRROMdata *romdata, const char dir[]) {
363  // set up ROM
364  int ncx = 159; // points in q
365  int ncy = 49; // points in chi
366  int N = ncx*ncy; // size of the data matrix for one SVD-mode
367 
368  int ret = XLAL_FAILURE;
369 
370  /* Create storage for structures */
371  if(romdata->setup)
372  {
373  XLALPrintError("WARNING: You tried to setup the SEOBNRv1ROMEffectiveSpin model that was already initialised. Ignoring\n");
374  return (XLAL_FAILURE);
375  }
376 
377  (romdata)->cvec_amp = gsl_vector_alloc(N*nk_amp);
378  (romdata)->cvec_phi = gsl_vector_alloc(N*nk_phi);
379  (romdata)->Bamp = gsl_matrix_alloc(nk_amp, nk_amp);
380  (romdata)->Bphi = gsl_matrix_alloc(nk_phi, 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);
383 
384  if(XLAL_SUCCESS==ret) romdata->setup=1;
385  else SEOBNRROMdata_Cleanup(romdata);
386 
387  return (ret);
388 }
389 
390 
391 /* Deallocate contents of the given SEOBNRROMdata structure */
392 static void SEOBNRROMdata_Cleanup(SEOBNRROMdata *romdata) {
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);
398  romdata->setup=0;
399 }
400 
401 /* Structure for internal use */
402 static void SEOBNRROMdata_coeff_Init(SEOBNRROMdata_coeff **romdatacoeff) {
403 
404  if(!romdatacoeff) exit(1);
405  /* Create storage for structures */
406  if(!*romdatacoeff)
407  *romdatacoeff=XLALCalloc(1,sizeof(SEOBNRROMdata_coeff));
408  else
409  SEOBNRROMdata_coeff_Cleanup(*romdatacoeff);
410 
411  (*romdatacoeff)->c_amp = gsl_vector_alloc(nk_amp);
412  (*romdatacoeff)->c_phi = gsl_vector_alloc(nk_phi);
413 }
414 
415 /* Deallocate contents of the given SEOBNRROMdata_coeff structure */
417  if(romdatacoeff->c_amp) gsl_vector_free(romdatacoeff->c_amp);
418  if(romdatacoeff->c_phi) gsl_vector_free(romdatacoeff->c_phi);
419  XLALFree(romdatacoeff);
420 }
421 
422 /* Return the closest higher power of 2 */
423 static size_t NextPow2(const size_t n) {
424  return 1 << (size_t) ceil(log2(n));
425 }
426 
427 /**
428  * Core function for computing the ROM waveform.
429  * Interpolate projection coefficient data and evaluate coefficients at desired (q, chi).
430  * Construct 1D splines for amplitude and phase.
431  * Compute strain waveform from amplitude and phase.
432 */
434  COMPLEX16FrequencySeries **hptilde,
435  COMPLEX16FrequencySeries **hctilde,
436  double phiRef,
437  double fRef,
438  double distance,
439  double inclination,
440  double Mtot_sec,
441  double q,
442  double chi,
443  const REAL8Sequence *freqs_in, /* Frequency points at which to evaluate the waveform (Hz) */
444  double deltaF
445  /* If deltaF > 0, the frequency points given in freqs are uniformly spaced with
446  * spacing deltaF. Otherwise, the frequency points are spaced non-uniformly.
447  * Then we will use deltaF = 0 to create the frequency series we return. */
448  )
449 {
450  /* Check output arrays */
451  if(!hptilde || !hctilde)
453  SEOBNRROMdata *romdata=&__lalsim_SEOBNRv1ROMSS_data;
454  if(*hptilde || *hctilde) {
455  XLALPrintError("(*hptilde) and (*hctilde) are supposed to be NULL, but got %p and %p",(*hptilde),(*hctilde));
457  }
458  int retcode=0;
459 
460  // 'Nudge' parameter values to allowed boundary values if close by
461  if (q < 1.0) nudge(&q, 1.0, 1e-6);
462  if (q > 100.0) nudge(&q, 100.0, 1e-6);
463  if (chi < -1.0) nudge(&chi, -1.0, 1e-6);
464  if (chi > 0.6) nudge(&chi, 0.6, 1e-6);
465 
466  /* If either spin > 0.6, model not available, exit */
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__);
470  }
471 
472  if (q > 100) {
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);
475  }
476 
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__);
480  }
481 
482  /* Find frequency bounds */
483  if (!freqs_in) XLAL_ERROR(XLAL_EFAULT);
484  double fLow = freqs_in->data[0];
485  double fHigh = freqs_in->data[freqs_in->length - 1];
486 
487  if(fRef==0.0)
488  fRef=fLow;
489 
490  /* Convert to geometric units for frequency */
491  double Mf_ROM_min = fmax(gA[0], gPhi[0]); // lowest allowed geometric frequency for ROM
492  double Mf_ROM_max = fmin(gA[nk_amp-1], gPhi[nk_phi-1]); // highest allowed geometric frequency for ROM
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;
497 
498  // Enforce allowed geometric frequency range
499  if (fLow_geom < Mf_ROM_min)
500  XLAL_ERROR(XLAL_EDOM, "Starting frequency Mflow=%g is smaller than lowest frequency in ROM Mf=%g.\n", fLow_geom, Mf_ROM_min);
501  if (fHigh_geom == 0)
502  fHigh_geom = Mf_ROM_max;
503  else if (fHigh_geom > Mf_ROM_max) {
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);
505  fHigh_geom = Mf_ROM_max;
506  }
507  else if (fHigh_geom < Mf_ROM_min)
508  XLAL_ERROR(XLAL_EDOM, "End frequency %g is smaller than starting frequency %g!\n", fHigh_geom, fLow_geom);
509  if (fRef_geom > Mf_ROM_max) {
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);
511  fRef_geom = Mf_ROM_max; // If fref > fhigh we reset fref to default value of cutoff frequency.
512  }
513  if (fRef_geom < Mf_ROM_min) {
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);
515  fRef_geom = Mf_ROM_min;
516  }
517 
518  /* Internal storage for w.f. coefficiencts */
519  SEOBNRROMdata_coeff *romdata_coeff=NULL;
520  SEOBNRROMdata_coeff_Init(&romdata_coeff);
521  REAL8 amp_pre;
522 
523  /* Interpolate projection coefficients and evaluate them at (q,chi) */
525  q, // Input: q-value for which projection coefficients should be evaluated
526  chi, // Input: chi-value for which projection coefficients should be evaluated
527  romdata->cvec_amp, // Input: data for spline coefficients for amplitude
528  romdata->cvec_phi, // Input: data for spline coefficients for phase
529  romdata->cvec_amp_pre, // Input: data for spline coefficients for amplitude prefactor
530  romdata_coeff->c_amp, // Output: interpolated projection coefficients for amplitude
531  romdata_coeff->c_phi, // Output: interpolated projection coefficients for phase
532  &amp_pre // Output: interpolated amplitude prefactor
533  );
534 
535  if(retcode!=0) {
536  SEOBNRROMdata_coeff_Cleanup(romdata_coeff);
537  XLAL_ERROR(retcode, "Parameter-space interpolation failed.");
538  }
539 
540  // Compute function values of amplitude an phase on sparse frequency points by evaluating matrix vector products
541  // amp_pts = B_A^T . c_A
542  // phi_pts = B_phi^T . c_phi
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);
547 
548  // Setup 1d splines in frequency
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);
552 
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);
556 
557 
558  size_t npts = 0;
559  LIGOTimeGPS tC = {0, 0};
560  UINT4 offset = 0; // Index shift between freqs and the frequency series
561  REAL8Sequence *freqs = NULL;
562  if (deltaF > 0) { // freqs contains uniform frequency grid with spacing deltaF; we start at frequency 0
563  /* Set up output array with size closest power of 2 */
564  npts = NextPow2(fHigh_geom / deltaF_geom) + 1;
565  if (fHigh_geom < fHigh * Mtot_sec) /* Resize waveform if user wants f_max larger than cutoff frequency */
566  npts = NextPow2(fHigh * Mtot_sec / deltaF_geom) + 1;
567 
568  XLALGPSAdd(&tC, -1. / deltaF); /* coalesce at t=0 */
569  *hptilde = XLALCreateCOMPLEX16FrequencySeries("hptilde: FD waveform", &tC, 0.0, deltaF, &lalStrainUnit, npts);
570  *hctilde = XLALCreateCOMPLEX16FrequencySeries("hctilde: FD waveform", &tC, 0.0, deltaF, &lalStrainUnit, npts);
571 
572  // Recreate freqs using only the lower and upper bounds
573  // Use fLow, fHigh and deltaF rather than geometric frequencies for numerical accuracy
574  double fHigh_temp = fHigh_geom / Mtot_sec;
575  UINT4 iStart = (UINT4) ceil(fLow / deltaF);
576  UINT4 iStop = (UINT4) ceil(fHigh_temp / deltaF);
577  freqs = XLALCreateREAL8Sequence(iStop - iStart);
578  if (!freqs) {
579  XLAL_ERROR(XLAL_EFUNC, "Frequency array allocation failed.");
580  }
581  for (UINT4 i=iStart; i<iStop; i++)
582  freqs->data[i-iStart] = i*deltaF_geom;
583 
584  offset = iStart;
585  } else { // freqs contains frequencies with non-uniform spacing; we start at lowest given frequency
586  npts = freqs_in->length;
587  *hptilde = XLALCreateCOMPLEX16FrequencySeries("hptilde: FD waveform", &tC, fLow, 0, &lalStrainUnit, npts);
588  *hctilde = XLALCreateCOMPLEX16FrequencySeries("hctilde: FD waveform", &tC, fLow, 0, &lalStrainUnit, npts);
589  offset = 0;
590 
591  freqs = XLALCreateREAL8Sequence(freqs_in->length);
592  if (!freqs) {
593  XLAL_ERROR(XLAL_EFUNC, "Frequency array allocation failed.");
594  }
595  for (UINT4 i=0; i<freqs_in->length; i++)
596  freqs->data[i] = freqs_in->data[i] * Mtot_sec;
597  }
598 
599 
600  if (!(*hptilde) || !(*hctilde))
601  {
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);
609  SEOBNRROMdata_coeff_Cleanup(romdata_coeff);
610  XLAL_ERROR(XLAL_EFUNC, "Waveform allocation failed.");
611  }
612  memset((*hptilde)->data->data, 0, npts * sizeof(COMPLEX16));
613  memset((*hctilde)->data->data, 0, npts * sizeof(COMPLEX16));
614 
615  XLALUnitMultiply(&(*hptilde)->sampleUnits, &(*hptilde)->sampleUnits, &lalSecondUnit);
616  XLALUnitMultiply(&(*hctilde)->sampleUnits, &(*hctilde)->sampleUnits, &lalSecondUnit);
617 
618  COMPLEX16 *pdata=(*hptilde)->data->data;
619  COMPLEX16 *cdata=(*hctilde)->data->data;
620 
621  REAL8 cosi = cos(inclination);
622  REAL8 pcoef = 0.5*(1.0 + cosi*cosi);
623  REAL8 ccoef = cosi;
624 
625  REAL8 s = 1.0/sqrt(2.0); // Scale polarization amplitude so that strain agrees with FFT of SEOBNRv1
626  double Mtot = Mtot_sec / LAL_MTSUN_SI;
627  double amp0 = Mtot * amp_pre * Mtot_sec * LAL_MRSUN_SI / (distance); // Correct overall amplitude to undo mass-dependent scaling used in single-spin ROM
628 
629  // Evaluate reference phase for setting phiRef correctly
630  double phase_change = gsl_spline_eval(spline_phi, fRef_geom, acc_phi) - 2*phiRef;
631 
632  // Assemble waveform from aplitude and phase
633  for (UINT4 i=0; i<freqs->length; i++) { // loop over frequency points in sequence
634  double f = freqs->data[i];
635  if (f > Mf_ROM_max) continue; // We're beyond the highest allowed frequency; since freqs may not be ordered, we'll just skip the current frequency and leave zero in the buffer
636  int j = i + offset; // shift index for frequency series if needed
637  double A = gsl_spline_eval(spline_amp, f, acc_amp);
638  double phase = gsl_spline_eval(spline_phi, f, acc_phi) - phase_change;
639  COMPLEX16 htilde = s*amp0*A * cexp(I*phase);
640  pdata[j] = pcoef * htilde;
641  cdata[j] = -I * ccoef * htilde;
642  }
643 
644  /* Correct phasing so we coalesce at t=0 (with the definition of the epoch=-1/deltaF above) */
645 
646  // Get SEOBNRv1 ringdown frequency for 22 mode
647  double Mf_final = SEOBNRROM_Ringdown_Mf_From_Mtot_q(Mtot_sec, q, chi, chi, SEOBNRv1);
648 
649  UINT4 L = freqs->length;
650  // prevent gsl interpolation errors
651  if (Mf_final > freqs->data[L-1])
652  Mf_final = freqs->data[L-1];
653  if (Mf_final < freqs->data[0])
654  {
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);
662  SEOBNRROMdata_coeff_Cleanup(romdata_coeff);
663  XLAL_ERROR(XLAL_EDOM, "f_ringdown < f_min");
664  }
665 
666  // Time correction is t(f_final) = 1/(2pi) dphi/df (f_final)
667  // We compute the dimensionless time correction t/M since we use geometric units.
668  REAL8 t_corr = gsl_spline_eval_deriv(spline_phi, Mf_final, acc_phi) / (2*LAL_PI);
669 
670  // Now correct phase
671  for (UINT4 i=0; i<freqs->length; i++) { // loop over frequency points in sequence
672  double f = freqs->data[i] - fRef_geom;
673  int j = i + offset; // shift index for frequency series if needed
674  pdata[j] *= cexp(-2*LAL_PI * I * f * t_corr);
675  cdata[j] *= cexp(-2*LAL_PI * I * f * t_corr);
676  }
677 
679 
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);
686  SEOBNRROMdata_coeff_Cleanup(romdata_coeff);
687 
688  return(XLAL_SUCCESS);
689 }
690 
691 /**
692  * @addtogroup LALSimIMRSEOBNRROM_c
693  *
694  * @brief Functions for producing SEOBNRv1 and v2 waveforms
695  * using reduced order models.
696  *
697  * @review SEOBNRv1/2_ROM_(Effective/Double)Spin reviewed by Frank Ohme, Sarah Caudill, Michael Puerrer, Ian Harry, John Veitch, Gareth Thomas. Review concluded with git hash 9dc5e84583bfe2707ac20638e7b89bf988d4d482 (July 2015).
698  *
699  * @{
700  *
701  * @name SEOBNRv1 Reduced Order Model (Effective Spin)
702  *
703  * @author Michael Puerrer, John Veitch
704  *
705  * @brief C code for SEOBNRv1 reduced order model (equal spin version).
706  * See CQG 31 195010, 2014, arXiv:1402.4146 for details.
707  *
708  * This is a frequency domain model that approximates the time domain SEOBNRv1 model with equal spins.
709  * Note that SEOBNRv2 supersedes SEOBNRv1.
710  *
711  * The binary data files are available at https://dcc.ligo.org/T1400701-v1.
712  * Put the untared data into a location in your LAL_DATA_PATH.
713  *
714  * @note Note that due to its construction the iFFT of the ROM has a small (~ 20 M) offset
715  * in the peak time that scales with total mass as compared to the time-domain SEOBNRv1 model.
716  *
717  * @note Due to non-smoothness in SEOBNRv1 at chi1=chi2 ~ -0.8 and 20 <= q <= 40 the
718  * ROM deviates from the SEOBNRv1 behavior there. See arXiv:1402.4146, Fig 7,
719  * Fig 11, and Fig 13 for details.
720  *
721  * @note Parameter ranges:
722  * * 1 <= q <= 100
723  * * -1 <= chi <= 0.6
724  * * Mtot >= 1.4Msun
725  *
726  * Equal spin chi = chi1 = chi2.
727  * Asymmetric mass-ratio q = max(m1/m2, m2/m1).
728  * Total mass Mtot.
729  *
730  * @{
731  */
732 
733 /**
734  * Compute waveform in LAL format at specified frequencies for the SEOBNRv1_ROM_EffectiveSpin model.
735  *
736  * XLALSimIMRSEOBNRv1ROMEffectiveSpin() returns the plus and cross polarizations as a complex
737  * frequency series with equal spacing deltaF and contains zeros from zero frequency
738  * to the starting frequency and zeros beyond the cutoff frequency in the ringdown.
739  *
740  * In contrast, XLALSimIMRSEOBNRv1ROMEffectiveSpinFrequencySequence() returns a
741  * complex frequency series with entries exactly at the frequencies specified in
742  * the sequence freqs (which can be unequally spaced). No zeros are added.
743  *
744  * If XLALSimIMRSEOBNRv1ROMEffectiveSpinFrequencySequence() is called with frequencies that
745  * are beyond the maxium allowed geometric frequency for the ROM, zero strain is returned.
746  * It is not assumed that the frequency sequence is ordered.
747  *
748  * This function is designed as an entry point for reduced order quadratures.
749  */
751  struct tagCOMPLEX16FrequencySeries **hptilde, /**< Output: Frequency-domain waveform h+ */
752  struct tagCOMPLEX16FrequencySeries **hctilde, /**< Output: Frequency-domain waveform hx */
753  const REAL8Sequence *freqs, /**< Frequency points at which to evaluate the waveform (Hz), need to be strictly monotonically increasing */
754  REAL8 phiRef, /**< Orbital phase at reference time */
755  REAL8 fRef, /**< Reference frequency (Hz); 0 defaults to fLow */
756  REAL8 distance, /**< Distance of source (m) */
757  REAL8 inclination, /**< Inclination of source (rad) */
758  REAL8 m1SI, /**< Mass of companion 1 (kg) */
759  REAL8 m2SI, /**< Mass of companion 2 (kg) */
760  REAL8 chi) /**< Effective aligned spin */
761 {
762 
763  /* Get masses in terms of solar mass */
764  double mass1 = m1SI / LAL_MSUN_SI;
765  double mass2 = m2SI / LAL_MSUN_SI;
766  double Mtot = mass1+mass2;
767  double q = mass2 / mass1;
768  if(q<1.0) q=1./q;
769  /* Total mass in seconds */
770  double Mtot_sec = Mtot * LAL_MTSUN_SI;
771 
772  // Load ROM data if not already loaded
773 #ifdef LAL_PTHREAD_LOCK
774  (void) pthread_once(&SEOBNRv1ROMEffectiveSpin_is_initialized, SEOBNRv1ROMEffectiveSpin_Init_LALDATA);
775 #else
777 #endif
778 
779  if(!SEOBNRv1ROMEffectiveSpin_IsSetup()) XLAL_ERROR(XLAL_EFAILED,"Error setting up SEOBNRv1ROMEffectiveSpin model - check your $LAL_DATA_PATH\n");
780 
781  // Call the internal core function with deltaF = 0 to indicate that freqs is non-uniformly
782  // spaced and we want the strain only at these frequencies
783  int retcode = SEOBNRv1ROMEffectiveSpinCore(hptilde, hctilde,
784  phiRef, fRef, distance, inclination, Mtot_sec, q, chi, freqs, 0);
785 
786  return(retcode);
787 }
788 
789 /**
790  * Compute waveform in LAL format for the SEOBNRv1_ROM_EffectiveSpin model.
791  *
792  * Returns the plus and cross polarizations as a complex frequency series with
793  * equal spacing deltaF and contains zeros from zero frequency to the starting
794  * frequency fLow and zeros beyond the cutoff frequency fHigh to the next power of 2 in
795  * the size of the frequency series.
796  */
798  struct tagCOMPLEX16FrequencySeries **hptilde, /**< Output: Frequency-domain waveform h+ */
799  struct tagCOMPLEX16FrequencySeries **hctilde, /**< Output: Frequency-domain waveform hx */
800  REAL8 phiRef, /**< Orbital phase at reference time */
801  REAL8 deltaF, /**< Sampling frequency (Hz) */
802  REAL8 fLow, /**< Starting GW frequency (Hz) */
803  REAL8 fHigh, /**< End frequency; 0 defaults to Mf=0.14 */
804  REAL8 fRef, /**< Reference frequency (Hz); 0 defaults to fLow */
805  REAL8 distance, /**< Distance of source (m) */
806  REAL8 inclination, /**< Inclination of source (rad) */
807  REAL8 m1SI, /**< Mass of companion 1 (kg) */
808  REAL8 m2SI, /**< Mass of companion 2 (kg) */
809  REAL8 chi) /**< Effective aligned spin */
810 {
811 
812  /* Get masses in terms of solar mass */
813  double mass1 = m1SI / LAL_MSUN_SI;
814  double mass2 = m2SI / LAL_MSUN_SI;
815  double Mtot = mass1+mass2;
816  double q = mass2 / mass1;
817  if(q<1.0) q=1./q;
818  /* Total mass in seconds */
819  double Mtot_sec = Mtot * LAL_MTSUN_SI;
820 
821  if(fRef==0.0)
822  fRef=fLow;
823 
824  // Load ROM data if not already loaded
825 #ifdef LAL_PTHREAD_LOCK
826  (void) pthread_once(&SEOBNRv1ROMEffectiveSpin_is_initialized, SEOBNRv1ROMEffectiveSpin_Init_LALDATA);
827 #else
829 #endif
830 
831  if(!SEOBNRv1ROMEffectiveSpin_IsSetup()) XLAL_ERROR(XLAL_EFAILED,"Error setting up SEOBNRv1ROMEffectiveSpin data - check your $LAL_DATA_PATH\n");
832 
833  // Use fLow, fHigh, deltaF to compute freqs sequence
834  // Instead of building a full sequency we only transfer the boundaries and let
835  // the internal core function do the rest (and properly take care of corner cases).
837  freqs->data[0] = fLow;
838  freqs->data[1] = fHigh;
839 
840  int retcode = SEOBNRv1ROMEffectiveSpinCore(hptilde, hctilde,
841  phiRef, fRef, distance, inclination, Mtot_sec, q, chi, freqs, deltaF);
842 
844 
845  return(retcode);
846 }
847 
848 /** @} */
849 /** @} */
850 
851 /** Setup SEOBNRv1ROMEffectiveSpin model using data files installed in $LAL_DATA_PATH
852  */
854 {
855  if (SEOBNRv1ROMEffectiveSpin_IsSetup()) return;
856 
857  // If we find one ROM datafile in a directory listed in LAL_DATA_PATH,
858  // then we expect the remaining datafiles to also be there.
859  char datafile[] = "SEOBNRv1ROM_SS_Phase_ciall.dat";
860 
861  char *path = XLAL_FILE_RESOLVE_PATH(datafile);
862  if (path==NULL)
863  XLAL_ERROR_VOID(XLAL_EIO, "Unable to resolve data file %s in $LAL_DATA_PATH\n", datafile);
864  char *dir = dirname(path);
865  int ret = SEOBNRv1ROMEffectiveSpin_Init(dir);
866  XLALFree(path);
867 
868  if(ret!=XLAL_SUCCESS)
869  XLAL_ERROR_VOID(XLAL_FAILURE, "Unable to find SEOBNRv1ROMEffectiveSpin data files in $LAL_DATA_PATH\n");
870 }
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 const double gA[]
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)
int s
Definition: bh_qnmode.c:137
double i
Definition: bh_ringdown.c:118
double e
Definition: bh_ringdown.c:117
sigmaKerr data[0]
#define XLAL_FILE_RESOLVE_PATH(fname)
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_MTSUN_SI
#define LAL_MRSUN_SI
double complex COMPLEX16
double REAL8
uint32_t UINT4
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
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.
static const INT4 m
static const INT4 q
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(...)
#define XLAL_ERROR(...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
int int XLALPrintWarning(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
XLAL_EIO
XLAL_EFAILED
XLAL_FAILURE
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
path
int N
int pthread_once(pthread_once_t *once_control, void(*init_routine)(void))
REAL8 * data
gsl_bspline_workspace * bwx
gsl_bspline_workspace * bwy