LALSimulation  5.4.0.1-b72065a

Detailed Description

General routines for generating binary inspiral waveforms.

Prototypes

int XLALSimInspiralPNPolarizationWaveformsFromModes (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8TimeSeries *v, REAL8TimeSeries *phi, REAL8 v0, REAL8 m1, REAL8 m2, REAL8 r, REAL8 i, int O)
 Given time series for a binary's orbital dynamical variables, construct the waveform polarizations h+ and hx as a sum of -2 spin-weighted spherical harmonic modes, h_lm. More...
 
int XLALSimInspiralPolarizationsFromSphHarmTimeSeries (REAL8TimeSeries **hp, REAL8TimeSeries **hc, SphHarmTimeSeries *hlms, REAL8 iota, REAL8 phiRef)
 Compute the polarizations from all the -2 spin-weighted spherical harmonic modes stored in 'hlms'. More...
 
int XLALSimInspiralPolarizationsFromChooseFDModes (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 UNUSED longAscNodes, const REAL8 UNUSED eccentricity, const REAL8 UNUSED meanPerAno, const REAL8 deltaF, const REAL8 f_min, const REAL8 f_max, REAL8 f_ref, LALDict *LALparams, const Approximant approximant)
 Function returning the Fourier domain polarizations for positive frequencies built from the individual modes computed with ChooseFDModes. More...
 
int XLALSimInspiralPolarizationsFromSphHarmFrequencySeries (COMPLEX16FrequencySeries **hp, COMPLEX16FrequencySeries **hc, SphHarmFrequencySeries *hlms, REAL8 theta, REAL8 phi)
 Return polarizations for positive frequencies built by summing the individual modes present in the input array SphHarmFrequencySeries *hlms computed with ChooseFDModes. More...
 
int XLALSimInspiralPNPolarizationWaveformsEccentric (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8TimeSeries *V, REAL8TimeSeries *Ecc, REAL8TimeSeries *U, REAL8TimeSeries *Phi, REAL8 m1, REAL8 m2, REAL8 r, REAL8 i, int ampO, int ph_O)
 Given time series for a binary's orbital dynamical variables, computes the radial and angular orbital motion; then constructs the waveform polarizations h+ and hx in terms of the relative separation r, the true anomaly phi and their time derivatives. More...
 
int XLALSimInspiralPrecessingPolarizationWaveforms (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8TimeSeries *V, REAL8TimeSeries *Phi, REAL8TimeSeries *S1x, REAL8TimeSeries *S1y, REAL8TimeSeries *S1z, REAL8TimeSeries *S2x, REAL8TimeSeries *S2y, REAL8TimeSeries *S2z, REAL8TimeSeries *LNhatx, REAL8TimeSeries *LNhaty, REAL8TimeSeries *LNhatz, REAL8TimeSeries *E1x, REAL8TimeSeries *E1y, REAL8TimeSeries *E1z, REAL8 m1, REAL8 m2, REAL8 r, INT4 ampO)
 Computes polarizations h+ and hx for a spinning, precessing binary when provided time series of all the dynamical quantities. More...
 
int XLALSimInspiralPrecessingPolarizationWaveformHarmonic (COMPLEX16 *hplus, COMPLEX16 *hcross, REAL8 v, REAL8 s1x, REAL8 s1y, REAL8 s1z, REAL8 s2x, REAL8 s2y, REAL8 s2z, REAL8 lnhx, REAL8 lnhy, REAL8 lnhz, REAL8 e1x, REAL8 e1y, REAL8 e1z, REAL8 dm, REAL8 eta, REAL8 v0, INT4 n, INT4 ampO)
 Computes polarizations h+ and hx for a spinning, precessing binary when provided a single value of all the dynamical quantities. More...
 

New Interface Generator Routines

void XLALDestroySimInspiralGenerator (LALSimInspiralGenerator *generator)
 Destroy LALSimInspiralGenerator object. More...
 
LALSimInspiralGenerator * XLALCreateSimInspiralGenerator (const LALSimInspiralGenerator *generator, LALDict *params)
 Create LALSimInspiralGenerator object. More...
 
LALSimInspiralGenerator * XLALSimInspiralChooseGenerator (Approximant approx, LALDict *params)
 Returns LALSimInspiralGenerator object from approximant. More...
 
const charXLALSimInspiralGeneratorName (LALSimInspiralGenerator *generator)
 Return approximant name from generator object. More...
 

New Interface Waveform Routines

int XLALSimInspiralGenerateTDWaveform (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, LALDict *params, LALSimInspiralGenerator *generator)
 Returns time-domain polarizations for a specific approximant. More...
 
int XLALSimInspiralGenerateTDModes (SphHarmTimeSeries **hlm, LALDict *params, LALSimInspiralGenerator *generator)
 Compute time-domain modes for a specific approximant. More...
 
int XLALSimInspiralGenerateFDWaveform (COMPLEX16FrequencySeries **hplus, COMPLEX16FrequencySeries **hcross, LALDict *params, LALSimInspiralGenerator *generator)
 Returns frequency-domain polarizations for a specific approximant. More...
 
int XLALSimInspiralGenerateFDModes (SphHarmFrequencySeries **hlm, LALDict *params, LALSimInspiralGenerator *generator)
 Compute frequency-domain modes for a specific approximant. More...
 

New Interface parse parameters Routines

void XLALSimInspiralParseDictionaryToChooseTDWaveform (REAL8 *m1, REAL8 *m2, REAL8 *S1x, REAL8 *S1y, REAL8 *S1z, REAL8 *S2x, REAL8 *S2y, REAL8 *S2z, REAL8 *distance, REAL8 *inclination, REAL8 *phiRef, REAL8 *longAscNodes, REAL8 *eccentricity, REAL8 *meanPerAno, REAL8 *deltaT, REAL8 *f_min, REAL8 *f_ref, LALDict *params)
 Insert all the input arguments needed by XALSimInspiralChooseTDWaveform() into a laldictionary. More...
 
void XLALSimInspiralParseDictionaryToChooseTDModes (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 *distance, INT4 *lmax, LALDict *params)
 Insert all the input arguments needed by XLALSimInspiralChooseTDModes() into a laldictionary. More...
 
void XLALSimInspiralParseDictionaryToChooseFDWaveform (REAL8 *m1, REAL8 *m2, REAL8 *S1x, REAL8 *S1y, REAL8 *S1z, REAL8 *S2x, REAL8 *S2y, REAL8 *S2z, REAL8 *distance, REAL8 *inclination, REAL8 *phiRef, REAL8 *longAscNodes, REAL8 *eccentricity, REAL8 *meanPerAno, REAL8 *deltaF, REAL8 *f_min, REAL8 *f_max, REAL8 *f_ref, LALDict *params)
 Insert all the input arguments needed by XLALSimInspiralChooseFDWaveform() into a laldictionary. More...
 
