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