Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
HeterodynedPulsarModel.c
Go to the documentation of this file.
1/*
2 * Copyright (C) 2018 Matthew Pitkin
3 *
4 * This program is free software; you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation; either version 2 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with with program; see the file COPYING. If not, write to the
16 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17 * MA 02110-1301 USA
18 */
19
20#include <lal/HeterodynedPulsarModel.h>
21
22/** Macro to round a value to the nearest integer. */
23#define ROUND(x) (floor(x+0.5))
24
25/** Macro to square a value. */
26#define SQUARE(x) ( (x) * (x) )
27
28/**
29* \brief The phase evolution difference compared to a heterodyned phase (for a pulsar)
30*
31* This function will calculate the difference in the phase evolution of a
32* source at a particular sky location as observed in a detector on Earth
33* compared with that used to heterodyne (or SpectralInterpolate) the data. If
34* the phase evolution is described by a Taylor expansion:
35* \f[
36* \phi(T) = \sum_{k=1}^n \frac{f^{(k-1)}}{k!} T^k,
37* \f]
38* where \f$ f^{(x)} \f$ is the xth time derivative of the gravitational-wave
39* frequency, and \f$ T \f$ is the pulsar proper time, then the phase difference
40* is given by
41* \f[
42* \Delta\phi(t) = \sum_{k=1}^n \left( \frac{\Delta f^{(k-1)}}{k!}(t+\delta t_1)^k + \frac{f^{(k-1)}_2}{k!} \sum_{i=0}^{i<k} \left(\begin{array}{c}k \\ i\end{array}\right) (\Delta t)^{k-i} (t+\delta t_1)^i \right),
43* \f]
44* where \f$ t \f$ is the signal arrival time at the detector minus the given
45* pulsar period epoch, \f$ \delta t_1 \f$ is the barycentring time delay (from
46* both solar system and binary orbital effects) calculated at the heterodyned
47* values, \f$ \Delta f^{(x)} = f_2^{(x)}-f1^{(x)} \f$ is the diffence in
48* frequency (derivative) between the current value ( \f$ f_2^{(x)} \f$ ) and the
49* heterodyne value ( \f$ f_1^{(x)} \f$ ), and
50* \f$ \Delta t = \delta t_2 - \delta t_1 \f$ is the difference between the
51* barycentring time delay calculated at the current values ( \f$ \delta t_1 \f$ )
52* and the heterodyned values. Frequency time derivatives up to any order are
53* allowed. The pulsar proper time is calculated by correcting the time of
54* arrival at Earth, \f$ t \f$ to the solar system barycentre and if necessary
55* the binary system barycenter, so
56* \f$ T = t + \delta{}t_{\rm SSB} + \delta{}t_{\rm BSB} \f$ .
57*
58* In this function the time delay needed to correct to the solar system
59* barycenter is only calculated if required, i.e., if an update is required
60* due to a change in the sky position. The same is true for the binary system
61* time delay, which is only calculated if it needs updating due to a change in
62* the binary system parameters. The function will also calculate any phase
63* evolution due to pulsar glitch parameters and/or \c FITWAVES parameters
64* (parameters that define sinusoidal terms used to represent strong red noise
65* in the signal) if required.
66*
67* This function can just return the phase evolution (and not phase evolution
68* difference) if \c ssbdts, \c bsbdts, and \c origparams are NULL, which will
69* be based on the parameters set in the \c params structure.
70*
71* \param params [in] A set of pulsar parameters
72* \param origparams [in] The original parameters used to heterodyne the data
73* \param datatimes [in] A vector of GPS times at which to calculate the phase difference
74* \param freqfactor [in] The multiplicative factor on the pulsar frequency for a particular model
75* \param ssbdts [in] The vector of SSB time delays used for the original heterodyne. If this is
76* \c NULL then the SSB delay for the given position will be calculated.
77* \param calcSSBDelay [in] Set to a non-zero value if the SSB delay needs to be calculated. This
78* would be required if using an updated sky position compared to that used for the heterodyne, or
79* if calculating the full phase evolution rather than a phase difference.
80* \param bsbdts [in] The vector of BSB time delays used for the original heterodyne. If this is
81* \c NULL then the BSB delay for the given binary parameters will be calculated.
82* \param calcBSBDelay [in] Set to a non-zero value if the BSB delay needs to be calculated. This
83* would be required if using updated binary parameters compared to that used for the heterodyne, or
84* if calculating the full phase evolution rather than a phase difference.
85* \param glphase [in] The vector of the glitch phase evolution (in EM cycles) used in the original
86* heterodyne. If this is \c NULL then the glitch phase for the given glitch parameters will be
87* calculated.
88* \param calcglphase [in] Set to a non-zero value if the glitch phase needs to be calculated. This
89* would be required if using a set of updated glitch parameters compared to that used for the
90* heterodyne, or if calculating the full phase evolution rather than a phase difference.
91* \param fitwavesphase [in] The vector of the phase evolution (in EM cycles) for any FITWAVES
92* (sinusoids used to fit red noise) parameters used in the original heterodyne. If this is \c NULL
93* then the FITWAVES phase for the given parameters will be calculated.
94* \param calcfitwaves [in] Set to a non-zero value if the FITWAVES phase needs to be calculated.
95* This would be required if using an updated set of FITWAVES parameters, or if calculating the
96* full phase evolution rather than a phase difference.
97* \param detector [in] A pointer to a \c LALDetector structure for a particular detector
98* \param ephem [in] A pointer to an \c EphemerisData structure containing solar system ephemeris information
99* \param tdat [in] A pointer to a \c TimeCorrectionData structure containing time system correction information
100* \param ttype [in] The \c TimeCorrectionType value
101*
102* \return A vector of rotational phase difference values (in cycles NOT radians)
103*
104* \sa XLALHeterodynedPulsarGetSSBDelay
105* \sa XLALHeterodynedPulsarGetBSBDelay
106*/
108 PulsarParameters *origparams,
109 const LIGOTimeGPSVector *datatimes,
110 REAL8 freqfactor,
111 REAL8Vector *ssbdts,
112 UINT4 calcSSBDelay,
113 REAL8Vector *bsbdts,
114 UINT4 calcBSBDelay,
115 REAL8Vector *glphase,
116 UINT4 calcglphase,
117 REAL8Vector *fitwavesphase,
118 UINT4 calcfitwaves,
119 const LALDetector *detector,
120 const EphemerisData *ephem,
122 TimeCorrectionType ttype )
123{
124 /* check inputs */
125 XLAL_CHECK_NULL( params != NULL, XLAL_EFUNC, "PulsarParameters must not be NULL" );
126 XLAL_CHECK_NULL( datatimes != NULL, XLAL_EFUNC, "datatimes must not be NULL" );
127 XLAL_CHECK_NULL( freqfactor > 0., XLAL_EFUNC, "freqfactor must be greater than zero" );
128 XLAL_CHECK_NULL( detector != NULL, XLAL_EFUNC, "LALDetector must not be NULL" );
129 XLAL_CHECK_NULL( ephem != NULL, XLAL_EFUNC, "EphemerisData must not be NULL" );
130
131 UINT4 i = 0, j = 0, k = 0, length = 0, nfreqs = 0;
132
133 REAL8 DT = 0., deltat = 0., deltatpow = 0., deltatpowinner = 1., taylorcoeff = 1., Ddelay = 0., Ddelaypow = 0.;
134
135 REAL8Vector *phis = NULL, *dts = NULL, *fixdts = NULL, *bdts = NULL, *fixbdts = NULL, *fixglph = NULL, *glph = NULL, *fixfitwavesph = NULL, *fitwavesph = NULL;
136
137 REAL8 pepoch = PulsarGetREAL8ParamOrZero( params, "PEPOCH" ); /* time of ephem info */
138 REAL8 T0 = pepoch, pepochorig = pepoch;
139 REAL8 *deltafs = NULL, *frequpdate = NULL;
140
141 length = datatimes->length;
142
143 /* allocate memory for phases */
144 phis = XLALCreateREAL8Vector( length );
145
146 /* get solar system barycentring time delays */
147 fixdts = ssbdts;
148 if ( calcSSBDelay ) {
149 if ( origparams != NULL && ssbdts != NULL ) {
150 dts = XLALHeterodynedPulsarGetSSBDelay( params, datatimes, detector, ephem, tdat, ttype );
151 } else {
152 fixdts = XLALHeterodynedPulsarGetSSBDelay( params, datatimes, detector, ephem, tdat, ttype );
153 }
154 XLAL_CHECK_NULL( length == fixdts->length, XLAL_EFUNC, "Lengths of time stamp vector and SSB delay vector are not the same" );
155 }
156
157 /* get the binary system barycentring time delays */
158 fixbdts = bsbdts;
159 if ( calcBSBDelay ) {
160 REAL8Vector *whichssb = dts != NULL ? dts : fixdts;
161 if ( origparams != NULL && bsbdts != NULL ) {
162 bdts = XLALHeterodynedPulsarGetBSBDelay( params, datatimes, whichssb, ephem );
163 } else {
164 fixbdts = XLALHeterodynedPulsarGetBSBDelay( params, datatimes, whichssb, ephem );
165 }
166 XLAL_CHECK_NULL( length == fixbdts->length, XLAL_EFUNC, "Lengths of time stamp vector and BSB delay vector are not the same" );
167 }
168
169 /* get the glitch parameter phase evolution */
170 fixglph = glphase;
171 if ( calcglphase ) {
172 const REAL8Vector *whichssb = dts != NULL ? dts : fixdts;
173 const REAL8Vector *whichbsb = bdts != NULL ? bdts : fixbdts;
174 if ( origparams != NULL && glphase != NULL ) {
175 glph = XLALHeterodynedPulsarGetGlitchPhase( params, datatimes, whichssb, whichbsb );
176 } else {
177 fixglph = XLALHeterodynedPulsarGetGlitchPhase( params, datatimes, whichssb, whichbsb );
178 }
179 }
180
181 /* set the "new" frequencies, with zeros where required */
182 if ( origparams != NULL ) {
183 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( origparams, "F" );
184 const REAL8Vector *freqs = PulsarGetREAL8VectorParam( params, "F" );
185 nfreqs = tmpvec->length >= freqs->length ? tmpvec->length : freqs->length;
186 pepochorig = PulsarGetREAL8ParamOrZero( origparams, "PEPOCH" );
187 } else {
188 const REAL8Vector *freqs = PulsarGetREAL8VectorParam( params, "F" );
189 nfreqs = freqs->length;
190 }
191
192 frequpdate = XLALCalloc( nfreqs, sizeof( REAL8 ) );
193 const REAL8Vector *freqsu = PulsarGetREAL8VectorParam( params, "F" );
194 /* initialise deltafs to zero */
195 deltafs = XLALCalloc( nfreqs, sizeof( REAL8 ) );
196 deltat = ( pepochorig - pepoch );
197
198 for ( i = 0; i < freqsu->length; i++ ) {
199 frequpdate[i] = freqsu->data[i];
200
201 /* if epochs are different update "new" frequency to heterodyne epoch */
202 if ( pepoch != pepochorig ) {
203 REAL8 deltatpowu = deltat;
204 for ( j = i + 1; j < freqsu->length; j++ ) {
205 frequpdate[i] += ( freqsu->data[j] * deltatpowu ) / gsl_sf_fact( j - i );
206
207 deltatpowu *= deltat;
208 }
209 }
210 }
211
212 if ( pepoch != pepochorig ) {
213 // correct the epoch value
214 T0 = pepochorig;
215 }
216
217 /* get vector of frequencies and frequency differences */
218 if ( origparams != NULL ) {
219 if ( PulsarCheckParam( origparams, "F" ) ) {
220 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( origparams, "F" );
221
222 for ( i = 0; i < nfreqs; i++ ) {
223 /* deal with cases where there are fewer frequency derivatives in the original heterodyne file */
224 if ( i >= tmpvec->length ) {
225 deltafs[i] = frequpdate[i];
226 } else {
227 deltafs[i] = frequpdate[i] - tmpvec->data[i];
228 }
229 }
230 } else {
231 XLAL_ERROR_NULL( XLAL_EINVAL, "No frequencies given in heterodyned parameter file!" );
232 }
233 } else {
234 /* set deltafs to (negative) frequencies */
235 for ( i = 0; i < nfreqs; i++ ) {
236 deltafs[i] = -frequpdate[i];
237 }
238 }
239
240 /* get FITWAVES phase */
241 fixfitwavesph = fitwavesphase;
242 if ( calcfitwaves ) {
243 const REAL8Vector *whichssb = dts != NULL ? dts : fixdts;
244 if ( origparams != NULL && fitwavesphase != NULL ) {
245 fitwavesph = XLALHeterodynedPulsarGetFITWAVESPhase( params, datatimes, whichssb, frequpdate[0] );
246 } else {
247 fixfitwavesph = XLALHeterodynedPulsarGetFITWAVESPhase( params, datatimes, whichssb, frequpdate[0] );
248 }
249 }
250
251 for ( i = 0; i < length; i++ ) {
252 REAL8 deltaphi = 0., innerphi = 0.; /* change in phase */
253 Ddelay = 0.; /* change in SSB/BSB delay */
254
255 REAL8 realT = XLALGPSGetREAL8( &datatimes->data[i] ); /* time of data */
256 DT = realT - T0; /* time diff between data and start of data */
257
258 /* get difference in solar system barycentring time delays */
259 deltat = DT;
260 if ( dts != NULL ) {
261 Ddelay += ( dts->data[i] - fixdts->data[i] );
262 }
263 if ( fixdts != NULL ) {
264 deltat += fixdts->data[i];
265 }
266
267 /* get difference in binary system barycentring time delays */
268 if ( bdts != NULL && fixbdts != NULL ) {
269 Ddelay += ( bdts->data[i] - fixbdts->data[i] );
270 }
271 if ( fixbdts != NULL ) {
272 deltat += fixbdts->data[i];
273 }
274
275 /* get the change in phase (compared to the heterodyned phase) */
276 deltatpow = deltat;
277 for ( j = 0; j < nfreqs; j++ ) {
278 taylorcoeff = gsl_sf_fact( j + 1 );
279 deltaphi += deltafs[j] * deltatpow / taylorcoeff;
280 if ( Ddelay != 0. ) {
281 innerphi = 0.;
282 deltatpowinner = 1.; /* this starts as one as it is first raised to the power of zero */
283 Ddelaypow = pow( Ddelay, j + 1 );
284 for ( k = 0; k < j + 1; k++ ) {
285 innerphi += gsl_sf_choose( j + 1, k ) * Ddelaypow * deltatpowinner;
286 deltatpowinner *= deltat; /* raise power */
287 Ddelaypow /= Ddelay; /* reduce power */
288 }
289 deltaphi += innerphi * frequpdate[j] / taylorcoeff;
290 }
291 deltatpow *= deltat;
292 }
293
294 /* get change in phase from glitches */
295 if ( glph != NULL && fixglph != NULL ) {
296 deltaphi += ( glph->data[i] - fixglph->data[i] );
297 } else if ( fixglph != NULL ) {
298 deltaphi += fixglph->data[i];
299 }
300
301 /* get change in phase from FITWAVES parameters */
302 if ( fitwavesph != NULL && fixfitwavesph != NULL ) {
303 deltaphi += ( fitwavesph->data[i] - fixfitwavesph->data[i] );
304 } else if ( fixfitwavesph != NULL ) {
305 deltaphi += fixfitwavesph->data[i];
306 }
307
308 deltaphi *= freqfactor; /* multiply by frequency factor */
309 phis->data[i] = deltaphi - floor( deltaphi ); /* only need to keep the fractional part of the phase */
310 }
311
312 /* free memory */
313 if ( dts != NULL ) {
315 }
316 if ( bdts != NULL ) {
318 }
319 if ( ssbdts == NULL ) {
320 XLALDestroyREAL8Vector( fixdts );
321 }
322 if ( bsbdts == NULL ) {
323 XLALDestroyREAL8Vector( fixbdts );
324 }
325 if ( glph != NULL ) {
327 }
328 if ( glphase == NULL ) {
329 XLALDestroyREAL8Vector( fixglph );
330 }
331 if ( fitwavesph != NULL ) {
332 XLALDestroyREAL8Vector( fitwavesph );
333 }
334 if ( fitwavesphase == NULL ) {
335 XLALDestroyREAL8Vector( fixfitwavesph );
336 }
337 XLALFree( deltafs );
338 XLALFree( frequpdate );
339
340 return phis;
341}
342
343
344/**
345 * \brief Computes the delay between a GPS time at Earth and the solar system barycentre
346 *
347 * This function calculates the time delay between a GPS time at a specific
348 * location (e.g., a gravitational-wave detector) on Earth and the solar system
349 * barycentre. The delay consists of three components: the geometric time delay
350 * (Roemer delay) \f$ t_R = \mathbf{r}(t)\hat{n}/c \f$ (where \f$ \mathbf{r}(t) \f$
351 * is the detector's position vector at time \f$ t \f$ ), the special relativistic
352 * Einstein delay \f$ t_E \f$ , and the general relativistic Shapiro delay
353 * \f$ t_S \f$ .
354 *
355 * \param pars [in] A set of pulsar parameters
356 * \param datatimes [in] A vector of GPS times at Earth
357 * \param detector [in] Information on the detector position on the Earth
358 * \param ephem [in] Information on the solar system ephemeris
359 * \param tdat [in] Informaion on the time system corrections
360 * \param ttype [in] The type of time system corrections to perform
361 *
362 * \return A vector of time delays in seconds
363 *
364 * \sa XLALBarycenter
365 * \sa XLALBarycenterEarthNew
366 */
368 const LIGOTimeGPSVector *datatimes,
369 const LALDetector *detector,
370 const EphemerisData *ephem,
372 TimeCorrectionType ttype )
373{
374 /* check inputs */
375 XLAL_CHECK_NULL( pars != NULL, XLAL_EFUNC, "PulsarParameters must not be NULL" );
376 XLAL_CHECK_NULL( datatimes != NULL, XLAL_EFUNC, "datatimes must not be NULL" );
377 XLAL_CHECK_NULL( detector != NULL, XLAL_EFUNC, "LALDetector must not be NULL" );
378 XLAL_CHECK_NULL( ephem != NULL, XLAL_EFUNC, "EphemerisData must not be NULL" );
379
380 INT4 i = 0, length = 0;
381
382 BarycenterInput bary;
383
384 REAL8Vector *dts = NULL;
385
386 /* copy barycenter and ephemeris data */
387 bary.site.location[0] = detector->location[0] / LAL_C_SI;
388 bary.site.location[1] = detector->location[1] / LAL_C_SI;
389 bary.site.location[2] = detector->location[2] / LAL_C_SI;
390
391 REAL8 ra = 0.;
392 if ( PulsarCheckParam( pars, "RA" ) ) {
393 ra = PulsarGetREAL8Param( pars, "RA" );
394 } else if ( PulsarCheckParam( pars, "RAJ" ) ) {
395 ra = PulsarGetREAL8Param( pars, "RAJ" );
396 } else {
397 XLAL_ERROR_NULL( XLAL_EINVAL, "No source right ascension specified!" );
398 }
399 REAL8 dec = 0.;
400 if ( PulsarCheckParam( pars, "DEC" ) ) {
401 dec = PulsarGetREAL8Param( pars, "DEC" );
402 } else if ( PulsarCheckParam( pars, "DECJ" ) ) {
403 dec = PulsarGetREAL8Param( pars, "DECJ" );
404 } else {
405 XLAL_ERROR_NULL( XLAL_EINVAL, "No source declination specified!" );
406 }
407 REAL8 pmra = PulsarGetREAL8ParamOrZero( pars, "PMRA" );
408 REAL8 pmdec = PulsarGetREAL8ParamOrZero( pars, "PMDEC" );
409 REAL8 pepoch = PulsarGetREAL8ParamOrZero( pars, "PEPOCH" );
410 REAL8 posepoch = PulsarGetREAL8ParamOrZero( pars, "POSEPOCH" );
411
412 REAL8 cgw = PulsarGetREAL8ParamOrZero( pars, "CGW" );
413
414 REAL8 dist = 0.; /* distance in light seconds */
415 /* set distance (a DIST param takes precedence over PX) */
416 if ( PulsarCheckParam( pars, "DIST" ) ) {
417 dist = PulsarGetREAL8Param( pars, "DIST" ) / LAL_C_SI;
418 } else if ( PulsarCheckParam( pars, "PX" ) ) {
419 dist = ( LAL_AU_SI / LAL_C_SI ) / PulsarGetREAL8Param( pars, "PX" );
420 }
421
422 /* set the position and frequency epochs if not already set */
423 if ( pepoch == 0. && posepoch != 0. ) {
424 pepoch = posepoch;
425 } else if ( posepoch == 0. && pepoch != 0. ) {
426 posepoch = pepoch;
427 }
428
429 length = datatimes->length;
430
431 /* allocate memory for times delays */
432 dts = XLALCreateREAL8Vector( length );
433
434 /* set 1/distance if distance is given */
435 if ( dist != 0. ) {
436 bary.dInv = 1. / dist;
437 } else {
438 bary.dInv = 0.;
439 }
440
441 /* make sure ra and dec are wrapped within 0--2pi and -pi.2--pi/2 respectively */
442 ra = fmod( ra, LAL_TWOPI );
443 REAL8 absdec = fabs( dec );
444 if ( absdec > LAL_PI_2 ) {
445 UINT4 nwrap = floor( ( absdec + LAL_PI_2 ) / LAL_PI );
446 dec = ( dec > 0 ? 1. : -1. ) * ( nwrap % 2 == 1 ? -1. : 1. ) * ( fmod( absdec + LAL_PI_2, LAL_PI ) - LAL_PI_2 );
447 ra = fmod( ra + ( REAL8 )nwrap * LAL_PI, LAL_TWOPI ); /* move RA by pi */
448 }
449
450 EarthState earth;
451 EmissionTime emit;
452 for ( i = 0; i < length; i++ ) {
453 REAL8 realT = XLALGPSGetREAL8( &datatimes->data[i] );
454
455 bary.tgps = datatimes->data[i];
456 bary.delta = dec + ( realT - posepoch ) * pmdec;
457 bary.alpha = ra + ( realT - posepoch ) * pmra / cos( bary.delta );
458
459 /* call barycentring routines */
460 XLAL_CHECK_NULL( XLALBarycenterEarthNew( &earth, &bary.tgps, ephem, tdat, ttype ) == XLAL_SUCCESS, XLAL_EFUNC, "Barycentring routine failed" );
461 XLAL_CHECK_NULL( XLALBarycenter( &emit, &bary, &earth ) == XLAL_SUCCESS, XLAL_EFUNC, "Barycentring routine failed" );
462
463 if ( cgw > 0.0 ) {
464 /* account for the speed of GWs not being the same a the speed of light
465 NOTE: we are only accounting for this in the Roemer delay */
466 dts->data[i] = ( emit.deltaT - emit.roemer ) + ( emit.roemer / cgw );
467 } else {
468 dts->data[i] = emit.deltaT;
469 }
470 }
471
472 return dts;
473}
474
475
476/**
477 * \brief Computes the delay between a pulsar in a binary system and the barycentre of the system
478 *
479 * This function uses \c XLALBinaryPulsarDeltaTNew to calculate the time delay
480 * between for a pulsar in a binary system between the time at the pulsar and
481 * the time at the barycentre of the system. This includes Roemer delays and
482 * relativistic delays. The orbit may be described by different models and can
483 * be purely Keplarian or include various relativistic corrections.
484 *
485 * \param pars [in] A set of pulsar parameters
486 * \param datatimes [in] A vector of GPS times
487 * \param dts [in] A vector of solar system barycentre time delays
488 * \param edat [in] Solar system ephemeris information
489 *
490 * \return A vector of time delays in seconds
491 *
492 * \sa XLALBinaryPulsarDeltaTNew
493 */
495 const LIGOTimeGPSVector *datatimes,
496 const REAL8Vector *dts,
497 const EphemerisData *edat )
498{
499 /* check inputs */
500 XLAL_CHECK_NULL( pars != NULL, XLAL_EFUNC, "PulsarParameters must not be NULL" );
501 XLAL_CHECK_NULL( datatimes != NULL, XLAL_EFUNC, "datatimes must not be NULL" );
502 XLAL_CHECK_NULL( dts != NULL, XLAL_EFUNC, "SSB delay times must not be NULL" );
503 XLAL_CHECK_NULL( edat != NULL, XLAL_EFUNC, "EphemerisData must not be NULL" );
504
505 REAL8 cgw = PulsarGetREAL8ParamOrZero( pars, "CGW" );
506
507 REAL8Vector *bdts = NULL;
508 BinaryPulsarInput binput;
509 BinaryPulsarOutput boutput;
510 EarthState earth;
511
512 INT4 i = 0, length = datatimes->length;
513
514 bdts = XLALCreateREAL8Vector( length );
515 memset( bdts->data, 0, bdts->length * sizeof( REAL8 ) ); // set to zeros
516
517 /* check whether there's a binary model */
518 if ( PulsarCheckParam( pars, "BINARY" ) ) {
519 for ( i = 0; i < length; i++ ) {
520 binput.tb = XLALGPSGetREAL8( &datatimes->data[i] ) + dts->data[i];
521
522 XLALGetEarthPosVel( &earth, edat, &datatimes->data[i] );
523
524 binput.earth = earth; /* current Earth state */
525 XLALBinaryPulsarDeltaTNew( &boutput, &binput, pars );
526
527 if ( cgw > 0. ) {
528 /* account for the speed of GWs not being the same as the speed of light
529 * NOTE: here we have to assume that the Roemer delay, which is effected
530 * by the speed difference is dominating the delay term. */
531 bdts->data[i] = boutput.deltaT / cgw;
532 } else {
533 bdts->data[i] = boutput.deltaT;
534 }
535 }
536 }
537 return bdts;
538}
539
540
541/**
542 * \brief Get the position and velocity of the Earth at a given time
543 *
544 * This function will get the position and velocity of the Earth from the
545 * ephemeris data at the time t. It will be returned in an EarthState
546 * structure. This is based on the start of the \c XLALBarycenterEarth
547 * function.
548 *
549 * \param earth [out] The returned \c EarthState structure containing the Earth's position and velocity
550 * \param edat [in] The solar system ephemeris data
551 * \param tGPS [in] A \c LIGOTimeGPS structure containing a GPS time
552 *
553 * \sa XLALBarycenterEarth
554 */
556 const EphemerisData *edat,
557 const LIGOTimeGPS *tGPS )
558{
559 REAL8 tgps[2];
560
561 REAL8 t0e; /* time since first entry in Earth ephem. table */
562 INT4 ientryE; /* entry in look-up table closest to current time, tGPS */
563
564 REAL8 tinitE; /* time (GPS) of first entry in Earth ephem table */
565 REAL8 tdiffE; /* current time tGPS minus time of nearest entry in Earth ephem look-up table */
566 REAL8 tdiff2E; /* tdiff2 = tdiffE * tdiffE */
567
568 INT4 j;
569
570 /* check input */
571 if ( !earth || !tGPS || !edat || !edat->ephemE || !edat->ephemS ) {
572 XLAL_ERROR_VOID( XLAL_EINVAL, "Invalid NULL input 'earth', 'tGPS', 'edat','edat->ephemE' or 'edat->ephemS'" );
573 }
574
575 tgps[0] = ( REAL8 )tGPS->gpsSeconds; /* convert from INT4 to REAL8 */
576 tgps[1] = ( REAL8 )tGPS->gpsNanoSeconds;
577
578 tinitE = edat->ephemE[0].gps;
579
580 t0e = tgps[0] - tinitE;
581 ientryE = ROUND( t0e / edat->dtEtable ); /* finding Earth table entry */
582
583 if ( ( ientryE < 0 ) || ( ientryE >= edat->nentriesE ) ) {
584 XLAL_ERROR_VOID( XLAL_EDOM, "Input GPS time %f outside of Earth ephem range [%f, %f]\n", tgps[0], tinitE, tinitE +
585 edat->nentriesE * edat->dtEtable );
586 }
587
588 /* tdiff is arrival time minus closest Earth table entry; tdiff can be pos. or neg. */
589 tdiffE = t0e - edat->dtEtable * ientryE + tgps[1] * 1.e-9;
590 tdiff2E = tdiffE * tdiffE;
591
592 REAL8 *pos = edat->ephemE[ientryE].pos;
593 REAL8 *vel = edat->ephemE[ientryE].vel;
594 REAL8 *acc = edat->ephemE[ientryE].acc;
595
596 for ( j = 0; j < 3; j++ ) {
597 earth->posNow[j] = pos[j] + vel[j] * tdiffE + 0.5 * acc[j] * tdiff2E;
598 earth->velNow[j] = vel[j] + acc[j] * tdiffE;
599 }
600}
601
602
603/**
604 * \brief Computes the phase evolution due to the presence of glitches
605 *
606 * This function calculates the phase offset due to the presence on pulsar
607 * glitches. This is based on the equations in the \c formResiduals.C file of
608 * TEMPO2 from Eqn. 1 of \cite Yu2013 .
609 *
610 * \param params [in] A set of pulsar parameters
611 * \param datatimes [in] A vector of GPS times at Earth
612 * \param ssbdts [in] A vector of solar system barycentre time delays
613 * \param bsbdts [in] A vector of binary system barycentre time delays
614 *
615 * \return A vector of (EM) phases in cycles
616 */
618 const LIGOTimeGPSVector *datatimes,
619 const REAL8Vector *ssbdts,
620 const REAL8Vector *bsbdts )
621{
622 /* check inputs */
623 XLAL_CHECK_NULL( params != NULL, XLAL_EFUNC, "PulsarParameters must not be NULL" );
624 XLAL_CHECK_NULL( datatimes != NULL, XLAL_EFUNC, "datatimes must not be NULL" );
625 XLAL_CHECK_NULL( ssbdts != NULL, XLAL_EFUNC, "ssbdts must not be NULL" );
626 XLAL_CHECK_NULL( ssbdts->length == datatimes->length, XLAL_EFUNC, "datatimes and ssbdts must be the same length" );
627 if ( bsbdts != NULL ) {
628 XLAL_CHECK_NULL( ssbdts->length == bsbdts->length, XLAL_EFUNC, "bsbdts and ssbdts must be the same length" );
629 }
630
631 UINT4 i = 0, length = 0;
632 REAL8Vector *glphase = NULL;
633 REAL8 tdt = 0.;
634
635 length = datatimes->length;
636
637 /* allocate memory for phases */
638 glphase = XLALCreateREAL8Vector( length );
639 memset( glphase->data, 0, glphase->length * sizeof( REAL8 ) ); // set to zeros
640
641 /* see if pulsar has glitch parameters */
642 if ( PulsarCheckParam( params, "GLEP" ) ) {
643 /* in case not all parameters are set and they are not all the same length, copy into arrays
644 of the same length that are initialised to zero - all glitches must have an epoch */
645 const REAL8Vector *glep = PulsarGetREAL8VectorParam( params, "GLEP" );
646 UINT4 glnum = glep->length;
647
648 /* get phase offsets */
649 REAL8 *glph = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
650 if ( PulsarCheckParam( params, "GLPH" ) ) {
651 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLPH" );
652 memcpy( glph, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
653 }
654
655 /* get frequencies offsets */
656 REAL8 *glf0 = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
657 if ( PulsarCheckParam( params, "GLF0" ) ) {
658 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLF0" );
659 memcpy( glf0, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
660 }
661
662 /* get frequency derivative offsets */
663 REAL8 *glf1 = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
664 if ( PulsarCheckParam( params, "GLF1" ) ) {
665 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLF1" );
666 memcpy( glf1, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
667 }
668
669 /* get second frequency derivative offsets */
670 REAL8 *glf2 = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
671 if ( PulsarCheckParam( params, "GLF2" ) ) {
672 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLF2" );
673 memcpy( glf2, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
674 }
675
676 /* get decaying frequency component offset derivative */
677 REAL8 *glf0d = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
678 if ( PulsarCheckParam( params, "GLF0D" ) ) {
679 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLF0D" );
680 memcpy( glf0d, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
681 }
682
683 /* get decaying frequency component decay time constant */
684 REAL8 *gltd = XLALCalloc( glnum, sizeof( REAL8 ) ); /* initialise to zeros */
685 if ( PulsarCheckParam( params, "GLTD" ) ) {
686 const REAL8Vector *tmpvec = PulsarGetREAL8VectorParam( params, "GLTD" );
687 memcpy( gltd, tmpvec->data, tmpvec->length * sizeof( REAL8 ) );
688 }
689
690 for ( i = 0; i < length; i++ ) {
691 glphase->data[i] = 0.;
692 tdt = XLALGPSGetREAL8( &datatimes->data[i] ) + ssbdts->data[i];
693 if ( bsbdts != NULL ) {
694 tdt += bsbdts->data[i];
695 }
696 for ( UINT4 k = 0; k < glnum; k++ ) {
697 if ( tdt >= glep->data[k] ) {
698 REAL8 dtg = 0, expd = 1.;
699 dtg = tdt - glep->data[k]; /* time since glitch */
700 if ( gltd[k] != 0. ) {
701 expd = exp( -dtg / gltd[k] ); /* decaying part of glitch */
702 }
703 glphase->data[i] += glph[k] + glf0[k] * dtg + 0.5 * glf1[k] * dtg * dtg + ( 1. / 6. ) * glf2[k] * dtg * dtg * dtg + glf0d[k] * gltd[k] * ( 1. - expd );
704 }
705 }
706 }
707
708 XLALFree( glph );
709 XLALFree( glf0 );
710 XLALFree( glf1 );
711 XLALFree( glf2 );
712 XLALFree( glf0d );
713 XLALFree( gltd );
714 }
715
716 return glphase;
717}
718
719
720/**
721 * \brief Computes the phase evolution due to the presence of FITWAVES parameters
722 *
723 * This function calculates the phase due to the presence on FITWAVES parameters in
724 * the pulsar timing model. These are used to fit multiple sinusoids to account for
725 * strong red noise effects (e.g., timing noise) in the pulsar model.
726 *
727 * \param params [in] A set of pulsar parameters
728 * \param datatimes [in] A vector of GPS times at Earth
729 * \param ssbdts [in] A vector of solar system barycentre time delays
730 * \param freq [in] The EM frequency.
731 *
732 * \return A vector of (EM) phases in cycles
733 */
735 const LIGOTimeGPSVector *datatimes,
736 const REAL8Vector *ssbdts,
737 REAL8 freq )
738{
739 /* check inputs */
740 XLAL_CHECK_NULL( params != NULL, XLAL_EFUNC, "PulsarParameters must not be NULL" );
741 XLAL_CHECK_NULL( datatimes != NULL, XLAL_EFUNC, "datatimes must not be NULL" );
742 XLAL_CHECK_NULL( ssbdts != NULL, XLAL_EFUNC, "ssbdts must not be NULL" );
743 XLAL_CHECK_NULL( ssbdts->length == datatimes->length, XLAL_EFUNC, "datatimes and ssbdts must be the same length" );
744 XLAL_CHECK_NULL( freq > 0.0, XLAL_EFUNC, "Frequency must be positive" );
745
746 UINT4 i = 0, length = 0, j = 0;
747 REAL8Vector *phasewave = NULL;
748
749 length = datatimes->length;
750
751 /* allocate memory for phases */
752 phasewave = XLALCreateREAL8Vector( length );
753 memset( phasewave->data, 0, phasewave->length * sizeof( REAL8 ) ); // set to zeros
754
755 /* get FITWAVES parameters */
756 if ( PulsarCheckParam( params, "WAVESIN" ) && PulsarCheckParam( params, "WAVECOS" ) ) {
757 const REAL8Vector *wsin = PulsarGetREAL8VectorParam( params, "WAVESIN" );
758 const REAL8Vector *wcos = PulsarGetREAL8VectorParam( params, "WAVECOS" );
759
760 XLAL_CHECK_NULL( wsin->length == wcos->length, XLAL_EFUNC, "Lengths of FITWAVES sin and cos parameters are not the same" );
761 UINT4 nwaves = wsin->length;
762
763 REAL8 waveepoch = PulsarGetREAL8ParamOrZero( params, "WAVEEPOCH" );
764 REAL8 waveom = PulsarGetREAL8ParamOrZero( params, "WAVE_OM" );
765
766 REAL8 twave = 0., dtwave = 0.;
767 for ( i = 0; i < length; i++ ) {
768 twave = 0.0;
769 /* note: this does not include the binary times delays */
770 dtwave = ( XLALGPSGetREAL8( &datatimes->data[i] ) + ssbdts->data[i] - waveepoch ) / 86400.; /* in days */
771
772 for ( j = 0; j < nwaves; j++ ) {
773 twave += wsin->data[j] * sin( waveom * ( REAL8 )( j + 1. ) * dtwave );
774 twave += wcos->data[j] * cos( waveom * ( REAL8 )( j + 1. ) * dtwave );
775 }
776
777 phasewave->data[i] = freq * twave;
778 }
779 }
780
781 return phasewave;
782}
783
784
785/**
786 * \brief The amplitude model of a complex heterodyned signal from the
787 * \f$ l=2, m=1,2 \f$ harmonics of a rotating neutron star.
788 *
789 * This function calculates the complex heterodyned time series model for a
790 * rotating neutron star. It will currently calculate the model for emission
791 * from the \f$ l=m=2 \f$ harmonic (which gives emission at twice the rotation
792 * frequency) and/or the \f$ l=2 \f$ and \f$ m=1 \f$ harmonic (which gives emission
793 * at the rotation frequency). See LIGO T1200265-v3. Further harmonics can be
794 * added and are defined by the \c freqFactor value, which is the multiple of
795 * the spin-frequency at which emission is produced.
796 *
797 * The antenna pattern functions are contained in a 1D lookup table, so within
798 * this function the correct value for the given time is interpolated from this
799 * lookup table using linear interpolation.
800 *
801 * \param pars [in] A set of pulsar parameters
802 * \param freqfactor [in] The harmonic frequency of the signal in units of the
803 * pulsar rotation frequency
804 * \param varyphase [in] Set to a non-zero value is the signal phase is
805 * different to the heterodyne phase (or if wanting the signal output at all
806 * time stamps).
807 * \param useroq [in] Set to a non-zero value if a reduced order quadrature
808 * likelihood is being used
809 * \param nonGR [in] Set to a non-zero value to indicate a non-GR polarisation
810 * is used
811 * \param timestamps [in] A vector of GPS times at which to calculate the signal
812 * \param resp [in] A detector response function look-up table
813 *
814 * \return A \c COMPLEX16TimeSeries the amplitude time series
815 */
817 REAL8 freqfactor,
818 UINT4 varyphase,
819 UINT4 useroq,
820 UINT4 nonGR,
822 const DetResponseTimeLookupTable *resp )
823{
824 COMPLEX16TimeSeries *csignal = NULL;
825
826 /* create signal if not already allocated */
827 LIGOTimeGPS gpstime; /* dummy time stamp */
828 gpstime.gpsSeconds = 0;
829 gpstime.gpsNanoSeconds = 0;
830
831 /* check for transient signal */
832 INT4 rectWindow = 0, expWindow = 0;
833 REAL8 transientStartTime = PulsarGetREAL8ParamOrZero( pars, "TRANSIENTSTARTTIME" );
834 REAL8 transientTau = PulsarGetREAL8ParamOrZero( pars, "TRANSIENTTAU" );
835
836 if ( PulsarCheckParam( pars, "TRANSIENTWINDOWTYPE" ) ) {
837 const CHAR *transientWindow = PulsarGetStringParam( pars, "TRANSIENTWINDOWTYPE" );
838 rectWindow = !strcmp( transientWindow, "RECT" );
839 expWindow = !strcmp( transientWindow, "EXP" );
840 }
841
842 if ( varyphase || useroq || rectWindow || expWindow ) {
843 XLAL_CHECK_NULL( timestamps->length > 0, XLAL_EFUNC, "length must be a positive integer" );
844
845 /* in these cases use signal length to specify the length of the signal */
847 } else {
848 if ( nonGR ) {
849 /* for non-GR signals there are six values that need to be held */
850 csignal = XLALCreateCOMPLEX16TimeSeries( "", &gpstime, 0., 1., &lalSecondUnit, 6 );
851 } else {
852 /* for GR signals just two values are required */
853 csignal = XLALCreateCOMPLEX16TimeSeries( "", &gpstime, 0., 1., &lalSecondUnit, 2 );
854 }
855 }
856
857 /* switch to waveform parameters if required */
859
860 /* check inputs */
861 XLAL_CHECK_NULL( pars != NULL, XLAL_EFUNC, "PulsarParameters is NULL" );
862 XLAL_CHECK_NULL( freqfactor > 0., XLAL_EFUNC, "freqfactor must be a positive number" );
863 XLAL_CHECK_NULL( resp != NULL, XLAL_EFUNC, "Response look-up table is NULL" );
864 XLAL_CHECK_NULL( freqfactor != 1. || nonGR != 1, XLAL_EFUNC, "A non-GR model is not defined for a frequency harmonic of 1" );
865
866 UINT4 i = 0;
867
868 REAL8 T, twopsi, t0;
869 REAL8 cosiota = 0.;
870 if ( PulsarCheckParam( pars, "IOTA" ) ) {
871 /* use IOTA if present */
872 cosiota = cos( PulsarGetREAL8Param( pars, "IOTA" ) );
873 } else {
874 cosiota = PulsarGetREAL8ParamOrZero( pars, "COSIOTA" );
875 }
876 REAL8 siniota = sin( acos( cosiota ) );
877 REAL8 s2psi = 0., c2psi = 0., spsi = 0., cpsi = 0.;
878
879 twopsi = 2.*PulsarGetREAL8ParamOrZero( pars, "PSI" );
880 s2psi = sin( twopsi );
881 c2psi = cos( twopsi );
882
883 t0 = XLALGPSGetREAL8( &resp->t0 ); /* epoch of detector response look-up table */
884
885 if ( nonGR == 1 ) {
886 spsi = sin( PulsarGetREAL8ParamOrZero( pars, "PSI" ) );
887 cpsi = cos( PulsarGetREAL8ParamOrZero( pars, "PSI" ) );
888 }
889
890 COMPLEX16 expPhi;
891 COMPLEX16 Cplus = 0., Ccross = 0., Cx = 0., Cy = 0., Cl = 0., Cb = 0.;
892
893 /* get the amplitude and phase factors */
894 if ( freqfactor == 1. ) {
895 /* the l=2, m=1 harmonic at the rotation frequency. */
896 if ( nonGR ) { /* amplitude if nonGR is specified */
897 COMPLEX16 expPhiTensor, expPsiTensor, expPhiScalar, expPsiScalar, expPhiVector, expPsiVector;
898
899 expPhiTensor = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0TENSOR_F" ) );
900 expPsiTensor = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSITENSOR_F" ) );
901 expPhiScalar = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0SCALAR_F" ) );
902 expPsiScalar = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSISCALAR_F" ) );
903 expPhiVector = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0VECTOR_F" ) );
904 expPsiVector = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSIVECTOR_F" ) );
905
906 Cplus = 0.5 * expPhiTensor * PulsarGetREAL8ParamOrZero( pars, "HPLUS_F" );
907 Ccross = 0.5 * expPhiTensor * PulsarGetREAL8ParamOrZero( pars, "HCROSS_F" ) * expPsiTensor;
908 Cx = 0.5 * expPhiVector * PulsarGetREAL8ParamOrZero( pars, "HVECTORX_F" );
909 Cy = 0.5 * expPhiVector * PulsarGetREAL8ParamOrZero( pars, "HVECTORY_F" ) * expPsiVector;
910 Cb = 0.5 * expPhiScalar * PulsarGetREAL8ParamOrZero( pars, "HSCALARB_F" );
911 Cl = 0.5 * expPhiScalar * PulsarGetREAL8ParamOrZero( pars, "HSCALARL_F" ) * expPsiScalar;
912 } else {
913 expPhi = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI21" ) );
914 Cplus = -0.25 * PulsarGetREAL8ParamOrZero( pars, "C21" ) * siniota * cosiota * expPhi;
915 Ccross = 0.25 * I * PulsarGetREAL8ParamOrZero( pars, "C21" ) * siniota * expPhi;
916 }
917 } else if ( freqfactor == 2. ) {
918 /* the l=2, m=2 harmonic at twice the rotation frequency */
919 if ( nonGR ) { /* amplitude if nonGR is specified */
920 COMPLEX16 expPhiTensor, expPsiTensor, expPhiScalar, expPsiScalar, expPhiVector, expPsiVector;
921
922 expPhiTensor = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0TENSOR" ) );
923 expPsiTensor = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSITENSOR" ) );
924 expPhiScalar = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0SCALAR" ) );
925 expPsiScalar = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSISCALAR" ) );
926 expPhiVector = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI0VECTOR" ) );
927 expPsiVector = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PSIVECTOR" ) );
928
929 Cplus = 0.5 * expPhiTensor * PulsarGetREAL8ParamOrZero( pars, "HPLUS" );
930 Ccross = 0.5 * expPhiTensor * PulsarGetREAL8ParamOrZero( pars, "HCROSS" ) * expPsiTensor;
931 Cx = 0.5 * expPhiVector * PulsarGetREAL8ParamOrZero( pars, "HVECTORX" );
932 Cy = 0.5 * expPhiVector * PulsarGetREAL8ParamOrZero( pars, "HVECTORY" ) * expPsiVector;
933 Cb = 0.5 * expPhiScalar * PulsarGetREAL8ParamOrZero( pars, "HSCALARB" );
934 Cl = 0.5 * expPhiScalar * PulsarGetREAL8ParamOrZero( pars, "HSCALARL" ) * expPsiScalar;
935 } else { /* just GR tensor mode amplitudes */
936 expPhi = cexp( I * PulsarGetREAL8ParamOrZero( pars, "PHI22" ) );
937 Cplus = -0.5 * PulsarGetREAL8ParamOrZero( pars, "C22" ) * ( 1. + cosiota * cosiota ) * expPhi;
938 Ccross = I * PulsarGetREAL8ParamOrZero( pars, "C22" ) * cosiota * expPhi;
939 }
940 } else {
941 XLAL_ERROR_NULL( XLAL_EINVAL, "Error... currently unknown frequency factor (%.2lf) for models.", freqfactor );
942 }
943
944 if ( varyphase || useroq || rectWindow || expWindow ) { /* have to compute the full time domain signal */
945 REAL8 tsv;
946
947 /* set lookup table parameters */
948 tsv = LAL_DAYSID_SI / ( REAL8 )resp->ntimebins;
949
950 for ( i = 0; i < csignal->data->length; i++ ) {
951 REAL8 plus00, plus01, cross00, cross01, plus = 0., cross = 0.;
952 REAL8 x00, x01, y00, y01, b00, b01, l00, l01;
953 REAL8 timeGPS, timeScaled, timeMin, timeMax;
954 REAL8 plusT = 0., crossT = 0., x = 0., y = 0., xT = 0., yT = 0., b = 0., l = 0.;
955 INT4 timebinMin, timebinMax;
956
957 /* set the time bin for the lookup table */
958 /* sidereal day in secs*/
959 timeGPS = XLALGPSGetREAL8( &timestamps->data[i] );
960 T = fmod( timeGPS - t0, LAL_DAYSID_SI ); /* convert GPS time to sidereal day */
961
962 /* set signal to zero for transient signals as required */
963 if ( rectWindow && ( timeGPS < transientStartTime || timeGPS > transientStartTime + transientTau ) ) {
964 csignal->data->data[i] = 0.0 + I * 0.0;
965 continue;
966 } else if ( expWindow && timeGPS < transientStartTime ) {
967 csignal->data->data[i] = 0.0 + I * 0.0;
968 continue;
969 }
970
971 timebinMin = ( INT4 )fmod( floor( T / tsv ), resp->ntimebins );
972 timeMin = timebinMin * tsv;
973 timebinMax = ( INT4 )fmod( timebinMin + 1, resp->ntimebins );
974 timeMax = timeMin + tsv;
975
976 /* get values of matrix for linear interpolation */
977 plus00 = resp->fplus->data[timebinMin];
978 plus01 = resp->fplus->data[timebinMax];
979
980 cross00 = resp->fcross->data[timebinMin];
981 cross01 = resp->fcross->data[timebinMax];
982
983 /* rescale time for linear interpolation on a unit square */
984 timeScaled = ( T - timeMin ) / ( timeMax - timeMin );
985
986 plus = plus00 + ( plus01 - plus00 ) * timeScaled;
987 cross = cross00 + ( cross01 - cross00 ) * timeScaled;
988
989 plusT = plus * c2psi + cross * s2psi;
990 crossT = cross * c2psi - plus * s2psi;
991
992 if ( nonGR ) {
993 x00 = resp->fx->data[timebinMin];
994 x01 = resp->fx->data[timebinMax];
995 y00 = resp->fy->data[timebinMin];
996 y01 = resp->fy->data[timebinMax];
997 b00 = resp->fb->data[timebinMin];
998 b01 = resp->fb->data[timebinMax];
999 l00 = resp->fl->data[timebinMin];
1000 l01 = resp->fl->data[timebinMax];
1001
1002 x = x00 + ( x01 - x00 ) * timeScaled;
1003 y = y00 + ( y01 - y00 ) * timeScaled;
1004 b = b00 + ( b01 - b00 ) * timeScaled;
1005 l = l00 + ( l01 - l00 ) * timeScaled;
1006
1007 xT = x * cpsi + y * spsi;
1008 yT = y * cpsi - x * spsi;
1009 }
1010
1011 /* create the complex signal amplitude model appropriate for the harmonic */
1012 csignal->data->data[i] = ( Cplus * plusT ) + ( Ccross * crossT );
1013
1014 /* add non-GR components if required */
1015 if ( nonGR ) {
1016 csignal->data->data[i] += ( Cx * xT ) + ( Cy * yT ) + Cb * b + Cl * l;
1017 }
1018
1019 /* add exponential window if required */
1020 if ( expWindow ) {
1021 csignal->data->data[i] *= exp( -( timeGPS - transientStartTime ) / transientTau );
1022 }
1023 }
1024 } else { /* just have to calculate the values to multiply the pre-summed data */
1025 /* for tensor-only models (e.g. the default of GR) calculate the two components of
1026 * the single model value - things multiplied by a(t) and things multiplied by b(t)
1027 * (both these will have real and imaginary components). NOTE: the values input into
1028 * csignal->data->data are not supposed to be identical to the above
1029 * relationships between the amplitudes and polarisation angles, as these are the
1030 * multiplicative coefficients of the a(t) and b(t) summations.
1031 */
1032
1033 /* put multiples of a(t) in first value and b(t) in second */
1034 if ( !nonGR ) {
1035 csignal->data->data[0] = ( Cplus * c2psi - Ccross * s2psi );
1036 csignal->data->data[1] = ( Cplus * s2psi + Ccross * c2psi );
1037 } else {
1038 csignal->data->data[0] = ( Cplus * c2psi - Ccross * s2psi );
1039 csignal->data->data[1] = ( Cplus * s2psi + Ccross * c2psi );
1040 csignal->data->data[2] = ( Cx * cpsi - Cy * spsi );
1041 csignal->data->data[3] = ( Cx * spsi + Cy * cpsi );
1042 csignal->data->data[4] = Cb;
1043 csignal->data->data[5] = Cl;
1044 }
1045 }
1046
1047 return csignal;
1048}
1049
1050
1051/**
1052 * \brief Generate the model of the neutron star signal
1053 *
1054 * Firstly the time varying amplitude of the signal will be calculated based
1055 * on the antenna pattern and amplitude parameters. Then, if searching over
1056 * phase parameters, the phase evolution of the signal will be calculated. The
1057 * difference between the new phase model, \f$ \phi(t)_n \f$ , and that used to
1058 * heterodyne the data, \f$ \phi(t)_h \f$ , will be calculated and the complex
1059 * signal model, \f$ M \f$ , modified accordingly:
1060 * \f[
1061 * M'(t) = M(t)\exp{i((\phi(t)_n - \phi(t)_h))}.
1062 * \f]
1063 * This does not try to undo the signal modulation in the data, but instead
1064 * replicates the modulation in the model, hence the positive phase difference
1065 * rather than a negative phase in the exponential function.
1066 *
1067 * The model is described in Equations 7 and 8 of \cite Pitkin2017 .
1068 *
1069 * \param pars [in] A \c PulsarParameters structure containing the model parameters
1070 * \param origpars [in] A \c PulsarParameters structure containing the original heterodyne parameters
1071 * \param freqfactor [in] The harmonic frequency of the signal in units of the
1072 * pulsar rotation frequency
1073 * \param usephase [in] Set to a non-zero value is the signal phase is
1074 * different to the heterodyne phase (or if wanting the signal output at all
1075 * time stamps).
1076 * \param useroq [in] Set to a non-zero value if a reduced order quadrature
1077 * likelihood is being used
1078 * \param nonGR [in] Set to a non-zero value to indicate a non-GR polarisation
1079 * is used
1080 * \param timestamps [in] A vector of GPS times at which to calculate the signal
1081 * \param hetssbdelays [in] The vector of SSB time delays used for the original heterodyne. If this is
1082 * \c NULL then the SSB delay for the given position will be calculated.
1083 * \param calcSSBDelay [in] Set to a non-zero value if the SSB delay needs to be recalculated at
1084 * an updated sky position compared to that used for the heterodyne.
1085 * \param hetbsbdelays [in] The vector of BSB time delays used for the original heterodyne. If this is
1086 * \c NULL then the BSB delay for the given binary parameters will be calculated.
1087 * \param calcBSBDelay [in] Set to a non-zero value if the BSB delay needs to be calulated at
1088 * a set of updated binary system parameters.
1089 * \param glphase [in] The vector containing the glitch phase evolution used for the original heterodyne.
1090 * If this is \c NULL then the glitch phase for the given glitch parameters will be calculated.
1091 * \param calcglphase [in] Set to a non-zero value if the glitch phase needs to be calulated at
1092 * a set of updated glitch parameters.
1093 * \param fitwavesphase [in] The vector of FITWAVES phases used for the original heterodyne. If this is
1094 * \c NULL then the FITWAVES phase for the given FITWAVES parameters will be calculated.
1095 * \param calcfitwaves [in] Set to a non-zero value if the FITWAVES phase needs to be calulated at
1096 * a set of updated FITWAVES parameters.
1097 * \param resp [in] A loop-up table for the detector response function.
1098 * \param ephem [in] A pointer to an \c EphemerisData structure containing solar system ephemeris information
1099 * \param tdat [in] A pointer to a \c TimeCorrectionData structure containing time system correction information
1100 * \param ttype [in] The \c TimeCorrectionType value
1101 *
1102 * \sa XLALHeterodynedPulsarGetAmplitudeModel
1103 * \sa XLALHeterodynedPulsarGetPhaseModel
1104 */
1106 PulsarParameters *origpars,
1107 REAL8 freqfactor,
1108 UINT4 usephase,
1109 UINT4 useroq,
1110 UINT4 nonGR,
1112 REAL8Vector *hetssbdelays,
1113 UINT4 calcSSBDelay,
1114 REAL8Vector *hetbsbdelays,
1115 UINT4 calcBSBDelay,
1116 REAL8Vector *glphase,
1117 UINT4 calcglphase,
1118 REAL8Vector *fitwavesphase,
1119 UINT4 calcfitwaves,
1120 const DetResponseTimeLookupTable *resp,
1121 const EphemerisData *ephem,
1122 const TimeCorrectionData *tdat,
1123 TimeCorrectionType ttype )
1124{
1125 UINT4 i = 0;
1126 COMPLEX16TimeSeries *csignal = NULL;
1127
1128 /* get amplitude model */
1130 freqfactor,
1131 usephase,
1132 useroq,
1133 nonGR,
1134 timestamps,
1135 resp );
1136
1137 // include phase change if required
1138 if ( usephase ) {
1139 REAL8Vector *dphi = NULL;
1141 origpars,
1142 timestamps,
1143 freqfactor,
1144 hetssbdelays,
1145 calcSSBDelay,
1146 hetbsbdelays,
1147 calcBSBDelay,
1148 glphase,
1149 calcglphase,
1150 fitwavesphase,
1151 calcfitwaves,
1152 resp->det,
1153 ephem,
1154 tdat,
1155 ttype );
1156
1157 COMPLEX16 M = 0., expp = 0.;
1158
1159 for ( i = 0; i < dphi->length; i++ ) {
1160 /* phase factor by which to multiply the (almost) DC signal model. NOTE: this does not try to undo
1161 * the signal modulation in the data, but instead replicates it in the model, hence the positive
1162 * phase rather than a negative phase in the cexp function. */
1163 expp = cexp( LAL_TWOPI * I * dphi->data[i] );
1164
1165 M = csignal->data->data[i];
1166
1167 /* heterodyne */
1168 csignal->data->data[i] = M * expp;
1169 }
1170
1171 XLALDestroyREAL8Vector( dphi );
1172 }
1173
1174 return csignal;
1175}
1176
1177
1178/**
1179 * \brief Creates a lookup table of the detector antenna pattern
1180 *
1181 * This function creates a lookup table of the tensor, vector and scalar
1182 * antenna patterns for a given detector orientation and source sky position.
1183 * For the tensor modes these are the functions given by equations 10-13 in
1184 * \cite JKS98 , whilst for the vector and scalar modes they are the \f$ \psi \f$
1185 * independent parts of e.g. equations 5-8 of \cite Nishizawa2009 . We remove
1186 * the \f$ \psi \f$ dependance by setting \f$ \psi=0 \f$ .
1187 *
1188 * If \c avedt is a value over 60 seconds then the antenna pattern will
1189 * actually be the mean value from 60 second intervals within that timespan.
1190 * This accounts for the fact that each data point is actually an averaged
1191 * value over the given timespan.
1192 *
1193 * \param t0 [in] initial GPS time of the data
1194 * \param det [in] structure containing the detector
1195 * \param alpha [in] the source right ascension in radians
1196 * \param delta [in] the source declination in radians
1197 * \param timeSteps [in] the number of grid bins to use in time
1198 * \param avedt [in] average the antenna pattern over this timespan
1199 */
1201 const LALDetector *det,
1202 REAL8 alpha,
1203 REAL8 delta,
1204 UINT4 timeSteps,
1205 REAL8 avedt )
1206{
1207 LIGOTimeGPS gps;
1208
1209 // allocate structure
1211 resp = XLALMalloc( sizeof( *resp ) );
1212 resp->det = XLALMalloc( sizeof( LALDetector * ) );
1213 memcpy( &resp->det, &det, sizeof( det ) );
1214 resp->alpha = alpha;
1215 resp->delta = delta;
1216 resp->psi = 0.;
1217 resp->ntimebins = timeSteps;
1218 resp->fplus = XLALCreateREAL8Vector( timeSteps );
1219 resp->fcross = XLALCreateREAL8Vector( timeSteps );
1220 resp->fx = XLALCreateREAL8Vector( timeSteps );
1221 resp->fy = XLALCreateREAL8Vector( timeSteps );
1222 resp->fb = XLALCreateREAL8Vector( timeSteps );
1223 resp->fl = XLALCreateREAL8Vector( timeSteps );
1224 XLALGPSSetREAL8( &resp->t0, t0 ); // set the epoch
1225
1226 REAL8 T = 0., Tstart = 0., Tav = 0., psi = 0.;
1227
1228 REAL8 fplus = 0., fcross = 0., fx = 0., fy = 0., fb = 0., fl = 0.;
1229 REAL8 tsteps = ( REAL8 )timeSteps;
1230
1231 UINT4 j = 0, k = 0, nav = 0;
1232
1233 /* number of points to average */
1234 if ( avedt == 60. ) {
1235 nav = 1;
1236 } else {
1237 nav = floor( avedt / 60. ) + 1;
1238 }
1239
1240 for ( j = 0 ; j < timeSteps ; j++ ) {
1241 resp->fplus->data[j] = 0.;
1242 resp->fcross->data[j] = 0.;
1243 resp->fx->data[j] = 0.;
1244 resp->fy->data[j] = 0.;
1245 resp->fb->data[j] = 0.;
1246 resp->fl->data[j] = 0.;
1247
1248 /* central time of lookup table point */
1249 T = t0 + ( REAL8 )j * LAL_DAYSID_SI / tsteps;
1250
1251 if ( nav % 2 ) { /* is odd */
1252 Tstart = T - 0.5 * ( ( REAL8 )nav - 1. ) * 60.;
1253 } else { /* is even */
1254 Tstart = T - ( 0.5 * ( REAL8 )nav - 1. ) * 60. - 30.;
1255 }
1256
1257 for ( k = 0; k < nav; k++ ) {
1258 Tav = Tstart + 60.*( REAL8 )k;
1259
1260 XLALGPSSetREAL8( &gps, Tav );
1261 XLALComputeDetAMResponseExtraModes( &fplus, &fcross, &fb, &fl, &fx, &fy, det->response,
1262 alpha, delta, psi,
1264
1265 resp->fplus->data[j] += fplus;
1266 resp->fcross->data[j] += fcross;
1267 resp->fx->data[j] += fx;
1268 resp->fy->data[j] += fy;
1269 resp->fb->data[j] += fb;
1270 resp->fl->data[j] += fl;
1271 }
1272
1273 resp->fplus->data[j] /= ( REAL8 )nav;
1274 resp->fcross->data[j] /= ( REAL8 )nav;
1275 resp->fx->data[j] /= ( REAL8 )nav;
1276 resp->fy->data[j] /= ( REAL8 )nav;
1277 resp->fb->data[j] /= ( REAL8 )nav;
1278 resp->fl->data[j] /= ( REAL8 )nav;
1279 }
1280
1281 return resp;
1282}
1283
1284
1285/** \brief Free memory for antenna response look-up table structure
1286 *
1287 * \param resp [in] the look-up table structure to be freed
1288 */
1290{
1291 if ( resp == NULL ) {
1292 return;
1293 }
1294
1295 if ( resp->fplus ) {
1297 }
1298 if ( resp->fcross ) {
1300 }
1301 if ( resp->fx ) {
1302 XLALDestroyREAL8Vector( resp->fx );
1303 }
1304 if ( resp->fy ) {
1305 XLALDestroyREAL8Vector( resp->fy );
1306 }
1307 if ( resp->fb ) {
1308 XLALDestroyREAL8Vector( resp->fb );
1309 }
1310 if ( resp->fl ) {
1311 XLALDestroyREAL8Vector( resp->fl );
1312 }
1313
1314 XLALFree( resp );
1315}
1316
1317
1318/**
1319 * \brief Convert source parameters into amplitude and phase notation parameters
1320 *
1321 * Convert the physical source parameters into the amplitude and phase notation
1322 * given in Eqns 62-65 of \cite Jones2015 .
1323 *
1324 * Note that \c phi0 is essentially the rotational phase of the pulsar. Also,
1325 * note that if using \f$ h_0 \f$ , and therefore the convention for a signal as
1326 * defined in \cite JKS98 , the sign of the waveform model is the opposite of
1327 * that in \cite Jones2015 , and therefore a sign flip is required in the
1328 * amplitudes.
1329 */
1331{
1332 REAL8 sinlambda, coslambda, sinlambda2, coslambda2, sin2lambda;
1333 REAL8 theta, sintheta, costheta2, sintheta2, sin2theta;
1336 REAL8 I21 = PulsarGetREAL8ParamOrZero( params, "I21" );
1337 REAL8 I31 = PulsarGetREAL8ParamOrZero( params, "I31" );
1338 REAL8 C21 = PulsarGetREAL8ParamOrZero( params, "C21" );
1339 REAL8 C22 = PulsarGetREAL8ParamOrZero( params, "C22" );
1340 REAL8 phi21 = PulsarGetREAL8ParamOrZero( params, "PHI21" );
1341 REAL8 phi22 = PulsarGetREAL8ParamOrZero( params, "PHI22" );
1342 REAL8 lambda = PulsarGetREAL8ParamOrZero( params, "LAMBDAPIN" );
1344 REAL8 q22 = PulsarGetREAL8ParamOrZero( params, "Q22" );
1345
1346 REAL8 dist = 0.; /* distance in metres */
1347 /* set distance (a DIST param takes precedence over PX) */
1348 if ( PulsarCheckParam( params, "DIST" ) ) {
1349 dist = PulsarGetREAL8Param( params, "DIST" );
1350 } else if ( PulsarCheckParam( params, "PX" ) ) {
1352 }
1353
1354 if ( ( q22 != 0. && dist != 0. ) && ( I21 == 0. && I21 == 0. && C21 == 0. && C22 == 0. ) ) {
1355 /* convert mass quadrupole to C22 parameter */
1356 if ( !PulsarCheckParam( params, "F" ) ) {
1357 XLAL_ERROR_VOID( XLAL_EINVAL, "Error... not frequency values given" );
1358 }
1359
1361 if ( f0 == 0. ) {
1362 XLAL_ERROR_VOID( XLAL_EINVAL, "Error... frequency is zero" );
1363 }
1364
1365 /* convert q22 to h0 */
1366 h0 = q22 * sqrt( 8. * LAL_PI / 15. ) * 16. * LAL_PI * LAL_PI * LAL_G_SI * f0 * f0 / ( LAL_C_SI * LAL_C_SI * LAL_C_SI * LAL_C_SI * dist );
1367
1368 C22 = -0.5 * h0;
1369 PulsarAddParam( params, "C22", &C22, PULSARTYPE_REAL8_t );
1370
1371 if ( phi22 == 0. ) {
1372 phi22 = 2.*phi0;
1373 phi22 = phi22 - LAL_TWOPI * floor( phi22 / LAL_TWOPI );
1374 PulsarAddParam( params, "PHI22", &phi22, PULSARTYPE_REAL8_t );
1375 }
1376 } else if ( h0 != 0. ) {
1377 C22 = -0.5 * h0; /* note the change in sign so that the triaxial model conforms to the convertion in JKS98 */
1378 PulsarAddParam( params, "C22", &C22, PULSARTYPE_REAL8_t );
1379
1380 if ( phi22 == 0. ) {
1381 phi22 = 2.*phi0;
1382 phi22 = phi22 - LAL_TWOPI * floor( phi22 / LAL_TWOPI );
1383 PulsarAddParam( params, "PHI22", &phi22, PULSARTYPE_REAL8_t );
1384 }
1385 } else if ( ( I21 != 0. || I31 != 0. ) && ( C22 == 0. && C21 == 0. ) ) {
1386 sinlambda = sin( lambda );
1387 coslambda = cos( lambda );
1388 sin2lambda = sin( 2. * lambda );
1389 sinlambda2 = SQUARE( sinlambda );
1390 coslambda2 = SQUARE( coslambda );
1391
1392 theta = acos( costheta );
1393 sintheta = sin( theta );
1394 sin2theta = sin( 2. * theta );
1395 sintheta2 = SQUARE( sintheta );
1396 costheta2 = SQUARE( costheta );
1397
1398 REAL8 A22 = I21 * ( sinlambda2 - coslambda2 * costheta2 ) - I31 * sintheta2;
1399 REAL8 B22 = I21 * sin2lambda * costheta;
1400 REAL8 A222 = SQUARE( A22 );
1401 REAL8 B222 = SQUARE( B22 );
1402
1403 REAL8 A21 = I21 * sin2lambda * sintheta;
1404 REAL8 B21 = sin2theta * ( I21 * coslambda2 - I31 );
1405 REAL8 A212 = SQUARE( A21 );
1406 REAL8 B212 = SQUARE( B21 );
1407
1408 C22 = 2.*sqrt( A222 + B222 );
1409 C21 = 2.*sqrt( A212 + B212 );
1410
1411 PulsarAddParam( params, "C22", &C22, PULSARTYPE_REAL8_t );
1412 PulsarAddParam( params, "C21", &C21, PULSARTYPE_REAL8_t );
1413
1414 phi22 = fmod( 2.*phi0 - atan2( B22, A22 ), LAL_TWOPI );
1415 phi21 = fmod( phi0 - atan2( B21, A21 ), LAL_TWOPI );
1416
1417 PulsarAddParam( params, "PHI22", &phi22, PULSARTYPE_REAL8_t );
1418 PulsarAddParam( params, "PHI21", &phi21, PULSARTYPE_REAL8_t );
1419 }
1420}
void XLALBinaryPulsarDeltaTNew(BinaryPulsarOutput *output, BinaryPulsarInput *input, PulsarParameters *params)
function to calculate the binary system delay using new parameter structure
REAL8Vector * XLALHeterodynedPulsarPhaseDifference(PulsarParameters *params, PulsarParameters *origparams, const LIGOTimeGPSVector *datatimes, REAL8 freqfactor, REAL8Vector *ssbdts, UINT4 calcSSBDelay, REAL8Vector *bsbdts, UINT4 calcBSBDelay, REAL8Vector *glphase, UINT4 calcglphase, REAL8Vector *fitwavesphase, UINT4 calcfitwaves, const LALDetector *detector, const EphemerisData *ephem, const TimeCorrectionData *tdat, TimeCorrectionType ttype)
The phase evolution difference compared to a heterodyned phase (for a pulsar)
void XLALPulsarSourceToWaveformParams(PulsarParameters *params)
Convert source parameters into amplitude and phase notation parameters.
COMPLEX16TimeSeries * XLALHeterodynedPulsarGetAmplitudeModel(PulsarParameters *pars, REAL8 freqfactor, UINT4 varyphase, UINT4 useroq, UINT4 nonGR, const LIGOTimeGPSVector *timestamps, const DetResponseTimeLookupTable *resp)
The amplitude model of a complex heterodyned signal from the harmonics of a rotating neutron star.
COMPLEX16TimeSeries * XLALHeterodynedPulsarGetModel(PulsarParameters *pars, PulsarParameters *origpars, REAL8 freqfactor, UINT4 usephase, UINT4 useroq, UINT4 nonGR, const LIGOTimeGPSVector *timestamps, REAL8Vector *hetssbdelays, UINT4 calcSSBDelay, REAL8Vector *hetbsbdelays, UINT4 calcBSBDelay, REAL8Vector *glphase, UINT4 calcglphase, REAL8Vector *fitwavesphase, UINT4 calcfitwaves, const DetResponseTimeLookupTable *resp, const EphemerisData *ephem, const TimeCorrectionData *tdat, TimeCorrectionType ttype)
Generate the model of the neutron star signal.
void XLALDestroyDetResponseTimeLookupTable(DetResponseTimeLookupTable *resp)
Free memory for antenna response look-up table structure.
DetResponseTimeLookupTable * XLALDetResponseLookupTable(REAL8 t0, const LALDetector *det, REAL8 alpha, REAL8 delta, UINT4 timeSteps, REAL8 avedt)
Creates a lookup table of the detector antenna pattern.
#define SQUARE(x)
Macro to square a value.
REAL8Vector * XLALHeterodynedPulsarGetSSBDelay(PulsarParameters *pars, const LIGOTimeGPSVector *datatimes, const LALDetector *detector, const EphemerisData *ephem, const TimeCorrectionData *tdat, TimeCorrectionType ttype)
Computes the delay between a GPS time at Earth and the solar system barycentre.
#define ROUND(x)
Macro to round a value to the nearest integer.
void XLALGetEarthPosVel(EarthState *earth, const EphemerisData *edat, const LIGOTimeGPS *tGPS)
Get the position and velocity of the Earth at a given time.
REAL8Vector * XLALHeterodynedPulsarGetBSBDelay(PulsarParameters *pars, const LIGOTimeGPSVector *datatimes, const REAL8Vector *dts, const EphemerisData *edat)
Computes the delay between a pulsar in a binary system and the barycentre of the system.
REAL8Vector * XLALHeterodynedPulsarGetGlitchPhase(PulsarParameters *params, const LIGOTimeGPSVector *datatimes, const REAL8Vector *ssbdts, const REAL8Vector *bsbdts)
Computes the phase evolution due to the presence of glitches.
REAL8Vector * XLALHeterodynedPulsarGetFITWAVESPhase(PulsarParameters *params, const LIGOTimeGPSVector *datatimes, const REAL8Vector *ssbdts, REAL8 freq)
Computes the phase evolution due to the presence of FITWAVES parameters.
int j
int k
static double double delta
const REAL8Vector * PulsarGetREAL8VectorParam(const PulsarParameters *pars, const CHAR *name)
Return a REAL8Vector parameter.
REAL8 PulsarGetREAL8ParamOrZero(const PulsarParameters *pars, const CHAR *name)
Return a REAL8 parameter if it exists, otherwise return zero.
int PulsarCheckParam(const PulsarParameters *pars, const CHAR *name)
Check for the existence of the parameter name in the PulsarParameters structure.
void PulsarAddParam(PulsarParameters *pars, const CHAR *name, void *value, PulsarParamType type)
Add a parameter and value to the PulsarParameters structure.
REAL8 PulsarGetREAL8Param(const PulsarParameters *pars, const CHAR *name)
Return a REAL8 parameter.
const CHAR * PulsarGetStringParam(const PulsarParameters *pars, const CHAR *name)
Return a string parameter.
REAL8 PulsarGetREAL8VectorParamIndividual(const PulsarParameters *pars, const CHAR *name)
Return an individual REAL8 value from a REAL8Vector parameter.
@ PULSARTYPE_REAL8_t
LIGOTimeGPSVector * timestamps
int l
double theta
const double wcos
const double costheta
void XLALComputeDetAMResponseExtraModes(double *fplus, double *fcross, double *fb, double *fl, double *fx, double *fy, const REAL4 D[3][3], const double ra, const double dec, const double psi, const double gmst)
int XLALBarycenter(EmissionTime *emit, const BarycenterInput *baryinput, const EarthState *earth)
Transforms from detector arrival time in GPS (as specified in the LIGOTimeGPS structure) to pulse em...
int XLALBarycenterEarthNew(EarthState *earth, const LIGOTimeGPS *tGPS, const EphemerisData *edat, const TimeCorrectionData *tdat, TimeCorrectionType ttype)
Computes the position and orientation of the Earth, at some arrival time, but unlike XLALBarycenterEa...
TimeCorrectionType
Enumerated type denoting the time system type to be produced in the solar system barycentring routine...
Definition: LALBarycenter.h:72
#define LAL_PI_2
#define LAL_DAYSID_SI
#define LAL_PI
#define LAL_TWOPI
#define LAL_G_SI
double complex COMPLEX16
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
void * XLALMalloc(size_t n)
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
COMPLEX16TimeSeries * XLALCreateCOMPLEX16TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
const LALUnit lalSecondUnit
REAL8Vector * XLALCreateREAL8Vector(UINT4 length)
void XLALDestroyREAL8Vector(REAL8Vector *vector)
#define XLAL_ERROR_VOID(...)
#define XLAL_ERROR_NULL(...)
#define XLAL_CHECK_NULL(assertion,...)
XLAL_SUCCESS
XLAL_EFUNC
XLAL_EDOM
XLAL_EINVAL
REAL8 XLALGreenwichMeanSiderealTime(const LIGOTimeGPS *gpstime)
LIGOTimeGPS * XLALGPSSetREAL8(LIGOTimeGPS *epoch, REAL8 t)
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
int gpstime
Definition: hwinject.c:365
M
list y
T
pos
dist
int T0
double alpha
double t0
Basic input structure to LALBarycenter.c.
REAL8 alpha
Source right ascension in ICRS J2000 coords (radians).
REAL8 dInv
1/(distance to source), in 1/sec.
LIGOTimeGPS tgps
input GPS arrival time.
REAL8 delta
Source declination in ICRS J2000 coords (radians)
LALDetector site
detector site info.
structure containing the input parameters for the binary delay function
REAL8 tb
Time of arrival (TOA) at the SSB.
EarthState earth
The current Earth state (for e.g.
structure containing the output parameters for the binary delay function
REAL8 deltaT
deltaT to add to TDB in order to account for binary
COMPLEX16Sequence * data
COMPLEX16 * data
Pulsar signal-generation routines for heterodyned data.
REAL8 alpha
the source right ascension in radians
REAL8Vector * fb
scalar breathing mode polarisation response
REAL8 delta
the source declination in radians
REAL8Vector * fl
scalar longitudinal mode polarisation response
REAL8 psi
the polarisation angle in radians
UINT4 ntimebins
number of time bins in the look-up table
REAL8Vector * fy
vector "y" polarisation response
REAL8Vector * fx
vector "x" polarisation response
REAL8Vector * fcross
tensor cross polarisation response
LIGOTimeGPS t0
GPS time epoch for the look-up table.
REAL8Vector * fplus
tensor plus polarisation response
Basic output structure of LALBarycenterEarth.c.
REAL8 velNow[3]
dimensionless velocity of Earth's center at tgps, extrapolated from JPL DE405 ephemeris
REAL8 posNow[3]
Cartesian coords of Earth's center at tgps, extrapolated from JPL DE405 ephemeris; units= sec.
Basic output structure produced by LALBarycenter.c.
REAL8 deltaT
(TDB) - (GPS)
REAL8 roemer
the Roemer delay
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
REAL4 response[3][3]
REAL8 location[3]
INT4 gpsNanoSeconds
A vector of 'timestamps' of type LIGOTimeGPS.
Definition: SFTfileIO.h:188
LIGOTimeGPS * data
array of timestamps
Definition: SFTfileIO.h:193
UINT4 length
number of timestamps
Definition: SFTfileIO.h:192
The PulsarParameters structure to contain a set of pulsar parameters.
REAL8 * data
This structure will contain a vector of time corrections used during conversion from TT to TDB/TCB/Te...