void XLALSimInspiralParseDictionaryToChooseFDModes (REAL8 *m1, REAL8 *m2, REAL8 *S1x, REAL8 *S1y, REAL8 *S1z, REAL8 *S2x, REAL8 *S2y, REAL8 *S2z, REAL8 *deltaF, REAL8 *f_min, REAL8 *f_max, REAL8 *f_ref, REAL8 *phiRef, REAL8 *distance, REAL8 *inclination, LALDict *params)
 Insert all the input arguments needed by XLALSimInspiralChooseFDModes() into a laldictionary. More...
 

General Waveform Switching Generation Routines

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 waveforms, all known spin effects up to given PN order are included Returns the waveform in the time domain. More...
 
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 waveforms, all known spin effects up to given PN order are included Returns the waveform in the frequency domain. More...
 
SphHarmTimeSeriesXLALSimInspiralTDModesFromPolarizations (REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 phiRef, REAL8 longAscNodes, REAL8 eccentricity, REAL8 meanPerAno, REAL8 deltaT, REAL8 f_min, REAL8 f_ref, LALDict *LALparams, Approximant approximant)
 Generates an time domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned, suitable for injection into data, and decomposed into the (2, \(\pm\) 2), spin -2 weighted spherical harmonic modes. More...
 
int XLALSimInspiralTDFromTD (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 longAscNodes, REAL8 eccentricity, REAL8 meanPerAno, REAL8 deltaT, REAL8 f_min, REAL8 f_ref, LALDict *LALparams, Approximant approximant)
 Helper routines for XLALSimInspiralTD(): performs conditioning of a TD waveform. More...
 
int XLALSimInspiralTDFromFD (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 longAscNodes, REAL8 eccentricity, REAL8 meanPerAno, REAL8 deltaT, REAL8 f_min, REAL8 f_ref, LALDict *LALparams, Approximant approximant)
 Helper routines for XLALSimInspiralTD(): performs conditioning of a FD waveform and transforms it to TD. More...
 
int XLALSimInspiralTD (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 longAscNodes, REAL8 eccentricity, REAL8 meanPerAno, REAL8 deltaT, REAL8 f_min, REAL8 f_ref, LALDict *LALparams, Approximant approximant)
 Generates an time domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned and suitable for injection into data. More...
 
int XLALSimInspiralFD (COMPLEX16FrequencySeries **hptilde, COMPLEX16FrequencySeries **hctilde, REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 distance, REAL8 inclination, REAL8 phiRef, REAL8 longAscNodes, REAL8 eccentricity, REAL8 meanPerAno, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 f_ref, LALDict *LALparams, Approximant approximant)
 Generates a frequency domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned and suitable for injection into data. More...
 
int XLALSimInspiralChooseWaveform (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 meanPerAno, const REAL8 deltaT, const REAL8 f_min, const REAL8 f_ref, LALDict *LALpars, const Approximant approximant)
 

General Waveform Switching Mode Generation Routines

SphHarmTimeSeriesXLALSimInspiralChooseTDModes (UNUSED 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 distance, LALDict *params, int lmax, Approximant approximant)
 Interface to compute a set of -2 spin-weighted spherical harmonic modes for a binary inspiral for a given waveform approximant. More...
 
SphHarmFrequencySeriesXLALSimInspiralChooseFDModes (REAL8 m1, REAL8 m2, REAL8 S1x, REAL8 S1y, REAL8 S1z, REAL8 S2x, REAL8 S2y, REAL8 S2z, REAL8 deltaF, REAL8 f_min, REAL8 f_max, REAL8 f_ref, REAL8 phiRef, REAL8 distance, REAL8 inclination, LALDict *params, Approximant approximant)
 Interface to compute a set of -2 spin-weighted spherical harmonic modes for a binary merger for a given waveform approximant in the Fourier domain. More...
 
SphHarmTimeSeriesXLALSimInspiralModesTD (REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 f_min, REAL8 f_ref, REAL8 r, LALDict *LALpars, int lmax, Approximant approximant)
 Interface to compute a conditioned set of -2 spin-weighted spherical harmonic modes for a binary inspiral. More...
 
COMPLEX16TimeSeriesXLALSimInspiralChooseTDMode (REAL8 deltaT, REAL8 m1, REAL8 m2, REAL8 f_min, REAL8 f_ref, REAL8 r, REAL8 lambda1, REAL8 lambda2, LALSimInspiralWaveformFlags *waveFlags, LALSimInspiralTestGRParam *nonGRparams, int amplitudeO, int phaseO, int l, int m, Approximant approximant)
 Interface to compute a single -2 spin-weighted spherical harmonic mode for a binary inspiral of any available amplitude and phase PN order. More...
 

Routines for Generating Inspiral Waveforms from Orbital Data

int XLALSimInspiralPNPolarizationWaveforms (REAL8TimeSeries **hplus, REAL8TimeSeries **hcross, REAL8TimeSeries *V, REAL8TimeSeries *Phi, REAL8 v0, REAL8 m1, REAL8 m2, REAL8 r, REAL8 i, int ampO)
 Given time series for a binary's orbital dynamical variables, construct the waveform polarizations h+ and hx directly. More...
 

Function Documentation

◆ XLALDestroySimInspiralGenerator()

void XLALDestroySimInspiralGenerator ( LALSimInspiralGenerator *  generator)

Destroy LALSimInspiralGenerator object.

Definition at line 390 of file LALSimInspiral.c.

◆ XLALCreateSimInspiralGenerator()

LALSimInspiralGenerator* XLALCreateSimInspiralGenerator ( const LALSimInspiralGenerator *  generator,
LALDict *  params 
)

Create LALSimInspiralGenerator object.

The LALDict can be None for all the C approximants. For external python model, the LALDict must contain the file and class object where the generate_waveform methods are defined. Activating the option condition in the LALDict will generate the waveform with adequate conditioning for (inverse)Fourier transform.

Definition at line 411 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseGenerator()

LALSimInspiralGenerator* XLALSimInspiralChooseGenerator ( Approximant  approx,
LALDict *  params 
)

Returns LALSimInspiralGenerator object from approximant.

The LALDict can be None for all the C approximants. For external python model, the LALDict must contain the file and class object where the generate_waveform methods are defined. Activating the option condition in the LALDict will generate the waveform with adequate conditioning for (inverse)Fourier transform.

Definition at line 439 of file LALSimInspiral.c.

◆ XLALSimInspiralGeneratorName()

const char* XLALSimInspiralGeneratorName ( LALSimInspiralGenerator *  generator)

Return approximant name from generator object.

Definition at line 450 of file LALSimInspiral.c.

◆ XLALSimInspiralGenerateTDWaveform()

int XLALSimInspiralGenerateTDWaveform ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
LALDict *  params,
LALSimInspiralGenerator *  generator 
)

Returns time-domain polarizations for a specific approximant.

Equivalent to XLALSimInspiralChooseTDWaveform(). Equivalent to XLALSimInspiralTD() if the option condition is activated in the LALDict. The waveform arguments are inserted into the LALDict. The generator carries the info about the approximant and potentially extra data which could be recycled by the model to speed-up calculation.

