Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-8a6b96f
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Header LALSimInspiralPrecess.h

Detailed Description

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.

Author
Chris Pankow

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...
 

Function Documentation

◆ XLALSimInspiralPrecessionRotateModes()

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)\)

Parameters
h_lmspherical harmonic decomposed modes, modified in place
alphaalpha Euler angle time series
betabeta Euler angle time series
gamgamma Euler angle time series

Definition at line 36 of file LALSimInspiralPrecess.c.

◆ XLALSimInspiralConstantPrecessionConeWaveformModes()

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.

Parameters
h_lm_tmp(l,m) modes, modified in place
precess_freqPrecession frequency in Hz
aOpening angle of precession cone in rads
phi_precessinitial phase in cone of L around J
alpha_0azimuth btwn center of cone and line of sight
beta_0zenith btwn center of cone and line of sight

Definition at line 95 of file LALSimInspiralPrecess.c.

◆ XLALSimInspiralConstantPrecessionConeWaveform()

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

Parameters
hpOutput precessing plus polarization
hxOutput precessing cross polarization
h_lmInput non-precessing (l,m) modes
precess_freqPrecession frequency in Hz
aOpening angle of precession cone in rads
phi_precessinitial phase in cone of L around J
alpha_0azimuth btwn center of cone and line of sight
beta_0zenith btwn center of cone and line of sight

Definition at line 216 of file LALSimInspiralPrecess.c.

◆ XLALSimInspiralOrientationMatrixForL2()

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.

Parameters
mtx3x3 orientation matrix, modified in place
h22l=2, m=2 h sample
h2m2l=2, m=-2 h sample
h21l=2, m=1 h sample
h2m1l=2, m=-1 h sample
h20l=2, m=0 h sample

Definition at line 246 of file LALSimInspiralPrecess.c.

◆ XLALSimInspiralOrientationMatrixDirection()

int XLALSimInspiralOrientationMatrixDirection ( REAL8  vec[3],
REAL8  mtx[3][3] 
)

Determine a direction vector from the orientation matrix.

NOTE: vec should be initialized to be a unit vector – its value is used, then replaced

Parameters
vec3 dim unit vector, modified in place
mtx3x3 orientation matrix

Definition at line 301 of file LALSimInspiralPrecess.c.

◆ XLALSimInspiralPrecessionRotateModesOut()

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:

  • leaves unaltered the input SphericalHarmonicTimeSeries
  • The Wigner DMatrix rotation function XLALWignerDMatrix() is called with arguments (alpha, -beta, gamma), instead of (alpha, beta, gamma), to ensure that (h+ + i hx)(alpha,beta,gamma)=Sum_{lmm'} Y_lm(0,0) D_mm'(alpha,beta,gamma) h_lm'(0,0,0)
Parameters
hlm_outspherical harmonic decomposed modes, output
hlm_inspherical harmonic decomposed modes, input
alphaalpha Euler angle time series
betabeta Euler angle time series
gamgamma Euler angle time series

Definition at line 353 of file LALSimInspiralPrecess.c.