Functions for working with SSB times.
Prototypes | |
int | XLALAddBinaryTimes (SSBtimes **tSSBOut, const SSBtimes *tSSBIn, const PulsarDopplerParams *Doppler) |
Compute extra time-delays for a CW source in a (Keplerian) binary orbital system. More... | |
static double | gsl_E_solver (REAL8 E, void *par) |
Function implementing Eq. More... | |
int | XLALAddMultiBinaryTimes (MultiSSBtimes **multiSSBOut, const MultiSSBtimes *multiSSBIn, const PulsarDopplerParams *Doppler) |
Multi-IFO version of XLALAddBinaryTimes(). More... | |
SSBtimes * | XLALDuplicateSSBtimes (const SSBtimes *tSSB) |
Duplicate (ie allocate + copy) an input SSBtimes structure. More... | |
MultiSSBtimes * | XLALDuplicateMultiSSBtimes (const MultiSSBtimes *multiSSB) |
Duplicate (ie allocate + copy) an input MultiSSBtimes structure. More... | |
SSBtimes * | XLALGetSSBtimes (const DetectorStateSeries *DetectorStates, SkyPosition pos, LIGOTimeGPS refTime, SSBprecision precision) |
For a given DetectorStateSeries, calculate the time-differences \( \Delta T_\alpha\equiv T(t_\alpha) - T_0 \) , and their derivatives \( \dot{T}_\alpha \equiv d T / d t (t_\alpha) \) . More... | |
MultiSSBtimes * | XLALGetMultiSSBtimes (const MultiDetectorStateSeries *multiDetStates, SkyPosition skypos, LIGOTimeGPS refTime, SSBprecision precision) |
Multi-IFO version of XLALGetSSBtimes(). More... | |
int | XLALEarliestMultiSSBtime (LIGOTimeGPS *out, const MultiSSBtimes *multiSSB, const REAL8 Tsft) |
Find the earliest timestamp in a multi-SSB data structure. More... | |
int | XLALLatestMultiSSBtime (LIGOTimeGPS *out, const MultiSSBtimes *multiSSB, const REAL8 Tsft) |
Find the latest timestamp in a multi-SSB data structure. More... | |
void | XLALDestroySSBtimes (SSBtimes *tSSB) |
Destroy a SSBtimes structure. More... | |
void | XLALDestroyMultiSSBtimes (MultiSSBtimes *multiSSB) |
Destroy a MultiSSBtimes structure. More... | |
Data Structures | |
struct | SSBtimes |
Simple container for two REAL8-vectors, namely the SSB-timings DeltaT_alpha and Tdot_alpha, with one entry per SFT-timestamp. More... | |
struct | MultiSSBtimes |
Multi-IFO container for SSB timings. More... | |
Enumerations | |
enum | SSBprecision { SSBPREC_NEWTONIAN , SSBPREC_RELATIVISTIC , SSBPREC_RELATIVISTICOPT , SSBPREC_DMOFF , SSBPREC_LAST } |
The precision in calculating the barycentric transformation. More... | |
Variables | |
const UserChoices | SSBprecisionChoices |
Static array of all SSBprecision choices, for use by the UserInput module parsing routines. More... | |
int XLALAddBinaryTimes | ( | SSBtimes ** | tSSBOut, |
const SSBtimes * | tSSBIn, | ||
const PulsarDopplerParams * | Doppler | ||
) |
Compute extra time-delays for a CW source in a (Keplerian) binary orbital system.
\[ \newcommand{\Det}{{\mathrm{det}}} \newcommand{\tDet}{t_\Det} \newcommand{\tSSB}{t_{\mathrm{SSB}}} \newcommand{\tEM}{t_{\mathrm{em}}} \newcommand{\tRef}{t_{\mathrm{ref}}} \newcommand{\tPeri}{t_{\mathrm{p}}} \newcommand{\TPeri}{T_{\mathrm{p}}} \newcommand{\argp}{\omega} \newcommand{\sini}{\sin i} \newcommand{\fdot}{\dot{f}} \]
For given binary-orbital parameters in Doppler, compute the source-frame emission times \( \tEM(\tDet) \) of CW wavefronts arriving at the detector frame at times \( \tDet \) . The input to this function should contain the wavefront arrival times \( \tSSB(\tDet) \) in the solar-system barycenter (SSB) frame, the additional time differences due to the binary orbital motion are then added by this function.
NOTE: the output vector tSSBOut can be passed either as
(*tSSBOut)==NULL
): it gets allocated here, orNOTE2: it is OK to pass an identical in- and output-vector, i.e. (*tSSBOut) = tSSBIn
: in this case the binary-orbital time offsets will be added, modifying the input.
NOTE3: (FIXME!) the SSBtimes structure naming is very misleading: 'SSB' here really refers to the source emission time, not generally the solar-system-barycenter time, so this naming is only correct in the case of isolated NSs. The additional time-delay computed by this function is between the inertial-SSB time and the source emission time, so this is quite confusing right now and the corresponding variable-names and documentation should be fixed ASAP.
Notation: (all times are referring to a global Newtonian time axis):
\begin{align*} & \tDet \ldots \text{wave-front arrival time at detector}\\ & \tSSB \ldots \text{arrival time at SSB}\\ & \tEM \ldots \text{wave-front emission time at the CW source}\\ & \tRef \ldots \text{reference time at which} \{\phi_0, f, \fdot, \ldots\} \text{are defined} \end{align*}
The intrinsic CW phase in the source frame can be written as
\begin{equation} \label{eq:sourcePhase} \begin{split} \Phi(\tEM) &= \phi_0 + 2 \pi \left[ f \, \Delta\tau + \frac{1}{2} \fdot \, \Delta\tau^2 + \ldots \right]\,,\\ \Delta\tau &\equiv \tEM - \tRef\,. \end{split} \end{equation}
In order to relate this to the CW phase \( \Phi_\Det \) arriving at the detector, we need the timing relation \( \tEM(\tDet) \) , such that \( \Phi_\Det(\tDet) = \Phi\left(\tEM(\tDet)\right) \) . We also compute the explicit values of the time-derivative, i.e. \( \frac{d\tEM}{d \tDet}(\tDet) \) , which are required by XLALComputeFaFb(), for example.
XLALGetSSBtimes() computes the timing relation between detector time \( \tDet \) and SSB time \( \tSSB \) , namely
\begin{align} \Delta\tau_0(\tDet) &\equiv \tSSB(\tDet) - \tRef\,,\\ \dot{\tau}_0(\tDet) &\equiv \frac{d\tSSB}{d \tDet}(\tDet)\, \end{align}
which is passed to this function as an input in tSSBIn. We therefore need to compute the extra time delay due to the binary-orbital motion, namely
\begin{align} \Delta\tau(\tDet) &\equiv \tEM(\tDet) - \tRef = [\tEM(\tDet) - \tSSB(\tDet)] + [\tSSB(\tDet) - \tRef] \notag\\ &= [\tEM - \tSSB](\tDet) + \Delta\tau_0(\tDet) \,, \label{eq:7a}\\ \dot{\tau}(\tDet) &\equiv \frac{d\tEM}{d \tDet}(\tDet) = \frac{d\tEM}{d\tSSB}\left(\tSSB(\tDet)\right) \, \frac{d\tSSB}{d\tDet}(\tDet) \notag\\ &= \left[ \frac{d\tEM}{d\tSSB}\left(\tSSB(\tDet)\right) \right]\,\dot{\tau}_0(\tDet)\,, \label{eq:7b} \end{align}
The relation between \( \tEM \) and \( \tSSB \) contains an unknown offset due to the distance of the binary system from the SSB, which can be either constant, or changing at a constant or variable rate (eg due to accelerations in globular clusters). However, we can absorb this unknown effect into the definition of the pulsar parameters \( \{\phi_0,f,\fdot,\ldots\} \) , which are either unknown anyway, or affected by the exact same re-definitions if they are known from photon astronomy. We therefore ignore this effect and pretend that there is no time delay between the SSB and the binary-system barycenter (BSB), ie effectively \( \tSSB = t_{\mathrm{BSB}} \) .
The extra time-delay from the binary orbital motion can be written as
\begin{equation} \label{eq:1} \tSSB(\tEM) = \tEM + R(\tEM)\,, \end{equation}
where \( R \) is the radial distance from the BSB to the emitting NS along the line of sight, and here and in the following we are using units where \( c=1 \) . The sign convention is such that \( R>0 \) means that the NS is further away than the BSB, when \( R<0 \) it is closer.
In terms of orbital parameters, the radial distance \( R \) can be expressed as (see this figure).
\begin{equation} \label{eq:2} R = r\,\sini\,\sin(\argp + \upsilon)\,, \end{equation}
where \( r \) is the distance of the NS from the BSB (ie the focus of the ellips), \( i \) is the inclination angle between the orbital plane and the sky, \( \argp \) is the argument of periapse, and \( \upsilon \) is the true anomaly (ie the angle from the periapse to the current NS location around the BSB.
Using elementary trigonometry (cf. this figure and https://en.wikipedia.org/wiki/Eccentric_anomaly), one can see that the elliptical orbit can be described in terms of the true anomaly \( \upsilon \) as
\begin{equation} \label{eq:3} r(\upsilon) = \frac{a\,(1- e^2)}{1 + e\cos\upsilon}\,, \end{equation}
and in terms of the eccentric anomaly \( E \) as
\begin{equation} \label{eq:4} r(E) = a \, \left( 1 - e\,\cos E\right)\,. \end{equation}
From Eq. \eqref{eq:3} and Eq. \eqref{eq:4} we easily obtain the relations
\begin{equation} \begin{split} \cos\upsilon &= \frac{\cos E - e}{1 - e\cos E}\,,\\ \sin\upsilon &= \sin E \frac{\sqrt{1 - e^2}}{1 - e\cos E}\,, \end{split} \label{eq:9} \end{equation}
where in the second equation we have used the fact that \( \sin E \) and \( \sin\upsilon \) always have the same sign, as can be seen from this figure.
The (Keplerian) motion of the NS on this elliptical orbit is described by Kepler's equation:
\begin{equation} \label{eq:6} \tEM - \tPeri = \frac{P}{2\pi}\left( E - e\,\sin E\right)\,, \end{equation}
where we fixed the (discrete) gauge as \( E=0 \) at \( \tEM=\tPeri \) .
Following Teviet's strategy explained in LALGenerateEllipticSpinOrbitCW() we proceed by first expressing \( \tSSB(E) \) , solving it numerically for \( E(\tSSB) \) , from which we obtain by substitution \( \tEM(\tSSB) \) and \( d\tEM/d\tSSB \) required for Eq. \eqref{eq:7a}, Eq. \eqref{eq:7b}.
We start by writing Eq. \eqref{eq:1} using Eq. \eqref{eq:6} as
\begin{equation} \label{eq:8} \tSSB(E) = \tPeri + \frac{P}{2\pi}\left( E - e\,\sin E\right) + R(E)\,, \end{equation}
and using Eq. \eqref{eq:2}, Eq. \eqref{eq:3} with Eq. \eqref{eq:9}, we obtain
\begin{equation} \label{eq:R_E} R(E) = a\sini\left[ \sin\argp ( \cos E - e) + \cos\argp\,\sin E \sqrt{1 - e^2} \right]\,, \end{equation}
Before solving Eq. \eqref{eq:8} for \( E(\tSSB) \) , it is useful to rewrite it a bit further: First we note that at periapse (i.e. \( E=0 \) ),
\begin{equation} \TPeri \equiv \tSSB(E=0) = \tPeri + a\sini \,\sin\argp\,(1-e)\,, \label{eq:defTPeri} \end{equation}
which corresponds to the observed (SSB) time of periapse passage.
Furthermore, as pointed out in LALGenerateEllipticSpinOrbitCW(), the (Newtonian) binary NS timing is periodic, namely \( \tSSB(E + m\,2\pi) = \tSSB(E) + m\,P \) for integer \( m \) , and we can therefore simplify the numerical solution of Eq. \eqref{eq:8} for \( E(\tSSB) \) by restricting the solution to the interval \( E\in[0,2\pi) \) by mapping \( \tSSB \) into \( [\TPeri,\,\TPeri + P) \) , via
\begin{equation} x_0 \equiv \frac{2\pi}{P} (\tSSB - \TPeri) \mod 2\pi\,, \end{equation}
where we add \( 2\pi \) if \( x_0 < 0 \) to ensure that \( x_0 \in [0,\,2\pi) \) . We can therefore rewrite Eq. \eqref{eq:8} as
\begin{align} x_0 &= E + A\,\sin E + B \,( \cos E -1 )\,, \label{eq:E_tSSB}\\ A &\equiv \frac{2\pi}{P}\, a\sini \,\cos\argp\,\sqrt{1-e^2} - e\,, \label{eq:defA}\\ B &\equiv \frac{2\pi}{P}\, a\sini \,\sin\argp\,, \label{eq:defB} \end{align}
which (after some substitutions) can be seen to agree with Teviet's equation found in LALGenerateEllipticSpinOrbitCW().
This is solved numerically for \( E(\tSSB) \in [0,\,2\pi) \) , and plugging this back into Eq. \eqref{eq:R_E} and Eq. \eqref{eq:1} we obtain \( \tEM(E) \) , and from this the required timing relation Eq. \eqref{eq:7a}.
We start from Eq. \eqref{eq:1} and write,
\begin{equation} \label{eq:12} \frac{d\tEM}{d\tSSB} = \left[\frac{d\tSSB}{d\tEM}\right]^{-1} = \left[1 + \frac{dR}{d E}\,\frac{d E}{d\tEM}\right]^{-1}\,. \end{equation}
From Eq. \eqref{eq:6} we obtain
\begin{equation} \label{eq:11} \frac{d E}{d\tEM} = \frac{2\pi}{P} \frac{1}{1 - e\cos E}\,, \end{equation}
and Eq. \eqref{eq:R_E} yields
\begin{equation} \label{eq:13} \begin{split} \frac{dR}{d E} &= a\sini\left[ - \sin E \, \sin\argp + \cos E\,\cos\argp\,\sqrt{1- e^2}\right]\\ &= \frac{P}{2\pi}\left( (A + e)\,\cos E - B\,\sin E \right)\,, \end{split} \end{equation}
which results in the derivative of Eq. \eqref{eq:12} to be expressible as
\begin{equation} \label{eq:dtEM_dtSSB} \frac{d\tEM}{d\tSSB} = \frac{1 - e\cos E}{1 + A\,\cos E - B \,\sin E}\,. \end{equation}
[out] | tSSBOut | reference-time offsets in emission frame: \( \tEM(\tDet)-\tRef \) and \( d\tEM/d\tDet \) |
[in] | tSSBIn | reference-time offsets in SSB frame: \( \tSSB(\tDet)-\tRef \) and \( d\tSSB/d\tDet \) |
[in] | Doppler | pulsar Doppler parameters, includes binary orbit parameters */ |
Definition at line 259 of file SSBtimes.c.
|
static |
Function implementing Eq.
\eqref{eq:E_tSSB} to be solved via numerical root-finder for \( E(\tSSB) \)
Definition at line 378 of file SSBtimes.c.
int XLALAddMultiBinaryTimes | ( | MultiSSBtimes ** | multiSSBOut, |
const MultiSSBtimes * | multiSSBIn, | ||
const PulsarDopplerParams * | Doppler | ||
) |
Multi-IFO version of XLALAddBinaryTimes().
For given binary-orbital parameters in Doppler, compute the source-frame emission times \( \tEM(\tDet) \) of CW wavefronts arriving at the detector frame at times \( \tDet \) . The input to this function should contain the wavefront arrival times \( \tSSB(\tDet) \) in the solar-system barycenter (SSB) frame, the additional time differences due to the binary orbital motion are then added by this function.
NOTE: the output vector multiSSBOut can be passed either as
(*tSSBOut)==NULL
): it gets allocated here, orNOTE2: it is OK to pass an identical in- and output-vector, i.e. (*tSSBOut) = tSSBIn
: in this case the binary-orbital time offsets will be added, modifying the input.
NOTE3: (FIXME!) the SSBtimes structure naming is very misleading: 'SSB' here really refers to the source emission time, not generally the solar-system-barycenter time, so this naming is only correct in the case of isolated NSs. The additional time-delay computed by this function is between the inertial-SSB time and the source emission time, so this is quite confusing right now and the corresponding variable-names and documentation should be fixed ASAP.
[out] | multiSSBOut | output SSB times |
[in] | multiSSBIn | SSB-timings for all input detector-state series |
[in] | Doppler | pulsar Doppler parameters, includes binary orbit parameters |
Definition at line 413 of file SSBtimes.c.
Duplicate (ie allocate + copy) an input SSBtimes structure.
This can be useful for creating a copy before adding binary-orbital corrections in XLALAddBinaryTimes()
Definition at line 454 of file SSBtimes.c.
MultiSSBtimes * XLALDuplicateMultiSSBtimes | ( | const MultiSSBtimes * | multiSSB | ) |
Duplicate (ie allocate + copy) an input MultiSSBtimes structure.
Definition at line 487 of file SSBtimes.c.
SSBtimes * XLALGetSSBtimes | ( | const DetectorStateSeries * | DetectorStates, |
SkyPosition | pos, | ||
LIGOTimeGPS | refTime, | ||
SSBprecision | precision | ||
) |
For a given DetectorStateSeries, calculate the time-differences \( \Delta T_\alpha\equiv T(t_\alpha) - T_0 \) , and their derivatives \( \dot{T}_\alpha \equiv d T / d t (t_\alpha) \) .
[in] | DetectorStates | detector-states at timestamps t_i |
pos | source sky-location | |
refTime | SSB reference-time T_0 of pulsar-parameters | |
precision | relativistic or Newtonian SSB transformation? |
Definition at line 518 of file SSBtimes.c.
MultiSSBtimes * XLALGetMultiSSBtimes | ( | const MultiDetectorStateSeries * | multiDetStates, |
SkyPosition | skypos, | ||
LIGOTimeGPS | refTime, | ||
SSBprecision | precision | ||
) |
Multi-IFO version of XLALGetSSBtimes().
Get all SSB-timings for all input detector-series.
NOTE: this functions allocates the output-vector, use XLALDestroyMultiSSBtimes() to free this.
[in] | multiDetStates | detector-states at timestamps t_i |
skypos | source sky-position [in equatorial coords!] | |
refTime | SSB reference-time T_0 for SSB-timing | |
precision | use relativistic or Newtonian SSB timing? |
Definition at line 657 of file SSBtimes.c.
int XLALEarliestMultiSSBtime | ( | LIGOTimeGPS * | out, |
const MultiSSBtimes * | multiSSB, | ||
const REAL8 | Tsft | ||
) |
Find the earliest timestamp in a multi-SSB data structure.
out | output earliest GPS time |
multiSSB | input multi SSB SFT-midpoint timestamps |
Tsft | the length of an SFT |
Definition at line 691 of file SSBtimes.c.
int XLALLatestMultiSSBtime | ( | LIGOTimeGPS * | out, |
const MultiSSBtimes * | multiSSB, | ||
const REAL8 | Tsft | ||
) |
Find the latest timestamp in a multi-SSB data structure.
out | output latest GPS time |
multiSSB | input multi SSB SFT-midpoint timestamps |
Tsft | the length of an SFT |
Definition at line 754 of file SSBtimes.c.
void XLALDestroySSBtimes | ( | SSBtimes * | tSSB | ) |
Destroy a SSBtimes structure.
Note, this is "NULL-robust" in the sense that it will not crash on NULL-entries anywhere in this struct, so it can be used for failure-cleanup even on incomplete structs
Definition at line 822 of file SSBtimes.c.
void XLALDestroyMultiSSBtimes | ( | MultiSSBtimes * | multiSSB | ) |
Destroy a MultiSSBtimes structure.
Note, this is "NULL-robust" in the sense that it will not crash on NULL-entries anywhere in this struct, so it can be used for failure-cleanup even on incomplete structs
Definition at line 845 of file SSBtimes.c.
enum SSBprecision |
The precision in calculating the barycentric transformation.
Enumerator | |
---|---|
SSBPREC_NEWTONIAN | simple Newtonian: \( \tau = t + \vec{r}\cdot\vec{n}/c \) |
SSBPREC_RELATIVISTIC | detailed relativistic: \( \tau=\tau(t; \vec{n}, \vec{r}) \) |
SSBPREC_RELATIVISTICOPT | optimized relativistic, numerically equivalent to SSBPREC_RELATIVISTIC, but faster |
SSBPREC_DMOFF | switch off all demodulatoin terms |
SSBPREC_LAST | end marker |
Definition at line 45 of file SSBtimes.h.
|
extern |
Static array of all SSBprecision
choices, for use by the UserInput module parsing routines.
Definition at line 41 of file SSBtimes.c.