The parameters in the LALDict must be in SI units.

Definition at line 470 of file LALSimInspiral.c.

◆ XLALSimInspiralGenerateTDModes()

int XLALSimInspiralGenerateTDModes ( SphHarmTimeSeries **  hlm,
LALDict *  params,
LALSimInspiralGenerator *  generator 
)

Compute time-domain modes for a specific approximant.

Equivalent to XLALSimInspiralChooseTDModes(). The only difference is that the SphHarmSeries object needs to be passed as an argument to the function. The actual returned value is an integer which indicates success or error in the waveform evaluation (see https://lscsoft.docs.ligo.org/lalsuite/lal/group___x_l_a_l_error__h.html). The waveform arguments are inserted into the LALDict. The generator carries the info about the approximant and potentially extra data which could be recycled by the model to speed-up calculation.

The parameters in the LALDict must be in SI units.

Definition at line 493 of file LALSimInspiral.c.

◆ XLALSimInspiralGenerateFDWaveform()

int XLALSimInspiralGenerateFDWaveform ( COMPLEX16FrequencySeries **  hplus,
COMPLEX16FrequencySeries **  hcross,
LALDict *  params,
LALSimInspiralGenerator *  generator 
)

Returns frequency-domain polarizations for a specific approximant.

Equivalent to XLALSimInspiralChooseFDWaveform(). Equivalent to XLALSimInspiralFD() if the option condition is activated in the LALDict. The waveform arguments are inserted into the LALDict. The generator carries the info about the approximant and potentially extra data which could be recycled by the model to speed-up calculation.

The parameters in the LALDict must be in SI units.

Definition at line 515 of file LALSimInspiral.c.

◆ XLALSimInspiralGenerateFDModes()

int XLALSimInspiralGenerateFDModes ( SphHarmFrequencySeries **  hlm,
LALDict *  params,
LALSimInspiralGenerator *  generator 
)

Compute frequency-domain modes for a specific approximant.

Equivalent to XLALSimInspiralChooseFDModes. The only difference is that the SphHarmSeries object needs to be passed as an argument to the function. The actual returned value is an integer which indicates success or error in the waveform evaluation (see https://lscsoft.docs.ligo.org/lalsuite/lal/group___x_l_a_l_error__h.html). The waveform arguments are inserted into the LALDict. The generator carries the info about the approximant and potentially extra data which could be recycled by the model to speed-up calculation.

The parameters in the LALDict must be in SI units.

Definition at line 537 of file LALSimInspiral.c.

◆ XLALSimInspiralParseDictionaryToChooseTDWaveform()

void XLALSimInspiralParseDictionaryToChooseTDWaveform ( REAL8 m1,
REAL8 m2,
REAL8 S1x,
REAL8 S1y,
REAL8 S1z,
REAL8 S2x,
REAL8 S2y,
REAL8 S2z,
REAL8 distance,
REAL8 inclination,
REAL8 phiRef,
REAL8 longAscNodes,
REAL8 eccentricity,
REAL8 meanPerAno,
REAL8 deltaT,
REAL8 f_min,
REAL8 f_ref,
LALDict *  params 
)

Insert all the input arguments needed by XALSimInspiralChooseTDWaveform() into a laldictionary.

Parameters
[out]m1mass of companion 1 (kg)
[out]m2mass of companion 2 (kg)
[out]S1xx-component of the dimensionless spin of object 1
[out]S1yy-component of the dimensionless spin of object 1
[out]S1zz-component of the dimensionless spin of object 1
[out]S2xx-component of the dimensionless spin of object 2
[out]S2yy-component of the dimensionless spin of object 2
[out]S2zz-component of the dimensionless spin of object 2
[out]distancedistance of source (m)
[out]inclinationinclination of source (rad)
[out]phiRefreference orbital phase (rad)
[out]longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
[out]eccentricityeccentrocity at reference epoch
[out]meanPerAnomean anomaly of periastron
[out]deltaTsampling interval (s)
[out]f_minstarting GW frequency (Hz)
[out]f_refreference frequency (Hz)
paramsInput lal dictionary with waveform (and optional) parameters *

Definition at line 562 of file LALSimInspiral.c.

◆ XLALSimInspiralParseDictionaryToChooseTDModes()

void XLALSimInspiralParseDictionaryToChooseTDModes ( 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 distance,
INT4 lmax,
LALDict *  params 
)

Insert all the input arguments needed by XLALSimInspiralChooseTDModes() into a laldictionary.

Parameters
[out]phiRefreference orbital phase (rad)
[out]deltaTsampling interval (s)
[out]m1mass of companion 1 (kg)
[out]m2mass of companion 2 (kg)
[out]S1xx-component of the dimensionless spin of object 1
[out]S1yy-component of the dimensionless spin of object 1
[out]S1zz-component of the dimensionless spin of object 1
[out]S2xx-component of the dimensionless spin of object 2
[out]S2yy-component of the dimensionless spin of object 2
[out]S2zz-component of the dimensionless spin of object 2
[out]f_minstarting GW frequency (Hz)
[out]f_refreference frequency (Hz)
[out]distancedistance of source (m)
[out]lmaxgenerate all modes with l <= lmax
paramsInput lal dictionary with waveform (and optional) parameters *

Definition at line 607 of file LALSimInspiral.c.

◆ XLALSimInspiralParseDictionaryToChooseFDWaveform()

void XLALSimInspiralParseDictionaryToChooseFDWaveform ( REAL8 m1,
REAL8 m2,
REAL8 S1x,
REAL8 S1y,
REAL8 S1z,
REAL8 S2x,
REAL8 S2y,
REAL8 S2z,
REAL8 distance,
REAL8 inclination,
REAL8 phiRef,
REAL8 longAscNodes,
REAL8 eccentricity,
REAL8 meanPerAno,
REAL8 deltaF,
REAL8 f_min,
REAL8 f_max,
REAL8 f_ref,
LALDict *  params 
)

Insert all the input arguments needed by XLALSimInspiralChooseFDWaveform() into a laldictionary.

Parameters
[out]m1mass of companion 1 (kg)
[out]m2mass of companion 2 (kg)
[out]S1xx-component of the dimensionless spin of object 1
[out]S1yy-component of the dimensionless spin of object 1
[out]S1zz-component of the dimensionless spin of object 1
[out]S2xx-component of the dimensionless spin of object 2
[out]S2yy-component of the dimensionless spin of object 2
[out]S2zz-component of the dimensionless spin of object 2
[out]distancedistance of source (m)
[out]inclinationinclination of source (rad)
[out]phiRefreference orbital phase (rad)
[out]longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
[out]eccentricityeccentrocity at reference epoch
[out]meanPerAnomean anomaly of periastron
[out]deltaFfrequency interval (Hz)
[out]f_minstarting GW frequency (Hz)
[out]f_maxending GW frequency (Hz)
[out]f_refreference frequency (Hz)
paramsInput lal dictionary with waveform (and optional) parameters *

Definition at line 646 of file LALSimInspiral.c.

◆ XLALSimInspiralParseDictionaryToChooseFDModes()

void XLALSimInspiralParseDictionaryToChooseFDModes ( REAL8 m1,
REAL8 m2,
REAL8 S1x,
REAL8 S1y,
REAL8 S1z,
REAL8 S2x,
REAL8 S2y,
REAL8 S2z,
REAL8 deltaF,
REAL8 f_min,
REAL8 f_max,
REAL8 f_ref,
REAL8 phiRef,
REAL8 distance,
REAL8 inclination,
LALDict *  params 
)

Insert all the input arguments needed by XLALSimInspiralChooseFDModes() into a laldictionary.

Parameters
[out]m1mass of companion 1 (kg)
[out]m2mass of companion 2 (kg)
[out]S1xx-component of the dimensionless spin of object 1
[out]S1yy-component of the dimensionless spin of object 1
[out]S1zz-component of the dimensionless spin of object 1
[out]S2xx-component of the dimensionless spin of object 2
[out]S2yy-component of the dimensionless spin of object 2
[out]S2zz-component of the dimensionless spin of object 2
[out]deltaFsampling interval (s)
[out]f_minstarting GW frequency (Hz)
[out]f_maxending GW frequency (Hz)
[out]f_refreference GW frequency (Hz)
[out]phiRefreference phase (rad)
[out]distancedistance of source (m)
[out]inclinationinclination of source (rad)
paramsInput lal dictionary with waveform (and optional) parameters *

Definition at line 694 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseTDWaveform()

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 waveforms, all known spin effects up to given PN order are included Returns the waveform in the time domain.

The parameters passed must be in SI units.

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
paramsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 750 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseFDWaveform()

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 waveforms, all known spin effects up to given PN order are included Returns the waveform in the frequency domain.

Parameters
hptildeFD plus polarization
hctildeFD cross polarization
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentricity at reference epoch
meanPerAnomean anomaly of periastron
deltaFsampling interval (Hz)
f_minstarting GW frequency (Hz)
f_maxending GW frequency (Hz)
f_refReference frequency (Hz)
paramsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 818 of file LALSimInspiral.c.

◆ XLALSimInspiralTDModesFromPolarizations()

SphHarmTimeSeries* XLALSimInspiralTDModesFromPolarizations ( REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  phiRef,
REAL8  longAscNodes,
REAL8  eccentricity,
REAL8  meanPerAno,
REAL8  deltaT,
REAL8  f_min,
REAL8  f_ref,
LALDict *  LALparams,
Approximant  approximant 
)

Generates an time domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned, suitable for injection into data, and decomposed into the (2, \(\pm\) 2), spin -2 weighted spherical harmonic modes.

NOTE: This is an algebraic decomposition, and will only be correct for approximants which use only the dominant 2, \(\pm\) 2 mode.

For spinning waveforms, all known spin effects up to given PN order are included

This routine can generate FD approximants and transform them into the time domain. Waveforms are generated beginning at a slightly lower starting frequency and tapers are put in this early region so that the waveform smoothly turns on. Artifacts at the very end of the waveform are also tapered. The resulting waveform is high-pass filtered at frequency f_min so that it should have little content at lower frequencies.

This routine used to have one additional parameter relative to XLALSimInspiralChooseTDWaveform: the redshift, z, of the waveform, which is now stuffed into the LALDict structure. This should be set to zero (default value) for sources in the nearby universe (that are nearly at rest relative to the earth). For sources at cosmological distances, the mass parameters m1 and m2 should be interpreted as the physical (source frame) masses of the bodies and the distance parameter r is the comoving (transverse) distance. If the calling routine has already applied cosmological "corrections" to m1 and m2 and regards r as a luminosity distance then the redshift factor should again be set to zero.

Note
The parameters passed must be in SI units.
Parameters
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 917 of file LALSimInspiral.c.

◆ XLALSimInspiralTDFromTD()

int XLALSimInspiralTDFromTD ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  longAscNodes,
REAL8  eccentricity,
REAL8  meanPerAno,
REAL8  deltaT,
REAL8  f_min,
REAL8  f_ref,
LALDict *  LALparams,
Approximant  approximant 
)

Helper routines for XLALSimInspiralTD(): performs conditioning of a TD waveform.

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 978 of file LALSimInspiral.c.

◆ XLALSimInspiralTDFromFD()

int XLALSimInspiralTDFromFD ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  longAscNodes,
REAL8  eccentricity,
REAL8  meanPerAno,
REAL8  deltaT,
REAL8  f_min,
REAL8  f_ref,
LALDict *  LALparams,
Approximant  approximant 
)

