Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALSimInspiralWaveformCache.c
Go to the documentation of this file.
1/*
2 * Copyright (C) 2013 Evan Ochsner and Will M. Farr
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#include <math.h>
22#include <lal/LALSimInspiral.h>
23#include <lal/LALSimIMR.h>
24#include <lal/FrequencySeries.h>
25#include <lal/Sequence.h>
26#include <lal/LALConstants.h>
27#include <lal/LALSimInspiralEOS.h>
28
31
32/**
33 * Bitmask enumerating which parameters have changed, to determine
34 * if the requested waveform can be transformed from a cached waveform
35 * or if it must be generated from scratch.
36 */
37typedef enum {
42 INCLINATION = 8
44
47 REAL8 phiRef,
48 REAL8 deltaTF,
49 REAL8 m1,
50 REAL8 m2,
51 REAL8 S1x, REAL8 S1y, REAL8 S1z,
52 REAL8 S2x, REAL8 S2y, REAL8 S2z,
53 REAL8 f_min, REAL8 f_ref, REAL8 f_max,
54 REAL8 r,
55 REAL8 i,
56 LALDict *LALpars,
58 REAL8Sequence *frequencies);
59
61 REAL8Sequence *newFrequencies,
62 REAL8Sequence *cachedFrequencies);
63
65 REAL8TimeSeries *hplus,
66 REAL8TimeSeries *hcross,
67 REAL8 phiRef,
69 REAL8 m1, REAL8 m2,
70 REAL8 S1x, REAL8 S1y, REAL8 S1z,
71 REAL8 S2x, REAL8 S2y, REAL8 S2z,
72 REAL8 f_min, REAL8 f_ref,
73 REAL8 r,
74 REAL8 i,
75 LALDict *LALpars,
77
81 REAL8 phiRef,
83 REAL8 m1, REAL8 m2,
84 REAL8 S1x, REAL8 S1y, REAL8 S1z,
85 REAL8 S2x, REAL8 S2y, REAL8 S2z,
86 REAL8 f_min, REAL8 f_ref, REAL8 f_max,
87 REAL8 r,
88 REAL8 i,
89 LALDict *LALpars,
91 REAL8Sequence *frequencies);
92
93
94/**
95 * @addtogroup LALSimInspiralWaveformCache_h
96 * @{
97 */
98
99/**
100 * Chooses between different approximants when requesting a waveform to be generated
101 * Returns the waveform in the time domain.
102 * The parameters passed must be in SI units.
103 *
104 * This version allows caching of waveforms. The most recently generated
105 * waveform and its parameters are stored. If the next call requests a waveform
106 * that can be obtained by a simple transformation, then it is done.
107 * This bypasses the waveform generation and speeds up the code.
108 */
110 REAL8TimeSeries **hplus, /**< +-polarization waveform */
111 REAL8TimeSeries **hcross, /**< x-polarization waveform */
112 REAL8 phiRef, /**< reference orbital phase (rad) */
113 REAL8 deltaT, /**< sampling interval (s) */
114 REAL8 m1, /**< mass of companion 1 (kg) */
115 REAL8 m2, /**< mass of companion 2 (kg) */
116 REAL8 S1x, /**< x-component of the dimensionless spin of object 1 */
117 REAL8 S1y, /**< y-component of the dimensionless spin of object 1 */
118 REAL8 S1z, /**< z-component of the dimensionless spin of object 1 */
119 REAL8 S2x, /**< x-component of the dimensionless spin of object 2 */
120 REAL8 S2y, /**< y-component of the dimensionless spin of object 2 */
121 REAL8 S2z, /**< z-component of the dimensionless spin of object 2 */
122 REAL8 f_min, /**< starting GW frequency (Hz) */
123 REAL8 f_ref, /**< reference GW frequency (Hz) */
124 REAL8 r, /**< distance of source (m) */
125 REAL8 i, /**< inclination of source (rad) */
126 LALDict *LALpars, /**< LALDictionary containing non-mandatory variables/flags */
127 Approximant approximant, /**< post-Newtonian approximant to use for waveform production */
128 LALSimInspiralWaveformCache *cache /**< waveform cache structure; use NULL for no caching */
129 )
130{
131 int status;
132 size_t j;
133 REAL8 phasediff, dist_ratio, incl_ratio_plus, incl_ratio_cross;
134 REAL8 cosrot, sinrot;
135 CacheVariableDiffersBitmask changedParams;
136
137 // If nonGRparams are not NULL, don't even try to cache.
138 if ( !XLALSimInspiralWaveformParamsNonGRAreDefault(LALpars) || (!cache) ||
149
150
151 return XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2, S1x, S1y, S1z, S2x, S2y, S2z,
152 r, i, phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars,
154
155 // Check which parameters have changed
156 changedParams = CacheArgsDifferenceBitmask(cache, phiRef, deltaT,
157 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, 0., r, i,
158 LALpars, approximant, NULL);
159
160 // No parameters have changed! Copy the cached polarizations
161 if( changedParams == NO_DIFFERENCE ) {
162 *hplus = XLALCutREAL8TimeSeries(cache->hplus, 0,
163 cache->hplus->data->length);
164 if (*hplus == NULL) return XLAL_ENOMEM;
165 *hcross = XLALCutREAL8TimeSeries(cache->hcross, 0,
166 cache->hcross->data->length);
167 if (*hcross == NULL) {
169 *hplus = NULL;
170 return XLAL_ENOMEM;
171 }
172
173 return XLAL_SUCCESS;
174 }
175
176 // Intrinsic parameters have changed. We must generate a new waveform
177 if( (changedParams & INTRINSIC) != 0 ) {
178
179 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2, S1x, S1y, S1z, S2x, S2y, S2z,
180 r, i, phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars,
182 if (status == XLAL_FAILURE) return status;
183
184 // FIXME: Need to add hlms, dynamic variables, etc. in cache
185 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
186 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i, LALpars, approximant);
187 }
188
190 // case 1: Precessing waveforms
192 // If polarizations are not cached we must generate a fresh waveform
193 // FIXME: Will need to check hlms and/or dynamical variables as well
194 if( cache->hplus == NULL || cache->hcross == NULL) {
195 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
196 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
197 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars,
199 if (status == XLAL_FAILURE) return status;
200
201 // FIXME: Need to add hlms, dynamic variables, etc. in cache
202 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
203 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
204 LALpars, approximant);
205 }
206
207 if( changedParams & INCLINATION ) {
208 // FIXME: For now just treat as intrinsic parameter.
209 // Will come back and put in transformation
210 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
211 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
212 phiRef, 0., 0., 0., deltaT, f_min, f_ref,
213 LALpars, approximant);
214 if (status == XLAL_FAILURE) return status;
215
216 // FIXME: Need to add hlms, dynamic variables, etc. in cache
217 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
218 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
219 LALpars, approximant);
220 }
221 if( changedParams & PHI_REF ) {
222 // FIXME: For now just treat as intrinsic parameter.
223 // Will come back and put in transformation
224 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
225 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
226 phiRef, 0., 0., 0., deltaT, f_min, f_ref,
227 LALpars, approximant);
228 if (status == XLAL_FAILURE) return status;
229
230 // FIXME: Need to add hlms, dynamic variables, etc. in cache
231 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
232 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
233 LALpars, approximant);
234 }
235 if( (changedParams & DISTANCE) != 0 ) {
236 // Return rescaled copy of cached polarizations
237 dist_ratio = cache->r / r;
238 *hplus = XLALCreateREAL8TimeSeries(cache->hplus->name,
239 &(cache->hplus->epoch), cache->hplus->f0,
240 cache->hplus->deltaT, &(cache->hplus->sampleUnits),
241 cache->hplus->data->length);
242 if (*hplus == NULL) return XLAL_ENOMEM;
243
244 *hcross = XLALCreateREAL8TimeSeries(cache->hcross->name,
245 &(cache->hcross->epoch), cache->hcross->f0,
246 cache->hcross->deltaT, &(cache->hcross->sampleUnits),
247 cache->hcross->data->length);
248 if (*hcross == NULL) {
250 *hplus = NULL;
251 return XLAL_ENOMEM;
252 }
253
254 for (j = 0; j < cache->hplus->data->length; j++) {
255 (*hplus)->data->data[j] = cache->hplus->data->data[j]
256 * dist_ratio;
257 (*hcross)->data->data[j] = cache->hcross->data->data[j]
258 * dist_ratio;
259 }
260 }
261
262 return XLAL_SUCCESS;
263 }
264
265 // case 2: Non-precessing, ampO = 0
266 else if( ampO==0 && (approximant==TaylorT1 || approximant==TaylorT2
269 // If polarizations are not cached we must generate a fresh waveform
270 if( cache->hplus == NULL || cache->hcross == NULL) {
271 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
272 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
273 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars, approximant);
274 if (status == XLAL_FAILURE) return status;
275
276 // FIXME: Need to add hlms, dynamic variables, etc. in cache
277 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
278 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
279 LALpars, approximant);
280 }
281
282 // Set transformation coefficients for identity transformation.
283 // We'll adjust them depending on which extrinsic parameters changed.
284 dist_ratio = incl_ratio_plus = incl_ratio_cross = cosrot = 1.;
285 phasediff = sinrot = 0.;
286
287 if( changedParams & PHI_REF ) {
288 // Only 2nd harmonic present, so {h+,hx} rotates by 2*deltaphiRef
289 phasediff = 2.*(phiRef - cache->phiRef);
290 cosrot = cos(phasediff);
291 sinrot = sin(phasediff);
292 }
293 if( changedParams & INCLINATION) {
294 // Rescale h+, hx by ratio of new/old inclination dependence
295 incl_ratio_plus = (1.0 + cos(i)*cos(i))
296 / (1.0 + cos(cache->i)*cos(cache->i));
297 incl_ratio_cross = cos(i) / cos(cache->i);
298 }
299 if( changedParams & DISTANCE ) {
300 // Rescale h+, hx by ratio of (1/new_dist)/(1/old_dist) = old/new
301 dist_ratio = cache->r / r;
302 }
303
304 // Create the output polarizations
305 *hplus = XLALCreateREAL8TimeSeries(cache->hplus->name,
306 &(cache->hplus->epoch), cache->hplus->f0,
307 cache->hplus->deltaT, &(cache->hplus->sampleUnits),
308 cache->hplus->data->length);
309 if (*hplus == NULL) return XLAL_ENOMEM;
310 *hcross = XLALCreateREAL8TimeSeries(cache->hcross->name,
311 &(cache->hcross->epoch), cache->hcross->f0,
312 cache->hcross->deltaT, &(cache->hcross->sampleUnits),
313 cache->hcross->data->length);
314 if (*hcross == NULL) {
316 *hplus = NULL;
317 return XLAL_ENOMEM;
318 }
319
320 // Get new polarizations by transforming the old
321 incl_ratio_plus *= dist_ratio;
322 incl_ratio_cross *= dist_ratio;
323 // FIXME: Do changing phiRef and inclination commute?!?!
324 for (j = 0; j < cache->hplus->data->length; j++) {
325 (*hplus)->data->data[j] = incl_ratio_plus
326 * (cosrot*cache->hplus->data->data[j]
327 - sinrot*cache->hcross->data->data[j]);
328 (*hcross)->data->data[j] = incl_ratio_cross
329 * (sinrot*cache->hplus->data->data[j]
330 + cosrot*cache->hcross->data->data[j]);
331 }
332
333 return XLAL_SUCCESS;
334 }
335 // case 3: Non-precessing, ampO > 0
336 // FIXME: EOBNRv2HM and TEOBResumS actually ignores ampO. If it's given with ampO==0,
337 // it will fall to the catch-all and not be cached.
338 else if( (ampO==-1 || ampO>0) && (approximant==TaylorT1
341 || approximant==TEOBResumS) ) {
342 // If polarizations are not cached we must generate a fresh waveform
343 // FIXME: Add in check that hlms non-NULL
344 if( cache->hplus == NULL || cache->hcross == NULL) {
345 // FIXME: This will change to a code-path: inputs->hlms->{h+,hx}
346 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
347 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
348 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars, approximant);
349 if (status == XLAL_FAILURE) return status;
350
351 // FIXME: Need to add hlms, dynamic variables, etc. in cache
352 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
353 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
354 LALpars, approximant);
355 }
356
357 if( changedParams & INCLINATION) {
358 // FIXME: For now just treat as intrinsic parameter.
359 // Will come back and put in transformation
360 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
361 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
362 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars, approximant);
363 if (status == XLAL_FAILURE) return status;
364
365 // FIXME: Need to add hlms, dynamic variables, etc. in cache
366 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
367 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
368 LALpars, approximant);
369
370 }
371 if( changedParams & PHI_REF ) {
372 // FIXME: For now just treat as intrinsic parameter.
373 // Will come back and put in transformation
374 status = XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
375 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
376 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars,
378 if (status == XLAL_FAILURE) return status;
379
380 // FIXME: Need to add hlms, dynamic variables, etc. in cache
381 return StoreTDHCache(cache, *hplus, *hcross, phiRef, deltaT, m1, m2,
382 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, r, i,
383 LALpars, approximant);
384
385 }
386 if( changedParams & DISTANCE ) {
387 // Return rescaled copy of cached polarizations
388 dist_ratio = cache->r / r;
389 *hplus = XLALCreateREAL8TimeSeries(cache->hplus->name,
390 &(cache->hplus->epoch), cache->hplus->f0,
391 cache->hplus->deltaT, &(cache->hplus->sampleUnits),
392 cache->hplus->data->length);
393 if (*hplus == NULL) return XLAL_ENOMEM;
394
395 *hcross = XLALCreateREAL8TimeSeries(cache->hcross->name,
396 &(cache->hcross->epoch), cache->hcross->f0,
397 cache->hcross->deltaT, &(cache->hcross->sampleUnits),
398 cache->hcross->data->length);
399 if (*hcross == NULL) {
401 *hplus = NULL;
402 return XLAL_ENOMEM;
403 }
404
405 for (j = 0; j < cache->hplus->data->length; j++) {
406 (*hplus)->data->data[j] = cache->hplus->data->data[j]
407 * dist_ratio;
408 (*hcross)->data->data[j] = cache->hcross->data->data[j]
409 * dist_ratio;
410 }
411 }
412
413 return XLAL_SUCCESS;
414 }
415
416 // Catch-all. Unsure what to do, don't try to cache.
417 // Basically, you requested a waveform type which is not setup for caching
418 // b/c of lack of interest or it's unclear what/how to cache for that model
419 else {
420 return XLALSimInspiralChooseTDWaveform(hplus, hcross, m1, m2,
421 S1x, S1y, S1z, S2x, S2y, S2z, r, i,
422 phiRef, 0., 0., 0., deltaT, f_min, f_ref, LALpars, approximant);
423 }
424}
425
426/**
427 * Chooses between different approximants when requesting a waveform to be generated
428 * Returns the waveform in the frequency domain.
429 * The parameters passed must be in SI units.
430 *
431 * This version allows caching of waveforms. The most recently generated
432 * waveform and its parameters are stored. If the next call requests a waveform
433 * that can be obtained by a simple transformation, then it is done.
434 * This bypasses the waveform generation and speeds up the code.
435 */
437 COMPLEX16FrequencySeries **hptilde, /**< +-polarization waveform */
438 COMPLEX16FrequencySeries **hctilde, /**< x-polarization waveform */
439 REAL8 phiRef, /**< reference orbital phase (rad) */
440 REAL8 deltaF, /**< sampling interval (Hz) */
441 REAL8 m1, /**< mass of companion 1 (kg) */
442 REAL8 m2, /**< mass of companion 2 (kg) */
443 REAL8 S1x, /**< x-component of the dimensionless spin of object 1 */
444 REAL8 S1y, /**< y-component of the dimensionless spin of object 1 */
445 REAL8 S1z, /**< z-component of the dimensionless spin of object 1 */
446 REAL8 S2x, /**< x-component of the dimensionless spin of object 2 */
447 REAL8 S2y, /**< y-component of the dimensionless spin of object 2 */
448 REAL8 S2z, /**< z-component of the dimensionless spin of object 2 */
449 REAL8 f_min, /**< starting GW frequency (Hz) */
450 REAL8 f_max, /**< ending GW frequency (Hz) */
451 REAL8 f_ref, /**< Reference GW frequency (Hz) */
452 REAL8 r, /**< distance of source (m) */
453 REAL8 i, /**< inclination of source (rad) */
454 LALDict *LALpars, /**< LALDictionary containing non-mandatory variables/flags */
455 Approximant approximant, /**< post-Newtonian approximant to use for waveform production */
456 LALSimInspiralWaveformCache *cache, /**< waveform cache structure */
457 REAL8Sequence *frequencies /**< sequence of frequencies for which the waveform will be computed. Pass in NULL (or None in python) for standard f_min to f_max sequence. */
458 )
459{
460 int status;
461 size_t j;
462 REAL8 dist_ratio, incl_ratio_plus, incl_ratio_cross, phase_diff;
463 COMPLEX16 exp_dphi;
464 CacheVariableDiffersBitmask changedParams;
465
466
467 // If nonGRparams are not NULL, don't even try to cache.
468 if ( !XLALSimInspiralWaveformParamsNonGRAreDefault(LALpars) || (!cache) ) {
469 if (frequencies != NULL)
470 return XLALSimInspiralChooseFDWaveformSequence(hptilde, hctilde, phiRef,
471 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_ref,
472 r, i,
473 LALpars, approximant,frequencies);
474 else
475 return XLALSimInspiralChooseFDWaveform(hptilde, hctilde, m1, m2,
476 S1x, S1y, S1z, S2x, S2y, S2z, r, i, phiRef,
477 0., 0., 0., deltaF, f_min, f_max, f_ref,
478 LALpars,
480 }
481
482 // Check which parameters have changed
483 changedParams = CacheArgsDifferenceBitmask(cache, phiRef, deltaF,
484 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, f_max, r, i,
485 LALpars, approximant, frequencies);
486
487 // No parameters have changed! Copy the cached polarizations
488 if( changedParams == NO_DIFFERENCE ) {
489 *hptilde = XLALCutCOMPLEX16FrequencySeries(cache->hptilde, 0,
490 cache->hptilde->data->length);
491 if (*hptilde == NULL) return XLAL_ENOMEM;
492 *hctilde = XLALCutCOMPLEX16FrequencySeries(cache->hctilde, 0,
493 cache->hctilde->data->length);
494 if (*hctilde == NULL) {
496 *hptilde = NULL;
497 return XLAL_ENOMEM;
498 }
499
500 return XLAL_SUCCESS;
501 }
502
503 // Intrinsic parameters have changed. We must generate a new waveform
504 if( (changedParams & INTRINSIC) != 0 ) {
505 if ( frequencies != NULL ){
506 status = XLALSimInspiralChooseFDWaveformSequence(hptilde, hctilde, phiRef,
507 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_ref,
508 r, i, LALpars, approximant, frequencies);
509 }
510 else {
511 status = XLALSimInspiralChooseFDWaveform(hptilde, hctilde, m1, m2,
512 S1x, S1y, S1z, S2x, S2y, S2z,
513 r, i, phiRef, 0., 0., 0.,
514 deltaF, f_min, f_max, f_ref,
515 LALpars, approximant);
516 }
517 if (status == XLAL_FAILURE) return status;
518
519 return StoreFDHCache(cache, *hptilde, *hctilde, phiRef, deltaF, m1, m2,
520 S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, f_max, r, i, LALpars, approximant, frequencies);
521 }
522
523 // case 1: Non-precessing, 2nd harmonic only
527 || approximant == IMRPhenomC ) {
528 // If polarizations are not cached we must generate a fresh waveform
529 // FIXME: Will need to check hlms and/or dynamical variables as well
530 if( cache->hptilde == NULL || cache->hctilde == NULL) {
531 if ( frequencies != NULL ){
532 status = XLALSimInspiralChooseFDWaveformSequence(hptilde, hctilde, phiRef,
533 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_ref,
534 r, i, LALpars, approximant,frequencies);
535 }
536 else {
537 status = XLALSimInspiralChooseFDWaveform(hptilde, hctilde, m1, m2,
538 S1x, S1y, S1z, S2x, S2y, S2z, r, i, phiRef,
539 0., 0., 0., deltaF, f_min, f_max, f_ref,
540 LALpars, approximant);
541 }
542 if (status == XLAL_FAILURE) return status;
543
544 return StoreFDHCache(cache, *hptilde, *hctilde, phiRef, deltaF,
545 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_ref, f_max, r, i,
546 LALpars, approximant, frequencies);
547 }
548
549 // Set transformation coefficients for identity transformation.
550 // We'll adjust them depending on which extrinsic parameters changed.
551 dist_ratio = incl_ratio_plus = incl_ratio_cross = 1.;
552 phase_diff = 0.;
553 exp_dphi = 1.;
554
555 if( changedParams & PHI_REF ) {
556 // Only 2nd harmonic present, so {h+,hx} \propto e^(2 i phiRef)
557 phase_diff = 2.*(phiRef - cache->phiRef);
558 exp_dphi = cpolar(1., phase_diff);
559 }
560 if( changedParams & INCLINATION) {
561 // Rescale h+, hx by ratio of new/old inclination dependence
562 incl_ratio_plus = (1.0 + cos(i)*cos(i))
563 / (1.0 + cos(cache->i)*cos(cache->i));
564 incl_ratio_cross = cos(i) / cos(cache->i);
565 }
566 if( changedParams & DISTANCE ) {
567 // Rescale h+, hx by ratio of (1/new_dist)/(1/old_dist) = old/new
568 dist_ratio = cache->r / r;
569 }
570
571 // Create the output polarizations
573 &(cache->hptilde->epoch), cache->hptilde->f0,
574 cache->hptilde->deltaF, &(cache->hptilde->sampleUnits),
575 cache->hptilde->data->length);
576 if (*hptilde == NULL) return XLAL_ENOMEM;
577
579 &(cache->hctilde->epoch), cache->hctilde->f0,
580 cache->hctilde->deltaF, &(cache->hctilde->sampleUnits),
581 cache->hctilde->data->length);
582 if (*hctilde == NULL) {
584 *hptilde = NULL;
585 return XLAL_ENOMEM;
586 }
587
588 // Get new polarizations by transforming the old
589 incl_ratio_plus *= dist_ratio;
590 incl_ratio_cross *= dist_ratio;
591 for (j = 0; j < cache->hptilde->data->length; j++) {
592 (*hptilde)->data->data[j] = exp_dphi * incl_ratio_plus
593 * cache->hptilde->data->data[j];
594 (*hctilde)->data->data[j] = exp_dphi * incl_ratio_cross
595 * cache->hctilde->data->data[j];
596 }
597
598 return XLAL_SUCCESS;
599 }
600
601 // case 2: Precessing
602 /*else if( approximant == SpinTaylorF2 ) {
603
604 }*/
605
606 // Catch-all. Unsure what to do, don't try to cache.
607 // Basically, you requested a waveform type which is not setup for caching
608 // b/c of lack of interest or it's unclear what/how to cache for that model
609 else {
610 if ( frequencies != NULL ){
611 return XLALSimInspiralChooseFDWaveformSequence(hptilde, hctilde, phiRef,
612 m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_ref,
613 r, i, LALpars, approximant,frequencies);
614 }
615 else {
616 return XLALSimInspiralChooseFDWaveform(hptilde, hctilde, m1, m2,
617 S1x, S1y, S1z, S2x, S2y, S2z,
618 r, i, phiRef, 0., 0., 0.,
619 deltaF, f_min, f_max, f_ref,
620 LALpars, approximant);
621 }
622 }
623
624}
625
626/**
627 * Construct and initialize a waveform cache. Caches are used to
628 * avoid re-computation of waveforms that differ only by simple
629 * scaling relations in extrinsic parameters.
630 */
632{
635
636 return cache;
637}
638
639/**
640 * Destroy a waveform cache.
641 */
643{
644 if (cache != NULL) {
649 if(cache->LALpars) XLALDestroyDict(cache->LALpars);
650 XLALFree(cache);
651 }
652}
653
654/** @} */
655
656/**
657 * Function to compare the requested arguments to those stored in the cache,
658 * returns a bitmask which determines if a cached waveform can be recycled.
659 */
662 REAL8 phiRef,
663 REAL8 deltaTF,
664 REAL8 m1,
665 REAL8 m2,
666 REAL8 S1x, REAL8 S1y, REAL8 S1z,
667 REAL8 S2x, REAL8 S2y, REAL8 S2z,
668 REAL8 f_min, REAL8 f_ref, REAL8 f_max,
669 REAL8 r,
670 REAL8 i,
671 LALDict *LALpars,
673 REAL8Sequence *frequencies
674 )
675{
677 if (cache == NULL) return INTRINSIC;
678
679 if ( !XLALSimInspiralWaveformFlagsEqual(LALpars, cache->LALpars) )
680 return INTRINSIC;
681
682 if ( deltaTF != cache->deltaTF) return INTRINSIC;
683 if ( m1 != cache->m1) return INTRINSIC;
684 if ( m2 != cache->m2) return INTRINSIC;
685 if ( S1x != cache->S1x) return INTRINSIC;
686 if ( S1y != cache->S1y) return INTRINSIC;
687 if ( S1z != cache->S1z) return INTRINSIC;
688 if ( S2x != cache->S2x) return INTRINSIC;
689 if ( S2y != cache->S2y) return INTRINSIC;
690 if ( S2z != cache->S2z) return INTRINSIC;
691 if ( f_min != cache->f_min) return INTRINSIC;
692 if ( f_ref != cache->f_ref) return INTRINSIC;
693 if ( f_max != cache->f_max) return INTRINSIC;
702
705
706 if ( approximant != cache->approximant) return INTRINSIC;
707
708 if (r != cache->r) difference = difference | DISTANCE;
709 if (phiRef != cache->phiRef) difference = difference | PHI_REF;
710 if (i != cache->i) difference = difference | INCLINATION;
711
712 if (FrequenciesAreDifferent(frequencies,cache->frequencies)) return INTRINSIC;
713
714 return difference;
715}
716
717/**
718 * Function to compare two frequencies sequences.
719 * Returns 1 if different, 0 if the same sequences (including if NULL pointers)
720 */
722 REAL8Sequence *newFrequencies,
723 REAL8Sequence *cachedFrequencies
724 )
725{
726 size_t j;
727 if ( newFrequencies == NULL && cachedFrequencies == NULL) return 0;
728 if ( newFrequencies == NULL && cachedFrequencies != NULL) return 1;
729 if ( newFrequencies != NULL && cachedFrequencies == NULL) return 1;
730 if ( newFrequencies->length != cachedFrequencies->length) return 1;
731 for ( j = 0; j < newFrequencies->length; j++){
732 if ( newFrequencies->data[j] != cachedFrequencies->data[j]) return 1;
733 }
734 return 0;
735}
736
737/** Store the output TD hplus and hcross in the cache. */
739 REAL8TimeSeries *hplus,
740 REAL8TimeSeries *hcross,
741 REAL8 phiRef,
743 REAL8 m1, REAL8 m2,
744 REAL8 S1x, REAL8 S1y, REAL8 S1z,
745 REAL8 S2x, REAL8 S2y, REAL8 S2z,
746 REAL8 f_min, REAL8 f_ref,
747 REAL8 r,
748 REAL8 i,
749 LALDict *LALpars,
751 )
752{
753 /* Clear any frequency-domain data. */
754 if (cache->hptilde != NULL) {
756 cache->hptilde = NULL;
757 }
758
759 if (cache->hctilde != NULL) {
761 cache->hctilde = NULL;
762 }
763
764 /* Store params in cache */
765 cache->phiRef = phiRef;
766 cache->deltaTF = deltaT;
767 cache->m1 = m1;
768 cache->m2 = m2;
769 cache->S1x = S1x;
770 cache->S1y = S1y;
771 cache->S1z = S1z;
772 cache->S2x = S2x;
773 cache->S2y = S2y;
774 cache->S2z = S2z;
775 cache->f_min = f_min;
776 cache->f_ref = f_ref;
777 cache->r = r;
778 cache->i = i;
779 if(cache->LALpars) XLALDestroyDict(cache->LALpars);
780 cache->LALpars = XLALDictDuplicate(LALpars);
781 cache->approximant = approximant;
782 cache->frequencies = NULL;
783
784 // Copy over the waveforms
785 // NB: XLALCut... creates a new Series object and copies data and metadata
788 if (hplus == NULL || hcross == NULL || hplus->data == NULL || hcross->data == NULL){
789 XLALPrintError("We have null pointers for h+, hx in StoreTDHCache \n");
790 XLALPrintError("Houston-S, we've got a problem SOS, SOS, SOS, the waveform generator returns NULL!!!... m1 = %.18e, m2 = %.18e, fMin = %.18e, spin1 = {%.18e, %.18e, %.18e}, spin2 = {%.18e, %.18e, %.18e} \n",
791 m1, m2, (double)f_min, S1x, S1y, S1z, S2x, S2y, S2z);
792 return XLAL_ENOMEM;
793 }
794 cache->hplus = XLALCutREAL8TimeSeries(hplus, 0, hplus->data->length);
795 if (cache->hplus == NULL) return XLAL_ENOMEM;
796 cache->hcross = XLALCutREAL8TimeSeries(hcross, 0, hcross->data->length);
797 if (cache->hcross == NULL) {
799 cache->hplus = NULL;
800 return XLAL_ENOMEM;
801 }
802
803 return XLAL_SUCCESS;
804}
805
806/** Store the output FD hptilde and hctilde in cache. */
810 REAL8 phiRef,
812 REAL8 m1, REAL8 m2,
813 REAL8 S1x, REAL8 S1y, REAL8 S1z,
814 REAL8 S2x, REAL8 S2y, REAL8 S2z,
815 REAL8 f_min, REAL8 f_ref, REAL8 f_max,
816 REAL8 r,
817 REAL8 i,
818 LALDict *LALpars,
820 REAL8Sequence *frequencies
821 )
822{
823 /* Clear any time-domain data. */
824 if (cache->hplus != NULL) {
826 cache->hplus = NULL;
827 }
828
829 if (cache->hcross != NULL) {
831 cache->hcross = NULL;
832 }
833
834 /* Store params in cache */
835 cache->phiRef = phiRef;
836 cache->deltaTF = deltaT;
837 cache->m1 = m1;
838 cache->m2 = m2;
839 cache->S1x = S1x;
840 cache->S1y = S1y;
841 cache->S1z = S1z;
842 cache->S2x = S2x;
843 cache->S2y = S2y;
844 cache->S2z = S2z;
845 cache->f_min = f_min;
846 cache->f_ref = f_ref;
847 cache->f_max = f_max;
848 cache->r = r;
849 cache->i = i;
850 if(cache->LALpars) XLALDestroyDict(cache->LALpars);
851 cache->LALpars = XLALDictDuplicate(LALpars);
852 cache->approximant = approximant;
853
855 cache->frequencies = NULL;
856 if (frequencies != NULL){
857 cache->frequencies = XLALCopyREAL8Sequence(frequencies);
858 }
859
860 // Copy over the waveforms
861 // NB: XLALCut... creates a new Series object and copies data and metadata
864 cache->hptilde = XLALCutCOMPLEX16FrequencySeries(hptilde, 0,
865 hptilde->data->length);
866 if (cache->hptilde == NULL) return XLAL_ENOMEM;
867 cache->hctilde = XLALCutCOMPLEX16FrequencySeries(hctilde, 0,
868 hctilde->data->length);
869 if (cache->hctilde == NULL) {
871 cache->hptilde = NULL;
872 return XLAL_ENOMEM;
873 }
874
875 return XLAL_SUCCESS;
876}
877
878/**
879 * Wrapper similar to XLALSimInspiralChooseFDWaveform() for waveforms to be generated a specific freqencies.
880 * Returns the waveform in the frequency domain at the frequencies of the REAL8Sequence frequencies.
881 */
883 COMPLEX16FrequencySeries **hptilde, /**< FD plus polarization */
884 COMPLEX16FrequencySeries **hctilde, /**< FD cross polarization */
885 REAL8 phiRef, /**< reference orbital phase (rad) */
886 REAL8 m1, /**< mass of companion 1 (kg) */
887 REAL8 m2, /**< mass of companion 2 (kg) */
888 REAL8 S1x, /**< x-component of the dimensionless spin of object 1 */
889 REAL8 S1y, /**< y-component of the dimensionless spin of object 1 */
890 REAL8 S1z, /**< z-component of the dimensionless spin of object 1 */
891 REAL8 S2x, /**< x-component of the dimensionless spin of object 2 */
892 REAL8 S2y, /**< y-component of the dimensionless spin of object 2 */
893 REAL8 S2z, /**< z-component of the dimensionless spin of object 2 */
894 REAL8 f_ref, /**< Reference frequency (Hz) */
895 REAL8 distance, /**< distance of source (m) */
896 REAL8 inclination, /**< inclination of source (rad) */
897 LALDict *LALpars, /**< LALDictionary containing non-mandatory variables/flags */
898 Approximant approximant, /**< post-Newtonian approximant to use for waveform production */
899 REAL8Sequence *frequencies /**< sequence of frequencies for which the waveform will be computed. Pass in NULL (or None in python) for standard f_min to f_max sequence. */
900)
901{
902 int ret;
903 unsigned int j;
904 REAL8 pfac, cfac;
905 PNPhasingSeries pfa;
906
907 /* Support variables for precessing wfs*/
908 //REAL8 incl;
909 //REAL8 spin1x,spin1y,spin1z;
910 //REAL8 spin2x,spin2y,spin2z;
911
912 /* Variables for IMRPhenomP and IMRPhenomPv2 */
913 REAL8 chi1_l, chi2_l, chip, thetaJN, alpha0, phi_aligned, zeta_polariz, cos_2zeta, sin_2zeta;
914 COMPLEX16 PhPpolp, PhPpolc;
915
916 LALDict *LALparams_aux;
917
918 /* General sanity checks that will abort
919 *
920 * If non-GR approximants are added, include them in
921 * XLALSimInspiralApproximantAcceptTestGRParams()
922 */
924 XLALPrintError("XLAL Error - %s: Passed in non-NULL testGRparams for an approximant that does not use them\n", __func__);
926 }
927 if (!frequencies) XLAL_ERROR(XLAL_EFAULT);
928 REAL8 f_min = frequencies->data[0];
929
930 /* General sanity check the input parameters - only give warnings! */
931 if( m1 < 0.09 * LAL_MSUN_SI )
932 XLALPrintWarning("XLAL Warning - %s: Small value of m1 = %e (kg) = %e (Msun) requested...Perhaps you have a unit conversion error?\n", __func__, m1, m1/LAL_MSUN_SI);
933 if( m2 < 0.09 * LAL_MSUN_SI )
934 XLALPrintWarning("XLAL Warning - %s: Small value of m2 = %e (kg) = %e (Msun) requested...Perhaps you have a unit conversion error?\n", __func__, m2, m2/LAL_MSUN_SI);
935 if( m1 + m2 > 1000. * LAL_MSUN_SI )
936 XLALPrintWarning("XLAL Warning - %s: Large value of total mass m1+m2 = %e (kg) = %e (Msun) requested...Signal not likely to be in band of ground-based detectors.\n", __func__, m1+m2, (m1+m2)/LAL_MSUN_SI);
937 if( S1x*S1x + S1y*S1y + S1z*S1z > 1.000001 )
938 XLALPrintWarning("XLAL Warning - %s: S1 = (%e,%e,%e) with norm > 1 requested...Are you sure you want to violate the Kerr bound?\n", __func__, S1x, S1y, S1z);
939 if( S2x*S2x + S2y*S2y + S2z*S2z > 1.000001 )
940 XLALPrintWarning("XLAL Warning - %s: S2 = (%e,%e,%e) with norm > 1 requested...Are you sure you want to violate the Kerr bound?\n", __func__, S2x, S2y, S2z);
941 if( f_min < 1. )
942 XLALPrintWarning("XLAL Warning - %s: Small value of fmin = %e requested...Check for errors, this could create a very long waveform.\n", __func__, f_min);
943 if( f_min > 40.000001 )
944 XLALPrintWarning("XLAL Warning - %s: Large value of fmin = %e requested...Check for errors, the signal will start in band.\n", __func__, f_min);
945
946 /* The non-precessing waveforms return h(f) for optimal orientation
947 * (i=0, Fp=1, Fc=0; Lhat pointed toward the observer)
948 * To get generic polarizations we multiply by inclination dependence
949 * and note hc(f) \propto -I * hp(f)
950 * Non-precessing waveforms multiply hp by pfac, hc by -I*cfac
951 */
952 cfac = cos(inclination);
953 pfac = 0.5 * (1. + cfac*cfac);
954
957
958 switch (approximant)
959 {
960 /* inspiral-only models */
961 case TaylorF2:
962 /* Waveform-specific sanity checks */
964 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
966 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
967 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
968 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
969
970 /* Call the waveform driver routine */
972 XLAL_CHECK(ret == XLAL_SUCCESS, XLAL_EFUNC, "Failed to set quadparams from Universal relation.\n");
974 S1z, S2z, S1z*S1z, S2z*S2z,
975 S1z*S2z, LALpars);
976 ret = XLALSimInspiralTaylorF2Core(hptilde, frequencies, phiRef,
977 m1, m2, f_ref, 0., distance, LALpars, &pfa);
979 /* Produce both polarizations */
980 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
981 &((*hptilde)->epoch), (*hptilde)->f0, 0.0,
982 &((*hptilde)->sampleUnits), (*hptilde)->data->length);
983 for(j = 0; j < (*hptilde)->data->length; j++) {
984 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j];
985 (*hptilde)->data->data[j] *= pfac;
986 }
987 break;
988 /* inspiral-merger-ringdown models */
990 /* Waveform-specific sanity checks */
992 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
993 if (!checkAlignedSpinsEqual(S1z, S2z))
994 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
995 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
996 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
997 if( !checkTidesZero(lambda1, lambda2) )
998 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
999 if( !checkTidesZero(lambda1, lambda2) )
1000 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1001
1002 ret = XLALSimIMRSEOBNRv1ROMEffectiveSpinFrequencySequence(hptilde, hctilde, frequencies,
1003 phiRef, f_ref, distance, inclination, m1, m2, XLALSimIMRPhenomBComputeChi(m1, m2, S1z, S2z));
1004 break;
1005
1007 /* Waveform-specific sanity checks */
1009 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1010 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1011 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1012 if( !checkTidesZero(lambda1, lambda2) )
1013 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1014
1015 ret = XLALSimIMRSEOBNRv1ROMDoubleSpinFrequencySequence(hptilde, hctilde, frequencies,
1016 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z);
1017 break;
1018
1020 /* Waveform-specific sanity checks */
1022 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1023 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1024 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1025 if( !checkTidesZero(lambda1, lambda2) )
1026 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1027
1028 ret = XLALSimIMRSEOBNRv2ROMEffectiveSpinFrequencySequence(hptilde, hctilde, frequencies,
1029 phiRef, f_ref, distance, inclination, m1, m2, XLALSimIMRPhenomBComputeChi(m1, m2, S1z, S2z));
1030 break;
1031
1033 /* Waveform-specific sanity checks */
1035 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1036 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1037 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1038 if( !checkTidesZero(lambda1, lambda2) )
1039 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1040
1041 ret = XLALSimIMRSEOBNRv2ROMDoubleSpinFrequencySequence(hptilde, hctilde, frequencies,
1042 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z);
1043 break;
1044
1046 /* Waveform-specific sanity checks */
1048 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1049 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1050 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1051 if( !checkTidesZero(lambda1, lambda2) )
1052 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1053
1054 ret = XLALSimIMRSEOBNRv2ROMDoubleSpinHIFrequencySequence(hptilde, hctilde, frequencies,
1055 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, -1);
1056 break;
1057
1058 case SEOBNRv4_ROM:
1059 /* Waveform-specific sanity checks */
1061 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1062 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1063 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1064 if( !checkTidesZero(lambda1, lambda2) )
1065 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1066
1067 ret = XLALSimIMRSEOBNRv4ROMFrequencySequence(hptilde, hctilde, frequencies,
1068 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, -1, LALpars, NoNRT_V);
1069 break;
1070
1071 case SEOBNRv4HM_ROM:
1072 /* Waveform-specific sanity checks */
1074 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1075 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1076 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1077 if( !checkTidesZero(lambda1, lambda2) )
1078 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1079
1080 int nModes = 5;
1081 ret = XLALSimIMRSEOBNRv4HMROMFrequencySequence(hptilde, hctilde, frequencies,
1082 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, -1, nModes, LALpars);
1083 break;
1084
1085 case SEOBNRv5_ROM:
1086 /* Waveform-specific sanity checks */
1088 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1089 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1090 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1091 if( !checkTidesZero(lambda1, lambda2) )
1092 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1093
1094 int nModesv5 = 1;
1095 ret = XLALSimIMRSEOBNRv5HMROMFrequencySequence(hptilde, hctilde, frequencies,
1096 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, -1, nModesv5, LALpars, NoNRT_V);
1097 break;
1098
1099 case SEOBNRv5HM_ROM:
1100 /* Waveform-specific sanity checks */
1102 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1103 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1104 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1105 if( !checkTidesZero(lambda1, lambda2) )
1106 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1107
1108 int nModesv5hm = 7;
1109 ret = XLALSimIMRSEOBNRv5HMROMFrequencySequence(hptilde, hctilde, frequencies,
1110 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, -1, nModesv5hm, LALpars, NoNRT_V);
1111 break;
1112
1114 /* Waveform-specific sanity checks */
1116 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1117 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1118 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1119
1121 XLAL_CHECK(XLAL_SUCCESS == ret, ret, "Failed to set QuadMon from Lambdas for SEOBNRv4_ROM_NRTidal");
1122
1123 ret = XLALSimIMRSEOBNRv4ROMNRTidalFrequencySequence(hptilde, hctilde, frequencies,
1124 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidal_V);
1125 break;
1126
1128 /* Waveform-specific sanity checks */
1130 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1131 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1132 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1133
1135 XLAL_CHECK(XLAL_SUCCESS == ret, ret, "Failed to set QuadMon from Lambdas for SEOBNRv4_ROM_NRTidalv2");
1136
1137 ret = XLALSimIMRSEOBNRv4ROMNRTidalFrequencySequence(hptilde, hctilde, frequencies,
1138 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidalv2_V);
1139 break;
1140
1142 /* Waveform-specific sanity checks */
1144 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1145 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1146 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1147 if (m1 < m2)
1148 XLAL_ERROR(XLAL_EFUNC, "m1 = %e, m2=%e. m1 should be greater than or equal to m2 for SEOBNRv4_ROM_NRTidalv2_NSBH", m1,m2);
1149 if( lambda1 != 0 )
1150 XLAL_ERROR(XLAL_EFUNC, "lambda1 = %f. lambda1 should be zero for SEOBNRv4_ROM_NRTidalv2_NSBH", lambda1);
1151 if( lambda2 < 0 )
1152 XLAL_ERROR(XLAL_EFUNC, "lambda2 = %f. lambda2 should be nonnegative for SEOBNRv4_ROM_NRTidalv2_NSBH", lambda2);
1153 if( lambda2 > 5000 )
1154 XLAL_ERROR(XLAL_EDOM, "lambda2 = %f. lambda2 should be < 5000", lambda2);
1155 if( S2z !=0 )
1156 XLAL_PRINT_WARNING("WARNING: S2z = %f. SEOBNRv4_ROM_NRTidalv2_NSBH is calibrated to NR data for which the NS spin is zero.",S2z);
1157 if( m2 < 1 * LAL_MSUN_SI )
1158 XLAL_PRINT_WARNING("WARNING: m2=%e MSun. SEOBNRv4_ROM_NRTidalv2_NSBH is calibrated to NR data for which the NS mass is >=1 solar mass.",m2/LAL_MSUN_SI);
1159 if( m2 > 3 * LAL_MSUN_SI )
1160 XLAL_ERROR(XLAL_EDOM, "m2=%e Msun. NS Mass should be <=3 solar masses",m2/LAL_MSUN_SI);
1161 if (m1/m2 > 100)
1162 XLAL_ERROR(XLAL_EDOM, "m1/m2=%e mass ratio should be < 100",m1/m2);
1164 XLAL_CHECK(XLAL_SUCCESS == ret, ret, "Failed to set QuadMon from Lambdas for SEOBNRv4_ROM_NRTidalv2_NSBH");
1165
1166 ret = XLALSimIMRSEOBNRv4ROMNRTidalFrequencySequence(hptilde, hctilde, frequencies,
1167 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidalv2NSBH_V);
1168 break;
1169
1171 /* Waveform-specific sanity checks */
1173 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1174 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1175 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1176
1177 ret = XLALSimIMRSEOBNRv4TSurrogateFrequencySequence(hptilde, hctilde, frequencies,
1178 phiRef, f_ref, distance, inclination,
1179 m1, m2, S1z, S2z, lambda1, lambda2,
1181 break;
1182
1184 /* Waveform-specific sanity checks */
1186 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1187 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1188 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1189
1191 XLAL_CHECK(XLAL_SUCCESS == ret, ret, "Failed to set QuadMon from Lambdas for SEOBNRv4_ROM_NRTidalv2");
1192
1193 ret = XLALSimIMRSEOBNRv5ROMNRTidalFrequencySequence(hptilde, hctilde, frequencies,
1194 phiRef, f_ref, distance, inclination, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidalv3_V);
1195 break;
1196
1198 /* Waveform-specific sanity checks */
1200 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1201 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1202 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1203
1204 ret = XLALSimIMRLackeyTidal2013FrequencySequence(hptilde, hctilde, frequencies,
1205 phiRef, f_ref, distance, inclination, m1, m2, S1z, lambda2);
1206 break;
1207
1208 case IMRPhenomP:
1209 /* Waveform-specific sanity checks */
1211 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");/* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1213 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1214 /* Default is (2,2) or l=2 modes. */
1215 if( !checkTidesZero(lambda1, lambda2) )
1216 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1217 /* Tranform to model parameters */
1218 if(f_ref==0.0)
1219 f_ref = f_min; /* Default reference frequency is minimum frequency */
1221 &chi1_l, &chi2_l, &chip, &thetaJN, &alpha0, &phi_aligned, &zeta_polariz,
1222 m1, m2, f_ref, phiRef, inclination,
1223 S1x, S1y, S1z,
1224 S2x, S2y, S2z, IMRPhenomPv1_V);
1225 /* Call the waveform driver routine */
1226 ret = XLALSimIMRPhenomPFrequencySequence(hptilde, hctilde, frequencies,
1227 chi1_l, chi2_l, chip, thetaJN,
1228 m1, m2, distance, alpha0, phi_aligned, f_ref, IMRPhenomPv1_V, NoNRT_V, NULL);
1229 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1230 cos_2zeta = cos(2.0*zeta_polariz);
1231 sin_2zeta = sin(2.0*zeta_polariz);
1232 for (UINT4 idx=0;idx<(*hptilde)->data->length;idx++) {
1233 PhPpolp=(*hptilde)->data->data[idx];
1234 PhPpolc=(*hctilde)->data->data[idx];
1235 (*hptilde)->data->data[idx] = cos_2zeta*PhPpolp + sin_2zeta*PhPpolc;
1236 (*hctilde)->data->data[idx] = cos_2zeta*PhPpolc - sin_2zeta*PhPpolp;
1237 }
1238 break;
1239
1240 case IMRPhenomPv2:
1241 /* Waveform-specific sanity checks */
1243 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");/* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1245 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1246 /* Default is (2,2) or l=2 modes. */
1247 if( !checkTidesZero(lambda1, lambda2) )
1248 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1249 /* Tranform to model parameters */
1250 if(f_ref==0.0)
1251 f_ref = f_min; /* Default reference frequency is minimum frequency */
1253 &chi1_l, &chi2_l, &chip, &thetaJN, &alpha0, &phi_aligned, &zeta_polariz,
1254 m1, m2, f_ref, phiRef, inclination,
1255 S1x, S1y, S1z,
1256 S2x, S2y, S2z, IMRPhenomPv2_V);
1257 /* Call the waveform driver routine */
1258 ret = XLALSimIMRPhenomPFrequencySequence(hptilde, hctilde, frequencies,
1259 chi1_l, chi2_l, chip, thetaJN,
1260 m1, m2, distance, alpha0, phi_aligned, f_ref, IMRPhenomPv2_V, NoNRT_V, NULL);
1261 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1262 cos_2zeta = cos(2.0*zeta_polariz);
1263 sin_2zeta = sin(2.0*zeta_polariz);
1264 for (UINT4 idx=0;idx<(*hptilde)->data->length;idx++) {
1265 PhPpolp=(*hptilde)->data->data[idx];
1266 PhPpolc=(*hctilde)->data->data[idx];
1267 (*hptilde)->data->data[idx] = cos_2zeta*PhPpolp + sin_2zeta*PhPpolc;
1268 (*hctilde)->data->data[idx] = cos_2zeta*PhPpolc - sin_2zeta*PhPpolp;
1269 }
1270 break;
1271
1272 case IMRPhenomD:
1273 /* Waveform-specific sanity checks */
1275 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1276 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1277 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1278 if( !checkTidesZero(lambda1, lambda2) )
1279 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1280
1281 ret = XLALSimIMRPhenomDFrequencySequence(hptilde, frequencies,
1282 phiRef, f_ref, m1, m2, S1z, S2z, distance, LALpars, NoNRT_V);
1283 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1284 /* Produce both polarizations */
1285 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1286 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1287 &((*hptilde)->sampleUnits), (*hptilde)->data->length);
1288 for(j = 0; j < (*hptilde)->data->length; j++) {
1289 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j];
1290 (*hptilde)->data->data[j] *= pfac;
1291 }
1292 break;
1293
1294 case IMRPhenomD_NRTidal:
1295 /* Waveform-specific sanity checks */
1297 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1298 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1299 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1300
1301 ret = XLALSimIMRPhenomDNRTidalFrequencySequence(hptilde, frequencies,
1302 phiRef, f_ref, distance, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidal_V);
1303 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1304 /* Produce both polarizations */
1305 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1306 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1307 &((*hptilde)->sampleUnits), (*hptilde)->data->length);
1308 for(j = 0; j < (*hptilde)->data->length; j++) {
1309 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j];
1310 (*hptilde)->data->data[j] *= pfac;
1311 }
1312 break;
1313
1315 /* Waveform-specific sanity checks */
1317 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1318 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1319 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1320
1322 XLAL_CHECK(XLAL_SUCCESS == ret, ret, "Failed to set QuadMon from Lambdas for IMRPhenomD_NRTidalv2");
1323
1324 ret = XLALSimIMRPhenomDNRTidalFrequencySequence(hptilde, frequencies,
1325 phiRef, f_ref, distance, m1, m2, S1z, S2z, lambda1, lambda2, LALpars, NRTidalv2_V);
1326 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1327 /* Produce both polarizations */
1328 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1329 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1330 &((*hptilde)->sampleUnits), (*hptilde)->data->length);
1331 for(j = 0; j < (*hptilde)->data->length; j++) {
1332 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j];
1333 (*hptilde)->data->data[j] *= pfac;
1334 }
1335 break;
1336
1338 /* Waveform-specific sanity checks */
1340 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");/* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1342 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1343 /* Tranform to model parameters */
1344 if(f_ref==0.0)
1345 f_ref = f_min; /* Default reference frequency is minimum frequency */
1347 &chi1_l, &chi2_l, &chip, &thetaJN, &alpha0, &phi_aligned, &zeta_polariz,
1348 m1, m2, f_ref, phiRef, inclination,
1349 S1x, S1y, S1z,
1350 S2x, S2y, S2z, IMRPhenomPv2NRTidal_V);
1351 /* Call the waveform driver routine */
1352 ret = XLALSimIMRPhenomPFrequencySequence(hptilde, hctilde, frequencies,
1353 chi1_l, chi2_l, chip, thetaJN,
1354 m1, m2, distance, alpha0, phi_aligned, f_ref, IMRPhenomPv2NRTidal_V, NRTidal_V, LALpars);
1355 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1356 for (UINT4 idx=0;idx<(*hptilde)->data->length;idx++) {
1357 PhPpolp=(*hptilde)->data->data[idx];
1358 PhPpolc=(*hctilde)->data->data[idx];
1359 (*hptilde)->data->data[idx] =cos(2.*zeta_polariz)*PhPpolp+sin(2.*zeta_polariz)*PhPpolc;
1360 (*hctilde)->data->data[idx]=cos(2.*zeta_polariz)*PhPpolc-sin(2.*zeta_polariz)*PhPpolp;
1361 }
1362 break;
1363
1365 /* Waveform-specific sanity checks */
1367 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");/* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1369 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1370 /* Tranform to model parameters */
1371 if(f_ref==0.0)
1372 f_ref = f_min; /* Default reference frequency is minimum frequency */
1374 &chi1_l, &chi2_l, &chip, &thetaJN, &alpha0, &phi_aligned, &zeta_polariz,
1375 m1, m2, f_ref, phiRef, inclination,
1376 S1x, S1y, S1z,
1377 S2x, S2y, S2z, IMRPhenomPv2NRTidal_V);
1378 /* Call the waveform driver routine */
1379 ret = XLALSimIMRPhenomPFrequencySequence(hptilde, hctilde, frequencies,
1380 chi1_l, chi2_l, chip, thetaJN,
1381 m1, m2, distance, alpha0, phi_aligned, f_ref, IMRPhenomPv2NRTidal_V, NRTidalv2_V, LALpars);
1382 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1383 for (UINT4 idx=0;idx<(*hptilde)->data->length;idx++) {
1384 PhPpolp=(*hptilde)->data->data[idx];
1385 PhPpolc=(*hctilde)->data->data[idx];
1386 (*hptilde)->data->data[idx] =cos(2.*zeta_polariz)*PhPpolp+sin(2.*zeta_polariz)*PhPpolc;
1387 (*hctilde)->data->data[idx]=cos(2.*zeta_polariz)*PhPpolc-sin(2.*zeta_polariz)*PhPpolp;
1388 }
1389 break;
1390
1391 case IMRPhenomHM:
1392 if (!checkTransverseSpinsZero(S1x, S1y, S2x, S2y))
1393 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1394 if (!checkTidesZero(lambda1, lambda2))
1395 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1396 ret = XLALSimIMRPhenomHM(hptilde, hctilde, frequencies, m1, m2,
1397 S1z, S2z, distance, inclination, phiRef, 0., f_ref,
1398 LALpars);
1399 if (ret == XLAL_FAILURE)
1401 break;
1402
1403 case IMRPhenomXAS:
1404 {
1405 /* Waveform-specific sanity checks */
1407 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1408 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1409 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1410 if( !checkTidesZero(lambda1, lambda2) )
1411 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1412
1413 /*
1414 This is the factor that comes from Y_22star + (-1)^l * Y_2-2 without the dependence in inclination, that is included in pfac and cfac
1415 We add the azimuthal part exp^{i*m*beta} of the spherical harmonics Ylm(inclination, beta),
1416 with beta = PI/2 - phiRef, phiRef is included in the individual mode
1417 */
1418 COMPLEX16 Ylmfactor = 2.0*sqrt(5.0 / (64.0 * LAL_PI)) * cexp(-I*2*(LAL_PI/2 ));
1419 /* The factor for hc is the same but opposite sign */
1420
1421 ret = XLALSimIMRPhenomXASFrequencySequence(hptilde, frequencies,
1422 m1, m2, S1z, S2z, distance, phiRef, f_ref, LALpars);
1423 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1424
1425 /* Produce both polarizations */
1426 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1427 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1428 &((*hptilde)->sampleUnits), (*hptilde)->data->length
1429 );
1430 for(j = 0; j < (*hptilde)->data->length; j++)
1431 {
1432 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j] * Ylmfactor;
1433 (*hptilde)->data->data[j] *= pfac * Ylmfactor;
1434 }
1435 }
1436 break;
1437
1439 {
1440 /* Waveform-specific sanity checks */
1442 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1443 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1444 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1445
1446 if(LALpars)
1448
1449 /*
1450 This is the factor that comes from Y_22star + (-1)^l * Y_2-2 without the dependence in inclination, that is included in pfac and cfac
1451 We add the azimuthal part exp^{i*m*beta} of the spherical harmonics Ylm(inclination, beta),
1452 with beta = PI/2 - phiRef, phiRef is included in the individual mode
1453 */
1454 COMPLEX16 Ylmfactor = 2.0*sqrt(5.0 / (64.0 * LAL_PI)) * cexp(-I*2*(LAL_PI/2 ));
1455 /* The factor for hc is the same but opposite sign */
1456
1457 ret = XLALSimIMRPhenomXASFrequencySequence(hptilde, frequencies,
1458 m1, m2, S1z, S2z, distance, phiRef, f_ref, LALpars);
1459 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1460
1461 if(LALpars)
1463
1464 /* Produce both polarizations */
1465 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1466 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1467 &((*hptilde)->sampleUnits), (*hptilde)->data->length
1468 );
1469 for(j = 0; j < (*hptilde)->data->length; j++)
1470 {
1471 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j] * Ylmfactor;
1472 (*hptilde)->data->data[j] *= pfac * Ylmfactor;
1473 }
1474
1475 }
1476 break;
1477
1479 {
1480 /* Waveform-specific sanity checks */
1482 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1483 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1484 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1485
1486 if(LALpars)
1488
1489 /*
1490 This is the factor that comes from Y_22star + (-1)^l * Y_2-2 without the dependence in inclination, that is included in pfac and cfac
1491 We add the azimuthal part exp^{i*m*beta} of the spherical harmonics Ylm(inclination, beta),
1492 with beta = PI/2 - phiRef, phiRef is included in the individual mode
1493 */
1494 COMPLEX16 Ylmfactor = 2.0*sqrt(5.0 / (64.0 * LAL_PI)) * cexp(-I*2*(LAL_PI/2 ));
1495 /* The factor for hc is the same but opposite sign */
1496
1497 ret = XLALSimIMRPhenomXASFrequencySequence(hptilde, frequencies,
1498 m1, m2, S1z, S2z, distance, phiRef, f_ref, LALpars);
1499 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1500
1501 if(LALpars)
1503
1504 /* Produce both polarizations */
1505 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1506 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1507 &((*hptilde)->sampleUnits), (*hptilde)->data->length
1508 );
1509 for(j = 0; j < (*hptilde)->data->length; j++)
1510 {
1511 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j] * Ylmfactor;
1512 (*hptilde)->data->data[j] *= pfac * Ylmfactor;
1513 }
1514
1515 }
1516 break;
1517
1518 case IMRPhenomXHM:
1519 /* Waveform-specific sanity checks */
1521 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1522 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1523 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1524 if( !checkTidesZero(lambda1, lambda2) )
1525 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1526
1527
1528 ret = XLALSimIMRPhenomXHMFrequencySequence(hptilde, hctilde, frequencies,
1529 m1, m2, S1z, S2z, distance, inclination, phiRef, f_ref, LALpars);
1530
1531 if (ret == XLAL_FAILURE)
1533 break;
1534
1535 case IMRPhenomXP:
1536 /* Waveform-specific sanity checks */
1538 {
1539 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1540 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1541 }
1543 {
1544 /* Default is (2,2) or l=2 modes. */
1545 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1546 }
1547 if( !checkTidesZero(lambda1, lambda2) )
1548 {
1549 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1550 }
1551 if(f_ref==0.0)
1552 {
1553 /* Default reference frequency is minimum frequency */
1554 f_ref = f_min;
1555 }
1556
1557 /* Call the main waveform driver. Note that we pass the full spin vectors
1558 with XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame being
1559 effectively called in the initialization of the pPrec struct
1560 */
1562 hptilde, hctilde,
1563 frequencies,
1564 m1, m2,
1565 S1x, S1y, S1z,
1566 S2x, S2y, S2z,
1567 distance, inclination,
1568 phiRef, f_ref, LALpars
1569 );
1570 if (ret == XLAL_FAILURE)
1571 {
1573 }
1574 break;
1575
1577 /* Waveform-specific sanity checks */
1579 {
1580 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1581 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1582 }
1584 {
1585 /* Default is (2,2) or l=2 modes. */
1586 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1587 }
1588
1589 if(LALpars)
1591
1592 if(f_ref==0.0)
1593 {
1594 /* Default reference frequency is minimum frequency */
1595 f_ref = f_min;
1596 }
1597
1598 /* Call the main waveform driver. Note that we pass the full spin vectors
1599 with XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame being
1600 effectively called in the initialization of the pPrec struct
1601 */
1603 hptilde, hctilde,
1604 frequencies,
1605 m1, m2,
1606 S1x, S1y, S1z,
1607 S2x, S2y, S2z,
1608 distance, inclination,
1609 phiRef, f_ref, LALpars
1610 );
1611
1612 if(LALpars)
1614
1615 if (ret == XLAL_FAILURE)
1616 {
1618 }
1619 break;
1620
1622 /* Waveform-specific sanity checks */
1624 {
1625 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1626 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1627 }
1629 {
1630 /* Default is (2,2) or l=2 modes. */
1631 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1632 }
1633
1634 if(LALpars)
1636
1637 if(f_ref==0.0)
1638 {
1639 /* Default reference frequency is minimum frequency */
1640 f_ref = f_min;
1641 }
1642
1643 /* Call the main waveform driver. Note that we pass the full spin vectors
1644 with XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame being
1645 effectively called in the initialization of the pPrec struct
1646 */
1648 hptilde, hctilde,
1649 frequencies,
1650 m1, m2,
1651 S1x, S1y, S1z,
1652 S2x, S2y, S2z,
1653 distance, inclination,
1654 phiRef, f_ref, LALpars
1655 );
1656
1657 if(LALpars)
1659
1660 if (ret == XLAL_FAILURE)
1661 {
1663 }
1664 break;
1665
1666 case IMRPhenomXPHM:
1667 /* Waveform-specific sanity checks */
1669 {
1670 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1671 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1672 }
1674 {
1675 /* Default is (2,2) or l=2 modes. */
1676 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1677 }
1678 if( !checkTidesZero(lambda1, lambda2) )
1679 {
1680 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1681 }
1682 // if(f_ref==0.0)
1683 // {
1684 // /* Default reference frequency is minimum frequency */
1685 // f_ref = f_min;
1686 // }
1687
1689 hptilde, hctilde,
1690 frequencies,
1691 m1, m2,
1692 S1x, S1y, S1z,
1693 S2x, S2y, S2z,
1694 distance, inclination,
1695 phiRef, f_ref, LALpars
1696 );
1697
1698 if (ret == XLAL_FAILURE)
1699 {
1701 }
1702 break;
1703
1704 case IMRPhenomXPNR:
1705
1706 if (LALpars == NULL){
1707 LALparams_aux = XLALCreateDict();
1708 }
1709 else{
1710 LALparams_aux = XLALDictDuplicate(LALpars);
1711 }
1712 /* Waveform-specific sanity checks */
1714 {
1715 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1716 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1717 }
1719 {
1720 /* Default is (2,2) or l=2 modes. */
1721 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1722 }
1723 if( !checkTidesZero(lambda1, lambda2) )
1724 {
1725 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1726 }
1727 // if(f_ref==0.0)
1728 // {
1729 // /* Default reference frequency is minimum frequency */
1730 // f_ref = f_min;
1731 // }
1732
1733
1734 /* XO4 uses previous version of XHM */
1736
1737
1738 /* Toggle on PNR angles */
1739 if(!XLALDictContains(LALparams_aux, "PNRUseTunedAngles")){
1741 }
1742
1743 /* Toggle on tuned coprecessing strain */
1744 if(!XLALDictContains(LALparams_aux, "PNRUseTunedCoprec")){
1746 }
1747 if(!XLALDictContains(LALparams_aux, "PNRForceXHMAlignment")){
1749 }
1750
1751 /* Toggle on antisymmetric contributions */
1752 if(!XLALDictContains(LALparams_aux, "AntisymmetricWaveform")){
1754 }
1755
1756 /* Toggle on reviewed PrecVersion */
1757 if(!XLALDictContains(LALparams_aux, "PrecVersion")){
1759 }
1760
1761 if(!XLALDictContains(LALparams_aux, "FinalSpinMod")){
1763 }
1764
1765 /* Toggle the conditional precession multibanding (for consistency, even though this is turned off anyway)*/
1766 if(!XLALDictContains(LALparams_aux, "MBandPrecVersion")){
1768 }
1769
1771 hptilde, hctilde,
1772 frequencies,
1773 m1, m2,
1774 S1x, S1y, S1z,
1775 S2x, S2y, S2z,
1776 distance, inclination,
1777 phiRef, f_ref, LALparams_aux
1778 );
1779
1780 XLALDestroyDict(LALparams_aux);
1781
1782 if (ret == XLAL_FAILURE)
1783 {
1785 }
1786
1787 break;
1788
1789case IMRPhenomXO4a:
1790
1791 if (LALpars == NULL){
1792 LALparams_aux = XLALCreateDict();
1793 }
1794 else{
1795 LALparams_aux = XLALDictDuplicate(LALpars);
1796 }
1797 /* Waveform-specific sanity checks */
1799 {
1800 /* Default is LAL_SIM_INSPIRAL_FRAME_AXIS_ORBITAL_L : z-axis along direction of orbital angular momentum. */
1801 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralFrameAxis provided, but this approximant does not use that flag.");
1802 }
1804 {
1805 /* Default is (2,2) or l=2 modes. */
1806 XLAL_ERROR(XLAL_EINVAL, "Non-default LALSimInspiralModesChoice provided, but this approximant does not use that flag.");
1807 }
1808 if( !checkTidesZero(lambda1, lambda2) )
1809 {
1810 XLAL_ERROR(XLAL_EINVAL, "Non-zero tidal parameters were given, but this is approximant doe not have tidal corrections.");
1811 }
1812 // if(f_ref==0.0)
1813 // {
1814 // /* Default reference frequency is minimum frequency */
1815 // f_ref = f_min;
1816 // }
1817
1818
1819 /* XO4 uses previous version of XHM */
1821
1822
1823 /* Toggle on PNR angles */
1824 if(!XLALDictContains(LALparams_aux, "PNRUseTunedAngles")){
1826 }
1827
1828 /* Toggle on tuned coprecessing strain */
1829 if(!XLALDictContains(LALparams_aux, "PNRUseTunedCoprec")){
1831 }
1832 if(!XLALDictContains(LALparams_aux, "PNRForceXHMAlignment")){
1834 }
1835
1836 /* Toggle on antisymmetric contributions */
1837 if(!XLALDictContains(LALparams_aux, "AntisymmetricWaveform")){
1839 }
1840
1841 /* Toggle on reviewed PrecVersion */
1842 if(!XLALDictContains(LALparams_aux, "PrecVersion")){
1844 }
1845
1847 hptilde, hctilde,
1848 frequencies,
1849 m1, m2,
1850 S1x, S1y, S1z,
1851 S2x, S2y, S2z,
1852 distance, inclination,
1853 phiRef, f_ref, LALparams_aux
1854 );
1855
1856 XLALDestroyDict(LALparams_aux);
1857
1858 if (ret == XLAL_FAILURE)
1859 {
1861 }
1862
1863 break;
1864
1865 case IMRPhenomNSBH:
1866 /* Waveform-specific sanity checks */
1868 XLAL_ERROR(XLAL_EINVAL, "Non-default flags given, but this approximant does not support this case.");
1869 if( !checkTransverseSpinsZero(S1x, S1y, S2x, S2y) )
1870 XLAL_ERROR(XLAL_EINVAL, "Non-zero transverse spins were given, but this is a non-precessing approximant.");
1871 if(f_ref==0.0)
1872 f_ref = f_min;
1873
1874 ret = XLALSimIMRPhenomNSBHFrequencySequence(hptilde, frequencies,
1875 phiRef, f_ref, distance, m1, m2, S1z, S2z, LALpars);
1876 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1877 /* Produce both polarizations */
1878 *hctilde = XLALCreateCOMPLEX16FrequencySeries("FD hcross",
1879 &((*hptilde)->epoch), (*hptilde)->f0, (*hptilde)->deltaF,
1880 &((*hptilde)->sampleUnits), (*hptilde)->data->length);
1881 for(j = 0; j < (*hptilde)->data->length; j++) {
1882 (*hctilde)->data->data[j] = -I*cfac * (*hptilde)->data->data[j];
1883 (*hptilde)->data->data[j] *= pfac;
1884 }
1885 break;
1886
1887 default:
1888 XLALPrintError("FD version of approximant not implemented in lalsimulation\n");
1890 }
1891
1892 if (ret == XLAL_FAILURE) XLAL_ERROR(XLAL_EFUNC);
1893
1894 return ret;
1895}
int XLALDictContains(const LALDict *dict, const char *key)
void XLALDestroyDict(LALDict *dict)
LALDict * XLALDictDuplicate(const LALDict *orig)
LALDict * XLALCreateDict(void)
int XLALSimIMRPhenomHM(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1z, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, const REAL8 deltaF, REAL8 f_ref, LALDict *extraParams)
int XLALSimIMRSEOBNRv4HMROMFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2, INT4 nk_max, UINT4 nModes, LALDict *LALParams)
int XLALSimIMRSEOBNRv5HMROMFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2, INT4 nk_max, UINT4 nModes, LALDict *LALParams, NRTidal_version_type NRTidal_version)
#define nModes
int XLALSimInspiralApproximantAcceptTestGRParams(Approximant approx)
int XLALSimInspiralSetQuadMonParamsFromLambdas(LALDict *LALparams)
if you do NOT provide a quadparam[1,2] term and you DO provide lamdba[1,2] then we calculate quad-mon...
static void UNUSED XLALSimInspiralPNPhasing_F2(PNPhasingSeries *pfa, const REAL8 m1, const REAL8 m2, const REAL8 chi1L, const REAL8 chi2L, const REAL8 chi1sq, const REAL8 chi2sq, const REAL8 chi1dotchi2, LALDict *p)
CacheVariableDiffersBitmask
Bitmask enumerating which parameters have changed, to determine if the requested waveform can be tran...
int XLALSimInspiralChooseFDWaveformSequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 phiRef, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_ref, REAL8 distance, REAL8 inclination, LALDict *LALpars, Approximant approximant, REAL8Sequence *frequencies)
Wrapper similar to XLALSimInspiralChooseFDWaveform() for waveforms to be generated a specific freqenc...
static CacheVariableDiffersBitmask CacheArgsDifferenceBitmask(LALSimInspiralWaveformCache *cache, REAL8 phiRef, REAL8 deltaTF, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_min, REAL8 f_ref, REAL8 f_max, REAL8 r, REAL8 i, LALDict *LALpars, Approximant approximant, REAL8Sequence *frequencies)
Function to compare the requested arguments to those stored in the cache, returns a bitmask which det...
static int FrequenciesAreDifferent(REAL8Sequence *newFrequencies, REAL8Sequence *cachedFrequencies)
Function to compare two frequencies sequences.
static int StoreFDHCache(LALSimInspiralWaveformCache *cache, COMPLEX16FrequencySeries *hptilde, COMPLEX16FrequencySeries *hctilde, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_min, REAL8 f_ref, REAL8 f_max, REAL8 r, REAL8 i, LALDict *LALpars, Approximant approximant, REAL8Sequence *frequencies)
Store the output FD hptilde and hctilde in cache.
static int StoreTDHCache(LALSimInspiralWaveformCache *cache, REAL8TimeSeries *hplus, REAL8TimeSeries *hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_min, REAL8 f_ref, REAL8 r, REAL8 i, LALDict *LALpars, Approximant approximant)
Store the output TD hplus and hcross in the cache.
int XLALSimInspiralWaveformParamsDTau440IsDefault(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupTidalHexadecapolarLambda1(LALDict *params)
int XLALSimInspiralWaveformParamsDOmega550IsDefault(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupTidalHexadecapolarLambda2(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXPNRUseTunedAngles(LALDict *params, INT4 value)
REAL8 XLALSimInspiralWaveformParamsLookupdQuadMon1(LALDict *params)
int XLALSimInspiralWaveformParamsDOmega220IsDefault(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupTidalLambda2(LALDict *params)
int XLALSimInspiralWaveformParamsDTau330IsDefault(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXAntisymmetricWaveform(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsInsertPhenomXPNRForceXHMAlignment(LALDict *params, INT4 value)
REAL8 XLALSimInspiralWaveformParamsLookupTidalLambda1(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXPFinalSpinMod(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsInsertPhenomXHMReleaseVersion(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsInsertPhenomXPNRUseTunedCoprec(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsDOmega210IsDefault(LALDict *params)
int XLALSimInspiralWaveformParamsDTau220IsDefault(LALDict *params)
INT4 XLALSimInspiralWaveformParamsLookupPNAmplitudeOrder(LALDict *params)
int XLALSimInspiralWaveformParamsModesChoiceIsDefault(LALDict *params)
int XLALSimInspiralWaveformParamsDOmega330IsDefault(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupdQuadMon2(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXTidalFlag(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsDOmega440IsDefault(LALDict *params)
int XLALSimInspiralWaveformParamsInsertPhenomXPHMMBandVersion(LALDict *params, INT4 value)
int XLALSimInspiralWaveformParamsInsertPhenomXPrecVersion(LALDict *params, INT4 value)
REAL8 XLALSimInspiralWaveformParamsLookupTidalOctupolarLambda1(LALDict *params)
int XLALSimInspiralWaveformParamsFrameAxisIsDefault(LALDict *params)
int XLALSimInspiralWaveformParamsDTau550IsDefault(LALDict *params)
REAL8 XLALSimInspiralWaveformParamsLookupTidalOctupolarLambda2(LALDict *params)
int XLALSimInspiralWaveformParamsDTau210IsDefault(LALDict *params)
static REAL8 pfac(int n)
double i
Definition: bh_ringdown.c:118
#define checkTransverseSpinsZero(s1x, s1y, s2x, s2y)
#define checkAlignedSpinsEqual(s1z, s2z)
#define checkTidesZero(lambda1, lambda2)
COMPLEX16FrequencySeries * XLALCreateCOMPLEX16FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
COMPLEX16FrequencySeries * XLALCutCOMPLEX16FrequencySeries(const COMPLEX16FrequencySeries *series, size_t first, size_t length)
void XLALDestroyCOMPLEX16FrequencySeries(COMPLEX16FrequencySeries *series)
#define LAL_MSUN_SI
#define LAL_PI
#define cpolar(r, th)
double complex COMPLEX16
double REAL8
uint32_t UINT4
int32_t INT4
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
@ IMRPhenomPv1_V
version 1: based on IMRPhenomC
Definition: LALSimIMR.h:74
@ IMRPhenomPv2_V
version 2: based on IMRPhenomD
Definition: LALSimIMR.h:75
@ IMRPhenomPv2NRTidal_V
version Pv2_NRTidal: based on IMRPhenomPv2; NRTides added before precession; can be used with both NR...
Definition: LALSimIMR.h:76
@ SEOBNRv4TSurrogate_CUBIC
use cubic splines in frequency
Definition: LALSimIMR.h:91
@ NRTidalv3_V
version NRTidalv3
Definition: LALSimIMR.h:83
@ NRTidal_V
version NRTidal: based on https://arxiv.org/pdf/1706.02969.pdf
Definition: LALSimIMR.h:81
@ NoNRT_V
special case for PhenomPv2 BBH baseline
Definition: LALSimIMR.h:87
@ NRTidalv2NSBH_V
version NRTidalv2: https://arxiv.org/abs/1905.06011 with amplitude corrections for NSBH (used for SEO...
Definition: LALSimIMR.h:86
@ NRTidalv2_V
version NRTidalv2: https://arxiv.org/abs/1905.06011
Definition: LALSimIMR.h:82
int XLALSimIMRPhenomPCalculateModelParametersFromSourceFrame(REAL8 *chi1_l, REAL8 *chi2_l, REAL8 *chip, REAL8 *thetaJN, REAL8 *alpha0, REAL8 *phi_aligned, REAL8 *zeta_polariz, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 f_ref, const REAL8 phiRef, const REAL8 incl, const REAL8 s1x, const REAL8 s1y, const REAL8 s1z, const REAL8 s2x, const REAL8 s2y, const REAL8 s2z, IMRPhenomP_version_type IMRPhenomP_version)
Function to map LAL parameters (masses, 6 spin components, phiRef and inclination at f_ref) (assumed ...
int XLALSimIMRPhenomNSBHFrequencySequence(COMPLEX16FrequencySeries **htilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 mBH_SI, REAL8 mNS_SI, REAL8 chi_BH, REAL8 chi_NS, LALDict *extraParams)
Compute waveform in LAL format at specified frequencies for the IMRPhenomNSBH model.
int XLALSimIMRPhenomDFrequencySequence(COMPLEX16FrequencySeries **htilde, const REAL8Sequence *freqs, const REAL8 phi0, const REAL8 fRef_in, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 chi1, const REAL8 chi2, const REAL8 distance, LALDict *extraParams, NRTidal_version_type NRTidal_version)
Compute waveform in LAL format at specified frequencies for the IMRPhenomD model.
double XLALSimIMRPhenomBComputeChi(const REAL8 m1, const REAL8 m2, const REAL8 s1z, const REAL8 s2z)
Compute the dimensionless, spin-aligned parameter chi as used in the IMRPhenomB waveform.
int XLALSimIMRPhenomPFrequencySequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, const REAL8 chi1_l, const REAL8 chi2_l, const REAL8 chip, const REAL8 thetaJ, REAL8 m1_SI, const REAL8 m2_SI, const REAL8 distance, const REAL8 alpha0, const REAL8 phic, const REAL8 f_ref, IMRPhenomP_version_type IMRPhenomP_version, NRTidal_version_type NRTidal_version, LALDict *extraParams)
Driver routine to compute the precessing inspiral-merger-ringdown phenomenological waveform IMRPhenom...
int XLALSimIMRPhenomXHMFrequencySequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1z, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef, LALDict *lalParams)
Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies spe...
int XLALSimIMRPhenomXASFrequencySequence(COMPLEX16FrequencySeries **htilde22, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1L, REAL8 chi2L, REAL8 distance, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
Compute waveform in LAL format at specified frequencies for the IMRPhenomX model.
int XLALSimIMRPhenomXPHMFrequencySequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
Returns hptilde and hctilde as a complex frequency series with entries exactly at the frequencies spe...
int XLALSimIMRPhenomXPFrequencySequence(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1x, REAL8 chi1y, REAL8 chi1z, REAL8 chi2x, REAL8 chi2y, REAL8 chi2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, REAL8 fRef_In, LALDict *lalParams)
Compute waveform in LAL format at specified frequencies for the IMRPhenomXP model.
int XLALSimIMRSEOBNRv1ROMDoubleSpinFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2)
Compute waveform in LAL format at specified frequencies for the SEOBNRv1_ROM_DoubleSpin model.
int XLALSimIMRSEOBNRv4ROMFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2, INT4 nk_max, LALDict *LALparams, NRTidal_version_type NRTidal_version)
Compute waveform in LAL format at specified frequencies for the SEOBNRv4_ROM model.
int XLALSimIMRSEOBNRv2ROMEffectiveSpinFrequencySequence(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 SEOBNRv2_ROM_EffectiveSpin model.
int XLALSimIMRSEOBNRv2ROMDoubleSpinFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2)
Compute waveform in LAL format at specified frequencies for the SEOBNRv2_ROM_DoubleSpin 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.
int XLALSimIMRSEOBNRv2ROMDoubleSpinHIFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2, INT4 nk_max)
Compute waveform in LAL format at specified frequencies for the SEOBNRv2_ROM_DoubleSpin_HI model.
int XLALSimIMRSEOBNRv4TSurrogateFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1SI, REAL8 m2SI, REAL8 chi1, REAL8 chi2, REAL8 lambda1, REAL8 lambda2, SEOBNRv4TSurrogate_spline_order spline_order)
Compute waveform in LAL format at specified frequencies for the SEOBNRv4T_surrogate model.
int XLALSimIMRSEOBNRv4ROMNRTidalFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1, REAL8 chi2, REAL8 Lambda1, REAL8 Lambda2, LALDict *LALparams, NRTidal_version_type NRTidal_version)
Compute waveform in LAL format at specified frequencies for the SEOBNRv4_ROM_NRTidal tidal model base...
int XLALSimIMRPhenomDNRTidalFrequencySequence(COMPLEX16FrequencySeries **htilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1, REAL8 chi2, REAL8 lambda1, REAL8 lambda2, LALDict *extraParams, NRTidal_version_type NRTidal_version)
Compute waveform in LAL format at specified frequencies for the IMRPhenomD_NRTidal tidal model based ...
int XLALSimIMRSEOBNRv5ROMNRTidalFrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 m1_SI, REAL8 m2_SI, REAL8 chi1, REAL8 chi2, REAL8 Lambda1, REAL8 Lambda2, LALDict *LALparams, NRTidal_version_type NRTidal_version)
Compute waveform in LAL format at specified frequencies for the SEOBNRv5_ROM_NRTidal tidal model base...
int XLALSimIMRLackeyTidal2013FrequencySequence(struct tagCOMPLEX16FrequencySeries **hptilde, struct tagCOMPLEX16FrequencySeries **hctilde, const REAL8Sequence *freqs, REAL8 phiRef, REAL8 fRef, REAL8 distance, REAL8 inclination, REAL8 mBH_SI, REAL8 mNS_SI, REAL8 chi_BH, REAL8 Lambda)
Compute waveform in LAL format at specified frequencies for the Lackey et al (2013) tidal model based...
int XLALSimInspiralChooseFDWaveform(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, const REAL8 m1, const REAL8 m2, const REAL8 S1x, const REAL8 S1y, const REAL8 S1z, const REAL8 S2x, const REAL8 S2y, const REAL8 S2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, const REAL8 longAscNodes, const REAL8 eccentricity, const REAL8 UNUSED meanPerAno, const REAL8 deltaF, const REAL8 f_min, const REAL8 f_max, REAL8 f_ref, LALDict *params, const Approximant approximant)
Chooses between different approximants when requesting a waveform to be generated For spinning wavefo...
int XLALSimInspiralChooseTDWaveform(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, const REAL8 m1, const REAL8 m2, const REAL8 S1x, const REAL8 S1y, const REAL8 S1z, const REAL8 S2x, const REAL8 S2y, const REAL8 S2z, const REAL8 distance, const REAL8 inclination, const REAL8 phiRef, const REAL8 longAscNodes, const REAL8 eccentricity, const REAL8 UNUSED meanPerAno, const REAL8 deltaT, const REAL8 f_min, REAL8 f_ref, LALDict *params, const Approximant approximant)
Chooses between different approximants when requesting a waveform to be generated For spinning wavefo...
Approximant
Enum that specifies the PN approximant to be used in computing the waveform.
@ EOBNRv2HM
UNDOCUMENTED.
@ SpinTaylorT4
Spinning case T4 models (lalsimulation's equivalent of SpinTaylorFrameless).
@ TaylorF2RedSpinTidal
TaylorF2 waveforms for non-precessing spins, defined in terms of a single (reduced-spin) parameter [A...
@ SEOBNRv2_ROM_DoubleSpin_HI
High resolution low-mass double-spin frequency domain reduced order model of spin-aligned EOBNR model...
@ SEOBNRv4_ROM_NRTidal
Low-mass double-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv4 [Bohe ...
@ IMRPhenomXAS_NRTidalv2
Spin non-precessing EOBNR model v4 with higher modes post-adiabatic dynamics (time domain) and TGR ri...
@ IMRPhenomP
Frequency domain (generic spins) inspiral-merger-ringdown templates of Hannam et al....
@ SEOBNRv4_ROM_NRTidalv2
based on NRTidalv2; https://arxiv.org/abs/1905.06011.
@ IMRPhenomXP
Frequency domain, precessing phenomenological IMR waveform model.
@ IMRPhenomC
Frequency domain (non-precessing spins) inspiral-merger-ringdown templates of Santamaria et al [Santa...
@ SEOBNRv4HM_ROM
Low-mass double-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv4hm.
@ IMRPhenomPv2_NRTidal
Frequency domain tidal version of IMRPhenomPv2, using NRTidal framework from arXiv:1706....
@ IMRPhenomXP_NRTidalv3
Tidal extension of IMRPhenomXP based on NRTidalv3.
@ SEOBNRv1
Spin-aligned EOBNR model.
@ Lackey_Tidal_2013_SEOBNRv2_ROM
Frequency domain tidal model based on reduced order model of SEOBNRv2.
@ EOBNRv2
UNDOCUMENTED.
@ IMRPhenomXO4a
Frequency domain, precessing with subdominant modes phenomenological IMR waveform model with NR-tuned...
@ IMRPhenomNSBH
NSBH Tidal model.
@ IMRPhenomXPHM
Frequency domain, precessing with subdominant modes phenomenological IMR waveform model.
@ IMRPhenomD
Frequency domain (non-precessing spins) inspiral-merger-ringdown templates of Husa et al,...
@ IMRPhenomXHM
Frequency domain, non-precessing phenomenological IMR waveform model with subdominant modes ([arXiv:2...
@ SEOBNRv2_ROM_EffectiveSpin
Single-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv2.
@ IMRPhenomD_NRTidal
Uses arxiv:1706.02969 to upgrad IMRPhenomD to a tidal approximant.
@ IMRPhenomXAS_NRTidalv3
Tidal extension of IMRPhenomXAS based on NRTidalv3.
@ TaylorF2RedSpin
TaylorF2 waveforms for non-precessing spins, defined in terms of a single (reduced-spin) parameter [A...
@ IMRPhenomXPNR
Frequency domain, precessing with subdominant modes phenomenological IMR waveform model with SpinTayl...
@ SEOBNRv1_ROM_EffectiveSpin
Single-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv1 See [Purrer:201...
@ IMRPhenomD_NRTidalv2
NRTidalv2; https://arxiv.org/abs/1905.06011.
@ IMRPhenomA
Time domain (non-spinning) inspiral-merger-ringdown waveforms generated from the inverse FFT of IMRPh...
@ IMRPhenomHM
Frequency domain with higher modes (non-precessing spins) inspiral-merger-ringdown templates,...
@ IMRPhenomPv2
Frequency domain (generic spins) inspiral-merger-ringdown templates of Hannam et al....
@ IMRPhenomPv2_NRTidalv2
Frequency domain tidal version; based on https://arxiv.org/abs/1905.06011.
@ TaylorT3
Time domain Taylor approximant in which phase is explicitly given as a function of time; outputs a ti...
@ IMRPhenomXAS
Frequency domain, non-precessing phenomenological IMR waveform model ([arXiv:2001....
@ SpinTaylorT5
Spinning case T5 models, which is a variant of the spinning version of the original TaylorT2 (see ) d...
@ SEOBNRv5_ROM_NRTidalv3
based on NRTidalv3 (arXiv:2311.07456);
@ TaylorT4
UNDOCUMENTED.
@ SEOBNRv4_ROM_NRTidalv2_NSBH
NSBH model based on SEOBNRv4_ROM_NRTidalv2.
@ SEOBNRv1_ROM_DoubleSpin
Double-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv1 See [Purrer:201...
@ TaylorF2
The standard stationary phase approximation; Outputs a frequency-domain wave.
@ SEOBNRv4_ROM
Low-mass double-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv4.
@ SEOBNRv2_ROM_DoubleSpin
Double-spin frequency domain reduced order model of spin-aligned EOBNR model SEOBNRv2.
@ IMRPhenomB
Time domain (non-precessing spins) inspiral-merger-ringdown waveforms generated from the inverse FFT ...
@ TaylorT1
Time domain Taylor approximant in which the energy and flux are both kept as Taylor expansions and a ...
@ SEOBNRv5_ROM
Time domain, precessing phenomenological IMR waveform model with subdominant modes ([arXiv: 20XY....
@ IMRPhenomXP_NRTidalv2
Tidal extension of IMRPhenomXP based on [arXiv:1905.06011].
@ TEOBResumS
Resummed Spin-aligned Tidal EOB.
@ SEOBNRv5HM_ROM
External Python model.
@ TaylorT2
Time domain Taylor approximant in which the phase evolution is obtained by iteratively solving post-...
@ SEOBNRv4T_surrogate
Double-spin frequency domain surrogate model of spin-aligned tidal EOBNR model SEOBNRv4T.
@ LAL_SIM_INSPIRAL_TESTGR_PARAMS
These approximants cannot accept testGR params as input params.
int XLALSimInspiralTaylorF2Core(COMPLEX16FrequencySeries **htilde, const REAL8Sequence *freqs, const REAL8 phi_ref, const REAL8 m1_SI, const REAL8 m2_SI, const REAL8 f_ref, const REAL8 shft, const REAL8 r, LALDict *LALparams, PNPhasingSeries *pfaP)
int XLALSimInspiralWaveformParamsNonGRAreDefault(LALDict *params)
LALSimInspiralWaveformCache * XLALCreateSimInspiralWaveformCache(void)
Construct and initialize a waveform cache.
int XLALSimInspiralChooseTDWaveformFromCache(REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 phiRef, REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_min, REAL8 f_ref, REAL8 r, REAL8 i, LALDict *LALpars, Approximant approximant, LALSimInspiralWaveformCache *cache)
Chooses between different approximants when requesting a waveform to be generated Returns the wavefor...
void XLALDestroySimInspiralWaveformCache(LALSimInspiralWaveformCache *cache)
Destroy a waveform cache.
int XLALSimInspiralChooseFDWaveformFromCache(COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 phiRef, REAL8 deltaF, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 f_min, REAL8 f_max, REAL8 f_ref, REAL8 r, REAL8 i, LALDict *LALpars, Approximant approximant, LALSimInspiralWaveformCache *cache, REAL8Sequence *frequencies)
Chooses between different approximants when requesting a waveform to be generated Returns the wavefor...
bool XLALSimInspiralWaveformParamsFlagsAreDefault(LALDict *params)
Returns true if waveFlags is non-NULL and all of its fields have default value; returns false otherwi...
bool XLALSimInspiralWaveformFlagsEqual(LALDict *LALpars1, LALDict *LALpars2)
Checks if all flags in two LALSimInspiralWaveformFlags structs are equal.
static const INT4 r
void XLALDestroyREAL8Sequence(REAL8Sequence *sequence)
REAL8Sequence * XLALCopyREAL8Sequence(REAL8Sequence *sequence)
void XLALDestroyREAL8TimeSeries(REAL8TimeSeries *series)
REAL8TimeSeries * XLALCutREAL8TimeSeries(const REAL8TimeSeries *series, size_t first, size_t length)
REAL8TimeSeries * XLALCreateREAL8TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
#define XLAL_PRINT_WARNING(...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
int int XLALPrintWarning(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
XLAL_EINVAL
XLAL_FAILURE
string status
string approximant
CHAR name[LALNameLength]
COMPLEX16Sequence * data
COMPLEX16 * data
COMPLEX16FrequencySeries * hctilde
COMPLEX16FrequencySeries * hptilde
REAL8Sequence * data
LALUnit sampleUnits
LIGOTimeGPS epoch
CHAR name[LALNameLength]
REAL8 * data
double f_min
Definition: unicorn.c:22
double deltaT
Definition: unicorn.c:24
double f_max
Definition: unicorn.c:23