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.
Prototypes | |
int | XLALSimInspiralPrecessionRotateModes (SphHarmTimeSeries *h_lm, REAL8TimeSeries *alpha, REAL8TimeSeries *beta, REAL8TimeSeries *gam) |
Takes in the h_lm spherical harmonic decomposed modes and rotates the modes by Euler angles alpha, beta, and gamma using the Wigner D matrices. More... | |
int | XLALSimInspiralConstantPrecessionConeWaveformModes (SphHarmTimeSeries **h_lm_tmp, double precess_freq, double a, double phi_precess, double alpha_0, double beta_0) |
Takes in the l=2, abs(m)=2 decomposed modes as a strain time series and imposes the effect of a constant cone of precession. More... | |
int | XLALSimInspiralConstantPrecessionConeWaveform (REAL8TimeSeries **hp, REAL8TimeSeries **hx, SphHarmTimeSeries *h_lm, double precess_freq, double a, double phi_precess, double alpha_0, double beta_0) |
Takes in the spherical harmonic decomposed modes as a strain time series and imposes the effect of a constant cone of precession. More... | |
int | XLALSimInspiralOrientationMatrixForL2 (REAL8 mtx[3][3], COMPLEX16 h22, COMPLEX16 h2m2, COMPLEX16 h21, COMPLEX16 h2m1, COMPLEX16 h20) |
Determine the `‘preferred direction’' matrix for the l=2 subspace from the relevant mode decomposed h. More... | |
int | XLALSimInspiralOrientationMatrixDirection (REAL8 vec[3], REAL8 mtx[3][3]) |
Determine a direction vector from the orientation matrix. More... | |
int | XLALSimInspiralPrecessionRotateModesOut (SphHarmTimeSeries **hlm_out, SphHarmTimeSeries *hlm_in, const REAL8TimeSeries *alpha, const REAL8TimeSeries *beta, const REAL8TimeSeries *gam) |
Takes in the h_lm spherical harmonic decomposed modes and rotates the modes by Euler angles alpha, beta, and gamma using the Wigner D matrices, analog to int XLALSimInspiralPrecessionRotateModes but with 2 crucial differences: More... | |
int XLALSimInspiralPrecessionRotateModes | ( | SphHarmTimeSeries * | h_lm, |
REAL8TimeSeries * | alpha, | ||
REAL8TimeSeries * | beta, | ||
REAL8TimeSeries * | gam | ||
) |
Takes in the h_lm spherical harmonic decomposed modes and rotates the modes by Euler angles alpha, beta, and gamma using the Wigner D matrices.
e.g.
\(\tilde{h}_{l,m}(t) = D^l_{m',m} h_{l,m'}(t)\)
h_lm | spherical harmonic decomposed modes, modified in place |
alpha | alpha Euler angle time series |
beta | beta Euler angle time series |
gam | gamma Euler angle time series |
Definition at line 36 of file LALSimInspiralPrecess.c.
int XLALSimInspiralConstantPrecessionConeWaveformModes | ( | SphHarmTimeSeries ** | h_lm_tmp, |
double | precess_freq, | ||
double | a, | ||
double | phi_precess, | ||
double | alpha_0, | ||
double | beta_0 | ||
) |
Takes in the l=2, abs(m)=2 decomposed modes as a strain time series and imposes the effect of a constant cone of precession.
This is accomplished by taking the axis of the binary rotational plane and rotating the Y_lm such that it appears to "precess" around a fixed J direction. Note that h_2_2 and h_22 are modified in place.
Future revisions will change the first two pointers to this: COMPLEX16TimeSeries** h_lm
and add unsigned int l
Thus the input h_lm will be considered a list of pointers to the h_lm for a given l and the appropriate action will be taken for all of the submodes.
h_lm_tmp | (l,m) modes, modified in place |
precess_freq | Precession frequency in Hz |
a | Opening angle of precession cone in rads |
phi_precess | initial phase in cone of L around J |
alpha_0 | azimuth btwn center of cone and line of sight |
beta_0 | zenith btwn center of cone and line of sight |
Definition at line 95 of file LALSimInspiralPrecess.c.
int XLALSimInspiralConstantPrecessionConeWaveform | ( | REAL8TimeSeries ** | hp, |
REAL8TimeSeries ** | hx, | ||
SphHarmTimeSeries * | h_lm, | ||
double | precess_freq, | ||
double | a, | ||
double | phi_precess, | ||
double | alpha_0, | ||
double | beta_0 | ||
) |
Takes in the spherical harmonic decomposed modes as a strain time series and imposes the effect of a constant cone of precession.
The result is returned in the physical waveforms hp, hx, after they have been resummed from the modified h_lm waveforms.
NOTE: the h_lm modes will be modified in place
hp | Output precessing plus polarization |
hx | Output precessing cross polarization |
h_lm | Input non-precessing (l,m) modes |
precess_freq | Precession frequency in Hz |
a | Opening angle of precession cone in rads |
phi_precess | initial phase in cone of L around J |
alpha_0 | azimuth btwn center of cone and line of sight |
beta_0 | zenith btwn center of cone and line of sight |
Definition at line 216 of file LALSimInspiralPrecess.c.
int XLALSimInspiralOrientationMatrixForL2 | ( | REAL8 | mtx[3][3], |
COMPLEX16 | h22, | ||
COMPLEX16 | h2m2, | ||
COMPLEX16 | h21, | ||
COMPLEX16 | h2m1, | ||
COMPLEX16 | h20 | ||
) |
Determine the `‘preferred direction’' matrix for the l=2 subspace from the relevant mode decomposed h.
Modifies mtx in place.
mtx | 3x3 orientation matrix, modified in place |
h22 | l=2, m=2 h sample |
h2m2 | l=2, m=-2 h sample |
h21 | l=2, m=1 h sample |
h2m1 | l=2, m=-1 h sample |
h20 | l=2, m=0 h sample |
Definition at line 246 of file LALSimInspiralPrecess.c.
Determine a direction vector from the orientation matrix.
NOTE: vec should be initialized to be a unit vector – its value is used, then replaced
vec | 3 dim unit vector, modified in place |
mtx | 3x3 orientation matrix |
Definition at line 301 of file LALSimInspiralPrecess.c.
int XLALSimInspiralPrecessionRotateModesOut | ( | SphHarmTimeSeries ** | hlm_out, |
SphHarmTimeSeries * | hlm_in, | ||
const REAL8TimeSeries * | alpha, | ||
const REAL8TimeSeries * | beta, | ||
const REAL8TimeSeries * | gam | ||
) |
Takes in the h_lm spherical harmonic decomposed modes and rotates the modes by Euler angles alpha, beta, and gamma using the Wigner D matrices, analog to int XLALSimInspiralPrecessionRotateModes but with 2 crucial differences:
hlm_out | spherical harmonic decomposed modes, output |
hlm_in | spherical harmonic decomposed modes, input |
alpha | alpha Euler angle time series |
beta | beta Euler angle time series |
gam | gamma Euler angle time series |
Definition at line 353 of file LALSimInspiralPrecess.c.