Helper routines for XLALSimInspiralTD(): performs conditioning of a FD waveform and transforms it to TD.

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1077 of file LALSimInspiral.c.

◆ XLALSimInspiralTD()

int XLALSimInspiralTD ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  longAscNodes,
REAL8  eccentricity,
REAL8  meanPerAno,
REAL8  deltaT,
REAL8  f_min,
REAL8  f_ref,
LALDict *  LALparams,
Approximant  approximant 
)

Generates an time domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned and suitable for injection into data.

For spinning waveforms, all known spin effects up to given PN order are included

This routine can generate FD approximants and transform them into the time domain. Waveforms are generated beginning at a slightly lower starting frequency and tapers are put in this early region so that the waveform smoothly turns on. Artifacts at the very end of the waveform are also tapered. The resulting waveform is high-pass filtered at frequency f_min so that it should have little content at lower frequencies.

If calling with precessing time-domain approximants for which the reference frequency is the starting frequency, or if calling with NR_hdf5 approximant, the starting frequency is not altered. Uses XLALSimInspiralGetSpinFreqFromApproximant to determine appropriate behaviour. Similarly, if calling time-domain models for which a starting frequency of zero is allowed (as set in XLALSimInspiralGetAllowZeroMinFreqFromApproximant), the starting frequency is never altered, independent of f_min.

This routine used to have one additional parameter relative to XLALSimInspiralChooseTDWaveform: the redshift, z, of the waveform, which is now stuffed into the LALDict structure. This should be set to zero (default value) for sources in the nearby universe (that are nearly at rest relative to the earth). For sources at cosmological distances, the mass parameters m1 and m2 should be interpreted as the physical (source frame) masses of the bodies and the distance parameter r is the comoving (transverse) distance. If the calling routine has already applied cosmological "corrections" to m1 and m2 and regards r as a luminosity distance then the redshift factor should again be set to zero.

Note
The parameters passed must be in SI units.
Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1260 of file LALSimInspiral.c.

◆ XLALSimInspiralFD()

int XLALSimInspiralFD ( COMPLEX16FrequencySeries **  hptilde,
COMPLEX16FrequencySeries **  hctilde,
REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  distance,
REAL8  inclination,
REAL8  phiRef,
REAL8  longAscNodes,
REAL8  eccentricity,
REAL8  meanPerAno,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  f_ref,
LALDict *  LALparams,
Approximant  approximant 
)

