LALSimulation  5.4.0.1-89842e6
Inspiral Simulation Packages

Detailed Description

This package provides routines for generating simulated gravitational waveforms from binary inspiral.

Various types of inspiral approximants are supported for producing waveforms in the time- or frequency-domains. The high-level routines for generating simulated inspiral waveforms are XLALSimInspiralChooseTDWaveform() (for time-domain waveforms) and XLALSimInspiralChooseFDWaveform() (for frequency-domain waveforms). The following examples show their basic usage.

Generate a time-domain waveform:

#include <lal/LALSimInspiral.h>
#include <lal/LALConstants.h>
...
double m1 = 1.4 * LAL_MSUN_SI; // mass of body 1
double m2 = 1.4 * LAL_MSUN_SI; // mass of body 2
double S1x = 0.0; // x-component of dimensionless spin of body 1
double S1y = 0.0; // y-component of dimensionless spin of body 1
double S1z = 0.0; // z-component of dimensionless spin of body 1
double S2x = 0.0; // x-component of dimensionless spin of body 2
double S2y = 0.0; // y-component of dimensionless spin of body 2
double S2z = 0.0; // z-component of dimensionless spin of body 2
double r = 1e6 * LAL_PC_SI; // distance
double inclination = 0.0; // angle between L and view direction, \iota in @image
double phiRef = 0.0; // orbital phase at reference, helf of main GW phase at reference
double longAscNodes=0.0; // longitude of ascending nodes, degenerate with the polarization angle, related to Omega in @image by Omega=longAscNodes+pi/2
double eccentricity=0.0; // eccentricity at reference epoch
double meanPerAno=0.0; // mean anomaly at reference epoch, i.e. the ratio of time passed since last periastron passage to the time interval between two periastron passages, times 2pi. Note: This is not a geometric angle that can be visualized in @image
double deltaT = 1.0/16384.0; // series sampling interval
double f_min = 40.0; // start frequency of inspiral
double f_ref = 0.0; // reference frequency: 0 means waveform end
LALDict *LALpars = NULL; // structure containing variable with default values
Approximant approximant = TaylorT2; // post-Newtonian approximant
REAL8TimeSeries *hplus = NULL; // plus polarization to be returned
REAL8TimeSeries *hcross = NULL; // cross polarization to be returned
...
XLALSimInspiralChooseTDWaveform( &hplus, &hcross, m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, r, i, phiRef, longAscNodes, eccentricity, meanPerAno, deltaT, f_min, f_ref, LALpars, approximant);
double i
Definition: bh_ringdown.c:118
#define LAL_MSUN_SI
#define LAL_PC_SI
Approximant
Enum that specifies the PN approximant to be used in computing the waveform.
@ TaylorT2
Time domain Taylor approximant in which the phase evolution is obtained by iteratively solving post-...
static const INT4 r
double f_min
Definition: unicorn.c:22
double deltaT
Definition: unicorn.c:24

Generate a frequency-domain waveform:

