LALPulsar 7.1.2.1-bf6a62b
ReadPulsarParFile.h
Go to the documentation of this file.
1/*
2* Copyright (C) 2013 Bernd Machenschalk, Jolien Creighton, Matt 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/**
21 * \defgroup ReadPulsarParFile_h Header ReadPulsarParFile.h
22 * \ingroup lalpulsar_general
23 * \author Matt Pitkin
24 * \date 2013
25 * \brief Functions to read <tt>TEMPO(2)</tt> pulsar parameter files
26 *
27 * Radio astronomers fit pulsar parameters using <tt>TEMPO(2)</tt> which will output
28 * the parameters in a <tt>.par</tt> file. The values allowed in this file can be
29 * found in the TEMPO documentation.
30 *
31 * The function XLALReadTEMPOParFile() reads the parameters into a linked list
32 * #PulsarParameters structure, from which the parameters can be accessed
33 * using the appropriate access function. These use a hash table to quick
34 * look-up. The parameters are assigned names, which are used as the hash table
35 * keys, which are fully uppercase versions of the TEMPO parameter names.
36 *
37 * All parameters read in are converted into SI units.
38 *
39 * Functions are is also included which converts a string containing the right ascension or
40 * declination in the format <tt>ddd/hh:mm:ss.s</tt> or <tt>ddd/hhmmss.s</tt>
41 * (as is given in the <tt>.par</tt> file) into a \c REAL8 value in
42 * radians.
43 *
44*/
45/** @{ */
46
47#ifndef _READPULSARPARFILE_H
48#define _READPULSARPARFILE_H
49
50#include <ctype.h>
51
52#include <lal/LALStdlib.h>
53#include <lal/StringVector.h>
54#include <lal/LALBarycenter.h>
55#include <lal/Date.h>
56#include <lal/LALHashTbl.h>
57
58#ifdef __cplusplus
59extern "C" {
60#endif
61
62/**\name Error Codes */ /** @{ */
63#define READPULSARPARFILEH_ENULLOUTPUT 1
64
65#define READPULSARPARFILEH_MSGENULLOUTPUT "Output was Null"
66/** @} */
67
68#define PULSAR_HASHTABLE_SIZE 512
69#define PULSAR_PARNAME_MAX 128
70
71#define GPS0MJD 44244.0 /* start of GPS time (MJD 44244) */
72/* the difference between TDT (or TT) and TAI (see e.g. Eqn 15 of T.-Y. Huang et
73 * al, A&A, 220, p. 329, 1989) */
74#define TDT_TAI 32.184
75/* the difference between TDT/TT and the GPS epoch */
76#define GPS_TDT (TDT_TAI + XLAL_EPOCH_GPS_TAI_UTC)
77
78/* to do conversions consistent with adopted constants in TEMPO2 */
79/* see: https://bitbucket.org/psrsoft/tempo2/src/fb2279f50513573e1c6cd9dacb582c4b6016fc92/DDSmodel.C#lines-44 */
80#define LALPULSAR_TEMPO2_MTSUN_SI 4.925490947e-6 /* value of GMsun/c^3 used in TEMPO2 */
81#define LALPULSAR_TEMPO2_MSUN_SI (LALPULSAR_TEMPO2_MTSUN_SI * LAL_MPL_SI / LAL_TPL_SI) /* derived value of Msun */
82
83/** An enumerated type for denoting the type of a variable. Several LAL types are supported. */
84typedef enum tagPulsarParamType {
91
92extern size_t PulsarTypeSize[5];
93
94
95/** \brief The PulsarParam list node structure
96 *
97 * This structure contains a pulsar parameter defined by the name.
98 *
99 * This should only be accessed using the accessor functions below.
100*/
101typedef struct tagPulsarParam {
102 CHAR name[PULSAR_PARNAME_MAX]; /**< Parameter name */
103 void *value; /**< Parameter value */
104 void *err; /**< Parameter error/uncertainty */
105 UINT4Vector *fitFlag; /**< Set to 1 if the parameter has been fit in the par file */
106 PulsarParamType type; /**< Parameter type e.g. REAL8, CHAR, INT4 */
107 struct tagPulsarParam *next;
109
110
111/** \brief The PulsarParameters structure to contain a set of pulsar parameters
112 *
113 * This is implemented as a linked list of PulsarParam structures and should
114 * only be accessed using the accessor functions below. This is intended as a replacement to the
115 * #BinaryPulsarParams structure.
116 */
117typedef struct tagPulsarParameters {
118 PulsarParam *head; /**< A linked list of #PulsarParam structures */
119 INT4 nparams; /**< The total number of parameters in the structure */
120 LALHashTbl *hash_table; /**< Hash table of parameters */
122
123
124/** A structure to contain all pulsar parameters and associated errors.
125 The structure does not have to be used for a binary pulsar, but can just contain the parameters for
126 an isolated pulsar. All parameters are in the same units as given by TEMPO.
127 */
128typedef struct tagBinaryPulsarParams {
129 CHAR *name; /**< pulsar name */
130 CHAR *jname; /**< pulsar J name */
131 CHAR *bname; /**< pulsar B name */
132
133 CHAR *model; /**< TEMPO binary model e.g. BT, DD, ELL1 */
134
135 CHAR *units; /**< The time system used e.g. TDB */
136 CHAR *ephem; /**< The JPL solar system ephemeris used e.g. DE405 */
137
138 REAL8 f0; /**< spin frequency (Hz) */
139 REAL8 f1; /**< frequency first derivative (Hz/s) */
140 REAL8 f2; /**< frequency second derivative (Hz/s^2) */
141 REAL8 f3; /**< frequency third derivative (Hz/s^3) */
142 REAL8 f4; /**< frequency fourth derivative (Hz/s^4) */
143 REAL8 f5; /**< frequency fifth derivative (Hz/s^5) */
144 REAL8 f6; /**< frequency sixth derivative (Hz/s^6) */
145 REAL8 f7; /**< frequency seventh derivative (Hz/s^7) */
146 REAL8 f8; /**< frequency eighth derivative (Hz/s^8) */
147 REAL8 f9; /**< frequency ninth derivative (Hz/s^9) */
148 REAL8 f10; /**< frequency tenth derivative (Hz/s^10) */
149
150 REAL8 ra; /**< right ascension (rads) */
151 REAL8 dec; /**< declination (rads) */
152 REAL8 pmra; /**< proper motion in RA (rads/s) */
153 REAL8 pmdec; /**< proper motion in dec (rads/s) */
154
155 REAL8 posepoch; /**< position epoch */
156 REAL8 pepoch; /**< period/frequency epoch */
157
158 REAL8 startTime; /**< start of parfile applicable time */
159 REAL8 finishTime; /**< finish of parfile applicable time */
160
161 /* all parameters will be in the same units as used in TEMPO */
162
163 /* Keplerian parameters */
164 REAL8 e; /**< orbital eccentricity */
165 REAL8 Pb; /**< orbital period (seconds) */
166 REAL8 w0; /**< logitude of periastron (radians) */
167 REAL8 x; /**< projected semi-major axis/speed of light (light secs) */
168 REAL8 T0; /**< time of orbital perisastron as measured in TDB (MJD) */
169
170 /* add extra parameters for the BT1P and BT2P models which contain two and
171 three orbits respectively (only the first one can be relativistic, so the
172 others only have the Keplerian parameters) */
174 REAL8 Pb2, Pb3;
175 REAL8 w02, w03;
177 REAL8 T02, T03;
178
180
181 /* for low eccentricity orbits (ELL1 model) use Laplace parameters */
182 /* (eps1 = e*sin(w), eps2 = e*cos(w)) instead of e, w. */
183 /* see Appendix A, Ch. Lange etal, MNRAS (2001) */
184 REAL8 eps1; /**< e*sin(w) */
185 REAL8 eps2; /**< e*cos(w) */
188 REAL8 Tasc; /**< time of the ascending node (used rather than T0) */
189 UINT4 nEll; /**< set to zero if have eps time derivs (default) set to 1 if have wdot */
190
191 /* Post-Keplarian parameters */
192 /* for Blandford-Teukolsky (BT) model */
193 REAL8 wdot; /**< precesion of longitude of periastron w = w0 + wdot(tb-T0) (degs/year) */
194 REAL8 gamma; /**< gravitational redshift and time dilation parameter (s)*/
195 REAL8 Pbdot; /**< rate of change of Pb (dimensionless) */
196 REAL8 xdot; /**< rate of change of x(=asini/c) - optional */
197 REAL8 edot; /**< rate of change of e */
198
199 /* for Epstein-Haugan (EH) model */
200 REAL8 s; /**< Shapiro 'shape' parameter sin i */
202
203 /* for Damour-Deruelle (DD) model */
204 /*REAL8 r; Shapiro 'range' parameter - defined internally as Gm2/c^3 */
206 REAL8 dth; /* (10^-6) */
207 REAL8 a0, b0; /**< abberation delay parameters */
208
209 /* for DD (General Relativity) (DDGR) - assumes GR is correct model */
210 /* we do not need wdot, gamma, Pbdot, s, r, xdot and edot */
211 REAL8 M; /**< M = m1 + m2 (m1 = pulsar mass, m2 = companion mass) (in kg) */
212 REAL8 m2; /**< companion mass (in kg) */
213
214 /* for the DDS model we need the shapmax parameter */
216
217 /* orbital frequency coefficients for BTX model (only for one orbit at the
218 moment i.e. a two body system) */
219 REAL8 *fb; /**< orbital frequency coefficients for BTX model */
220 INT4 nfb; /**< the number of fb coefficients */
221
222 REAL8 px; /**< pulsar parallax (converted from milliarcsecs to rads) */
223 REAL8 dist; /**< pulsar distance (in kiloparsecs) */
224
225 REAL8 DM; /**< dispersion measure */
226 REAL8 DM1; /**< first derivative of dispersion measure */
227
228 /* sinusoidal parameters for fitting quasi-periodic timing noise */
231 REAL8 wave_om; /**< fundamental frequency timing noise terms */
234
235 /* gravitational wave parameters */
236 REAL8 h0; /**< gravitational wave amplitude */
237 REAL8 Q22; /**< gravitational wave l=m=2 mass quadrupole moment */
238 REAL8 cosiota;/**< cosine of the pulsars inclination angle */
239 REAL8 iota; /**< inclination angle */
240 REAL8 psi; /**< polarisation angle */
241 REAL8 phi0; /**< initial phase */
242 REAL8 Aplus; /**< 0.5*h0*(1+cos^2iota) */
243 REAL8 Across; /**< h0*cosiota */
244
245 /* scalar-tensor-vector mode amplitude parameters */
246 REAL8 hPlus; /**< amplitude for the tensor plus polarisation */
247 REAL8 hCross; /**< amplitude for the tensor cross polarisation */
248 REAL8 psiTensor; /**< phase polarisation angle for tensor polarisation mode */
249 REAL8 phi0Tensor; /**< initial phase for tensor modes (to be used instead of phi0 or phi22) */
250 REAL8 hScalarB; /**< amplitude for scalar breathing polarisation mode */
251 REAL8 hScalarL; /**< amplitude for scalar longitudinal polarisation mode */
252 REAL8 psiScalar; /**< phase polarisation angle for scalar polarisation mode */
253 REAL8 phi0Scalar; /**< initial phase for the scalar modes */
254 REAL8 hVectorX; /**< amplitude for vector "x" polarisation mode */
255 REAL8 hVectorY; /**< amplitude for vector "y" polarisation mode */
256 REAL8 psiVector; /**< phase polarisation angle for vector polarisation mode */
257 REAL8 phi0Vector; /**< initial phase for vector modes */
258
259 /* pinned superfluid gw parameters*/
260 REAL8 I21; /**< parameter for pinsf model.**/
261 REAL8 I31; /**< parameter for pinsf model.**/
262 REAL8 lambdapin; /**< this is a longitude like angle between pinning axis and line of sight */
263 REAL8 costheta; /**< cosine of angle between rotation axis and pinning axis */
265
266 /* complex amplitude and phase parameters for l=2, m=1 and 2 harmonics */
271
272 /* parameters for Kopeikin terms */
273 REAL8 daop; /**< parameter for the Kopeikin annual orbital parallax */
274 INT4 daopset; /**< set if daop is set from the par file */
275 REAL8 kin; /**< binary inclination angle */
279
280 /******** errors read in from a .par file **********/
292
295
300
301 REAL8 eErr, e2Err, e3Err;
302 REAL8 PbErr, Pb2Err, Pb3Err;
303 REAL8 w0Err, w02Err, w03Err;
304 REAL8 xErr, x2Err, x3Err;
305 REAL8 T0Err, T02Err, T03Err;
306
308
314
320
323
324 /*REAL8 rErr; Shapiro 'range' parameter - defined internally as Gm2/c^3 */
327 REAL8 a0Err, b0Err;
328
331
333
336
339
340 /* gravitational wave parameters */
358
371
372 /* timing noise fitting parameters */
374
375 REAL8 cgw; /**< The speed of gravitational waves as a fraction of the speed of light <tt>LAL_C_SI</tt> */
378
379
380/* DEFINE FUNCTIONS */
381/** \brief Get the required parameter value from the #PulsarParameters structure
382 *
383 * This function will return a void pointer to the parameter value. This should be cast into the required
384 * variable type once returned.
385 */
386void *PulsarGetParam( const PulsarParameters *pars, const CHAR *name );
387
388/** \brief Get the required parameter error value from the #PulsarParameters structure
389 *
390 * This function will return a void pointer to the parameter error value. The parameter error will be either
391 * a \c REAL8 or a \c REAL8Vector and should be cast accordingly once returned.
392 */
393void *PulsarGetParamErr( const PulsarParameters *pars, const CHAR *name );
394
395/** \brief Get the fit flag array for a given parameter from the #PulsarParameters structure
396 *
397 * This function will return a \c UINT4 array to the parameter fit flag.
398 */
399const UINT4 *PulsarGetParamFitFlag( const PulsarParameters *pars, const CHAR *name );
400
401/** \brief Get the fit flag array for a given parameter from the #PulsarParameters structure
402 *
403 * This function will return a \c UINT4Vector array to the parameter fit flag.
404 */
406
407
408/** \brief Return a \c REAL8 parameter error value
409 *
410 * This function will call PulsarGetParamErr() for a \c REAL8 parameter and properly cast it for returning.
411 */
413
414/** \brief Return a \c REAL8Vector parameter error value
415 *
416 * This function will call PulsarGetParamErr() for a \c REAL8Vector parameter and properly cast it for returning.
417 */
419
420/** \brief Return an individual \c REAL8 value from a \c REAL8Vector parameter
421 *
422 * This function will call PulsarGetParam() for a \c REAL8Vector parameter's errors and then return a specific
423 * value from within the vector. The input \c name must be of the form e.g. \c FB1, which will get the \c REAL8Vector
424 * for the \c FB parameter and return the value from the \c 1 index.
425 */
427
428/** \brief Get the required parameter's type
429 */
430PulsarParamType PulsarGetParamType( const PulsarParameters *pars, const char *name );
431
432/** \brief Return a \c REAL8 parameter
433 *
434 * This function will call PulsarGetParam() for a \c REAL8 parameter and properly cast it for returning.
435 */
436REAL8 PulsarGetREAL8Param( const PulsarParameters *pars, const CHAR *name );
437
438/** \brief Return a \c REAL8 parameter if it exists, otherwise return zero
439 */
441
442/** \brief Return a \c UINT4 parameter
443 *
444 * This function will call PulsarGetParam() for a \c UINT4 parameter and properly cast it for returning.
445 */
446UINT4 PulsarGetUINT4Param( const PulsarParameters *pars, const CHAR *name );
447
448/** \brief Return a \c UINT4 parameter if it exists, otherwise return zero
449 */
451
452#ifdef SWIG /* SWIG interface directives */
453SWIGLAL( OWNS_THIS_ARG( const CHAR *, value ) );
454#endif
455
456/** \brief Return a string parameter
457 *
458 * This function will call PulsarGetParam() for a string parameter and properly cast it for returning.
459 * The return value should be copied e.g. with
460 * CHAR *str = XLALStringDuplicate( PulsarGetStringParam(pars, "NAME") );
461 * It also needs to be freed afterwards.
462 */
463const CHAR *PulsarGetStringParam( const PulsarParameters *pars, const CHAR *name );
464
465/** \brief Add a string parameter to the #PulsarParameters structure */
466void PulsarAddStringParam( PulsarParameters *pars, const CHAR *name, const CHAR *value );
467
468#ifdef SWIG /* SWIG interface directives */
469SWIGLAL_CLEAR( OWNS_THIS_ARG( const CHAR *, value ) );
470#endif
471
472#ifdef SWIG /* SWIG interface directives */
473SWIGLAL( OWNS_THIS_ARG( const REAL8Vector *, value ) );
474#endif
475
476/** \brief Return a \c REAL8Vector parameter
477 *
478 * This function will call PulsarGetParam() for a \c REAL8Vector parameter and properly cast it for returning.
479 */
481
482/** \brief Add a \c REAL8Vector parameter to the #PulsarParameters structure */
483void PulsarAddREAL8VectorParam( PulsarParameters *pars, const CHAR *name, const REAL8Vector *value );
484
485#ifdef SWIG /* SWIG interface directives */
486SWIGLAL_CLEAR( OWNS_THIS_ARG( const REAL8Vector *, value ) );
487#endif
488
489/** \brief Return an individual \c REAL8 value from a \c REAL8Vector parameter
490 *
491 * This function will call PulsarGetParam() for a \c REAL8Vector parameter and then return a specific value
492 * from within the vector. The input \c name must be of the form e.g. \c FB1, which will get the \c REAL8Vector
493 * for the \c FB parameter and return the value from the \c 1 index.
494 */
496
497/** \brief Add a parameter and value to the #PulsarParameters structure
498 *
499 * This function adds a new parameter, and associated value to the #PulsarParameters structure. If the parameter
500 * already exists then the old value is replaced with the new value.
501 */
502void PulsarAddParam( PulsarParameters *pars, const CHAR *name, void *value, PulsarParamType type );
503
504/** \brief Add a \c REAL8 parameter to the #PulsarParameters structure */
505void PulsarAddREAL8Param( PulsarParameters *pars, const CHAR *name, REAL8 value );
506
507/** \brief Add a \c UINT4 parameter to the #PulsarParameters structure */
508void PulsarAddUINT4Param( PulsarParameters *pars, const CHAR *name, UINT4 value );
509
510/** \brief Free all the parameters from a #PulsarParameters structure */
512
513/** \brief Remove a given parameter from a #PulsarParameters structure */
514void PulsarRemoveParam( PulsarParameters *pars, const CHAR *name );
515
516/** \brief Set the value of a parameter in the #PulsarParameters structure
517 *
518 * Set the value of the parameter given by \c name in the #PulsarParameters structure. The parameter must already
519 * exist in the structure, otherwise it should be added using PulsarAddParam().
520 */
521void PulsarSetParam( PulsarParameters *pars, const CHAR *name, const void *value );
522
523/** \brief Set the value of the error of a parameter in the #PulsarParameters structure
524 *
525 * Set the value of the error on the parameter given by \c name in the #PulsarParameters structure.
526 * The parameter must already exist in the structure, otherwise it should be added using PulsarAddParam().
527 * If the .par file contains a 1 before the error value (because that parameter has been included in the
528 * <tt>TEMPO(2)</tt> fitting procedure) then that must be input as the fit flag (this can be a vector for e.g. FB
529 * values with multiple parameters, in which case \c nfits will be the number of values in that vector).
530 */
531void PulsarSetParamErr( PulsarParameters *pars, const CHAR *name, void *value, const UINT4 *fitFlag, UINT4 len );
532
533/** \brief Set the error value for a \c REAL8 parameter */
534void PulsarSetREAL8ParamErr( PulsarParameters *pars, const CHAR *name, REAL8 value, UINT4 fitFlag );
535
536/** \brief Set the error values for a \c REAL8Vector parameter */
537void PulsarSetREAL8VectorParamErr( PulsarParameters *pars, const CHAR *name, const REAL8Vector *value, const UINT4 *fitFlag );
538
539/** \brief Check for the existence of the parameter \c name in the #PulsarParameters structure */
540int PulsarCheckParam( const PulsarParameters *pars, const CHAR *name );
541
542/** \brief Function to free memory from pulsar parameters */
544
545/** \brief Function to copy a #PulsarParameters structure */
546void PulsarCopyParams( PulsarParameters *origin, PulsarParameters *target );
547
548/** Conversion functions from units used in TEMPO parameter files into SI units */
549
550/** Convert the input string into a double precision floating point number */
551void ParConvToFloat( const CHAR *in, void *out );
552/** Convert the input string into a unsigned integer number */
553void ParConvToInt( const CHAR *in, void *out );
554/** Copy the input string into the output pointer */
555void ParConvToString( const CHAR *in, void *out );
556/** Convert the input string from degrees to radians */
557void ParConvDegsToRads( const CHAR *in, void *out );
558/** Convert the input string from milliarcsecs per year to radians per second */
559void ParConvMasPerYrToRadPerSec( const CHAR *in, void *out );
560/** Convert the input string from seconds to radians */
561void ParConvSecsToRads( const CHAR *in, void *out );
562/** Convert the input string from arcseconds to radians */
563void ParConvArcsecsToRads( const CHAR *in, void *out );
564/** Convert the input string from milliarcsecs to radians */
565void ParConvMasToRads( const CHAR *in, void *out );
566/** Convert the input string from 1/acrseconds to 1/radians */
567void ParConvInvArcsecsToInvRads( const CHAR *in, void *out );
568/** Convert the input string from days to seconds */
569void ParConvDaysToSecs( const CHAR *in, void *out );
570/** Convert the input string from kiloparsecs to metres */
571void ParConvKpcToMetres( const CHAR *in, void *out );
572
573/** Convert the binary system parameter from a string to a double, but make the check (as performed by TEMPO2)
574 * that this is > 1e-7 then it's in units of 1e-12, so needs converting by that factor. It also checks if the
575 * number is too large (> 10000) and if so sets to to zero.
576 */
577void ParConvBinaryUnits( const CHAR *in, void *out );
578/** Convert the input string from a TT MJD value into a GPS time */
579void ParConvMJDToGPS( const CHAR *in, void *out );
580/** Convert the input string from degrees per year to radians per second */
581void ParConvDegPerYrToRadPerSec( const CHAR *in, void *out );
582/** Convert the input string from solar masses to kilograms */
583void ParConvSolarMassToKg( const CHAR *in, void *out );
584/** Convert a right ascension input string in the format "hh:mm:ss.s" into radians */
585void ParConvRAToRads( const CHAR *in, void *out );
586/** Convert a declination input string in the format "dd:mm:ss.s" into radians */
587void ParConvDecToRads( const CHAR *in, void *out );
588/** Convert an input string from microseconds into seconds */
589void ParConvMicrosecToSec( const CHAR *in, void *out );
590
591
592/** \brief Read in the parameters from a <tt>TEMPO(2)</tt> parameter file into a #PulsarParameters structure
593 *
594 * This function will read in a <tt>TEMPO(2)</tt> parameter file into a #PulsarParameters structure. The structure of this
595 * function is similar to those in the <tt>TEMPO2</tt> code \c readParfile.C.
596 *
597 * \param pulsarAndPath [in] The path to the pulsar parameter file
598 */
599PulsarParameters *XLALReadTEMPOParFile( const CHAR *pulsarAndPath );
600
601/** \brief This function will read in a TEMPO-style parameter correlation matrix
602 *
603 * This function will read in a TEMPO-style parameter correlation matrix file,
604 * which contains the correlations between parameters as fit by TEMPO. An
605 * example the format would be:
606\verbatim
607 RA DEC F0
608RA 1.000
609DEC 0.954 1.000
610F0 -0.007 0.124 1.000
611\endverbatim
612 *
613 * In the output all parameter names will sometimes be
614 * converted to a more convenient naming convention. If non-diagonal parameter
615 * correlation values are +/-1 then they will be converted to be +/-0.99999 to
616 * avoid some problems of the matrix becoming singular. The output matrix will
617 * have both the upper and lower triangle completed.
618 *
619 * \param cormat [out] A REAL8 array into which the correlation matrix will be output
620 * \param corfile [in] A string containing the path and filename of the
621 * TEMPO-style correlation matrix file
622 */
624
625/** function to print out all the pulsar parameters read in from a par file */
627
628
629/** Functions for converting times given in Terrestrial time TT, TDB, or TCB in
630 * MJD to times in GPS - this is important for epochs given in <tt>.par</tt>
631 * files which are in TDB(MJD). TT and GPS are different by a factor of 51.184
632 * secs, this is just the historical factor of 32.184 secs between TT and TAI
633 * (International Atomic Time) and the other 19 seconds come from the leap
634 * seonds added between the TAI and UTC up to the point of definition of GPS
635 * time at UTC 01/01/1980.
636 */
637/** @{ */
638REAL8
639XLALTTMJDtoGPS( REAL8 MJD );
640
641REAL8
643
644REAL8
646/** @} */
647
648/** @} */
649
650#ifdef __cplusplus
651}
652#endif
653
654#endif /* _READPULSARPARFILE_H */
const double b0
const char * name
Definition: SearchTiming.c:93
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
void PulsarSetParam(PulsarParameters *pars, const CHAR *name, const void *value)
Set the value of a parameter in the PulsarParameters structure.
void ParConvKpcToMetres(const CHAR *in, void *out)
Convert the input string from kiloparsecs to metres.
const REAL8Vector * PulsarGetREAL8VectorParam(const PulsarParameters *pars, const CHAR *name)
Return a REAL8Vector parameter.
void ParConvToString(const CHAR *in, void *out)
Copy the input string into the output pointer.
void ParConvDegsToRads(const CHAR *in, void *out)
Convert the input string from degrees to radians.
void ParConvInvArcsecsToInvRads(const CHAR *in, void *out)
Convert the input string from 1/acrseconds to 1/radians.
void ParConvDegPerYrToRadPerSec(const CHAR *in, void *out)
Convert the input string from degrees per year to radians per second.
void PulsarSetREAL8ParamErr(PulsarParameters *pars, const CHAR *name, REAL8 value, UINT4 fitFlag)
Set the error value for a REAL8 parameter.
void PulsarAddREAL8Param(PulsarParameters *pars, const CHAR *name, REAL8 value)
Add a REAL8 parameter to the PulsarParameters structure.
void PulsarCopyParams(PulsarParameters *origin, PulsarParameters *target)
Function to copy a PulsarParameters structure.
const REAL8Vector * PulsarGetREAL8VectorParamErr(const PulsarParameters *pars, const CHAR *name)
Return a REAL8Vector parameter error value.
void ParConvSecsToRads(const CHAR *in, void *out)
Convert the input string from seconds to radians.
REAL8 PulsarGetREAL8VectorParamErrIndividual(const PulsarParameters *pars, const CHAR *name)
Return an individual REAL8 value from a REAL8Vector parameter.
void ParConvBinaryUnits(const CHAR *in, void *out)
Convert the binary system parameter from a string to a double, but make the check (as performed by TE...
REAL8 PulsarGetREAL8ParamErr(const PulsarParameters *pars, const CHAR *name)
Return a REAL8 parameter error value.
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 PulsarSetREAL8VectorParamErr(PulsarParameters *pars, const CHAR *name, const REAL8Vector *value, const UINT4 *fitFlag)
Set the error values for a REAL8Vector parameter.
REAL8 XLALTCBMJDtoGPS(REAL8 MJD)
If you have an MJD arrival time on the Earth then this will convert it to the equivalent GPS time in ...
void * PulsarGetParam(const PulsarParameters *pars, const CHAR *name)
Get the required parameter value from the PulsarParameters structure.
PulsarParameters * XLALReadTEMPOParFile(const CHAR *pulsarAndPath)
Read in the parameters from a TEMPO(2) parameter file into a PulsarParameters structure.
void PulsarClearParams(PulsarParameters *pars)
Free all the parameters from a PulsarParameters structure.
void PulsarSetParamErr(PulsarParameters *pars, const CHAR *name, void *value, const UINT4 *fitFlag, UINT4 len)
Set the value of the error of a parameter in the PulsarParameters structure.
LALStringVector * XLALReadTEMPOCorFile(REAL8Array *cormat, CHAR *corfile)
This function will read in a TEMPO-style parameter correlation matrix.
void ParConvDaysToSecs(const CHAR *in, void *out)
Convert the input string from days to seconds.
void PrintPulsarParameters(BinaryPulsarParams params)
function to print out all the pulsar parameters read in from a par file
void PulsarAddParam(PulsarParameters *pars, const CHAR *name, void *value, PulsarParamType type)
Add a parameter and value to the PulsarParameters structure.
void ParConvSolarMassToKg(const CHAR *in, void *out)
Convert the input string from solar masses to kilograms.
REAL8 PulsarGetREAL8Param(const PulsarParameters *pars, const CHAR *name)
Return a REAL8 parameter.
void PulsarRemoveParam(PulsarParameters *pars, const CHAR *name)
Remove a given parameter from a PulsarParameters structure.
REAL8 XLALTDBMJDtoGPS(REAL8 MJD)
If you have an MJD arrival time on the Earth then this will convert it to the equivalent GPS time in ...
UINT4 PulsarGetUINT4Param(const PulsarParameters *pars, const CHAR *name)
Return a UINT4 parameter.
const CHAR * PulsarGetStringParam(const PulsarParameters *pars, const CHAR *name)
Return a string parameter.
REAL8 XLALTTMJDtoGPS(REAL8 MJD)
This function converts a MJD format time corrected to Terrestrial Time (TT) into an equivalent GPS ti...
void PulsarAddStringParam(PulsarParameters *pars, const CHAR *name, const CHAR *value)
Add a string parameter to the PulsarParameters structure.
size_t PulsarTypeSize[5]
PulsarParamType
An enumerated type for denoting the type of a variable.
void ParConvMicrosecToSec(const CHAR *in, void *out)
Convert an input string from microseconds into seconds.
void PulsarAddREAL8VectorParam(PulsarParameters *pars, const CHAR *name, const REAL8Vector *value)
Add a REAL8Vector parameter to the PulsarParameters structure.
void * PulsarGetParamErr(const PulsarParameters *pars, const CHAR *name)
Get the required parameter error value from the PulsarParameters structure.
REAL8 PulsarGetREAL8VectorParamIndividual(const PulsarParameters *pars, const CHAR *name)
Return an individual REAL8 value from a REAL8Vector parameter.
void ParConvMasPerYrToRadPerSec(const CHAR *in, void *out)
Convert the input string from milliarcsecs per year to radians per second.
const UINT4 * PulsarGetParamFitFlag(const PulsarParameters *pars, const CHAR *name)
Get the fit flag array for a given parameter from the PulsarParameters structure.
void ParConvMasToRads(const CHAR *in, void *out)
Convert the input string from milliarcsecs to radians.
const UINT4Vector * PulsarGetParamFitFlagAsVector(const PulsarParameters *pars, const CHAR *name)
Get the fit flag array for a given parameter from the PulsarParameters structure.
void PulsarFreeParams(PulsarParameters *pars)
Function to free memory from pulsar parameters.
void ParConvToInt(const CHAR *in, void *out)
Convert the input string into a unsigned integer number.
void ParConvMJDToGPS(const CHAR *in, void *out)
Convert the input string from a TT MJD value into a GPS time.
void ParConvRAToRads(const CHAR *in, void *out)
Convert a right ascension input string in the format "hh:mm:ss.s" into radians.
void ParConvToFloat(const CHAR *in, void *out)
Conversion functions from units used in TEMPO parameter files into SI units.
PulsarParamType PulsarGetParamType(const PulsarParameters *pars, const char *name)
Get the required parameter's type.
void ParConvDecToRads(const CHAR *in, void *out)
Convert a declination input string in the format "dd:mm:ss.s" into radians.
#define PULSAR_PARNAME_MAX
void ParConvArcsecsToRads(const CHAR *in, void *out)
Convert the input string from arcseconds to radians.
UINT4 PulsarGetUINT4ParamOrZero(const PulsarParameters *pars, const CHAR *name)
Return a UINT4 parameter if it exists, otherwise return zero.
void PulsarAddUINT4Param(PulsarParameters *pars, const CHAR *name, UINT4 value)
Add a UINT4 parameter to the PulsarParameters structure.
@ PULSARTYPE_REAL8Vector_t
@ PULSARTYPE_void_ptr_t
@ PULSARTYPE_UINT4_t
@ PULSARTYPE_string_t
@ PULSARTYPE_REAL8_t
A structure to contain all pulsar parameters and associated errors.
REAL8 phi0Vector
initial phase for vector modes
REAL8 iota
inclination angle
REAL8 w0
logitude of periastron (radians)
REAL8 hPlus
amplitude for the tensor plus polarisation
REAL8 hCross
amplitude for the tensor cross polarisation
REAL8 px
pulsar parallax (converted from milliarcsecs to rads)
REAL8 m2
companion mass (in kg)
REAL8 DM
dispersion measure
REAL8 f7
frequency seventh derivative (Hz/s^7)
REAL8 psiVector
phase polarisation angle for vector polarisation mode
REAL8 f6
frequency sixth derivative (Hz/s^6)
REAL8 pepoch
period/frequency epoch
REAL8 daop
parameter for the Kopeikin annual orbital parallax
REAL8 edot
rate of change of e
REAL8 f5
frequency fifth derivative (Hz/s^5)
REAL8 Q22
gravitational wave l=m=2 mass quadrupole moment
REAL8 x
projected semi-major axis/speed of light (light secs)
CHAR * name
pulsar name
REAL8 f9
frequency ninth derivative (Hz/s^9)
REAL8 e
orbital eccentricity
REAL8 lambdapin
this is a longitude like angle between pinning axis and line of sight
REAL8 f2
frequency second derivative (Hz/s^2)
REAL8 hScalarL
amplitude for scalar longitudinal polarisation mode
CHAR * ephem
The JPL solar system ephemeris used e.g.
INT4 nfb
the number of fb coefficients
REAL8 Aplus
0.5*h0*(1+cos^2iota)
UINT4 nEll
set to zero if have eps time derivs (default) set to 1 if have wdot
REAL8 M
M = m1 + m2 (m1 = pulsar mass, m2 = companion mass) (in kg)
REAL8 Across
h0*cosiota
REAL8 ra
right ascension (rads)
REAL8 s
Shapiro 'shape' parameter sin i.
REAL8 I21
parameter for pinsf model.
REAL8 xdot
rate of change of x(=asini/c) - optional
REAL8 phi0
initial phase
REAL8 f3
frequency third derivative (Hz/s^3)
REAL8 f8
frequency eighth derivative (Hz/s^8)
REAL8 I31
parameter for pinsf model.
REAL8 hVectorY
amplitude for vector "y" polarisation mode
REAL8 hScalarB
amplitude for scalar breathing polarisation mode
REAL8 psiTensor
phase polarisation angle for tensor polarisation mode
REAL8 T0
time of orbital perisastron as measured in TDB (MJD)
REAL8 * fb
orbital frequency coefficients for BTX model
REAL8 wave_om
fundamental frequency timing noise terms
REAL8 cosiota
cosine of the pulsars inclination angle
INT4 daopset
set if daop is set from the par file
REAL8 h0
gravitational wave amplitude
REAL8 f4
frequency fourth derivative (Hz/s^4)
REAL8 hVectorX
amplitude for vector "x" polarisation mode
REAL8 Pb
orbital period (seconds)
REAL8 psiScalar
phase polarisation angle for scalar polarisation mode
REAL8 phi0Tensor
initial phase for tensor modes (to be used instead of phi0 or phi22)
CHAR * model
TEMPO binary model e.g.
CHAR * units
The time system used e.g.
REAL8 cgw
The speed of gravitational waves as a fraction of the speed of light LAL_C_SI
REAL8 pmdec
proper motion in dec (rads/s)
REAL8 pmra
proper motion in RA (rads/s)
REAL8 Tasc
time of the ascending node (used rather than T0)
REAL8 Pbdot
rate of change of Pb (dimensionless)
CHAR * jname
pulsar J name
REAL8 DM1
first derivative of dispersion measure
REAL8 startTime
start of parfile applicable time
REAL8 posepoch
position epoch
REAL8 f10
frequency tenth derivative (Hz/s^10)
REAL8 dec
declination (rads)
REAL8 wdot
precesion of longitude of periastron w = w0 + wdot(tb-T0) (degs/year)
REAL8 finishTime
finish of parfile applicable time
REAL8 f0
spin frequency (Hz)
REAL8 f1
frequency first derivative (Hz/s)
REAL8 phi0Scalar
initial phase for the scalar modes
REAL8 gamma
gravitational redshift and time dilation parameter (s)
REAL8 kin
binary inclination angle
REAL8 costheta
cosine of angle between rotation axis and pinning axis
REAL8 psi
polarisation angle
REAL8 dist
pulsar distance (in kiloparsecs)
CHAR * bname
pulsar B name
The PulsarParam list node structure.
void * value
Parameter value.
struct tagPulsarParam * next
void * err
Parameter error/uncertainty.
UINT4Vector * fitFlag
Set to 1 if the parameter has been fit in the par file.
PulsarParamType type
Parameter type e.g.
The PulsarParameters structure to contain a set of pulsar parameters.
PulsarParam * head
A linked list of PulsarParam structures.
LALHashTbl * hash_table
Hash table of parameters.
INT4 nparams
The total number of parameters in the structure.