Generates a frequency domain inspiral waveform using the specified approximant; the resulting waveform is appropriately conditioned and suitable for injection into data.

For spinning waveforms, all known spin effects up to given PN order are included.

This routine can generate TD approximants and transform them into the frequency domain. Waveforms are generated beginning at a slightly lower starting frequency and tapers are put in this early region so that the waveform smoothly turns on.

If an FD approximant is used, this routine applies tapers in the frequency domain between the slightly-lower frequency and the requested f_min. Also, the phase of the waveform is adjusted to introduce a time shift. This time shift should allow the resulting waveform to be Fourier transformed into the time domain without wrapping the end of the waveform to the beginning.

This routine assumes that f_max is the Nyquist frequency of a corresponding time-domain waveform, so that deltaT = 0.5 / f_max. If deltaF is set to 0 then this routine computes a deltaF that is small enough to represent the Fourier transform of a time-domain waveform while ensuring the length of the time-domain signal if a power of 2. If deltaF is specified but f_max / deltaF is not a power of 2, and the waveform approximant is a time-domain approximant, then f_max is increased so that f_max / deltaF is the next power of 2. (If the user wishes to discard the extra high frequency content, this must be done separately.)

The user should take care to ensure that deltaF is sufficiently small to contain the full signal (time series duration = 1 / deltaF). If the provided deltaF is too large the signal will be abruptly truncated for time-domain waveform generators. For frequency-domain generators the signal will be aliased in the frequency domain.

Similarly, if the provided f_max is less than the ringdown frequency the underlying waveform generator may raise an error. If not, the frequency domain signal will be aliased in the frequency domain.

Some waveform approximants have built in checks for the maximum frequency and signal length

This routine used to have one additional parameter relative to XLALSimInspiralChooseTDWaveform: the redshift, z, of the waveform, which is now stuffed into the LALDict. This should be set to zero (default value) for sources in the nearby universe (that are nearly at rest relative to the earth). For sources at cosmological distances, the mass parameters m1 and m2 should be interpreted as the physical (source frame) masses of the bodies and the distance parameter r is the comoving (transverse) distance. If the calling routine has already applied cosmological "corrections" to m1 and m2 and regards r as a luminosity distance then the redshift factor should again be set to zero.

Note
The parameters passed must be in SI units.
Parameters
hptildeFD plus polarization
hctildeFD cross polarization
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentricity at reference epoch
meanPerAnomean anomaly of periastron
deltaFsampling interval (Hz)
f_minstarting GW frequency (Hz)
f_maxending GW frequency (Hz)
f_refReference frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1350 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseWaveform()

int XLALSimInspiralChooseWaveform ( 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  meanPerAno,
const REAL8  deltaT,
const REAL8  f_min,
const REAL8  f_ref,
LALDict *  LALpars,
const Approximant  approximant 
)
Deprecated:
Use XLALSimInspiralChooseTDWaveform() instead

Chooses between different approximants when requesting a waveform to be generated For spinning waveforms, all known spin effects up to given PN order are included

The parameters passed must be in SI units.

Parameters
hplus+-polarization waveform
hcrossx-polarization waveform
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentrocity at reference epoch
meanPerAnomean anomaly of periastron
deltaTsampling interval (s)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
LALparsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1401 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseTDModes()

SphHarmTimeSeries* XLALSimInspiralChooseTDModes ( UNUSED 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  distance,
LALDict *  params,
int  lmax,
Approximant  approximant 
)

Interface to compute a set of -2 spin-weighted spherical harmonic modes for a binary inspiral for a given waveform approximant.

PN Approximants (TaylorT1 - T4), EOBNRv2 (EOBNRv2HM), NRSur7dq2, NRSur7dq4 NRHybSur3dq8 and spin-precessing SpintaylorT1, T5, T4 are implemented.

The EOBNRv2 model returns the (2,2), (2,1), (3,3), (4,4), and (5,5) modes. Note that the inclination parameter is not passed to create hlm modes, hence to recover the correct h+,x one has to combine the hlm modes with Euler angles alpha=0, iota=inclination, psi=0,Pi/2 (according to the approximat) i.e. (h+ + I hx) (psi,iota,alpha)= e^(-2Ialpha) Sum_{l,m} Y_lm(-iota,-psi) h_lm, or equivalently rotate h_lm -> h'lm=DWigner(psi,iota,alpha) h_lm and then obtain (h+ + I hx) = Sum{l,m} Y_lm(0,0) h'_lm,

Parameters
phiRefreference orbital phase (rad). This variable is not used and only kept here for backwards compatibility
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
distancedistance of source (m)
paramsLAL dictionary containing accessory parameters
lmaxgenerate all modes with l <= lmax
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1455 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseFDModes()

SphHarmFrequencySeries* XLALSimInspiralChooseFDModes ( REAL8  m1,
REAL8  m2,
REAL8  S1x,
REAL8  S1y,
REAL8  S1z,
REAL8  S2x,
REAL8  S2y,
REAL8  S2z,
REAL8  deltaF,
REAL8  f_min,
REAL8  f_max,
REAL8  f_ref,
REAL8  phiRef,
REAL8  distance,
REAL8  inclination,
LALDict *  params,
Approximant  approximant 
)

Interface to compute a set of -2 spin-weighted spherical harmonic modes for a binary merger for a given waveform approximant in the Fourier domain.

Non-precessing models IMRPhenomXHM, SEOBNRv4HM_ROM, SEOBNRv5(HM)_ROM and IMRPhenomHM and the precessing IMRPhenomXPHM are implemented. By default, all the modes available in the model are returned, although the list can be specified through the ModeArray option in the LAL dictionary.

In the Fourier domain the modes span over the whole frequency regime (positive and negative frequencies). However, in the aligned spin case, the modes have support only in one half of the frequency regime. The LAL conventions establish that the negative modes (m<0) have support for positive frequencies while the positive modes (m>0) have support for negative frequencies (this is based in the right hand rule and Fourier transform definition in LAL). Due to the equatorial symmetry of non-precessng systems, there exist a relation between them: h_{lm}(f) = (-1)^l h*_{l-m}(-f). In the precessing case this symmetry is broken and all the modes have support for both positive and negative frequencies.

For this reason, the ouput SphHarmFrequencySeries object will return the modes in the whole frequency range. The frequencies of this object will be sorted as: -f_max,...,-f_min,...0,...,f_min,...,f_max. The values of the waveform will be sorted correspondingly. Consequently, in the aligned spin case, half of the frequency spectrum consists of zeros.

