LALPulsar 7.1.2.1-bf6a62b
ReadPulsarParFile.c
Go to the documentation of this file.
1/*
2* Copyright (C) 2013 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/// \addtogroup ReadPulsarParFile_h
21/// @{
22
23#include <config.h>
24#include <string.h>
25#include <math.h>
26#include <locale.h>
27#include <ctype.h>
28#ifdef HAVE_UNISTD_H
29#include <unistd.h>
30#endif
31#include <lal/LALConstants.h>
32#include <lal/LALStdlib.h>
33#include <lal/LALString.h>
34#include <lal/ComputeFstat.h>
35#include <lal/TranslateAngles.h>
36#include <lal/TranslateMJD.h>
37#include <lal/LALHashFunc.h>
38#include <lal/ReadPulsarParFile.h>
39
40#ifdef __GNUC__
41#define UNUSED __attribute__ ((unused))
42#else
43#define UNUSED
44#endif
45
46
47#define DAYSTOSECS 86400.0 /* number of seconds in an SI day */
48
49size_t PulsarTypeSize[5] = {
50 sizeof( UINT4 ),
51 sizeof( REAL8 ),
52 sizeof( REAL8Vector * ),
53 sizeof( CHAR ) *PULSAR_PARNAME_MAX,
54 sizeof( void * )
55};
56
57
58/** Array for conversion from lowercase to uppercase */
59static const CHAR a2A[256] = {
60 ['a'] = 'A', ['b'] = 'B', ['c'] = 'C', ['d'] = 'D', ['e'] = 'E', ['f'] = 'F', ['g'] = 'G', ['h'] = 'H',
61 ['i'] = 'I', ['j'] = 'J', ['k'] = 'K', ['l'] = 'L', ['m'] = 'M', ['n'] = 'N', ['o'] = 'O', ['p'] = 'P',
62 ['q'] = 'Q', ['r'] = 'R', ['s'] = 'S', ['t'] = 'T', ['u'] = 'U', ['v'] = 'V', ['w'] = 'W', ['x'] = 'X',
63 ['y'] = 'Y', ['z'] = 'Z'
64};
65
66
67/** \brief Convert string to uppercase */
68static void strtoupper( CHAR *s )
69{
70 /* convert everything to uppercase */
71 CHAR c;
72 for ( ; *s; ++s ) {
73 if ( ( c = a2A[( int )( *s )] ) ) {
74 *s = c;
75 }
76 }
77}
78
79/* hash functions copied from LALInference.c */
80typedef struct taghash_elem {
81 const char *name;
82 PulsarParam *itemPtr;
83} hash_elem;
84
85static void *new_elem( const char *name, PulsarParam *itemPtr )
86{
87 hash_elem e = { .name = name, .itemPtr = itemPtr };
88 return memcpy( XLALMalloc( sizeof( e ) ), &e, sizeof( e ) );
89};
90
91static void del_elem( void *elem )
92{
93 XLALFree( elem );
94}
95
96
97/* Compute a hash value based on element */
98static UINT8 PulsarHash( const void *elem );
99static UINT8 PulsarHash( const void *elem )
100{
101 if ( !elem ) {
103 }
104 size_t len = strnlen( ( ( const hash_elem * )elem )->name, PULSAR_PARNAME_MAX );
105 return ( XLALCityHash64( ( ( const hash_elem * )elem )->name, len ) );
106}
107
108
109static int PulsarHashElemCmp( const void *elem1, const void *elem2 );
110static int PulsarHashElemCmp( const void *elem1, const void *elem2 )
111{
112 if ( !elem1 || !elem2 ) {
114 }
115 return ( strncmp( ( ( const hash_elem * )elem1 )->name, ( ( const hash_elem * )elem2 )->name, PULSAR_PARNAME_MAX ) );
116}
117
118
119/** \brief Get a pointer to a parameter of a given name from a \c PulsarParameters structure
120 *
121 * Note this function can only be used internally.
122 */
124{
125 /* (this function is only to be used internally) */
126 /* Returns pointer to item for given item name. */
127 if ( pars == NULL ) {
128 return NULL;
129 }
130 if ( pars->nparams == 0 ) {
131 return NULL;
132 }
133
134 PulsarParam *this = pars->head;
135
136 /* loop through all values checking the name */
137 while ( this != NULL ) {
138 if ( !strcmp( this->name, name ) ) {
139 break;
140 } else {
141 this = this->next;
142 }
143 }
144
145 return ( this );
146}
147
148
149/** \brief Get a pointer to a parameter of a given name from a \c PulsarParameters structure
150 *
151 * This function will return a pointer to the parameter. It will initially try and use the parameter's hash name,
152 * otherwise it will use \c PulsarGetParamItemSlow to loop through all parameters.
153 *
154 * Note this function can only be used internally.
155 */
157{
158 /* (this function is only to be used internally) */
159 /* Returns pointer to item for given item name. */
160 hash_elem tmp;
161 const hash_elem *match = NULL;
162
163 CHAR upperName[PULSAR_PARNAME_MAX];
165 strtoupper( upperName );
166
167 if ( pars == NULL ) {
168 return NULL;
169 }
170 if ( pars->nparams == 0 ) {
171 return NULL;
172 }
173 if ( !pars->hash_table ) {
174 return PulsarGetParamItemSlow( pars, upperName );
175 }
176
177 tmp.name = upperName;
178 XLALHashTblFind( pars->hash_table, &tmp, ( const void ** )&match );
179 if ( !match ) {
180 return NULL;
181 }
182
183 return match->itemPtr;
184}
185
186
187void *PulsarGetParam( const PulsarParameters *pars, const CHAR *name )
188{
189 /* Return the value of variable name from the pars structure */
190 PulsarParam *item = PulsarGetParamItem( pars, name );
191 if ( !item ) {
192 XLAL_ERROR_NULL( XLAL_EFAILED, "Entry \"%s\" not found.", name );
193 }
194
195 return ( item->value );
196}
197
198
200{
201 return PulsarGetParamItem( pars, name )->type;
202}
203
204
206{
207 /* check type is a REAL8 */
208 if ( PulsarGetParamType( pars, name ) == PULSARTYPE_REAL8_t ) {
209 return *( REAL8 * )PulsarGetParam( pars, name );
210 } else {
211 XLAL_ERROR_REAL8( XLAL_EINVAL, "Used wrong type for required parameter" );
212 }
213}
214
215
217{
218 /* if parameter is there return it otherwise return zero */
219 return PulsarCheckParam( pars, name ) ? PulsarGetREAL8Param( pars, name ) : 0.;
220}
221
222
224{
225 /* check type is a UINT4 */
226 if ( PulsarGetParamType( pars, name ) == PULSARTYPE_UINT4_t ) {
227 return *( UINT4 * )PulsarGetParam( pars, name );
228 } else {
229 XLAL_ERROR( XLAL_EINVAL, "Used wrong type for required parameter" );
230 }
231}
232
233
235{
236 /* if parameter is there return it otherwise return zero */
237 return PulsarCheckParam( pars, name ) ? PulsarGetUINT4Param( pars, name ) : 0;
238}
239
240
242{
243 /* check type is a string */
244 if ( PulsarGetParamType( pars, name ) == PULSARTYPE_string_t ) {
245 CHAR *rvalue = *( CHAR ** )PulsarGetParam( pars, name );
246 return rvalue;
247 } else {
248 XLAL_ERROR_NULL( XLAL_EINVAL, "Used wrong type for required parameter" );
249 }
250}
251
252
254{
255 /* check type is a REAL8Vector */
257 return *( REAL8Vector ** )PulsarGetParam( pars, name );
258 } else {
259 XLAL_ERROR_NULL( XLAL_EINVAL, "Used wrong type for required parameter" );
260 }
261}
262
263
265{
266 /* split the input name into the parameter name and an index e.g. FB0 into FB and 0, or GLEP_1 */
267 CHAR *namecopy = NULL, *token;
268 namecopy = XLALStringDuplicate( name );
269 REAL8 val = 0.;
270
271 if ( strchr( name, '_' ) ) { /* check for underscore delimiting numbers */
272 token = strtok( namecopy, "_" );
273 } else { /* name delimited by just numbers */
274 token = strtok( namecopy, "0123456789" );
275 }
276
277 const REAL8Vector *vpars = PulsarGetREAL8VectorParam( pars, token );
278
279 /* get the index value of the parameter name */
280 INT4 idx = -1;
281 INT4 loc = strlen( token );
282 if ( name[loc] == '_' ) {
283 loc++;
284 }
285 if ( sscanf( name + loc, "%d", &idx ) != 1 ) {
286 XLAL_ERROR_REAL8( XLAL_EINVAL, "Input parameter (%s) problem", name );
287 }
288
289 /* for glitch parameters GL, or WAVE parameters the index starts at 1, so shift it */
290 if ( !strncmp( name, "GL", 2 ) || !strncmp( name, "WAVE", 4 ) ) {
291 idx--;
292 }
293 if ( idx < 0 || ( UINT4 )idx > vpars->length - 1 ) {
294 XLAL_ERROR_REAL8( XLAL_EINVAL, "Input parameter index %d is wrong", idx );
295 }
296
297 val = vpars->data[idx];
298
299 XLALFree( namecopy );
300
301 return val;
302}
303
304
305void *PulsarGetParamErr( const PulsarParameters *pars, const CHAR *name )
306{
307 /* Return the error value of variable name from the pars structure */
308 PulsarParam *item = PulsarGetParamItem( pars, name );
309 if ( !item ) {
310 XLAL_ERROR_NULL( XLAL_EFAILED, "Entry \"%s\" not found.", name );
311 }
312
313 return ( item->err );
314}
315
316
318{
319 /* return the fit flag value */
320 PulsarParam *item = PulsarGetParamItem( pars, name );
321 if ( !item ) {
322 XLAL_ERROR_NULL( XLAL_EFAILED, "Entry \"%s\" not found.", name );
323 }
324
325 if ( item->fitFlag == NULL ) {
326 return NULL;
327 } else {
328 if ( item->fitFlag->data == NULL ) {
329 return NULL;
330 } else {
331 return item->fitFlag->data;
332 }
333 }
334}
335
336
338{
339 /* return the fit flag value(s) */
340 PulsarParam *item = PulsarGetParamItem( pars, name );
341 if ( !item ) {
342 XLAL_ERROR_NULL( XLAL_EFAILED, "Entry \"%s\" not found.", name );
343 }
344
345 if ( item->fitFlag == NULL ) {
346 return NULL;
347 } else {
348 return item->fitFlag;
349 }
350}
351
352
354{
355 /* Return the error value of variable name from the pars structure */
356 /* check type is a REAL8 */
357 if ( PulsarGetParamType( pars, name ) == PULSARTYPE_REAL8_t ) {
358 void *val = PulsarGetParamErr( pars, name );
359 if ( val == NULL ) {
360 return 0.; /* if no error is present return 0 */
361 } else {
362 return *( REAL8 * )val;
363 }
364 } else {
365 XLAL_ERROR_REAL8( XLAL_EINVAL, "Used wrong type for required parameter" );
366 }
367}
368
369
371{
372 /* Return the error value of variable name from the pars structure */
373 /* check type is a REAL8 */
375 void *val = PulsarGetParamErr( pars, name );
376 if ( val == NULL ) {
377 return NULL; /* if no error is present return NULL */
378 } else {
379 return *( REAL8Vector ** )val;
380 }
381 } else {
382 XLAL_ERROR_NULL( XLAL_EINVAL, "Used wrong type for required parameter" );
383 }
384}
385
386
388{
389 /* split the input name into the parameter name and an index e.g. FB0 into FB and 0*/
390 CHAR *namecopy = NULL, *token;
391 namecopy = XLALStringDuplicate( name );
392 REAL8 val = 0.;
393
394 if ( strchr( name, '_' ) ) { /* check for underscore delimiting numbers */
395 token = strtok( namecopy, "_" );
396 } else { /* name delimited by just numbers */
397 token = strtok( namecopy, "0123456789" );
398 }
399
400 const REAL8Vector *vpars = PulsarGetREAL8VectorParamErr( pars, token );
401
402 /* get the index value of the parameter name */
403 INT4 idx = -1;
404 INT4 loc = strlen( token );
405 if ( name[loc] == '_' ) {
406 loc++;
407 }
408 if ( sscanf( name + loc, "%d", &idx ) != 1 ) {
409 XLAL_ERROR_REAL8( XLAL_EINVAL, "Input parameter (%s) problem", name );
410 }
411
412 /* for glitch parameters GL, or WAVE parameters the index starts at 1, so shift it */
413 if ( !strncmp( name, "GL", 2 ) || !strncmp( name, "WAVE", 4 ) ) {
414 idx--;
415 }
416 if ( idx < 0 || ( UINT4 )idx > vpars->length - 1 ) {
417 XLAL_ERROR_REAL8( XLAL_EINVAL, "Input parameter index %d is wrong", idx );
418 }
419
420 val = vpars->data[idx];
421
422 XLALFree( namecopy );
423
424 return val;
425}
426
427
428void PulsarAddParam( PulsarParameters *pars, const CHAR *name, void *value, PulsarParamType type )
429{
430 /* Add the variable name with type type and value value to pars */
431 /* If variable already exists, it will over-write the current value if type compatible*/
432 PulsarParam *old = NULL;
433
434 /* create the hash table if it does not exist */
435 if ( !pars->hash_table ) {
437 }
438
439 /* Check input value is accessible */
440 if ( !value ) {
441 XLAL_ERROR_VOID( XLAL_EFAULT, "Unable to access value through null pointer; trying to add \"%s\".", name );
442 }
443
444 /* Check the name doesn't already exist */
445 if ( PulsarCheckParam( pars, name ) ) {
446 old = PulsarGetParamItem( pars, name );
447
448 /* If the type is incompatible, it should be removed */
449 if ( old->type != type || old->type == PULSARTYPE_REAL8Vector_t ) {
450 PulsarRemoveParam( pars, name );
451 } else {
452 PulsarSetParam( pars, name, value );
453 return;
454 }
455 }
456
457 /* If we get this far it is safe to create a new node for this variable */
458 PulsarParam *new = XLALMalloc( sizeof( PulsarParam ) );
459 memset( new, 0, sizeof( PulsarParam ) );
460 if ( new ) {
461 new->value = ( void * )XLALMalloc( PulsarTypeSize[type] );
462 new->fitFlag = NULL;
463
464 /* for REAL8 and REAL8Vector parameters initialise errors and fit flags to zero */
466 REAL8Vector *tmpvector = *( REAL8Vector ** )value;
467 UINT4 dlength = tmpvector->length;
468
469 new->fitFlag = XLALCreateUINT4Vector( dlength );
470 memset( new->fitFlag->data, 0, sizeof( UINT4 )*dlength ); // initialise to zero
471
473 memset( err->data, 0, sizeof( REAL8 )*dlength ); // initialise to zero
474
475 new->err = ( void * )XLALMalloc( PulsarTypeSize[type] );
476 memcpy( new->err, ( void * )&err, PulsarTypeSize[type] );
477 } else if ( type == PULSARTYPE_REAL8_t ) {
478 new->fitFlag = XLALCreateUINT4Vector( 1 );
479 memset( new->fitFlag->data, 0, sizeof( UINT4 ) ); // initialise to zero
480
481 new->err = ( void * )XLALMalloc( sizeof( UINT4Vector * ) );
482 REAL8 err = 0.;
483 memcpy( new->err, ( void * )&err, PulsarTypeSize[type] );
484 }
485 }
486
487 if ( new == NULL || new->value == NULL ) {
488 XLAL_ERROR_VOID( XLAL_ENOMEM, "Unable to allocate memory for list item." );
489 }
490
491 CHAR upperName[PULSAR_PARNAME_MAX];
493 strtoupper( upperName );
494
495 XLALStringCopy( new->name, upperName, PULSAR_PARNAME_MAX );
496 new->type = type;
497 memcpy( new->value, value, PulsarTypeSize[type] );
498 new->next = pars->head;
499 pars->head = new;
500 hash_elem *elem = new_elem( new->name, new );
501 XLALHashTblAdd( pars->hash_table, ( void * )elem );
502 pars->nparams++;
503}
504
505
507/* Typed version of PulsarAddParam for REAL8 values.*/
508{
509 PulsarAddParam( pars, name, ( void * )&value, PULSARTYPE_REAL8_t );
510}
511
512
514/* Typed version of PulsarAddParam for UINT4 values.*/
515{
516 PulsarAddParam( pars, name, ( void * )&value, PULSARTYPE_UINT4_t );
517}
518
519
521/* Typed version of PulsarAddParam for REAL8Vector values.*/
522{
523 REAL8Vector *fval = NULL;
524 fval = XLALCreateREAL8Vector( value->length );
525 memcpy( fval->data, value->data, sizeof( REAL8 )*value->length );
526 PulsarAddParam( pars, name, ( void * )&fval, PULSARTYPE_REAL8Vector_t );
527}
528
529
530void PulsarAddStringParam( PulsarParameters *pars, const CHAR *name, const CHAR *value )
531/* Typed version of PulsarAddParam for string values.*/
532{
533 CHAR *sval = NULL;
536 PulsarAddParam( pars, name, ( void * )&sval, PULSARTYPE_string_t );
537}
538
539
540/* Check for existence of name */
541int PulsarCheckParam( const PulsarParameters *pars, const CHAR *name )
542{
543 /* convert name to uppercase */
544 CHAR upperName[PULSAR_PARNAME_MAX];
546 strtoupper( upperName );
547
548 if ( PulsarGetParamItem( pars, upperName ) ) {
549 return 1;
550 } else {
551 return 0;
552 }
553}
554
555
557{
558 /* Free all variables inside the linked list, leaving only the head struct */
559 PulsarParam *this, *next = NULL;
560
561 if ( !pars ) {
562 return;
563 }
564
565 this = pars->head;
566
567 if ( this ) {
568 next = this->next;
569 }
570
571 while ( this ) {
572 if ( this->type == PULSARTYPE_REAL8Vector_t ) {
573 if ( this->value ) {
574 XLALDestroyREAL8Vector( *( REAL8Vector ** )this->value );
575 }
576 if ( this->err ) {
577 XLALDestroyREAL8Vector( *( REAL8Vector ** )this->err );
578 }
579 } else if ( this->type == PULSARTYPE_string_t ) {
580 if ( this->value ) {
581 XLALFree( *( CHAR ** )this->value );
582 }
583 }
584 if ( this->fitFlag ) {
585 XLALDestroyUINT4Vector( this->fitFlag );
586 }
587 XLALFree( this->value );
588 XLALFree( this->err );
589 XLALFree( this );
590 this = next;
591 if ( this ) {
592 next = this->next;
593 }
594 }
595 pars->head = NULL;
596 pars->nparams = 0;
597
598 if ( pars->hash_table ) {
600 }
601 pars->hash_table = NULL;
602}
603
604
605/* remove a given parameter */
607{
608 PulsarParam *this;
609
610 if ( !pars ) {
612 }
613
614 CHAR upperName[PULSAR_PARNAME_MAX];
616 strtoupper( upperName );
617
618 this = pars->head;
619 PulsarParam *parent = NULL;
620
621 while ( this ) {
622 if ( !strcmp( this->name, upperName ) ) {
623 break;
624 } else {
625 parent = this;
626 this = this->next;
627 }
628 }
629 if ( !this ) {
630 XLAL_PRINT_WARNING( "Entry \"%s\" not found.", upperName );
631 return;
632 }
633
634 if ( !parent ) {
635 pars->head = this->next;
636 } else {
637 parent->next = this->next;
638 }
639 /* Remove from hash table */
641 elem.name = this->name;
642 XLALHashTblRemove( pars->hash_table, ( void * )&elem );
643 if ( this->type == PULSARTYPE_REAL8Vector_t ) {
644 if ( this->value ) {
645 XLALDestroyREAL8Vector( *( REAL8Vector ** )this->value );
646 }
647 if ( this->err ) {
648 XLALDestroyREAL8Vector( *( REAL8Vector ** )this->err );
649 }
650 } else if ( this->type == PULSARTYPE_string_t ) {
651 if ( this->value ) {
652 XLALFree( *( CHAR ** )this->value );
653 }
654 }
655 if ( this->fitFlag ) {
656 XLALDestroyUINT4Vector( this->fitFlag );
657 }
658 XLALFree( this->value );
659 XLALFree( this->err );
660 this->value = NULL;
661 this->err = NULL;
662 this->fitFlag = NULL;
663 XLALFree( this );
664
665 this = NULL;
666 pars->nparams--;
667}
668
669
670/* Set the value of parameter name in the pars structure to value */
671void PulsarSetParam( PulsarParameters *pars, const CHAR *name, const void *value )
672{
673 PulsarParam *item;
674
675 /* convert name to uppercase */
676 CHAR upperName[PULSAR_PARNAME_MAX];
678 strtoupper( upperName );
679
680 item = PulsarGetParamItem( pars, upperName );
681 if ( !item ) {
682 XLAL_ERROR_VOID( XLAL_EINVAL, "Entry \"%s\" not found.", upperName );
683 }
684 memcpy( item->value, value, PulsarTypeSize[item->type] );
685}
686
687
688/* Set the value of parameter name error in the pars structure to value */
689void PulsarSetParamErr( PulsarParameters *pars, const CHAR *name, void *value, const UINT4 *fitFlag, UINT4 len )
690{
691 PulsarParam *item;
692
693 /* convert name to uppercase */
694 CHAR upperName[PULSAR_PARNAME_MAX];
696 strtoupper( upperName );
697
698 item = PulsarGetParamItem( pars, upperName );
699 if ( !item ) {
700 XLAL_ERROR_VOID( XLAL_EINVAL, "Entry \"%s\" not found.", name );
701 }
702
703 if ( item->err != NULL ) {
704 // Remove previous error values if present
705 if ( item->type == PULSARTYPE_REAL8Vector_t ) {
706 XLALDestroyREAL8Vector( *( REAL8Vector ** )item->err );
707 }
708 XLALFree( item->err );
709 }
710
711 if ( ( item->err = ( void * )XLALMalloc( PulsarTypeSize[PulsarGetParamType( pars, name )] ) ) == NULL ) {
712 XLAL_ERROR_VOID( XLAL_ENOMEM, "Unable to allocate memory for list item." );
713 }
714
715 // Remove previous fit flag if present
716 if ( item->fitFlag != NULL ) {
718
719 // set the fit flag values
720 item->fitFlag = NULL;
721 item->fitFlag = XLALCreateUINT4Vector( len );
722 memcpy( item->fitFlag->data, fitFlag, len * sizeof( UINT4 ) );
723 }
724
725 // set error values
726 memcpy( item->err, value, PulsarTypeSize[item->type] );
727}
728
729
730/* Set the value of parameter name error in the pars structure to value */
731void PulsarSetREAL8ParamErr( PulsarParameters *pars, const CHAR *name, REAL8 value, UINT4 fitFlag )
732{
733 PulsarSetParamErr( pars, name, ( void * )&value, &fitFlag, 1 );
734}
735
736
737/* Set the value of parameter name error in the pars structure to value */
738void PulsarSetREAL8VectorParamErr( PulsarParameters *pars, const CHAR *name, const REAL8Vector *value, const UINT4 *fitFlag )
739{
740 REAL8Vector *evals = NULL;
741 evals = XLALCreateREAL8Vector( value->length );
742 memcpy( evals->data, value->data, sizeof( REAL8 )*value->length );
743
744 PulsarSetParamErr( pars, name, ( void * )&evals, fitFlag, value->length );
745}
746
747
748/* free memory */
750{
751 if ( !pars ) {
752 return;
753 }
754
755 PulsarClearParams( pars );
756 if ( pars->hash_table ) {
758 }
759 XLALFree( pars );
760}
761
762
763/* create a copy of the parameters */
765{
766 /* Check that the source and origin differ */
767 if ( origin == target ) {
768 return;
769 }
770
771 if ( !origin ) {
772 XLAL_ERROR_VOID( XLAL_EFAULT, "Unable to access origin pointer." );
773 }
774
775 /* Make sure the structure is initialised */
776 if ( !target ) {
777 XLAL_ERROR_VOID( XLAL_EFAULT, "Unable to copy to uninitialised PulsarParameters structure." );
778 }
779
780 PulsarParam *this, *next = NULL;
781
782 this = origin->head;
783
784 if ( this ) {
785 next = this->next;
786 }
787
788 /* first clear target variables */
789 PulsarClearParams( target );
790
791 while ( this ) {
792 if ( this->type == PULSARTYPE_REAL8Vector_t ) {
793 const REAL8Vector *old = *( REAL8Vector ** )this->value;
795 if ( new ) {
796 memcpy( new->data, old->data, new->length * sizeof( REAL8 ) );
797 } else {
798 XLAL_ERROR_VOID( XLAL_ENOMEM, "Unable to copy vector!\n" );
799 }
800
801 PulsarAddREAL8VectorParam( target, this->name, ( const REAL8Vector * )new );
802
803 if ( this->err ) {
804 const REAL8Vector *olderr = *( REAL8Vector ** )this->err;
805 REAL8Vector *newerr = XLALCreateREAL8Vector( olderr->length );
806 if ( newerr ) {
807 memcpy( newerr->data, olderr->data, newerr->length * sizeof( REAL8 ) );
808 } else {
809 XLAL_ERROR_VOID( XLAL_ENOMEM, "Unable to copy vector!\n" );
810 }
811
812 const UINT4 *oldfit = PulsarGetParamFitFlag( origin, this->name );
813 UINT4 *newfit = ( UINT4 * )XLALCalloc( newerr->length, sizeof( UINT4 ) );
814
815 if ( newfit ) {
816 memcpy( newfit, oldfit, newerr->length * sizeof( UINT4 ) );
817 } else {
818 XLAL_ERROR_VOID( XLAL_ENOMEM, "Unable to copy vector!\n" );
819 }
820
821 PulsarSetREAL8VectorParamErr( target, this->name, ( const REAL8Vector * )newerr, ( const UINT4 * )newfit );
822 }
823 } else if ( this->type == PULSARTYPE_string_t ) {
824 CHAR *parstr = XLALStringDuplicate( *( CHAR ** )this->value );
825 PulsarAddStringParam( target, this->name, ( const CHAR * )parstr );
826 } else {
827 PulsarAddParam( target, this->name, this->value, this->type );
828
829 if ( ( this->type == PULSARTYPE_REAL8_t ) && this->err ) {
830 REAL8 err = *( REAL8 * )this->err;
831
832 const UINT4 *oldfit = PulsarGetParamFitFlag( origin, this->name );
833 UINT4 newfit = oldfit[0];
834
835 PulsarSetREAL8ParamErr( target, this->name, err, newfit );
836 }
837 }
838
839 this = next;
840 if ( this ) {
841 next = this->next;
842 }
843 }
844
845 return;
846}
847
848
849/* create functions for converting par file input e.g. from strings to floats, converting times, converting units */
850enum {
859
860
861#define DEFINE_CONV_FACTOR_FUNCTION( name, convfactor, type ) \
862 void ParConv ## name ( const CHAR *inval, void *out ){ \
863 CHAR *in = XLALStringDuplicate( inval ); \
864 if ( type != CONVSTRING && type != CONVINT ) { \
865 REAL8 x; \
866 if ( type == CONVFLOAT || type == CONVBINUNITS ){ \
867 /* check for exponents expressed as 'D/d' rather than 'E/e' and switch */ \
868 for ( INT4 i = 0; i < (INT4)strlen(in); i++ ) { \
869 if ( in[i] == 'D' || in[i] == 'd' ) { in[i] = 'e'; } \
870 } \
871 x = atof(in); \
872 if ( type == CONVBINUNITS ){ \
873 /* TEMPO2 checks if this is > 1e-7 then it's in units of 1e-12, so needs converting */ \
874 if ( fabs( x ) > 1e-7 ) { x *= 1.e-12; } \
875 /* some of the parameter files in the ATNF catalogue have values of EDOT that are stupidly large e.g. O(1e33). \
876 * These can cause the time delay routines to fail, so if values of EDOT are greater than 10000 ignore them and \
877 * set it to zero */ \
878 if ( x > 10000. ) { x = 0.; } \
879 } \
880 else{ x *= convfactor; } \
881 } \
882 else if ( type == CONVHMS ) { /* convert to radians from hh:mm:ss.s format */ \
883 XLALTranslateHMStoRAD( &x, in ); \
884 } \
885 else if ( type == CONVDMS ) { /* convert to radians from hh:mm:ss.s format */ \
886 XLALTranslateDMStoRAD( &x, in ); \
887 } \
888 else if ( type == CONVMJD ) { /* convert an MJD to a GPS time */ \
889 x = XLALTTMJDtoGPS( atof(in) ); \
890 } \
891 memcpy(out, &x, sizeof(REAL8) ); \
892 } \
893 else if ( type == CONVSTRING ){ memcpy( out, in, strlen(in)+1 ); } \
894 else if ( type == CONVINT ) { /* convert to an integer */ \
895 UINT4 x; \
896 x = (UINT4)atoi(in); \
897 memcpy(out, &x, sizeof(UINT4) ); \
898 } \
899 XLALFree( in ); \
900 }
901
904DEFINE_CONV_FACTOR_FUNCTION( ToFloat, 1., CONVFLOAT ) /* just convert string to float */
905DEFINE_CONV_FACTOR_FUNCTION( ToInt, 0, CONVINT ) /* just convert string to float */
906DEFINE_CONV_FACTOR_FUNCTION( DegsToRads, LAL_PI_180, CONVFLOAT ) /* convert degrees to radians */
907DEFINE_CONV_FACTOR_FUNCTION( MasPerYrToRadPerSec, LAL_PI_180 / ( 3600.e3 * 365.25 * 86400. ), CONVFLOAT ) /* convert milliarcseconds/yr to radians/s */
908DEFINE_CONV_FACTOR_FUNCTION( SecsToRads, LAL_TWOPI / ( 24.*3600. ), CONVFLOAT ) /* convert seconds to radians */
909DEFINE_CONV_FACTOR_FUNCTION( ArcsecsToRads, LAL_PI_180 / 3600., CONVFLOAT ) /* convert arcseconds to radians */
910DEFINE_CONV_FACTOR_FUNCTION( MasToRads, LAL_PI_180 / 3600.0e3, CONVFLOAT ) /* convert milliarcseconds to radians */
911DEFINE_CONV_FACTOR_FUNCTION( InvArcsecsToInvRads, 3600. / LAL_PI_180, CONVFLOAT ) /* convert 1/arcsec to 1/rads */
912DEFINE_CONV_FACTOR_FUNCTION( DaysToSecs, DAYSTOSECS, CONVFLOAT ) /* convert days to seconds */
913DEFINE_CONV_FACTOR_FUNCTION( KpcToMetres, LAL_PC_SI * 1.e3, CONVFLOAT ) /* convert kiloparsecs to metres */
914DEFINE_CONV_FACTOR_FUNCTION( BinaryUnits, 1.e-12, CONVBINUNITS ) /* convert certain binary units as defined in TEMPO2 with factor */
915DEFINE_CONV_FACTOR_FUNCTION( MJDToGPS, 0, CONVMJD ) /* convert from MJD to GPS time */
916DEFINE_CONV_FACTOR_FUNCTION( DegPerYrToRadPerSec, LAL_PI_180 / ( 365.25 * DAYSTOSECS ), CONVFLOAT ) /* convert degs/year to rads/s */
917DEFINE_CONV_FACTOR_FUNCTION( SolarMassToKg, LALPULSAR_TEMPO2_MSUN_SI, CONVFLOAT ) /* convert solar masses to kg */
918DEFINE_CONV_FACTOR_FUNCTION( RAToRads, 0, CONVHMS ) /* convert right ascension to radians */
919DEFINE_CONV_FACTOR_FUNCTION( DecToRads, 0, CONVDMS ) /* convert declination to radians */
920DEFINE_CONV_FACTOR_FUNCTION( MicrosecToSec, 1.e-6, CONVFLOAT ) /* convert microseconds to seconds */
921
922/** Function type definition for a conversion function */
923typedef void ( *ParamConversionFunc )( const CHAR *in, void *out );
924
925/** A strcuture to contain all possible pulsar parameters that can be read in from a par file, and define
926 * the conversion function and type used for each */
927typedef struct tagParConversion {
928 CHAR name[PULSAR_PARNAME_MAX]; /** Parameter name */
929 ParamConversionFunc convfunc; /** Conversion function from string to required value */
930 ParamConversionFunc converrfunc; /** Conversion function for error from string to value */
931 PulsarParamType ptype; /** Parameter type */
934
935#define NUM_PARS 132 /* number of allowed parameters */
936
937/** Initialise conversion structure with most allowed TEMPO2 parameter names and conversion functions
938 * (convert all read in parameters to SI units where necessary). See http://arxiv.org/abs/astro-ph/0603381 and
939 * http://arxiv.org/abs/astro-ph/0603381 for parameter definitions.
940 *
941 * If requiring a new parameter to be read in in should be added to this structure and \c NUM_PARS should be
942 * incremented.
943 */
945 { .name = "F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* vector containing frequency and derivatives (Hz, Hz/s, Hz/s^2 ...) */
946 { .name = "P", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* vector containing period and derivatives (s, s/s, s/s^2 ...) */
947 { .name = "DIST", .convfunc = ParConvKpcToMetres, .converrfunc = ParConvKpcToMetres, .ptype = PULSARTYPE_REAL8_t }, /* distance to pulsar in metres */
948 { .name = "PX", .convfunc = ParConvMasToRads, .converrfunc = ParConvMasToRads, .ptype = PULSARTYPE_REAL8_t }, /* parallax (converted to radians) */
949 { .name = "DM", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* dispersion measure */
950 { .name = "DM1", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /*first derivative of the dispersion measure */
951
952 /* position parameters */
953 { .name = "RA", .convfunc = ParConvRAToRads, .converrfunc = ParConvSecsToRads, .ptype = PULSARTYPE_REAL8_t }, /* right ascension (converted to radians) */
954 { .name = "RAJ", .convfunc = ParConvRAToRads, .converrfunc = ParConvSecsToRads, .ptype = PULSARTYPE_REAL8_t }, /* right ascension (converted to radians) */
955 { .name = "DEC", .convfunc = ParConvDecToRads, .converrfunc = ParConvArcsecsToRads, .ptype = PULSARTYPE_REAL8_t }, /* declination (converted to radians) */
956 { .name = "DECJ", .convfunc = ParConvDecToRads, .converrfunc = ParConvArcsecsToRads, .ptype = PULSARTYPE_REAL8_t }, /* declination (converted to radians) */
957 { .name = "PMRA", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in right ascension (converted to radians/s) */
958 { .name = "PMDEC", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in declination (converted to radians/s) */
959 { .name = "ELONG", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t }, /* ecliptic longitude (converted from degs to rads) */
960 { .name = "ELAT", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t }, /* ecliptic latitude (converted from degs to rads) */
961 { .name = "PMELONG", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in ecliptic longitude (converted to radians/s) */
962 { .name = "PMELAT", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in ecliptic latitude (converted to radians/s) */
963 { .name = "BETA", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t }, /* an alias for ELONG also giving ecliptic latitude (converted from degs to rads) */
964 { .name = "LAMBDA", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t }, /* an alias for ELAT also giving ecliptic longitude (converted from degs to rads) */
965 { .name = "PMBETA", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in ecliptic latitude (converted to radians/s) */
966 { .name = "PMLAMBDA", .convfunc = ParConvMasPerYrToRadPerSec, .converrfunc = ParConvMasPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* proper motion in ecliptic longitude (converted to radians/s) */
967
968 /* epoch parameters */
969 { .name = "PEPOCH", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* period epoch (saved as GPS time) */
970 { .name = "POSEPOCH", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* position epoch (saved as GPS time) */
971 { .name = "DMEPOCH", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* dispersion measure epoch (saved as GPS time) */
972
973 /* glitch parameters */
974 { .name = "GLEP", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch epochs (saved as GPS times) */
975 { .name = "GLPH", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch phase offsets (rotational phase) */
976 { .name = "GLF0", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch frequency offsets (Hz) */
977 { .name = "GLF1", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch frequency derivative offsets (Hz/s) */
978 { .name = "GLF2", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch second frequency derivative offsets (Hz/s^2) */
979 { .name = "GLF0D", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch decaying frequency component (Hz) */
980 { .name = "GLTD", .convfunc = ParConvDaysToSecs, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8Vector_t }, /* glitch decay time (seconds) */
981
982 /* string parameters */
983 { .name = "NAME", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* pulsar name */
984 { .name = "PSR", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* pulsar name */
985 { .name = "PSRJ", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* pulsar J name */
986 { .name = "PSRB", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* pulsar B name */
987
988 /* binary parameters */
989 { .name = "BINARY", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* binary model type */
990 { .name = "A1", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* projected semi-major axis (light seconds) */
991 { .name = "OM", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t }, /* angle of periastron (convert degrees to radians) */
992 { .name = "ECC", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* eccentricity */
993 { .name = "PB", .convfunc = ParConvDaysToSecs, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* orbital period (convert days to seconds) */
994 { .name = "T0", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* time of perisatron (GPS) */
995 { .name = "TASC", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t }, /* time ascending noise for ELL1 model (GPS) */
996 { .name = "EPS1", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* e*sin(w0) for ELL1 model */
997 { .name = "EPS2", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* e*cos(w0) for ELL1 model */
998 { .name = "GAMMA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* relativistic parameter */
999 { .name = "OMDOT", .convfunc = ParConvDegPerYrToRadPerSec, .converrfunc = ParConvDegPerYrToRadPerSec, .ptype = PULSARTYPE_REAL8_t }, /* angle of periastron time derivative (degs/year converted to rad/s) */
1000 { .name = "XDOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* project semi-major axis time derivative (light sec/sec) */
1001 { .name = "PBDOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* period time derivative */
1002 { .name = "EDOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* eccentricity time derivative (1/s) */
1003 { .name = "EPS1DOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1004 { .name = "EPS2DOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1005 { .name = "XPBDOT", .convfunc = ParConvBinaryUnits, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1006 { .name = "SINI", .convfunc = ParConvToString, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_string_t }, /* sin of inclination angle */
1007 { .name = "MTOT", .convfunc = ParConvSolarMassToKg, .converrfunc = ParConvSolarMassToKg, .ptype = PULSARTYPE_REAL8_t }, /* total system mass (convert solar masses to kg) */
1008 { .name = "M2", .convfunc = ParConvSolarMassToKg, .converrfunc = ParConvSolarMassToKg, .ptype = PULSARTYPE_REAL8_t }, /* binary companion mass (convert solar masses to kg) */
1009 { .name = "DR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1010 { .name = "DTHETA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1011 { .name = "SHAPMAX", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* used for DDS model */
1012
1013 /* multiple system terms for BT1P and BT2P models */
1014 { .name = "A1_2", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1015 { .name = "A1_3", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1016 { .name = "OM_2", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t },
1017 { .name = "OM_3", .convfunc = ParConvDegsToRads, .converrfunc = ParConvDegsToRads, .ptype = PULSARTYPE_REAL8_t },
1018 { .name = "PB_2", .convfunc = ParConvDaysToSecs, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t },
1019 { .name = "PB_3", .convfunc = ParConvDaysToSecs, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t },
1020 { .name = "T0_2", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t },
1021 { .name = "T0_3", .convfunc = ParConvMJDToGPS, .converrfunc = ParConvDaysToSecs, .ptype = PULSARTYPE_REAL8_t },
1022 { .name = "ECC_2", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1023 { .name = "ECC_3", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1024 { .name = "FB", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* orbital frequency components for the BTX model */
1025
1026 /* Aberration delay parameters */
1027 { .name = "A0", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1028 { .name = "B0", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1029
1030 /* Kopeikin terms */
1031 { .name = "D_AOP", .convfunc = ParConvInvArcsecsToInvRads, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* inverse arsecs converted to inverse radians */
1032 { .name = "KIN", .convfunc = ParConvDegsToRads, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t },
1033 { .name = "KOM", .convfunc = ParConvDegsToRads, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t },
1034
1035 /* FITWAVES parameters */
1036 { .name = "WAVE_OM", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* fundamental frequency (Hz) of sinusoids used for whitening */
1037 { .name = "WAVEEPOCH", .convfunc = ParConvMJDToGPS, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* frequency epoch (converted from MJD to GPS) */
1038 { .name = "WAVESIN", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* amplitudes of the sine terms for the kth sinusoids */
1039 { .name = "WAVECOS", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8Vector_t }, /* amplitudes of the cosine terms for the kth sinusoids */
1040
1041 /* TEMPO parameters */
1042 { .name = "EPHEM", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* ephemeris type e.g. DE405 */
1043 { .name = "UNITS", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* TEMPO2 units e.g. TDB */
1044 { .name = "START", .convfunc = ParConvMJDToGPS, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* start of observations */
1045 { .name = "FINISH", .convfunc = ParConvMJDToGPS, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* end of observations */
1046 { .name = "NTOA", .convfunc = ParConvToInt, .converrfunc = NULL, .ptype = PULSARTYPE_UINT4_t }, /* number of TOAs in observation */
1047 { .name = "TRES", .convfunc = ParConvMicrosecToSec, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* timing residual (convert microseconds to seconds) */
1048 { .name = "CLK", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* The observatory clock */
1049 { .name = "T2CMETHOD", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* Method for transforming from terrestrial to celestial frame */
1050 { .name = "TIMEEPH", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* Which time ephemeris to use */
1051
1052 /* GW parameters */
1053 { .name = "H0", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* gravitational wave amplitude */
1054 { .name = "APLUS", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* plus polarisation component of GW amplitude */
1055 { .name = "ACROSS", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* cross polarisation component of GW amplitude */
1056 { .name = "PHI0", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* gravitational wave initial phase (radians) */
1057 { .name = "PSI", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* gravitational wave polarisation angle (radians) */
1058 { .name = "COSIOTA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* cosine of source inclination angle */
1059 { .name = "IOTA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* the source inclination angle in radians */
1060 { .name = "C22", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW amplitude of C22 component */
1061 { .name = "C21", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW amplitude of C21 component */
1062 { .name = "PHI22", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase of C22 component (radians) */
1063 { .name = "PHI21", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase of C21 component (radians) */
1064 { .name = "CGW", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* speed of gravitational waves as a fraction of the speed of light */
1065 { .name = "LAMBDAPIN", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* parameters from http://uk.arxiv.org/abs/0909.4035 */
1066 { .name = "COSTHETA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1067 { .name = "THETA", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1068 { .name = "I21", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1069 { .name = "I31", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t },
1070 { .name = "Q22", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* the l=m=2 mass quadrupole in kg m^2 */
1071
1072 /* GW non-GR polarisation mode amplitude parameters (assuming emission at twice the rotation frequency) */
1073 { .name = "HPLUS", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor plus polarisation amplitude */
1074 { .name = "HCROSS", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor cross polarisation amplitude */
1075 { .name = "PSITENSOR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor angle polarisation */
1076 { .name = "PHI0TENSOR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase for tensor modes */
1077 { .name = "HSCALARB", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar breathing mode polarisation amplitude */
1078 { .name = "HSCALARL", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar longitudinal polarisation amplitude */
1079 { .name = "PSISCALAR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar angle polarisation */
1080 { .name = "PHI0SCALAR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase for scalar modes */
1081 { .name = "HVECTORX", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector x-mode amplitude */
1082 { .name = "HVECTORY", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector y-mode amplitude */
1083 { .name = "PSIVECTOR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector angle polarisation */
1084 { .name = "PHI0VECTOR", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector polarisation initial phase */
1085
1086 /* GW non-GR polarisation mode amplitude parameters (assuming emission at the rotation frequency) */
1087 { .name = "H0_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GR 21 polarisation amplitude */
1088 { .name = "HPLUS_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor plus polarisation amplitude */
1089 { .name = "HCROSS_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor cross polarisation amplitude */
1090 { .name = "PSITENSOR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW tensor angle polarisation */
1091 { .name = "PHI0TENSOR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase for tensor modes */
1092 { .name = "HSCALARB_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar breathing mode polarisation amplitude */
1093 { .name = "HSCALARL_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar longitudinal polarisation amplitude */
1094 { .name = "PSISCALAR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW scalar angle polarisation */
1095 { .name = "PHI0SCALAR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* initial phase for scalar modes */
1096 { .name = "HVECTORX_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector x-mode amplitude */
1097 { .name = "HVECTORY_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector y-mode amplitude */
1098 { .name = "PSIVECTOR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector angle polarisation */
1099 { .name = "PHI0VECTOR_F", .convfunc = ParConvToFloat, .converrfunc = ParConvToFloat, .ptype = PULSARTYPE_REAL8_t }, /* GW vector polarisation initial phase */
1100
1101 /* transient signal parameters */
1102 { .name = "TRANSIENTWINDOWTYPE", .convfunc = ParConvToString, .converrfunc = NULL, .ptype = PULSARTYPE_string_t }, /* window type, currently "RECT" (rectangular window) or "EXP" (exponentially decaying window) */
1103 { .name = "TRANSIENTSTARTTIME", .convfunc = ParConvMJDToGPS, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t }, /* transient start time in MJD */
1104 { .name = "TRANSIENTTAU", .convfunc = ParConvDaysToSecs, .converrfunc = NULL, .ptype = PULSARTYPE_REAL8_t } /* transient decay time scale (duration for "RECT" and decay constant for "EXP") */
1105};
1106
1107/** \brief Parse a single line from a pulsar parameter file
1108 *
1109 * This will parse a line from the TEMPO-style pulsar parameter file containing the
1110 * parameter given by \c name. The parameter will be added to the \c par structure.
1111 *
1112 * This function can only be used internally.
1113 */
1114static INT4 ParseParLine( PulsarParameters *par, const CHAR *name, FILE *fp )
1115{
1116 INT4 nread = 0; /* number of values read from line */
1118 /* three potential values on the line */
1120 INT4 i = 0, tmpnum = 0;
1121 CHAR *nname = NULL; /* duplicate of name */
1122
1123 if ( par == NULL ) {
1124 XLAL_ERROR( XLAL_EINVAL, "Error... PulsarParameter structure is not initialised!\n" );
1125 }
1126 if ( name == NULL ) {
1127 XLAL_ERROR( XLAL_EINVAL, "Error... parameter name is not set!\n" );
1128 }
1129 if ( fp == NULL ) {
1130 XLAL_ERROR( XLAL_EINVAL, "Error... parameter file pointer is NULL!\n" );
1131 }
1132
1133 /* parse the line from the fp pointer */
1134 if ( fgets( str, PULSAR_PARNAME_MAX, fp ) == NULL ) {
1135 XLAL_PRINT_WARNING( "No value for parameter %s", name );
1136 return XLAL_SUCCESS;
1137 }
1138
1139 /* scan the line for values */
1140 nread = sscanf( str, "%s %s %s", str1, str2, str3 );
1141
1142 if ( !nread ) {
1143 XLAL_PRINT_WARNING( "No value for parameter %s", name );
1144 return XLAL_SUCCESS;
1145 }
1146
1147 /* check for parameters with more than one possible name */
1148 if ( !strcmp( name, "E" ) ) {
1149 nname = XLALStringDuplicate( "ECC" );
1150 } else if ( !strcmp( name, "E_1" ) ) {
1151 nname = XLALStringDuplicate( "ECC_1" );
1152 } else if ( !strcmp( name, "E_2" ) ) {
1153 nname = XLALStringDuplicate( "ECC_2" );
1154 } else if ( !strcmp( name, "E_3" ) ) {
1155 nname = XLALStringDuplicate( "ECC_3" );
1156 } else {
1157 nname = XLALStringDuplicate( name );
1158 }
1159
1160 /* perform parameter dependent inputs */
1161 for ( i = 0; i < NUM_PARS; i++ ) {
1162 /* this has to be hard-coded for the F, WAVE, FB and glitch (GL) vector parameters */
1163
1164 /* check for parameters requiring placement in vectors */
1165 INT4 isfreq = ( !strncmp( nname, "F", 1 ) && !strcmp( "F", pc[i].name ) && sscanf( nname + strlen( "F" ), "%d", &tmpnum ) == 1 );
1166 INT4 isperiod = ( !strncmp( nname, "P", 1 ) && !strcmp( "P", pc[i].name ) && sscanf( nname + strlen( "P" ), "%d", &tmpnum ) == 1 );
1167 INT4 iswave = ( ( !strncmp( nname, "WAVE", 4 ) && ( !strcmp( "WAVESIN", pc[i].name ) || !strcmp( "WAVECOS", pc[i].name ) ) ) && strcmp( "WAVE_OM", nname ) && strcmp( "WAVEEPOCH", nname ) );
1168 INT4 isfb = ( !strncmp( nname, "FB", 2 ) && !strcmp( "FB", pc[i].name ) );
1169 INT4 isgl = ( !strncmp( nname, "GL", 2 ) && strstr( nname, pc[i].name ) && !strncmp( nname, pc[i].name, strcspn( nname, "_" ) ) );
1170
1171 if ( !strcmp( nname, pc[i].name ) || isfreq || isperiod || iswave || isfb || isgl ) {
1172 UINT4 num = 0, nsize = 0;
1173
1174 if ( pc[i].convfunc == NULL ) {
1175 XLAL_PRINT_WARNING( "No conversion function for parameter %s. Skipping parameter.\n", pc[i].name );
1176 return XLAL_SUCCESS;
1177 }
1178
1179 /* add parameter */
1180 if ( isfreq || isperiod || isfb ) { /* add frequency/period and frequency/period derivative values, or FB values */
1181 REAL8Vector *ptr = NULL, *eptr = NULL;
1182 UINT4 *fitFlag = NULL;
1183
1184 if ( isfb ) {
1185 if ( strlen( nname ) > strlen( "FB" ) ) {
1186 if ( sscanf( nname + strlen( "FB" ), "%d", &num ) != 1 ) {
1187 XLAL_ERROR( XLAL_EINVAL, "Error...problem reading %s number from par file.\n", nname );
1188 }
1189 }
1190 } else if ( isfreq ) {
1191 if ( strlen( nname ) > strlen( "F" ) ) {
1192 if ( sscanf( nname + strlen( "F" ), "%d", &num ) != 1 ) {
1193 XLAL_ERROR( XLAL_EINVAL, "Error...problem reading %s number from par file.\n", nname );
1194 }
1195 }
1196 } else {
1197 if ( strlen( nname ) > strlen( "P" ) ) {
1198 if ( sscanf( nname + strlen( "P" ), "%d", &num ) != 1 ) {
1199 XLAL_ERROR( XLAL_EINVAL, "Error...problem reading %s number from par file.\n", nname );
1200 }
1201 }
1202 }
1203
1204 void *val = ( void * )XLALMalloc( PulsarTypeSize[PULSARTYPE_REAL8_t] );
1205
1206 pc[i].convfunc( str1, val );
1207
1208 nsize = num + 1;
1209 if ( PulsarCheckParam( par, pc[i].name ) ) {
1210 // count any previous values
1211 const REAL8Vector *cptr = PulsarGetREAL8VectorParam( par, pc[i].name );
1212 nsize = ( cptr->length > num + 1 ) ? cptr->length : ( num + 1 );
1213 }
1214
1215 // pointer for parameter values
1216 ptr = XLALCreateREAL8Vector( nsize );
1217 memset( ptr->data, 0, sizeof( REAL8 )*nsize );
1218
1219 // pointers for error and fit flags
1220 eptr = XLALCreateREAL8Vector( nsize );
1221 memset( eptr->data, 0, sizeof( REAL8 )*nsize ); // initialise with zeros
1222 fitFlag = ( UINT4 * )XLALCalloc( nsize, sizeof( UINT4 ) ); // set fit values to zero (no-fit) by default
1223
1224 if ( PulsarCheckParam( par, pc[i].name ) ) {
1225 // copy any previous values
1226 const REAL8Vector *cptr = PulsarGetREAL8VectorParam( par, pc[i].name );
1227 memcpy( ptr->data, cptr->data, cptr->length * sizeof( REAL8 ) );
1228
1229 /* copy any error values, as updating REAL8Vector parameters destroys the original PulsarParam structure */
1231 const UINT4 *ff = PulsarGetParamFitFlag( par, pc[i].name );
1232
1233 // copy any previous error values
1234 memcpy( eptr->data, ceptr->data, sizeof( REAL8 )*cptr->length );
1235 memcpy( fitFlag, ff, sizeof( UINT4 )*cptr->length );
1236 }
1237
1238 ptr->data[num] = *( REAL8 * )val;
1239
1240 // NOTE: this will destroy the previous held PulsarParam structure for this parameter
1242
1243 // (re-)add errors
1244 PulsarSetREAL8VectorParamErr( par, pc[i].name, ( const REAL8Vector * )eptr, ( const UINT4 * )fitFlag );
1245
1246 XLALFree( val );
1248 XLALDestroyREAL8Vector( eptr );
1249 XLALFree( fitFlag );
1250 } else if ( iswave && nread == 2 ) { /* add WAVE values */
1251 REAL8Vector *ptr1 = NULL, *ptr2 = NULL;
1252
1253 if ( strlen( nname ) > strlen( "WAVE" ) ) {
1254 if ( sscanf( nname + strlen( "WAVE" ), "%d", &num ) != 1 ) {
1255 XLAL_ERROR( XLAL_EINVAL, "Error...problem reading %s number from par file.\n", nname );
1256 }
1257 } else {
1258 break;
1259 }
1260
1261 num--; /* WAVE values start from WAVE1, so subtract 1 from the num for the vector index */
1262
1263 void *val1 = ( void * )XLALMalloc( PulsarTypeSize[PULSARTYPE_REAL8_t] );
1264 void *val2 = ( void * )XLALMalloc( PulsarTypeSize[PULSARTYPE_REAL8_t] );
1265
1266 pc[i].convfunc( str1, val1 );
1267 pc[i].convfunc( str2, val2 );
1268
1269 nsize = num + 1;
1270 if ( PulsarCheckParam( par, "WAVESIN" ) && PulsarCheckParam( par, "WAVECOS" ) ) {
1271 // count number of values
1272 const REAL8Vector *cptr = PulsarGetREAL8VectorParam( par, pc[i].name );
1273
1274 nsize = ( cptr->length > num + 1 ) ? cptr->length : ( num + 1 );
1275 }
1276
1277 // pointers for parameter values
1278 ptr1 = XLALCreateREAL8Vector( nsize );
1279 ptr2 = XLALCreateREAL8Vector( nsize );
1280 memset( ptr1->data, 0, sizeof( REAL8 )*nsize );
1281 memset( ptr2->data, 0, sizeof( REAL8 )*nsize );
1282
1283 if ( PulsarCheckParam( par, "WAVESIN" ) && PulsarCheckParam( par, "WAVECOS" ) ) {
1284 const REAL8Vector *cptr1 = PulsarGetREAL8VectorParam( par, "WAVESIN" );
1285 const REAL8Vector *cptr2 = PulsarGetREAL8VectorParam( par, "WAVECOS" );
1286
1287 // copy previous values
1288 memcpy( ptr1->data, cptr1->data, ptr1->length * sizeof( REAL8 ) );
1289 memcpy( ptr2->data, cptr2->data, ptr2->length * sizeof( REAL8 ) );
1290 }
1291
1292 ptr1->data[num] = *( REAL8 * )val1;
1293 ptr2->data[num] = *( REAL8 * )val2;
1294
1295 PulsarAddREAL8VectorParam( par, "WAVESIN", ( const REAL8Vector * )ptr1 );
1296 PulsarAddREAL8VectorParam( par, "WAVECOS", ( const REAL8Vector * )ptr2 );
1297
1298 XLALFree( val1 );
1299 XLALFree( val2 );
1300 XLALDestroyREAL8Vector( ptr1 );
1301 XLALDestroyREAL8Vector( ptr2 );
1302
1303 /* there are no errors on the wave parameters, so set defaults of zero and then break */
1304 UINT4 *fitFlag = ( UINT4 * )XLALCalloc( nsize, sizeof( UINT4 ) );
1305 REAL8Vector *eptr = XLALCreateREAL8Vector( nsize );
1306 memset( eptr->data, 0, sizeof( REAL8 )*nsize );
1307 PulsarSetREAL8VectorParamErr( par, "WAVESIN", ( const REAL8Vector * )eptr, ( const UINT4 * )fitFlag );
1308 PulsarSetREAL8VectorParamErr( par, "WAVECOS", ( const REAL8Vector * )eptr, ( const UINT4 * )fitFlag );
1309 XLALDestroyREAL8Vector( eptr );
1310 XLALFree( fitFlag );
1311
1312 break;
1313 } else if ( isgl ) { /* get glitch parameters */
1314 REAL8Vector *ptr = NULL, *eptr = NULL;
1315 UINT4 *fitFlag = NULL;
1316
1317 if ( sscanf( strstr( nname, "_" ) + 1, "%d", &num ) != 1 ) {
1318 XLAL_ERROR( XLAL_EINVAL, "Error...problem reading %s number from par file.\n", nname );
1319 }
1320
1321 num--; /* GL values start from e.g. GLF0_1, so subtract 1 from the num for the vector index */
1322
1323 void *val = ( void * )XLALMalloc( PulsarTypeSize[PULSARTYPE_REAL8_t] );
1324
1325 pc[i].convfunc( str1, val );
1326
1327 nsize = num + 1;
1328 if ( PulsarCheckParam( par, pc[i].name ) ) {
1329 // count previous values
1330 const REAL8Vector *cptr = PulsarGetREAL8VectorParam( par, pc[i].name );
1331 nsize = ( cptr->length > num + 1 ) ? cptr->length : ( num + 1 );
1332 }
1333
1334 // pointer for parameter values
1335 ptr = XLALCreateREAL8Vector( nsize );
1336 memset( ptr->data, 0, sizeof( REAL8 )*nsize );
1337
1338 // pointers for error and fit flags
1339 eptr = XLALCreateREAL8Vector( nsize );
1340 memset( eptr->data, 0, sizeof( REAL8 )*nsize ); // initialise with zeros
1341 fitFlag = ( UINT4 * )XLALCalloc( nsize, sizeof( UINT4 ) ); // set fit values to zero (no-fit) by default
1342
1343 if ( PulsarCheckParam( par, pc[i].name ) ) {
1344 // copy any previous values
1345 const REAL8Vector *cptr = PulsarGetREAL8VectorParam( par, pc[i].name );
1346 memcpy( ptr->data, cptr->data, cptr->length * sizeof( REAL8 ) );
1347
1348 /* copy any error values, as updating REAL8Vector parameters destroys the original PulsarParam structure */
1350 const UINT4 *ff = PulsarGetParamFitFlag( par, pc[i].name );
1351
1352 // copy any previous error values
1353 memcpy( eptr->data, ceptr->data, sizeof( REAL8 )*cptr->length );
1354 memcpy( fitFlag, ff, sizeof( UINT4 )*cptr->length );
1355 }
1356
1357 ptr->data[num] = *( REAL8 * )val;
1358
1360
1361 // (re-)add errors
1362 PulsarSetREAL8VectorParamErr( par, pc[i].name, ( const REAL8Vector * )eptr, ( const UINT4 * )fitFlag );
1363
1364 XLALFree( val );
1366 XLALDestroyREAL8Vector( eptr );
1367 XLALFree( fitFlag );
1368 } else {
1369 if ( pc[i].ptype != PULSARTYPE_string_t ) {
1370 void *val = ( void * )XLALMalloc( PulsarTypeSize[pc[i].ptype] );
1371 pc[i].convfunc( str1, val );
1372 PulsarAddParam( par, pc[i].name, val, pc[i].ptype );
1373 XLALFree( val );
1374 } else {
1375 CHAR *val = XLALStringDuplicate( str1 );
1376 PulsarAddStringParam( par, pc[i].name, ( const CHAR * )val );
1377 XLALFree( val );
1378 }
1379 }
1380
1381 /* check for error values */
1382 if ( ( nread == 2 && strcmp( str2, "1" ) ) || ( nread == 3 && !strcmp( str2, "1" ) ) ) {
1383 if ( pc[i].converrfunc == NULL ) {
1384 XLAL_PRINT_WARNING( "No conversion function for parameter %s error. No error being set.\n", pc[i].name );
1385 return XLAL_SUCCESS;
1386 }
1387
1388 void *val = ( void * )XLALMalloc( PulsarTypeSize[PULSARTYPE_REAL8_t] );
1389 UINT4 isFit = 0;
1390
1391 /* get the fit flag */
1392 if ( isfreq || isperiod || isfb || isgl ) { /* do this for REAL8Vector parameters */
1394 const UINT4 *fitFlag = PulsarGetParamFitFlag( par, pc[i].name );
1395
1396 XLAL_CHECK( ptr != NULL, XLAL_EFUNC, "No default errors set for parameter %s", pc[i].name );
1397
1398 // copy of parameter errors and fit flags
1399 REAL8Vector *eptr = XLALCreateREAL8Vector( ptr->length );
1400 UINT4 *ff = ( UINT4 * )XLALMalloc( sizeof( UINT4 ) * ptr->length );
1401
1402 memcpy( eptr->data, ptr->data, sizeof( REAL8 )*ptr->length );
1403 memcpy( ff, fitFlag, sizeof( UINT4 )*ptr->length );
1404
1405 if ( nread == 2 ) {
1406 pc[i].converrfunc( str2, val );
1407 } else {
1408 if ( !strcmp( str2, "1" ) ) {
1409 isFit = 1; /* a fit flag is set to one */
1410 }
1411 pc[i].converrfunc( str3, val );
1412 }
1413
1414 eptr->data[num] = *( REAL8 * )val;
1415 ff[num] = isFit;
1416
1417 PulsarSetREAL8VectorParamErr( par, pc[i].name, ( const REAL8Vector * )eptr, ( const UINT4 * )ff );
1418
1419 XLALDestroyREAL8Vector( eptr );
1420 XLALFree( ff );
1421 } else {
1422 if ( nread == 2 ) {
1423 pc[i].converrfunc( str2, val );
1424 } else {
1425 if ( !strcmp( str2, "1" ) ) {
1426 isFit = 1;
1427 }
1428 pc[i].converrfunc( str3, val );
1429 }
1430
1431 PulsarSetREAL8ParamErr( par, pc[i].name, *( REAL8 * )val, isFit );
1432 }
1433
1434 XLALFree( val );
1435 }
1436
1437 break;
1438 }
1439 }
1440
1441 XLALFree( nname );
1442
1443 return XLAL_SUCCESS;
1444}
1445
1447/* read in the pulsar parameter file */
1448PulsarParameters *XLALReadTEMPOParFile( const CHAR *pulsarAndPath )
1449{
1450 FILE *fp = NULL;
1451 CHAR str[PULSAR_PARNAME_MAX]; /* string to contain first value on line */
1452
1453 PulsarParameters *par = XLALCalloc( sizeof( *par ), 1 );
1454
1455 /* open file */
1456 if ( ( fp = fopen( pulsarAndPath, "r" ) ) == NULL ) {
1457 XLAL_PRINT_ERROR( "Error... Cannot open .par file %s\n", pulsarAndPath );
1458 XLALFree( par );
1460 }
1461
1462 int nread = 0; /* number of parameters read from a line in the par file */
1463 int UNUSED c;
1464
1465 /* read in the par file */
1466 while ( !feof( fp ) ) {
1467 /* Read in a line from the parameter file */
1468 nread = fscanf( fp, "%s", str );
1469 if ( nread == 1 ) {
1470 /* check for comment line and skip to end of line */
1471 if ( str[0] == '#' ) {
1472 c = fscanf( fp, "%*[^\n]" );
1473 continue;
1474 }
1475
1476 CHAR upperName[PULSAR_PARNAME_MAX];
1477 XLALStringCopy( upperName, str, PULSAR_PARNAME_MAX );
1478 strtoupper( upperName );
1479
1480 if ( XLAL_SUCCESS != ParseParLine( par, upperName, fp ) ) {
1481 XLAL_PRINT_WARNING( "Parameter \"%s\" could not be successfully parsed from par file", str );
1482 }
1483 }
1484 }
1485
1486 /* check for linked parameters SINI and KIN */
1487 if ( PulsarCheckParam( par, "SINI" ) ) {
1488 CHAR *sini = XLALStringDuplicate( PulsarGetStringParam( par, "SINI" ) );
1489 strtoupper( sini );
1490
1491 REAL8 sinid;
1492
1493 PulsarRemoveParam( par, "SINI" );
1494
1495 if ( !strcmp( sini, "KIN" ) ) {
1496 if ( PulsarCheckParam( par, "KIN" ) ) {
1497 sinid = sin( PulsarGetREAL8Param( par, "KIN" ) );
1498 PulsarAddParam( par, "SINI", &sinid, PULSARTYPE_REAL8_t );
1499 } else {
1500 XLAL_PRINT_ERROR( "Error... KIN not set in .par file %s\n", pulsarAndPath );
1502 XLALFree( par );
1504 }
1505 } else {
1506 sinid = atof( sini );
1507 PulsarAddParam( par, "SINI", &sinid, PULSARTYPE_REAL8_t );
1508 }
1509
1510 XLALFree( sini );
1511 }
1512
1513 fclose( fp );
1514
1515 return par;
1516}
1517
1519/* function to print out to screen all the pulsar parameters and there associated errors */
1521{
1522 fprintf( stderr, "xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx\n" );
1523 fprintf( stderr, "PULSAR %s :\n", params.name );
1524 fprintf( stderr, "sky position:\tra %.7lf +/- %.3le rads, dec %.7lf +/- %.3le rads\n", params.ra,
1525 params.raErr, params.dec, params.decErr );
1526 if ( params.pmra != 0. || params.pmdec != 0. )
1527 fprintf( stderr, "proper motion:\tra %.4le +/- %.1le rads/s, dec %.4le +/- %.1le rads/s\n",
1528 params.pmra, params.pmraErr, params.pmdec, params.pmdecErr );
1529 if ( params.pepoch != 0. || params.posepoch != 0. )
1530 fprintf( stderr, "epochs:\tperiod %lf (GPS), position %lf (GPS)\n", params.pepoch,
1531 params.posepoch );
1532 fprintf( stderr, "xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx\n\n" );
1533
1534 fprintf( stderr, "Frequency parameters\n" );
1535 if ( params.f0 != 0. ) {
1536 fprintf( stderr, "\tf0 = %.10lf +/- %.3le (Hz)\n", params.f0, params.f0Err );
1537 }
1538 if ( params.f1 != 0. ) {
1539 fprintf( stderr, "\tf1 = %.5le +/- %.3le (Hz/s)\n", params.f1, params.f1Err );
1540 }
1541 if ( params.f2 != 0. ) {
1542 fprintf( stderr, "\tf1 = %.5le +/- %.3le (Hz/s^2)\n", params.f2, params.f2Err );
1543 }
1544 /* print binary parameters */
1545 if ( params.model != NULL ) {
1546 fprintf( stderr, "\nBinary parameters:\tmodel %s\n", params.model );
1547
1548 fprintf( stderr, "Keplarian parameters:-\n" );
1549 if ( params.Pb != 0. ) {
1550 fprintf( stderr, "\tperiod = %lf +/- %.3le (days)\n", params.Pb, params.PbErr );
1551 }
1552 if ( params.x != 0. )
1553 fprintf( stderr, "\tprojected semi-major axis = %lf +/- %.3le (light sec)\n", params.x,
1554 params.xErr );
1555 if ( params.e != 0. ) {
1556 fprintf( stderr, "\teccentricity = %lf +/- %.3le\n", params.e, params.eErr );
1557 }
1558 if ( params.w0 != 0. )
1559 fprintf( stderr, "\tlongitude of periastron = %lf +/- %.3lf (degs)\n", params.w0,
1560 params.w0Err );
1561 if ( params.T0 != 0. ) {
1562 fprintf( stderr, "\ttime of periastron = %lf +/- %.3lf (GPS)\n", params.T0, params.T0Err );
1563 }
1564 if ( params.Tasc != 0. )
1565 fprintf( stderr, "\ttime of ascending node = %lf +/- %.3lf (GPS)\n", params.Tasc,
1566 params.TascErr );
1567 if ( params.eps1 != 0. )
1568 fprintf( stderr, "\tfirst Laplace-Lagrange parameter (eps1) = %le +/- %.3le\n", params.eps1,
1569 params.eps1Err );
1570 if ( params.eps2 != 0. )
1571 fprintf( stderr, "\tsecond Laplace-Lagrange parameter (eps1) = %le +/- %.3le\n", params.eps2,
1572 params.eps2Err );
1573 if ( params.eps2 != 0. )
1574 fprintf( stderr, "\tsecond Laplace-Lagrange parameter (eps1) = %le +/- %.3le\n", params.eps2,
1575 params.eps2Err );
1576
1577 /*fprintf(stderr, "Post-Newtonian parameters:-\n");
1578 if(params.gamma != 0.)
1579 fprintf(stderr, "\tGravitational redshift parameter = %le +/- %.3le\n", params.gamma,
1580 params.gammaErr);*/
1581
1582 }
1583}
1585
1587{
1588 FILE *fp = NULL;
1589 CHAR *firstline = XLALStringDuplicate( "" );
1590 CHAR onechar[2];
1591 INT4 i = 0, numPars = 0, c = 1, sl = 0;
1592 LALStringVector *tmpparams = NULL; /* temporary parameter names */
1593 LALStringVector *params = NULL;
1594 UINT4Vector *dims = NULL;
1595
1596 /* check the file exists */
1597 if ( access( corfile, F_OK ) != 0 ) {
1598 XLAL_PRINT_ERROR( "Error... correlation matrix file does not exist!\n" );
1600 }
1601
1602 /* open file */
1603 if ( ( fp = fopen( corfile, "r" ) ) == NULL ) {
1604 XLAL_PRINT_ERROR( "Error... cannot open correlation matrix file!\n" );
1606 }
1607
1608 /* read in first line of the file */
1609 while ( !strchr( fgets( onechar, 2, fp ), '\n' ) ) {
1610 firstline = XLALStringAppend( firstline, onechar );
1611 }
1612
1613 sl = strlen( firstline );
1614
1615 /* count the number of parameters */
1616 for ( i = 0; i < sl; i++ ) {
1617 /* use isspace as delimiters could be unknown generic whitespace */
1618 if ( !isspace( firstline[i] ) ) {
1619 if ( c ) {
1620 numPars++;
1621 c = 0;
1622 }
1623 } else {
1624 c = 1;
1625 }
1626 }
1627
1628 /* parse the line and put into the params vector */
1629 rewind( fp ); /* rewind to start of the file */
1630 for ( i = 0; i < numPars; i++ ) {
1631 CHAR tmpStr[128];
1632
1633 if ( fscanf( fp, "%s", tmpStr ) == EOF ) {
1634 XLAL_PRINT_ERROR( "Error... Problem reading first line of correlation\
1635 matrix!\n" );
1637 }
1638
1639 tmpparams = XLALAppendString2Vector( tmpparams, tmpStr );
1640
1641 /* convert some parameter names to a more common convention */
1642 if ( !XLALStringCaseCompare( tmpStr, "RAJ" ) ) { /* convert RAJ to ra */
1644 } else if ( !XLALStringCaseCompare( tmpStr, "DECJ" ) ) { /* convert DECJ to dec */
1646 } else {
1648 }
1649 }
1650
1651 dims = XLALCreateUINT4Vector( 2 );
1652 dims->data[0] = numPars;
1653 dims->data[1] = numPars;
1654
1655 /* set the correlation matrix to the correct size */
1656 cormat = XLALResizeREAL8Array( cormat, dims );
1657
1658 /* read through covariance values */
1659 for ( i = 0; i < numPars; i++ ) {
1660 CHAR tmpStr[128];
1661 INT4 j = 0;
1662
1663 if ( fscanf( fp, "%s", tmpStr ) == EOF ) {
1664 XLAL_PRINT_ERROR( "Error... problem reading in correlation matrix!\n" );
1666 }
1667
1668 if ( strcmp( tmpStr, tmpparams->data[i] ) ) {
1669 XLAL_PRINT_ERROR( "Error... problem reading in correlation matrix. \
1670Parameters not in consistent order!\n" );
1672 }
1673
1674 for ( j = 0; j < i + 1; j++ ) {
1675 REAL8 tmpval = 0.;
1676
1677 if ( fscanf( fp, "%lf", &tmpval ) == EOF ) {
1678 XLAL_PRINT_ERROR( "Error... problem reading in correlation matrix!\n" );
1680 }
1681
1682 /* if off diagonal values are +/-1 set to +/- 0.99999 */
1683 if ( j != i && fabs( tmpval ) == 1. ) {
1684 tmpval *= 0.99999;
1685 }
1686
1687 cormat->data[i * numPars + j] = tmpval;
1688
1689 /* set opposite elements */
1690 if ( j != i ) {
1691 cormat->data[j * numPars + i] = tmpval;
1692 }
1693 }
1694 }
1695
1696 XLALDestroyStringVector( tmpparams );
1697
1698 return params;
1699}
1700
1701
1702/* functions for converting times given in Terrestrial time TT or TDB in MJD to
1703times in GPS - this is important for epochs given in .par files which are in
1704TDB. TT and GPS are different by a factor of 51.184 secs, this is just the
1705historical factor of 32.184 secs between TT and TAI (International Atomic Time)
1706and the other 19 seconds come from the leap seonds added between the TAI and
1707UTC up to the point of definition of GPS time at UTC 01/01/1980 (see
1708http://www.stjarnhimlen.se/comp/time.html for details) */
1709
1710/* a very good paper describing the tranforms between different time systems
1711and why they are necessary can be found in Seidelmann and Fukushima, A&A 265
1712(1992) http://ukads.nottingham.ac.uk/abs/1992A%26A...265..833S */
1713
1714/** This function converts a MJD format time corrected to Terrestrial Time (TT)
1715 * into an equivalent GPS time */
1717{
1718
1719 LIGOTimeGPS GPS;
1720 REAL8 mjdInt, mjdFrac;
1721 mjdFrac = modf( MJD, &mjdInt );
1722 XLAL_CHECK_REAL8( XLALTranslateMJDTTtoGPS( &GPS, ( INT4 )mjdInt, mjdFrac ) != NULL, XLAL_EFUNC );
1723
1724 return XLALGPSGetREAL8( &GPS );
1725}
1726
1727/** If you have an MJD arrival time on the Earth then this will convert it to
1728 * the equivalent GPS time in TDB (see Table 1 of Seidelmann and Fukushima,
1729 * Astronomy & Astrophysics, 265, 833-838 (1992).
1730 *
1731 * Note that XLALBarycenter performs these TDBtoTT corrections (i.e. the
1732 * Einstein delay) when correcting a GPS time on the Earth to TDB. Also, for
1733 * TEMPO produced pulsar epochs given in MJD these are already in the TDB
1734 * system and an equivalent GPS time in the TDB can be calculated just using
1735 * \c XLALTTMJDtoGPS.
1736 */
1738{
1739 REAL8 GPS;
1740 REAL8 T, TDBtoTT;
1741
1742 /* Check not before the start of GPS time */
1743 XLAL_CHECK_REAL8( MJD >= GPS0MJD, XLAL_EDOM, "Input MJD time %.1f is not in range, must be > %.1f.\n", MJD, GPS0MJD );
1744
1745 /* use factors from Table 1 of Seidelmann and Fukushima, Astronomy &
1746 * Astrophysics, 265, 833-838 (1992) where TDB = TDT + P
1747 * and:
1748 * P = 0.0016568 sin(35999.37 degs x T + 357.5 degs) +
1749 0.0000224 sin(32964.5 degs x T + 246.0 degs) +
1750 0.0000138 sin(71998.7 degs x T + 355.0 degs) +
1751 0.0000048 sin(3034.9 degs x T + 25.0 degs) +
1752 0.0000047 sin(34777.3 degs x T + 230.0 degs)
1753 * and T is the elapsed time from J2000 (which has a Julian day date of
1754 * JD 2451545.0) in Julian centuries.*/
1755 T = MJD + ( XLAL_MJD_REF - XLAL_EPOCH_J2000_0_JD );
1756 T /= 36525.; /* covert days to Julian centuries */
1757
1758 /* time diff in seconds (the Einstein delay) */
1759 TDBtoTT = 0.0016568 * sin( ( 35999.37 * T + 357.5 ) * LAL_PI_180 ) +
1760 0.0000224 * sin( ( 32964.5 * T + 246.0 ) * LAL_PI_180 ) +
1761 0.0000138 * sin( ( 71998.7 * T + 355.0 ) * LAL_PI_180 ) +
1762 0.0000048 * sin( ( 3034.9 * T + 25.0 ) * LAL_PI_180 ) +
1763 0.0000047 * sin( ( 34777.3 * T + 230.0 ) * LAL_PI_180 );
1764
1765 /* convert TDB to TT (TDB-TDBtoTT) and then convert TT to GPS */
1766 /* there is the magical number factor of 32.184 + 19 leap seconds to the
1767 * start of GPS time */
1768 GPS = ( MJD - GPS0MJD ) * 86400. - GPS_TDT - TDBtoTT;
1769
1770 return GPS;
1771}
1772
1773
1774/** If you have an MJD arrival time on the Earth then this will convert it to
1775 * the equivalent GPS time in TCB (see Table 1 of Seidelmann and Fukushima,
1776 * Astronomy & Astrophysics, 265, 833-838, 1992).
1777 *
1778 * Note that for default TEMPO2 produced pulsar epochs given in MJD these are
1779 * already in the TCB system and an equivalent GPS time in the TCB can be
1780 * calculated just using \c XLALTTMJDtoGPS. */
1782{
1783 REAL8 GPS;
1784 REAL8 Tdiff;
1785 REAL8 TCBtoTDB;
1786
1787 /* Check not before the start of GPS time (MJD 44244) */
1788 XLAL_CHECK_REAL8( MJD >= GPS0MJD, XLAL_EDOM, "Input MJD time %.1f is not in\
1789 range, must be > %.1f.\n", MJD, GPS0MJD );
1790
1791 /* from Seidelmann and Fukushima we have a linear drift term:
1792 * TCB - TDB = 1.550506e-8 x (JD - 2443144.5) x 86400
1793 */
1794 Tdiff = ( MJD + XLAL_MJD_REF - 2443144.5 ) * 86400.;
1795 TCBtoTDB = 1.550506e-8 * Tdiff;
1796
1797 /* convert from TDB to GPS */
1798 GPS = XLALTDBMJDtoGPS( MJD );
1799
1800 /* add extra factor as the MJD was really in TCB not TDB) */
1801 GPS -= TCBtoTDB;
1802
1803 return GPS;
1804}
1805
1806/// @}
int j
ProcessParamsTable * ptr
#define c
const char * name
Definition: SearchTiming.c:93
#define fscanf
#define fprintf
double e
REAL8Array * XLALResizeREAL8Array(REAL8Array *, UINT4Vector *)
#define XLAL_MJD_REF
#define XLAL_EPOCH_J2000_0_JD
#define LAL_PI_180
uint64_t UINT8
double REAL8
char CHAR
uint32_t UINT4
int32_t INT4
UINT8 XLALCityHash64(const char *buf, size_t len)
int XLALHashTblFind(const LALHashTbl *ht, const void *x, const void **y)
int XLALHashTblAdd(LALHashTbl *ht, void *x)
int XLALHashTblRemove(LALHashTbl *ht, const void *x)
LALHashTbl * XLALHashTblCreate(LALHashTblDtorFcn dtor, LALHashTblHashFcn hash, LALHashTblCmpFcn cmp)
void XLALHashTblDestroy(LALHashTbl *ht)
void * XLALMalloc(size_t n)
void * XLALCalloc(size_t m, size_t n)
void XLALFree(void *p)
char char * XLALStringDuplicate(const char *s)
int XLALStringCaseCompare(const char *s1, const char *s2)
size_t XLALStringCopy(char *dst, const char *src, size_t size)
int char * XLALStringAppend(char *s, const char *append)
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.
static int PulsarHashElemCmp(const void *elem1, const void *elem2)
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.
static PulsarParam * PulsarGetParamItemSlow(const PulsarParameters *pars, const CHAR *name)
Get a pointer to a parameter of a given name from a PulsarParameters structure.
static void del_elem(void *elem)
void ParConvInvArcsecsToInvRads(const CHAR *in, void *out)
Convert the input string from 1/acrseconds to 1/radians.
static PulsarParam * PulsarGetParamItem(const PulsarParameters *pars, const CHAR *name)
Get a pointer to a parameter of a given name from a PulsarParameters structure.
#define DAYSTOSECS
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.
#define GPS_TDT
REAL8 PulsarGetREAL8ParamOrZero(const PulsarParameters *pars, const CHAR *name)
Return a REAL8 parameter if it exists, otherwise return zero.
static const CHAR a2A[256]
Array for conversion from lowercase to uppercase.
#define DEFINE_CONV_FACTOR_FUNCTION(name, convfactor, type)
ParConversion
static void strtoupper(CHAR *s)
Convert string to uppercase.
#define LALPULSAR_TEMPO2_MSUN_SI
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.
static UINT8 PulsarHash(const void *elem)
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
static INT4 ParseParLine(PulsarParameters *par, const CHAR *name, FILE *fp)
Parse a single line from a pulsar parameter 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.
static void * new_elem(const char *name, PulsarParam *itemPtr)
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.
#define NUM_PARS
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.
ParConversion pc[NUM_PARS]
Initialise conversion structure with most allowed TEMPO2 parameter names and conversion functions (co...
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.
#define GPS0MJD
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.
@ CONVSTRING
@ CONVINT
@ CONVDMS
@ CONVFLOAT
@ CONVMJD
@ CONVHMS
@ CONVBINUNITS
@ PULSARTYPE_REAL8Vector_t
@ PULSARTYPE_UINT4_t
@ PULSARTYPE_string_t
@ PULSARTYPE_REAL8_t
void XLALDestroyStringVector(LALStringVector *vect)
LALStringVector * XLALAppendString2Vector(LALStringVector *vect, const CHAR *string)
LIGOTimeGPS * XLALTranslateMJDTTtoGPS(LIGOTimeGPS *gps, INT4 mjdDays, REAL8 mjdFracDays)
void XLALDestroyUINT4Vector(UINT4Vector *vector)
REAL8Vector * XLALCreateREAL8Vector(UINT4 length)
void XLALDestroyREAL8Vector(REAL8Vector *vector)
UINT4Vector * XLALCreateUINT4Vector(UINT4 length)
#define XLAL_ERROR_VOID(...)
#define XLAL_ERROR_REAL8(...)
#define XLAL_ERROR_NULL(...)
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
#define XLAL_PRINT_WARNING(...)
#define XLAL_CHECK_REAL8(assertion,...)
#define XLAL_PRINT_ERROR(...)
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EDOM
XLAL_EIO
XLAL_EINVAL
XLAL_EFAILED
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
T
out
A structure to contain all pulsar parameters and associated errors.
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.
REAL8 * data
REAL8 * data
UINT4 * data
LALInferenceVariableItem * itemPtr
const char * name
double ra
double dec