#include <lal/LALSimInspiral.h>
#include <lal/LALConstants.h>
...
double m1 = 1.4 * LAL_MSUN_SI; // mass of body 1
double m2 = 1.4 * LAL_MSUN_SI; // mass of body 2
double S1x = 0.0; // x-component of dimensionless spin of body 1
double S1y = 0.0; // y-component of dimensionless spin of body 1
double S1z = 0.0; // z-component of dimensionless spin of body 1
double S2x = 0.0; // x-component of dimensionless spin of body 2
double S2y = 0.0; // y-component of dimensionless spin of body 2
double S2z = 0.0; // z-component of dimensionless spin of body 2
double distance = 1e6 * LAL_PC_SI; // distance
double inclination = 0.0; // angle between L and view direction, \iota in @image
double phiRef = 0; // orbital phase at reference, half of main GW phase at reference
double longAscNodes=0.0; // longitude of ascending nodes, degenerate with the polarization angle, related to Omega in @image by Omega=longAscNodes+pi/2
double eccentricity=0.0; // eccentricity at reference epoch
double meanPerAno=0.0; // mean anomaly at reference epoch, i.e. the ratio of time passed since last periastron passage to the time interval between two periastron passages, times 2pi. Note: This is not a geometric angle that can be visualized in @image
double deltaF = 1.; // frequency sampling interval
double f_min = 40.0; // start frequency of inspiral
double f_max = 0.0; // end frequency of inspiral: 0 means use default
double f_ref = 0.0; // reference frequency: 0 means waveform end
LALDict *LALpars = NULL; // structure containing variable with default values
Approximant approximant = TaylorF2; // post-Newtonian approximant
COMPLEX16FrequencySeries *hptilde = NULL; // plus polarization to be returned
COMPLEX16FrequencySeries *hctilde = NULL; // cross polarization to be returned
...
XLALSimInspiralChooseFDWaveform(&hptilde, &hctilde, m1, m2, S1x, S1y, S1z, S2x, S2y, S2z, f_min, f_max, f_ref, r, i, phiRef, longAscNodes, eccentricity, meanPerAno, deltaF, f_min, f_ref, LALpars, approximant);
@ TaylorF2
The standard stationary phase approximation; Outputs a frequency-domain wave.
double f_max
Definition: unicorn.c:23

Coordinate Systems

The diagram below illustrates how the source frame (x,y,z) of the binary is related to the wave frame (X,Y,Z) in which the gravitational waveform is defined.

Orbital Elements

The origin of the coordinate systems is the instantaneous center-of-mass of the binary system. The orbiting body shown in the diagram is body 1.

The binary's instantaneous orbital angular momentum L at the reference gravitational wave frequency f_ref defines the z-axis of the binary system, while the unit vector from body 2 to body 1 defines the x-axis of the binary system. The x-y-plane is therefore the orbital plane, at least at the moment the binary system is at the reference gravitational wave frequency.

The spin components for body 1, (S1x,S1y, S1z), and for body 2, (S2x,S2y, S2z), are defined in the source-frame. Therefore, when the spins are aligned with the orbital angular momentum, S1x = S1y = S2x = S2y = 0.

Note
The spin components transverse to the orbital angular momentum L at the reference gravitational wave frequency f_ref are given with respect to the triad x-y-z, with x-axis parallel to the vector pointing from body 2 to body 1.

The wave frame is defined by the Z-axis, which points toward the Earth, and some reference direction, defining the X-axis. The X-Y-plane is therefore the plane of the sky.

The plus- and cross-polarizations of the gravitational waveform are defined in this wave frame. Specifically, if \( h^{ij} \) is computed in the source frame, then

\[ h_+ = \frac12 ( \hat{P}_i \hat{P}_j - \hat{Q}_i \hat{Q}_j ) h^{ij} \]

and

\[ h_\times = \frac12 ( \hat{P}_i \hat{Q}_j + \hat{Q}_i \hat{P}_j ) h^{ij} \]

where \( \hat{P}_i \) are the components of the unit vector pointing along the X-axis and \( \hat{Q}_i \) are the components of the unit vector pointing along the Y-axis.

The orbital elements are:

See also
The coordinate systems used here follows those of

L. Blanchet, G. Faye, B. R. Iyer and S. 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, (2008) 165003 Erratum: [Class. Quant. Grav. 29, (2012) 239501, arXiv:0802.1249 [gr-qc].

Modules

 Header LALSimIMR.h
 Routines for generating inspiral-merger-ringdown waveforms.
 
 Header LALSimInspiral.h
 Routines for generating binary inspiral gravitational waveforms.
 
 Header LALSimInspiralPrecess.h
 Functions to take an arbitrary waveform time series and impose the effects of causing the viewing angle to precess about a cone of L around J. The cone currently has a constant opening angle.
 
 Header LALSimInspiralWaveformCache.h
 Routines for saving previously-computed waveforms for reuse.
 
 Header LALSimBlackHoleRingdown.h
 Routines to generate black hole ringdown waveforms.