It is relevant to mention why the arguments inclination and phiRef are needed for computing the h_lm. For AS models the argument inclination is irrelevant and will not be use since it only enters in the Ylm. However, for the precessing model, since the modes are returned in the J-frame, we need the inclination argument to carry out the Euler transformation from the co-precessing L-frame to the inertial J-frame. Regarding the argument phiRef, this affects the output of the precessing model due to the same reason as before, while for the AS models it would not affect the output of SEOBNRv4HM_ROM, SEOBNRv5(HM)_ROM but will change the output of IMRPhenomHM and IMRPhenomXHM (this is due to the internals workings of the models). If one wants to built the polarizations from the individual modes of ChooseFDModes must be aware of this behaviour. Ideally, one would call ChooseFDModes with phiRef=0 to obtain the h_lms, then build the Fourier domain polarizations as

h_+ (f) = 1/2 sum_{l=2} sum_{m=-l}^{m=l} ( h_lm(f) * Y_lm(theta, vphi) + h*_lm(-f) * Y*_lm(theta, vphi) ) h_x (f) = i/2 sum_{l=2} sum_{m=-l}^{m=l} ( h_lm(f) * Y_lm(theta, vphi) - h*_lm(-f) * Y*_lm(theta, vphi) )

where theta is the inclination and vphi is pi/2 - phiRef.

If one does this, one will find generally a very close result to ChooseFDWaveform with very small mismatches (~10^-9 for IMRPhenomXHM), and this is what is found if one uses the XLALSimInspiralPolarizationsFromSphHarmFrequencySeries function, however this is not close to machine precision. The reason is that IMRPhenomHM and IMRPhenomXHM use internally the phiRef argument to compute the h_lms because at that time phiRef was considered to be also a reference phase for the h_lms and not just the argument for the azimuthal angle in the Y_lms. To take this into account we provide also the function XLALSimInspiralPolarizationsFromChooseFDModes which build the polarizations in the proper way for each model returning a result close to machine precision with ChooseFDWaveform.

For the precessing model IMRPhenomXPHM, since the h_lms are returned in the J-frame one must build the polarizations using theta = theta_JN and vphi = 0. The parameter theta_JN is computed internally when using XLALSimInspiralPolarizationsFromChooseFDModes and again here the result is close to machine precision to ChooseFDWaveform. However, one would have to compute theta_JN personally when using XLALSimInspiralPolarizationsFromSphHarmFrequencySeries. For IMRPhenomXPHM this parameter can be computed using XLALSimIMRPhenomXPCalculateModelParametersFromSourceFrame. Eventhough one would still have to correct with the polarization angle.

By default all the modes available in the model will be returned, both positive and negative modes. The mode content of AS models can be adjusted through the ModeArray option in the LAL dictionary argument, and it accepts any set of modes e.g. (2,2),(2,-1),(3,3),... For IMRPhenomXPHM, we can specify both the modes in the co-precessing L-frame, which are used to do the twisting-up, and the ouput modes in the inertial J-frame. The set of modes in the L-frame are specified with the standard ModeArray option while the set of modes in the J-frame are specified in a new option called ModeArrayJframe. Notice that in IMRPhenomXPHM, ModeArray does not distinguish between positive or negative modes and it always twists-up both positive and negative modes, i.e. the sets {(2,2)}, {(2,-2)} or {(2,2),(2,-2)} would return the same result. For the modes in ModeArrayJframe, we can specify both positive or negative modes for example {(2,2),(2,-2),(2,-1),(3,3),...}.

Parameters
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
deltaFsampling interval (s)
f_minstarting GW frequency (Hz)
f_maxending GW frequency (Hz)
f_refreference GW frequency (Hz)
phiRefreference phase (rad)
distancedistance of source (m)
inclinationinclination of source (rad)
paramsLAL dictionary containing accessory parameters (optional mode array)
approximantapproximant to use for waveform production

Definition at line 1565 of file LALSimInspiral.c.

◆ XLALSimInspiralModesTD()

SphHarmTimeSeries* XLALSimInspiralModesTD ( REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  f_min,
REAL8  f_ref,
REAL8  r,
LALDict *  LALpars,
int  lmax,
Approximant  approximant 
)

Interface to compute a conditioned set of -2 spin-weighted spherical harmonic modes for a binary inspiral.

This routine is a wrapper for XLALSimInspiralChooseTDModes which applies waveform conditioning to the modes generated by that routine. The conditioning algorithm is analogous to that performed in XLALSimInspiralTD. Note that the modes are high-pass filtered at frequency f_min, which is specified for the m=2 mode, which means that the low-frequency part of the m=1 mode is removed by the filtering. The phasing is computed with any of the TaylorT1, T2, T3, T4 methods. It can also return the (2,2), (2,1), (3,3), (4,4), (5,5) modes of the EOBNRv2 model. Note that EOBNRv2 will ignore ampO, phaseO, lmax and f_ref arguments.

Parameters
deltaTSampling interval (s)
m1Mass of companion 1 (kg)
m2Mass of companion 2 (kg)
f_minStarting GW frequency (Hz)
f_refReference GW frequency (Hz)
rDistance of source (m)
LALparsLAL dictionary containing accesory parameters
lmaxGenerate all modes with l <= lmax
approximantPost-Newtonian approximant to use for waveform production
Returns
Linked list of SphHarmTimeSeries modes.

Definition at line 1650 of file LALSimInspiral.c.

◆ XLALSimInspiralChooseTDMode()

COMPLEX16TimeSeries* XLALSimInspiralChooseTDMode ( REAL8  deltaT,
REAL8  m1,
REAL8  m2,
REAL8  f_min,
REAL8  f_ref,
REAL8  r,
REAL8  lambda1,
REAL8  lambda2,
LALSimInspiralWaveformFlags *  waveFlags,
LALSimInspiralTestGRParam nonGRparams,
int  amplitudeO,
int  phaseO,
int  l,
int  m,
Approximant  approximant 
)

Interface to compute a single -2 spin-weighted spherical harmonic mode for a binary inspiral of any available amplitude and phase PN order.

The phasing is computed with any of the TaylorT1, T2, T3, T4 methods.

Parameters
deltaTsampling interval (s)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
f_minstarting GW frequency (Hz)
f_refreference GW frequency (Hz)
rdistance of source (m)
lambda1(tidal deformability of mass 1) / m1^5 (dimensionless)
lambda2(tidal deformability of mass 2) / m2^5 (dimensionless)
waveFlagsSet of flags to control special behavior of some waveform families. Pass in NULL (or None in python) for default flags
nonGRparamsLinked list of non-GR parameters. Pass in NULL (or None in python) for standard GR waveforms
amplitudeOtwice post-Newtonian amplitude order
phaseOtwice post-Newtonian order
ll index of mode
mm index of mode
approximantpost-Newtonian approximant to use for waveform production

Definition at line 1772 of file LALSimInspiral.c.

◆ XLALSimInspiralPNPolarizationWaveforms()

int XLALSimInspiralPNPolarizationWaveforms ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8TimeSeries V,
REAL8TimeSeries Phi,
REAL8  v0,
REAL8  m1,
REAL8  m2,
REAL8  r,
REAL8  i,
int  ampO 
)

