Functions to read TEMPO pulsar parameter files. More...
Prototypes | |
static void | strtoupper (CHAR *s) |
Convert string to uppercase. More... | |
static void * | new_elem (const char *name, PulsarParam *itemPtr) |
static void | del_elem (void *elem) |
static UINT8 | PulsarHash (const void *elem) |
static int | PulsarHashElemCmp (const void *elem1, const void *elem2) |
static PulsarParam * | PulsarGetParamItemSlow (const PulsarParameters *pars, const CHAR *name) |
Get a pointer to a parameter of a given name from a PulsarParameters structure. More... | |
static PulsarParam * | PulsarGetParamItem (const PulsarParameters *pars, const CHAR *name) |
Get a pointer to a parameter of a given name from a PulsarParameters structure. More... | |
void * | PulsarGetParam (const PulsarParameters *pars, const CHAR *name) |
Get the required parameter value from the PulsarParameters structure. More... | |
PulsarParamType | PulsarGetParamType (const PulsarParameters *pars, const char *name) |
Get the required parameter's type. More... | |
REAL8 | PulsarGetREAL8Param (const PulsarParameters *pars, const CHAR *name) |
Return a REAL8 parameter. More... | |
REAL8 | PulsarGetREAL8ParamOrZero (const PulsarParameters *pars, const CHAR *name) |
Return a REAL8 parameter if it exists, otherwise return zero. More... | |
UINT4 | PulsarGetUINT4Param (const PulsarParameters *pars, const CHAR *name) |
Return a UINT4 parameter. More... | |
UINT4 | PulsarGetUINT4ParamOrZero (const PulsarParameters *pars, const CHAR *name) |
Return a UINT4 parameter if it exists, otherwise return zero. More... | |
const CHAR * | PulsarGetStringParam (const PulsarParameters *pars, const CHAR *name) |
Return a string parameter. More... | |
const REAL8Vector * | PulsarGetREAL8VectorParam (const PulsarParameters *pars, const CHAR *name) |
Return a REAL8Vector parameter. More... | |
REAL8 | PulsarGetREAL8VectorParamIndividual (const PulsarParameters *pars, const CHAR *name) |
Return an individual REAL8 value from a REAL8Vector parameter. More... | |
void * | PulsarGetParamErr (const PulsarParameters *pars, const CHAR *name) |
Get the required parameter error value from the PulsarParameters structure. More... | |
const UINT4 * | PulsarGetParamFitFlag (const PulsarParameters *pars, const CHAR *name) |
Get the fit flag array for a given parameter from the PulsarParameters structure. More... | |
const UINT4Vector * | PulsarGetParamFitFlagAsVector (const PulsarParameters *pars, const CHAR *name) |
Get the fit flag array for a given parameter from the PulsarParameters structure. More... | |
REAL8 | PulsarGetREAL8ParamErr (const PulsarParameters *pars, const CHAR *name) |
Return a REAL8 parameter error value. More... | |
const REAL8Vector * | PulsarGetREAL8VectorParamErr (const PulsarParameters *pars, const CHAR *name) |
Return a REAL8Vector parameter error value. More... | |
REAL8 | PulsarGetREAL8VectorParamErrIndividual (const PulsarParameters *pars, const CHAR *name) |
Return an individual REAL8 value from a REAL8Vector parameter. More... | |
void | PulsarAddParam (PulsarParameters *pars, const CHAR *name, void *value, PulsarParamType type) |
Add a parameter and value to the PulsarParameters structure. More... | |
void | PulsarAddREAL8Param (PulsarParameters *pars, const CHAR *name, REAL8 value) |
Add a REAL8 parameter to the PulsarParameters structure. More... | |
void | PulsarAddUINT4Param (PulsarParameters *pars, const CHAR *name, UINT4 value) |
Add a UINT4 parameter to the PulsarParameters structure. More... | |
void | PulsarAddREAL8VectorParam (PulsarParameters *pars, const CHAR *name, const REAL8Vector *value) |
Add a REAL8Vector parameter to the PulsarParameters structure. More... | |
void | PulsarAddStringParam (PulsarParameters *pars, const CHAR *name, const CHAR *value) |
Add a string parameter to the PulsarParameters structure. More... | |
int | PulsarCheckParam (const PulsarParameters *pars, const CHAR *name) |
Check for the existence of the parameter name in the PulsarParameters structure. More... | |
void | PulsarClearParams (PulsarParameters *pars) |
Free all the parameters from a PulsarParameters structure. More... | |
void | PulsarRemoveParam (PulsarParameters *pars, const CHAR *name) |
Remove a given parameter from a PulsarParameters structure. More... | |
void | PulsarSetParam (PulsarParameters *pars, const CHAR *name, const void *value) |
Set the value of a parameter in the PulsarParameters structure. More... | |
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. More... | |
void | PulsarSetREAL8ParamErr (PulsarParameters *pars, const CHAR *name, REAL8 value, UINT4 fitFlag) |
Set the error value for a REAL8 parameter. More... | |
void | PulsarSetREAL8VectorParamErr (PulsarParameters *pars, const CHAR *name, const REAL8Vector *value, const UINT4 *fitFlag) |
Set the error values for a REAL8Vector parameter. More... | |
void | PulsarFreeParams (PulsarParameters *pars) |
Function to free memory from pulsar parameters. More... | |
void | PulsarCopyParams (PulsarParameters *origin, PulsarParameters *target) |
Function to copy a PulsarParameters structure. More... | |
DEFINE_CONV_FACTOR_FUNCTION (ToFloat, 1., CONVFLOAT) | |
A strcuture to contain all possible pulsar parameters that can be read in from a par file, and define the conversion function and type used for each. More... | |
static INT4 | ParseParLine (PulsarParameters *par, const CHAR *name, FILE *fp) |
Parse a single line from a pulsar parameter file. More... | |
PulsarParameters * | XLALReadTEMPOParFile (const CHAR *pulsarAndPath) |
Read in the parameters from a TEMPO(2) parameter file into a PulsarParameters structure. More... | |
void | XLALReadTEMPOParFileOrig (BinaryPulsarParams *output, CHAR *pulsarAndPath) |
function to read in a TEMPO parameter file More... | |
void | PrintPulsarParameters (BinaryPulsarParams params) |
function to print out all the pulsar parameters read in from a par file More... | |
LALStringVector * | XLALReadTEMPOCorFile (REAL8Array *cormat, CHAR *corfile) |
This function will read in a TEMPO-style parameter correlation matrix. More... | |
REAL8 | XLALTTMJDtoGPS (REAL8 MJD) |
This function converts a MJD format time corrected to Terrestrial Time (TT) into an equivalent GPS time. More... | |
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 TDB (see Table 1 of Seidelmann and Fukushima, Astronomy & Astrophysics, 265, 833-838 (1992). More... | |
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 TCB (see Table 1 of Seidelmann and Fukushima, Astronomy & Astrophysics, 265, 833-838, 1992). More... | |
Functions to read TEMPO pulsar parameter files.
Functions for reading pulsar parameters from TEMPO .par files.
Radio astronomers fit pulsar parameters using TEMPO(2) which will output the parameters in a .par
file. The values allowed in this file can be found in the TEMPO documentation. Two function are available to extract these parameters from the .par
files:
XLALReadTEMPOParFileOrig
- this will read parameters into a BinaryPulsarParams
structure and set any unused parameters to zero or NULL
. To use this you must know the correct parameter name within the structure. This is deprecated in favour of XLALReadTEMPOParFile
and is no longer maintained XLALReadTEMPOParFile
- this reads the parameters into a linked list PulsarParameters
structure, from which the parameters can be accessed using the appropriate access function. These use a hash table to quick look-up. The parameters are assigned names, which are used as the hash table keys, which are fully uppercase versions of the TEMPO parameter names. All parameters read in are converted into SI units.
Functions are is also included which converts a string containing the right ascension or declination in the format ddd/hh:mm:ss.s
or ddd/hhmmss.s
(as is given in the .par
file) into a REAL8
value in radians.
Definition in file ReadPulsarParFile.c.
Go to the source code of this file.
Data Structures | |
struct | hash_elem |
Macros | |
#define | DAYSTOSECS 86400.0 /* number of seconds in an SI day */ |
#define | DEFINE_CONV_FACTOR_FUNCTION(name, convfactor, type) |
#define | NUM_PARS 130 /* number of allowed parameters */ |
Enumerations | |
enum | { CONVFLOAT = 0 , CONVINT , CONVSTRING , CONVHMS , CONVDMS , CONVMJD , CONVBINUNITS } |
Variables | |
size_t | PulsarTypeSize [5] |
static const CHAR | a2A [256] |
Array for conversion from lowercase to uppercase. More... | |
ParConversion | |
ParConversion | pc [NUM_PARS] |
Initialise conversion structure with most allowed TEMPO2 parameter names and conversion functions (convert all read in parameters to SI units where necessary). More... | |
#define DAYSTOSECS 86400.0 /* number of seconds in an SI day */ |
Definition at line 83 of file ReadPulsarParFile.c.
Definition at line 897 of file ReadPulsarParFile.c.
#define NUM_PARS 130 /* number of allowed parameters */ |
Definition at line 969 of file ReadPulsarParFile.c.
anonymous enum |
Enumerator | |
---|---|
CONVFLOAT | |
CONVINT | |
CONVSTRING | |
CONVHMS | |
CONVDMS | |
CONVMJD | |
CONVBINUNITS |
Definition at line 886 of file ReadPulsarParFile.c.
|
static |
Convert string to uppercase.
Definition at line 104 of file ReadPulsarParFile.c.
|
static |
Definition at line 121 of file ReadPulsarParFile.c.
|
static |
Definition at line 127 of file ReadPulsarParFile.c.
Definition at line 135 of file ReadPulsarParFile.c.
Definition at line 146 of file ReadPulsarParFile.c.
|
static |
Get a pointer to a parameter of a given name from a PulsarParameters
structure.
Note this function can only be used internally.
Definition at line 159 of file ReadPulsarParFile.c.
|
static |
Get a pointer to a parameter of a given name from a PulsarParameters
structure.
This function will return a pointer to the parameter. It will initially try and use the parameter's hash name, otherwise it will use PulsarGetParamItemSlow
to loop through all parameters.
Note this function can only be used internally.
Definition at line 192 of file ReadPulsarParFile.c.
void* PulsarGetParam | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Get the required parameter value from the PulsarParameters
structure.
This function will return a void pointer to the parameter value. This should be cast into the required variable type once returned.
Definition at line 223 of file ReadPulsarParFile.c.
PulsarParamType PulsarGetParamType | ( | const PulsarParameters * | pars, |
const char * | name | ||
) |
Get the required parameter's type.
Definition at line 235 of file ReadPulsarParFile.c.
REAL8 PulsarGetREAL8Param | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a REAL8
parameter.
This function will call PulsarGetParam
for a REAL8
parameter and properly cast it for returning.
Definition at line 241 of file ReadPulsarParFile.c.
REAL8 PulsarGetREAL8ParamOrZero | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a REAL8
parameter if it exists, otherwise return zero.
Definition at line 252 of file ReadPulsarParFile.c.
UINT4 PulsarGetUINT4Param | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a UINT4
parameter.
This function will call PulsarGetParam
for a UINT4
parameter and properly cast it for returning.
Definition at line 259 of file ReadPulsarParFile.c.
UINT4 PulsarGetUINT4ParamOrZero | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a UINT4
parameter if it exists, otherwise return zero.
Definition at line 270 of file ReadPulsarParFile.c.
Return a string parameter.
This function will call PulsarGetParam
for a string parameter and properly cast it for returning. The return value should be copied e.g. with CHAR *str = XLALStringDuplicate( PulsarGetStringParam(pars, "NAME") ); It also needs to be freed afterwards.
Definition at line 277 of file ReadPulsarParFile.c.
const REAL8Vector* PulsarGetREAL8VectorParam | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a REAL8Vector
parameter.
This function will call PulsarGetParam
for a REAL8Vector
parameter and properly cast it for returning.
Definition at line 289 of file ReadPulsarParFile.c.
REAL8 PulsarGetREAL8VectorParamIndividual | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return an individual REAL8
value from a REAL8Vector
parameter.
This function will call PulsarGetParam
for a REAL8Vector
parameter and then return a specific value from within the vector. The input name
must be of the form e.g. FB1
, which will get the REAL8Vector
for the FB
parameter and return the value from the 1
index.
Definition at line 300 of file ReadPulsarParFile.c.
void* PulsarGetParamErr | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Get the required parameter error value from the PulsarParameters
structure.
This function will return a void pointer to the parameter error value. The parameter error will be either a REAL8
or a REAL8Vector
and should be cast accordingly once returned.
Definition at line 341 of file ReadPulsarParFile.c.
Get the fit flag array for a given parameter from the PulsarParameters
structure.
This function will return a UINT4
array to the parameter fit flag.
Definition at line 353 of file ReadPulsarParFile.c.
const UINT4Vector* PulsarGetParamFitFlagAsVector | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Get the fit flag array for a given parameter from the PulsarParameters
structure.
This function will return a UINT4Vector
array to the parameter fit flag.
Definition at line 373 of file ReadPulsarParFile.c.
REAL8 PulsarGetREAL8ParamErr | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a REAL8
parameter error value.
This function will call PulsarGetParamErr
for a REAL8
parameter and properly cast it for returning.
Definition at line 389 of file ReadPulsarParFile.c.
const REAL8Vector* PulsarGetREAL8VectorParamErr | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return a REAL8Vector
parameter error value.
This function will call PulsarGetParamErr
for a REAL8Vector
parameter and properly cast it for returning.
Definition at line 406 of file ReadPulsarParFile.c.
REAL8 PulsarGetREAL8VectorParamErrIndividual | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Return an individual REAL8
value from a REAL8Vector
parameter.
This function will call PulsarGetParam
for a REAL8Vector
parameter's errors and then return a specific value from within the vector. The input name
must be of the form e.g. FB1
, which will get the REAL8Vector
for the FB
parameter and return the value from the 1
index.
Definition at line 423 of file ReadPulsarParFile.c.
void PulsarAddParam | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
void * | value, | ||
PulsarParamType | type | ||
) |
Add a parameter and value to the PulsarParameters
structure.
This function adds a new parameter, and associated value to the PulsarParameters
structure. If the parameter already exists then the old value is replaced with the new value.
Definition at line 464 of file ReadPulsarParFile.c.
void PulsarAddREAL8Param | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
REAL8 | value | ||
) |
Add a REAL8
parameter to the PulsarParameters
structure.
Definition at line 542 of file ReadPulsarParFile.c.
void PulsarAddUINT4Param | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
UINT4 | value | ||
) |
Add a UINT4
parameter to the PulsarParameters
structure.
Definition at line 549 of file ReadPulsarParFile.c.
void PulsarAddREAL8VectorParam | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
const REAL8Vector * | value | ||
) |
Add a REAL8Vector
parameter to the PulsarParameters
structure.
Definition at line 556 of file ReadPulsarParFile.c.
void PulsarAddStringParam | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
const CHAR * | value | ||
) |
Add a string parameter to the PulsarParameters
structure.
Definition at line 566 of file ReadPulsarParFile.c.
int PulsarCheckParam | ( | const PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Check for the existence of the parameter name
in the PulsarParameters
structure.
Definition at line 577 of file ReadPulsarParFile.c.
void PulsarClearParams | ( | PulsarParameters * | pars | ) |
Free all the parameters from a PulsarParameters
structure.
Definition at line 592 of file ReadPulsarParFile.c.
void PulsarRemoveParam | ( | PulsarParameters * | pars, |
const CHAR * | name | ||
) |
Remove a given parameter from a PulsarParameters
structure.
Definition at line 642 of file ReadPulsarParFile.c.
void PulsarSetParam | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
const void * | value | ||
) |
Set the value of a parameter in the PulsarParameters
structure.
Set the value of the parameter given by name
in the PulsarParameters
structure. The parameter must already exist in the structure, otherwise it should be added using PulsarAddParam()
.
Definition at line 707 of file ReadPulsarParFile.c.
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.
Set the value of the error on the parameter given by name
in the PulsarParameters
structure. The parameter must already exist in the structure, otherwise it should be added using PulsarAddParam()
. If the .par file contains a 1 before the error value (because that parameter has been included in the TEMPO(2) fitting procedure) then that must be input as the fit flag (this can be a vector for e.g. FB values with multiple parameters, in which case nfits
will be the number of values in that vector).
Definition at line 725 of file ReadPulsarParFile.c.
void PulsarSetREAL8ParamErr | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
REAL8 | value, | ||
UINT4 | fitFlag | ||
) |
Set the error value for a REAL8
parameter.
Definition at line 767 of file ReadPulsarParFile.c.
void PulsarSetREAL8VectorParamErr | ( | PulsarParameters * | pars, |
const CHAR * | name, | ||
const REAL8Vector * | value, | ||
const UINT4 * | fitFlag | ||
) |
Set the error values for a REAL8Vector
parameter.
Definition at line 774 of file ReadPulsarParFile.c.
void PulsarFreeParams | ( | PulsarParameters * | pars | ) |
Function to free memory from pulsar parameters.
Definition at line 785 of file ReadPulsarParFile.c.
void PulsarCopyParams | ( | PulsarParameters * | origin, |
PulsarParameters * | target | ||
) |
Function to copy a PulsarParameters
structure.
Definition at line 800 of file ReadPulsarParFile.c.
DEFINE_CONV_FACTOR_FUNCTION | ( | ToFloat | , |
1. | , | ||
CONVFLOAT | |||
) |
A strcuture to contain all possible pulsar parameters that can be read in from a par file, and define the conversion function and type used for each.
Parameter name
Conversion function from string to required value
Conversion function for error from string to value
Parameter type
Definition at line 938 of file ReadPulsarParFile.c.
|
static |
Parse a single line from a pulsar parameter file.
This will parse a line from the TEMPO-style pulsar parameter file containing the parameter given by name
. The parameter will be added to the par
structure.
This function can only be used internally.
Definition at line 1146 of file ReadPulsarParFile.c.
PulsarParameters* XLALReadTEMPOParFile | ( | const CHAR * | pulsarAndPath | ) |
Read in the parameters from a TEMPO(2) parameter file into a PulsarParameters
structure.
This function will read in a TEMPO(2) parameter file into a PulsarParameters
structure. The structure of this function is similar to those in the TEMPO2 code readParfile.C
and this supersedes the XLALReadTEMPOParFileOrig
function.
pulsarAndPath | [in] The path to the pulsar parameter file |
Definition at line 1480 of file ReadPulsarParFile.c.
void XLALReadTEMPOParFileOrig | ( | BinaryPulsarParams * | output, |
CHAR * | pulsarAndPath | ||
) |
function to read in a TEMPO parameter file
Definition at line 1553 of file ReadPulsarParFile.c.
void PrintPulsarParameters | ( | BinaryPulsarParams | params | ) |
function to print out all the pulsar parameters read in from a par file
Definition at line 2849 of file ReadPulsarParFile.c.
LALStringVector* XLALReadTEMPOCorFile | ( | REAL8Array * | cormat, |
CHAR * | corfile | ||
) |
This function will read in a TEMPO-style parameter correlation matrix.
This function will read in a TEMPO-style parameter correlation matrix file, which contains the correlations between parameters as fit by TEMPO. An example the format would be:
RA DEC F0 RA 1.000 DEC 0.954 1.000 F0 -0.007 0.124 1.000
In the output all parameter names will sometimes be converted to a more convenient naming convention. If non-diagonal parameter correlation values are +/-1 then they will be converted to be +/-0.99999 to avoid some problems of the matrix becoming singular. The output matrix will have both the upper and lower triangle completed.
cormat | [out] A REAL8 array into which the correlation matrix will be output |
corfile | [in] A string containing the path and filename of the TEMPO-style correlation matrix file |
Definition at line 2915 of file ReadPulsarParFile.c.
This function converts a MJD format time corrected to Terrestrial Time (TT) into an equivalent GPS time.
Functions for converting times given in Terrestrial time TT, TDB, or TCB in MJD to times in GPS - this is important for epochs given in .par
files which are in TDB(MJD).
Definition at line 3045 of file ReadPulsarParFile.c.
If you have an MJD arrival time on the Earth then this will convert it to the equivalent GPS time in TDB (see Table 1 of Seidelmann and Fukushima, Astronomy & Astrophysics, 265, 833-838 (1992).
Note that XLALBarycenter performs these TDBtoTT corrections (i.e. the Einstein delay) when correcting a GPS time on the Earth to TDB. Also, for TEMPO produced pulsar epochs given in MJD these are already in the TDB system and an equivalent GPS time in the TDB can be calculated just using XLALTTMJDtoGPS
.
Definition at line 3066 of file ReadPulsarParFile.c.
If you have an MJD arrival time on the Earth then this will convert it to the equivalent GPS time in TCB (see Table 1 of Seidelmann and Fukushima, Astronomy & Astrophysics, 265, 833-838, 1992).
Note that for default TEMPO2 produced pulsar epochs given in MJD these are already in the TCB system and an equivalent GPS time in the TCB can be calculated just using XLALTTMJDtoGPS
.
Definition at line 3110 of file ReadPulsarParFile.c.
size_t PulsarTypeSize[5] |
Definition at line 85 of file ReadPulsarParFile.c.
Array for conversion from lowercase to uppercase.
Definition at line 95 of file ReadPulsarParFile.c.
ParConversion |
Definition at line 966 of file ReadPulsarParFile.c.
Initialise conversion structure with most allowed TEMPO2 parameter names and conversion functions (convert all read in parameters to SI units where necessary).
See http://arxiv.org/abs/astro-ph/0603381 and http://arxiv.org/abs/astro-ph/0603381 for parameter definitions.
If requiring a new parameter to be read in in should be added to this structure and NUM_PARS
should be incremented.
Definition at line 978 of file ReadPulsarParFile.c.