Given time series for a binary's orbital dynamical variables, construct the waveform polarizations h+ and hx directly.

NB: Valid only for non-precessing binaries!

Implements Equations (8.8) - (8.10) of: Luc Blanchet, Guillaume Faye, Bala R. Iyer and Siddhartha Sinha, "The third post-Newtonian gravitational wave polarisations and associated spherical harmonic modes for inspiralling compact binaries in quasi-circular orbits", Class. Quant. Grav. 25 165003 (2008); arXiv:0802.1249. NB: Be sure to check arXiv:0802.1249v3, which corrects a typo!

Note however, that we do not include the constant "memory" terms

Parameters
hplus+-polarization waveform [returned]
hcrossx-polarization waveform [returned]
Vpost-Newtonian (PN) parameter
Phiorbital phase
v0tail-term gauge choice (default = 1)
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
rdistance of source (m)
iinclination of source (rad)
ampOtwice PN order of the amplitude

Definition at line 1915 of file LALSimInspiral.c.

◆ XLALSimInspiralPNPolarizationWaveformsFromModes()

int XLALSimInspiralPNPolarizationWaveformsFromModes ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8TimeSeries v,
REAL8TimeSeries phi,
REAL8  v0,
REAL8  m1,
REAL8  m2,
REAL8  r,
REAL8  i,
int  O 
)

Given time series for a binary's orbital dynamical variables, construct the waveform polarizations h+ and hx as a sum of -2 spin-weighted spherical harmonic modes, h_lm.

NB: Valid only for non-precessing systems!

Implements Equation (11) of: Lawrence E. Kidder, "Using Full Information When Computing Modes of Post-Newtonian Waveforms From Inspiralling Compact Binaries in Circular Orbit", Physical Review D 77, 044016 (2008), arXiv:0710.0614v1 [gr-qc].

Parameters
hplus+-polarization waveform [returned]
hcrossx-polarization waveform [returned]
vpost-Newtonian parameter
phiorbital phase
v0tail-term gauge choice (default = 1)
m1mass of companion 1
m2mass of companion 2
rdistance of source
iinclination of source (rad)
Otwice post-Newtonian order

Definition at line 2244 of file LALSimInspiral.c.

◆ XLALSimInspiralPolarizationsFromSphHarmTimeSeries()

int XLALSimInspiralPolarizationsFromSphHarmTimeSeries ( REAL8TimeSeries **  hp,
REAL8TimeSeries **  hc,
SphHarmTimeSeries hlms,
REAL8  iota,
REAL8  phiRef 
)

Compute the polarizations from all the -2 spin-weighted spherical harmonic modes stored in 'hlms'.

Be sure that 'hlms' is the head of the linked list!

The computation done is: \(hp(t) - i hc(t) = \sum_l \sum_m h_lm(t) -2Y_lm(iota,psi)\)

iota and psi are the inclination and polarization angle of the observer relative to the source of GWs.

Parameters
hpPlus polarization time series [returned]
hcCross polarization time series [returned]
hlmsHead of linked list of waveform modes
iotainclination of viewer to source frame (rad)
phiRefreference phase (rad)

Definition at line 2291 of file LALSimInspiral.c.

◆ XLALSimInspiralPolarizationsFromChooseFDModes()

int XLALSimInspiralPolarizationsFromChooseFDModes ( 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 UNUSED  longAscNodes,
const REAL8 UNUSED  eccentricity,
const REAL8 UNUSED  meanPerAno,
const REAL8  deltaF,
const REAL8  f_min,
const REAL8  f_max,
REAL8  f_ref,
LALDict *  LALparams,
const Approximant  approximant 
)

Function returning the Fourier domain polarizations for positive frequencies built from the individual modes computed with ChooseFDModes.

The output should be equivalent to that from ChooseFDWaveform, close to machine precision. Some of the AS models use the argument phiRef internally to build the hlms and build the polarizations with an azimuthal angle different from Pi/2 - phiRef. In this function we take into account those differences among models. For the precessing model IMRPhenomXPHM, since the modes are returned in the J-frame, the polarizations are built using theta_JN instead of the inclination and azimuthal angle = 0.

Parameters
hptildeFD plus polarization
hctildeFD cross polarization
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
S1xx-component of the dimensionless spin of object 1
S1yy-component of the dimensionless spin of object 1
S1zz-component of the dimensionless spin of object 1
S2xx-component of the dimensionless spin of object 2
S2yy-component of the dimensionless spin of object 2
S2zz-component of the dimensionless spin of object 2
distancedistance of source (m)
inclinationinclination of source (rad)
phiRefreference orbital phase (rad)
longAscNodeslongitude of ascending nodes, degenerate with the polarization angle, Omega in documentation
eccentricityeccentricity at reference epoch
meanPerAnomean anomaly of periastron
deltaFsampling interval (Hz)
f_minstarting GW frequency (Hz)
f_maxending GW frequency (Hz)
f_refReference frequency (Hz)
LALparamsLAL dictionary containing accessory parameters
approximantpost-Newtonian approximant to use for waveform production

Definition at line 2329 of file LALSimInspiral.c.

◆ XLALSimInspiralPolarizationsFromSphHarmFrequencySeries()

int XLALSimInspiralPolarizationsFromSphHarmFrequencySeries ( COMPLEX16FrequencySeries **  hp,
COMPLEX16FrequencySeries **  hc,
SphHarmFrequencySeries hlms,
REAL8  theta,
REAL8  phi 
)

Return polarizations for positive frequencies built by summing the individual modes present in the input array SphHarmFrequencySeries *hlms computed with ChooseFDModes.

Notice that in general the output may not be close to machine precision with ChooseFDWaveform due to differences in computing the h_lms and in the use of the azimuthal angle. For IMRPhenomXPHM, the argument theta should not be the inclination but theta_JN and phiRef should be 0.

Parameters
hpPlus polarization frequency series [returned]
hcCross polarization frequency series [returned]
hlmsHead of linked list of waveform modes
thetaPolar angle for the Ylms (rad): Inclination for h_lms in L0-frame and theta_JN for J-frame.
phiAzimuthal angle for the Ylms (rad): pi/2 - phiRef for AS models and 0 for precessing.

Definition at line 2503 of file LALSimInspiral.c.

◆ XLALSimInspiralPNPolarizationWaveformsEccentric()

int XLALSimInspiralPNPolarizationWaveformsEccentric ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8TimeSeries V,
REAL8TimeSeries Ecc,
REAL8TimeSeries U,
REAL8TimeSeries Phi,
REAL8  m1,
REAL8  m2,
REAL8  r,
REAL8  i,
int  ampO,
int  ph_O 
)

Given time series for a binary's orbital dynamical variables, computes the radial and angular orbital motion; then constructs the waveform polarizations h+ and hx in terms of the relative separation r, the true anomaly phi and their time derivatives.

NB: Valid only for non-spinning binaries in inspiralling! and precessing eccentric orbits!

Implements Equations (3.7a) - (3.7c) and Equations (B2a) - (B2d), (B4a), (B4b) of: Sashwat Tanay, Maria Haney, and Achamveedu Gopakumar, "Frequency and time domain inspiral waveforms from comparable mass compact binaries in eccentric orbits", (2015); arXiv:TBD https://dcc.ligo.org/P1500148-v1

(note that above equations use x = v^2 as the PN expansion parameter)

as well as Equations (6a) and (6b) of: Thibault Damour, Achamveedu Gopakumar, and Bala R. Iyer, "Phasing of gravitational waves from inspiralling eccentric binaries", Phys. Rev. D 70 064028 (2004); arXiv:gr-qc/0404128.

Parameters
hplus+-polarization waveform [returned]
hcrossx-polarization waveform [returned]
Vpost-Newtonian (PN) parameter
Eccorbital eccentricity
Ueccentric anomaly of quasi-Keplerian orbit
Phiorbital phase
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
rdistance of source (m)
iinclination of source (rad)
ampOtwice PN order of the amplitude
ph_Otwice PN order of the phase

Definition at line 2573 of file LALSimInspiral.c.

◆ XLALSimInspiralPrecessingPolarizationWaveforms()

int XLALSimInspiralPrecessingPolarizationWaveforms ( REAL8TimeSeries **  hplus,
REAL8TimeSeries **  hcross,
REAL8TimeSeries V,
REAL8TimeSeries Phi,
REAL8TimeSeries S1x,
REAL8TimeSeries S1y,
REAL8TimeSeries S1z,
REAL8TimeSeries S2x,
REAL8TimeSeries S2y,
REAL8TimeSeries S2z,
REAL8TimeSeries LNhatx,
REAL8TimeSeries LNhaty,
REAL8TimeSeries LNhatz,
REAL8TimeSeries E1x,
REAL8TimeSeries E1y,
REAL8TimeSeries E1z,
REAL8  m1,
REAL8  m2,
REAL8  r,
INT4  ampO 
)

Computes polarizations h+ and hx for a spinning, precessing binary when provided time series of all the dynamical quantities.

Amplitude can be chosen between 1.5PN and Newtonian orders (inclusive).

Based on K.G. Arun, Alesssandra Buonanno, Guillaume Faye and Evan Ochsner "Higher-order spin effects in the amplitude and phase of gravitational waveforms emitted by inspiraling compact binaries: Ready-to-use gravitational waveforms", Phys Rev. D 79, 104023 (2009), arXiv:0810.5336

HOWEVER, the formulae have been adapted to use the output of the so-called "Frameless" convention for evolving precessing binary dynamics, which is not susceptible to hitting coordinate singularities.

FIXME: Clean up and commit Mathematica NB Showing correctness. Cite here.

NOTE: The vectors MUST be given in the so-called radiation frame where Z is the direction of propagation, X is the principal '+' axis and Y = Z x X For different convention (Z is the direction of initial total angular momentum, useful for GRB and comparison to NR, see XLALSimSpinInspiralGenerator())

Parameters
hplus+-polarization waveform [returned]
hcrossx-polarization waveform [returned]
Vpost-Newtonian parameter
Phiorbital phase
S1xSpin1 vector x component
S1ySpin1 vector y component
S1zSpin1 vector z component
S2xSpin2 vector x component
S2ySpin2 vector y component
S2zSpin2 vector z component
LNhatxunit orbital ang. mom. x comp.
LNhatyunit orbital ang. mom. y comp.
LNhatzunit orbital ang. mom. z comp.
E1xorbital plane basis vector x comp.
E1yorbital plane basis vector y comp.
E1zorbital plane basis vector z comp.
m1mass of companion 1 (kg)
m2mass of companion 2 (kg)
rdistance of source (m)
ampOtwice amp. post-Newtonian order

Definition at line 2889 of file LALSimInspiral.c.

◆ XLALSimInspiralPrecessingPolarizationWaveformHarmonic()

int XLALSimInspiralPrecessingPolarizationWaveformHarmonic ( COMPLEX16 hplus,
COMPLEX16 hcross,
REAL8  v,
REAL8  s1x,
REAL8  s1y,
REAL8  s1z,
REAL8  s2x,
REAL8  s2y,
REAL8  s2z,
REAL8  lnhx,
REAL8  lnhy,
REAL8  lnhz,
REAL8  e1x,
REAL8  e1y,
REAL8  e1z,
REAL8  dm,
REAL8  eta,
REAL8  v0,
INT4  n,
INT4  ampO 
)

Computes polarizations h+ and hx for a spinning, precessing binary when provided a single value of all the dynamical quantities.

Amplitude can be chosen between 1.5PN and Newtonian orders (inclusive).

Based on K.G. Arun, Alesssandra Buonanno, Guillaume Faye and Evan Ochsner "Higher-order spin effects in the amplitude and phase of gravitational waveforms emitted by inspiraling compact binaries: Ready-to-use gravitational waveforms", Phys Rev. D 79, 104023 (2009), arXiv:0810.5336

HOWEVER, the formulae have been adapted to use the output of the so-called "Frameless" convention for evolving precessing binary dynamics, which is not susceptible to hitting coordinate singularities.

This has been written to reproduce XLALSimInspiralPrecessingPolarizationWaveforms. If hplus and hcross are the output of XLALSimInspiralPrecessingPolarizationWaveforms, and hp(n) and hc(n) the output of this function for a given harmonic number, then

hplus = sum_{n=0}^5 hp(n)*exp(-i*n*Phi) + c.c. hcross = sum_{n=0}^5 hc(n)*exp(-i*n*Phi) + c.c.

NOTE: The vectors MUST be given in the so-called radiation frame where Z is the direction of propagation, X is the principal '+' axis and Y = Z x X For different convention (Z is the direction of initial total angular momentum, useful for GRB and comparison to NR, see XLALSimSpinInspiralGenerator()) FIXME take out v0 as it can be absorbed in a 4PN additional phase term see discussion in sec. VIII of Class.Quant.Grav. 25 (2008) 165003, arXiv 0802.1249.

Parameters
hplus+-polarization waveform [returned]
hcrossx-polarization waveform [returned]
vpost-Newtonian parameter
s1xSpin1 vector x component
s1ySpin1 vector y component
s1zSpin1 vector z component
s2xSpin2 vector x component
s2ySpin2 vector y component
s2zSpin2 vector z component
lnhxunit orbital ang. mom. x comp.
lnhyunit orbital ang. mom. y comp.
lnhzunit orbital ang. mom. z comp.
e1xorbital plane basis vector x comp.
e1yorbital plane basis vector y comp.
e1zorbital plane basis vector z comp.
dmdimensionless mass difference (m1 - m2)/(m1 + m2) > 0
etasymmetric mass ratio m1*m2/(m1 + m2)^2
v0tail-term gauge choice (default = 1)
nharmonic number
ampOtwice amp. post-Newtonian order

Definition at line 3182 of file LALSimInspiral.c.