Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
PulsarParametersWrapper.py
Go to the documentation of this file.
1# Copyright (C) 2018 Matthew Pitkin
2#
3# This program is free software; you can redistribute it and/or modify it
4# under the terms of the GNU General Public License as published by the
5# Free Software Foundation; either version 2 of the License, or (at your
6# option) any later version.
7#
8# This program is distributed in the hope that it will be useful, but
9# WITHOUT ANY WARRANTY; without even the implied warranty of
10# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
11# Public License for more details.
12#
13# You should have received a copy of the GNU General Public License along
14# with this program; if not, write to the Free Software Foundation, Inc.,
15# 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
16"""
17
18Wrapper class for the PulsarParameters structure, so that it can be accessed
19in a dictionary-like way.
20
21"""
22
23from __future__ import division, absolute_import, print_function
24
25import os
26import re
27
28import numpy as np
29
30from astropy import units as u
31
32import lal
33import lalpulsar
34
35# set units of parameters in the PulsarParameters structure
36PPUNITS = {
37 "F": u.Hz, # Hz
38 "P": u.s, # seconds
39 "DIST": u.m, # metres
40 "PX": u.rad, # radians
41 "DM": u.pc / (u.cm) ** 3, # cm^-3 pc
42 "DM1": u.pc / (u.cm**3 * u.yr), # pc cm^-3 yr^-1
43 "RA": u.rad, # radians
44 "RAJ": u.rad, # radians
45 "DEC": u.rad, # radians
46 "DECJ": u.rad, # radians
47 "PMRA": u.rad / u.s, # rad/s
48 "PMDEC": u.rad / u.s, # rad/s
49 "ELONG": u.rad, # rad
50 "ELAT": u.rad, # rad
51 "PMELONG": u.rad, # rad/s
52 "PMELAT": u.rad, # rad/s
53 "BETA": u.rad, # rad
54 "LAMBDA": u.rad, # rad
55 "PMBETA": u.rad, # rad/s
56 "PMLAMBDA": u.rad, # rad/s
57 "PEPOCH": u.s, # GPS seconds
58 "POSEPOCH": u.s, # GPS seconds
59 "DMEPOCH": u.s, # GPS seconds
60 "GLEP": u.s, # GPS seconds
61 "GLPH": u.rad, # rad
62 "GLF0": u.Hz, # Hz
63 "GLF1": u.Hz / u.s, # Hz/s
64 "GLF2": u.Hz / u.s**2, # Hz s^-2
65 "GLF0D": u.Hz, # Hz
66 "GLTD": u.s, # sec
67 "A1": u.s, # light seconds
68 "OM": u.rad, # rad
69 "PB": u.s, # seconds
70 "T0": u.s, # GPS seconds
71 "TASC": u.s, # GPS seconds
72 "EPS1": u.dimensionless_unscaled,
73 "EPS2": u.dimensionless_unscaled,
74 "GAMMA": u.s, # seconds
75 "OMDOT": u.rad / u.s, # rad/s
76 "XDOT": u.s / u.s, # light seconds/sec
77 "PBDOT": u.s / u.s, # s/s
78 "EDOT": 1.0 / u.s, # 1/sec
79 "EPSDOT1": 1.0 / u.s, # 1/sec
80 "EPSDOT2": 1.0 / u.s, # 1/sec
81 "XPBDOT": u.s / u.s, # s/s
82 "SINI": u.dimensionless_unscaled,
83 "MTOT": u.kg, # kg
84 "M2": u.kg, # kg
85 "DR": u.dimensionless_unscaled,
86 "DTHETA": u.dimensionless_unscaled,
87 "SHAPMAX": u.dimensionless_unscaled,
88 "A1_2": u.s, # light seconds
89 "A1_3": u.s, # light seconds
90 "OM_2": u.rad, # radians
91 "OM_3": u.rad, # radians
92 "PB_2": u.s, # seconds
93 "PB_3": u.s, # seconds
94 "T0_2": u.s, # GPS seconds
95 "T0_3": u.s, # GPS seconds
96 "FB": u.Hz, # Hz
97 "A0": u.s, # seconds
98 "B0": u.s, # seconds
99 "D_AOP": 1.0 / u.rad, # 1/rad
100 "KIN": u.rad, # radians
101 "KOM": u.rad, # radians
102 "WAVE_OM": u.Hz, # Hz
103 "WAVEEPOCH": u.s, # GPS seconds
104 "WAVESIN": u.s, # seconds
105 "WAVECOS": u.s, # seconds
106 "START": u.s, # GPS seconds
107 "FINISH": u.s, # GPS seconds
108 "TRES": u.s, # seconds
109 "H0": u.dimensionless_unscaled,
110 "APLUS": u.dimensionless_unscaled,
111 "ACROSS": u.dimensionless_unscaled,
112 "PHI0": u.rad, # radians
113 "PSI": u.rad, # radians
114 "COSIOTA": u.dimensionless_unscaled,
115 "IOTA": u.rad, # radians
116 "C22": u.dimensionless_unscaled,
117 "C21": u.dimensionless_unscaled,
118 "PHI22": u.rad, # radians
119 "PHI21": u.rad, # radians
120 "CGW": u.dimensionless_unscaled,
121 "LAMBDAPIN": u.rad, # radians
122 "COSTHETA": u.dimensionless_unscaled,
123 "THETA": u.rad,
124 "I21": u.dimensionless_unscaled,
125 "I31": u.dimensionless_unscaled,
126 "Q22": u.kg * u.m**2, # kg m^2
127 "H0_F": u.dimensionless_unscaled,
128 "HPLUS": u.dimensionless_unscaled,
129 "HCROSS": u.dimensionless_unscaled,
130 "PSITENSOR": u.rad, # radians
131 "PHI0TENSOR": u.rad, # radians
132 "HSCALARB": u.dimensionless_unscaled,
133 "HSCALARL": u.dimensionless_unscaled,
134 "PSISCALAR": u.rad, # radians
135 "PHI0SCALAR": u.rad, # radians
136 "HVECTORX": u.dimensionless_unscaled,
137 "HVECTORY": u.dimensionless_unscaled,
138 "PSIVECTOR": u.rad, # radians
139 "PHI0VECTOR": u.rad, # radians
140 "HPLUS_F": u.dimensionless_unscaled,
141 "HCROSS_F": u.dimensionless_unscaled,
142 "PSITENSOR_F": u.rad, # radians
143 "PHI0TENSOR_F": u.rad, # radians
144 "HSCALARB_F": u.dimensionless_unscaled,
145 "HSCALARL_F": u.dimensionless_unscaled,
146 "PSISCALAR_F": u.rad, # radians
147 "PHI0SCALAR_F": u.rad, # radians
148 "HVECTORX_F": u.dimensionless_unscaled,
149 "HVECTORY_F": u.dimensionless_unscaled,
150 "PSIVECTOR_F": u.rad, # radians
151 "PHI0VECTOR_F": u.rad, # radians
152 "TRANSIENTSTARTTIME": u.s, # GPS seconds
153 "TRANSIENTTAU": u.s, # seconds
154}
155
156# set units of parameters in a TEMPO-style parameter file if different from above
157TEMPOUNITS = {
158 "DIST": u.kpc, # kpc
159 "PX": u.mas, # milliarcsecs
160 "RA": u.hourangle, # hh:mm:ss.s
161 "RAJ": u.hourangle, # hh:mm:ss.s
162 "DEC": u.deg, # hh:mm:ss.s
163 "DECJ": u.deg, # hh:mm:ss.s
164 "PMRA": u.mas / u.yr, # milliarcsecs/year
165 "PMDEC": u.mas / u.yr, # milliarcsecs/year
166 "ELONG": u.deg, # degrees
167 "ELAT": u.deg, # degrees
168 "PMELONG": u.mas / u.yr, # milliarcsecs/year
169 "PMELAT": u.mas / u.yr, # milliarcsecs/year
170 "BETA": u.deg,
171 "LAMBDA": u.deg,
172 "PMBETA": u.mas / u.yr,
173 "PMLAMBDA": u.mas / u.yr,
174 "PEPOCH": u.d, # MJD(TT) (day)
175 "POSEPOCH": u.d, # MJD(TT) (day)
176 "DMEPOCH": u.d, # MJD(TT) (day)
177 "GLEP": u.d, # MJD(TT) (day)
178 "GLTD": u.d, # days
179 "OM": u.deg, # degs
180 "PB": u.d, # day
181 "T0": u.d, # MJD(TT) (day)
182 "TASC": u.d, # MJD(TT) (day)
183 "OMDOT": u.deg / u.yr, # deg/yr
184 "MTOT": u.solMass, # M_sun
185 "M2": u.solMass, # M_sun
186 "OM_2": u.deg, # degrees
187 "OM_3": u.deg, # degrees
188 "PB_2": u.d, # days
189 "PB_3": u.d, # days
190 "T0_2": u.d, # MJD(TT) (days)
191 "T0_3": u.d, # MJD(TT) (days)
192 "D_AOP": 1.0 / u.arcsec, # 1/arcsec
193 "KIN": u.deg, # degrees
194 "KOM": u.deg, # degrees
195 "WAVEEPOCH": u.d, # MJD(TT) (days)
196 "START": u.d, # MJD(TT) (days)
197 "FINISH": u.d, # MJD(TT) (days)
198 "TRES": u.us, # microsecs
199 "TRANSIENTSTARTTIME": u.d, # MJD(TT) (day)
200 "TRANSIENTTAU": u.d, # days
201}
202
203# set units of error values in tempo if different from above
204TEMPOERRUNITS = {
205 "RA": u.s, # second
206 "RAJ": u.s, # second
207 "DEC": u.arcsec, # arcsecond
208 "DECJ": u.arcsec, # arcsecond
209}
210
211
212class PulsarParametersPy(object):
213 """
214 A class to wrap the SWIG-wrapped lalpulsar.PulsarParameters structure.
215
216 This class lets you access the structure in a more Pythonic way, as well as providing
217 a nice format for holding pulsar (`.par`) parameter files.
218
219 The class can be used to set numerical values (double precision, unsigned integers), strings,
220 or vectors of floating point values, e.g.:
221
222 >>> from lalpulsar.PulsarParametersWrapper import PulsarParametersPy
223 >>> pppy = PulsarParametersPy() # an empty structure
224 >>> pppy['DECJ'] = 0.23 # set a numerical value
225 >>> pppy['BINARY'] = 'BT' # set a string value
226 >>> pppy['F'] = [10.2, 1.4e-11] # set a vector of float values
227
228 Args:
229 pp (PulsarParameters, str): a lalpulsar.PulsarParameters structure, or a string giving the
230 path to a TEMPO-style (`.par`) pulsar parameter file. If nothing is given then an empty
231 lalpulsar.PulsarParameters structure is created. The `read()` method can subsequently
232 be used to read in a `.par` file, or parameters can be added.
233
234 Examples:
235 An example of initialising the class with a previously created `lalpulsar.PulsarParameters`
236 structure is:
237
238 >>> import lalpulsar
239 >>> from lalpulsar.PulsarParametersWrapper import PulsarParametersPy
240 >>> # read in a pulsar parameter file
241 >>> pp = lalpulsar.ReadTEMPOParFile('apulsar.par')
242 >>> # view as a PulsarParametersPy object
243 >>> pppy = PulsarParametersPy(pp)
244
245 The same thing could be achieved more directly using:
246
247 >>> pppy = PulsarParametersPy('apulsar.par')
248
249 or, equivalently with:
250
251 >>> pppy = PulsarParametersPy()
252 >>> pppy.read('apulsar.par')
253
254 """
255
256 keynames = [] # parameter names in PulsarParameters structure
257 length = 0 # number of parameters
258 _pulsarparameters = None
259
260 def __init__(self, pp=None):
261 # if pp is None create empty PulsarParameters structure
262 if pp is None:
263 self._pulsarparameters = lalpulsar.PulsarParameters()
264 else:
265 # check if pp is a pulsar parameters type or a (par file)
266 if not isinstance(pp, lalpulsar.PulsarParameters) and isinstance(pp, str):
267 if os.path.isfile(pp):
268 # try reading in file
269 self.read(pp)
270 else:
271 raise ValueError("Input string does not point to a file")
272 elif isinstance(pp, lalpulsar.PulsarParameters):
273 self._pulsarparameters = pp
274 else:
275 raise ValueError(
276 "Expected 'lalpulsar.PulsarParameters' type, string, or None"
277 )
278
279 def __len__(self):
280 _ = self.keys() # length is counted in the keys() method
281
282 return self.lengthlength
283
284 def __str__(self):
285 return self.pp_to_str()
286
287 def __getitem__(self, key):
288 """
289 Get value from pulsar parameters
290 """
291
292 if self._pulsarparameters is None:
293 return None
294
295 # check if key finishes with "_ERR", in which case check for error value
296 geterr = False
297 tkey = key
298 if key[-4:].upper() == "_ERR":
299 geterr = True
300 tkey = key[:-4] # get the actual parameter key name
301
302 # check if the key is asking for an individual parameter from a vector parameter
303 # (e.g. 'F0' gets the first value from the 'F' vector.
304 # NOTE: this is problematic for glitch parameters, e.g., GLF0, which could provide
305 # values to multiple glitches, so this cannot be used to get individual glitch
306 # parameters).
307 sname = re.sub(r"_\d", "", tkey) if "_" in tkey else re.sub(r"\d", "", tkey)
308 sidx = None
309 indkey = None
310 if sname != tkey and tkey[0:2] != "GL":
311 # check additional index is an integer
312 try:
313 sidx = (
314 int(tkey.split("_")[-1]) if "_" in tkey else int(tkey[len(sname) :])
315 )
316 except ValueError:
317 pass
318
319 # change tkey for checking parameter exists
320 if sidx is not None:
321 indkey = tkey # key with index
322 tkey = sname
323
324 # check if parameter key is present otherwise switch back to original name
325 # (needed for, e.g., 'H0', 'PHI0', ...)
326 if not lalpulsar.PulsarCheckParam(self._pulsarparameters, tkey):
327 tkey = indkey
328
329 # check if parameter given by the key is present
330 if not lalpulsar.PulsarCheckParam(self._pulsarparameters, tkey):
331 return None
332
333 # get type of parameter
334 ptype = lalpulsar.PulsarGetParamType(self._pulsarparameters, tkey)
335
336 if ptype == lalpulsar.PULSARTYPE_REAL8_t:
337 if not geterr:
338 value = lalpulsar.PulsarGetREAL8Param(self._pulsarparameters, tkey)
339 else:
340 value = lalpulsar.PulsarGetREAL8ParamErr(self._pulsarparameters, tkey)
341 elif ptype == lalpulsar.PULSARTYPE_REAL8Vector_t:
342 if not geterr:
343 if sidx is None:
344 tmpvalue = lalpulsar.PulsarGetREAL8VectorParam(
345 self._pulsarparameters, tkey
346 )
347 value = (
348 tmpvalue.data
349 ) # 'data' in a REAL8Vector gets returned as a numpy array
350 else:
351 value = lalpulsar.PulsarGetREAL8VectorParamIndividual(
352 self._pulsarparameters, indkey
353 )
354 else:
355 if sidx is None:
356 tmpvalue = lalpulsar.PulsarGetREAL8VectorParamErr(
357 self._pulsarparameters, tkey
358 )
359 value = tmpvalue.data
360 else:
361 value = lalpulsar.PulsarGetREAL8VectorParamErrIndividual(
362 self._pulsarparameters, indkey
363 )
364 elif ptype == lalpulsar.PULSARTYPE_string_t:
365 if not geterr:
366 value = lalpulsar.PulsarGetStringParam(self._pulsarparameters, tkey)
367 else:
368 raise ValueError("String-type cannot have an error")
369 elif ptype == lalpulsar.PULSARTYPE_UINT4_t:
370 if not geterr:
371 value = lalpulsar.PulsarGetUINT4Param(self._pulsarparameters, tkey)
372 else:
373 raise ValueError("UINT4-type cannot have an error")
374 else:
375 raise ValueError("Unrecognised type")
376
377 return value
378
379 def __setitem__(self, key, value):
380 """
381 Set the value of a key
382 """
383
384 # if parameter exists remove it
385 if lalpulsar.PulsarCheckParam(self._pulsarparameters, key):
386 lalpulsar.PulsarRemoveParam(self._pulsarparameters, key)
387
388 if isinstance(value, float):
389 lalpulsar.PulsarAddREAL8Param(self._pulsarparameters, key, value)
390 elif isinstance(value, str):
391 lalpulsar.PulsarAddStringParam(self._pulsarparameters, key, value)
392 elif isinstance(value, int):
393 if value < 0.0: # store negative integers as floats
394 lalpulsar.PulsarAddREAL8Param(self._pulsarparameters, key, float(value))
395 else:
396 lalpulsar.PulsarAddUINT4Param(self._pulsarparameters, key, value)
397 elif isinstance(value, list) or isinstance(value, np.ndarray):
398 tarray = lal.CreateREAL8Vector(len(value))
399 for i, tv in enumerate(value):
400 if isinstance(tv, float):
401 tarray.data[i] = tv
402 else:
403 raise ValueError("Non-float value in list or array")
404 lalpulsar.PulsarAddREAL8VectorParam(self._pulsarparameters, key, tarray)
405 else:
406 raise ValueError("Data-type not one of know types")
407
408 def convert_to_units(self, name, value):
409 """
410 Convert parameter values to equivalent dimensional versions.
411
412 Args:
413 name (str): the name of the parameter to convert
414 value: the value to the parameter to convert
415
416 Returns:
417 :class:`astropy.unit.Unit`: a unit class with dimensions for a float parameters, or a
418 list containing unit classes for a list or :class:`numpy.ndarray`
419 """
420
421 if isinstance(value, str):
422 # don't make any changes for a string type
423 return value
424
425 uname = name.upper()
426
427 ppunit = u.dimensionless_unscaled
428 if uname in PPUNITS:
429 ppunit = PPUNITS[uname]
430
431 if isinstance(value, np.ndarray) or isinstance(value, list):
432 cvalue = []
433 bunit = ppunit
434 for v in value:
435 cvalue.append(v * ppunit)
436 if uname in ["F", "FB", "P"]: # frequency/period values
437 ppunit *= bunit # increment unit (e.g. Hz -> Hz/s, Hz/s -> Hz/s^2)
438 else:
439 cvalue = value * ppunit
440
441 return cvalue
442
443 def convert_to_tempo_units(self, name, value, iserr=False):
444 """
445 Convert from PulsarParameter units to TEMPO-par-file units
446
447 Args:
448 name (str): a parameter name
449 value (float, list, :class:`numpy.ndarray`): the value of the parameter
450 iserr (bool): state whether where converting the parameter's error value
451
452 Returns:
453 :class:`astropy.unit.Unit`: a unit class with dimensions for a float parameters, or a
454 list containing unit classes for a list or :class:`numpy.ndarray`
455 """
456
457 from astropy.coordinates import ICRS, Angle
458 from astropy.time import Time
459
460 # for certain binary parameters there is an anomoly that their value
461 # may have been rescaled (I think this is a hangover from a TEMPO definition compared to
462 # the TEMPO2 definition)
463 binaryunits = ["XDOT", "PBDOT", "EPS1DOT", "EPS2DOT", "XPBDOT"]
464
465 # the names of epoch parameters that are held as GPS times, but must be converted back to
466 # MJD for a TEMPO-style par file
467 epochpars = [
468 "POSEPOCH",
469 "PEPOCH",
470 "WAVEEPOCH",
471 "T0",
472 "TASC",
473 "T0_2",
474 "T0_3",
475 "START",
476 "FINISH",
477 "DMEPOCH",
478 "GLEP",
479 "TRANSIENTSTARTTIME",
480 ]
481
482 uname = name.upper()
483
484 ppunit = None
485 if uname in PPUNITS:
486 ppunit = PPUNITS[uname]
487
488 tempounit = None
489 if uname in TEMPOUNITS:
490 if not iserr:
491 tempounit = TEMPOUNITS[uname]
492 else:
493 if uname not in TEMPOERRUNITS:
494 tempounit = TEMPOUNITS[uname]
495 else:
496 tempounit = TEMPOERRUNITS[uname]
497
498 if ppunit is None:
499 if uname in binaryunits:
500 # for these binary parameters there is a TEMPO2 oddity that has to be corrected for
501 if abs(value) / 1e-12 > 1e-7:
502 value /= 1e-12
503
504 if isinstance(value, str):
505 return value
506 else:
507 return value * u.dimensionless_unscaled
508
509 # convert to dimensionful value
510 pvalue = self.convert_to_units(uname, value)
511
512 if ppunit == tempounit or tempounit is None:
513 tempounit = ppunit
514
515 if uname in ["RA", "RAJ"]:
516 if not iserr:
517 # convert right ascension in radians into a string output format
518 c = ICRS(pvalue, 0.0 * u.rad)
519 cvalue = c.ra.to_string(unit=tempounit, sep=":", precision=12, pad=True)
520 else:
521 angle = Angle(pvalue)
522 cvalue = (
523 angle.hms[0] * (60.0**2) + angle.hms[1] * 60.0 + angle.hms[2]
524 ) * tempounit
525 elif uname in ["DEC", "DECJ"] and not iserr:
526 c = ICRS(0.0 * u.rad, pvalue)
527 cvalue = c.dec.to_string(unit=tempounit, sep=":", precision=12, pad=True)
528 elif uname in epochpars and not iserr:
529 if isinstance(pvalue, list):
530 cvalue = []
531 for i, pv in enumerate(pvalue):
532 t = Time(pv, format="gps", scale="tt")
533 cvalue.append(t.mjd * tempounit)
534 else:
535 t = Time(pvalue, format="gps", scale="tt")
536 cvalue = t.mjd * tempounit
537 elif uname in binaryunits and not iserr:
538 # for these binary parameters there is a TEMPO2 oddity that has to be corrected for
539 if abs(pvalue.value) / 1e-12 > 1e-7:
540 pvalue.value /= 1e-12
541 cvalue = pvalue.to(tempounit)
542 else:
543 # perform conversion
544 if isinstance(pvalue, list):
545 cvalue = []
546 bunit = tempounit
547 for i, pv in enumerate(pvalue):
548 cvalue.append(pv.to(tempounit))
549 if uname in ["F", "FB", "P"]: # frequency/period values
550 tempounit *= (
551 bunit # increment unit (e.g. Hz -> Hz/s, Hz/s -> Hz/s^2)
552 )
553 else:
554 cvalue = pvalue.to(tempounit)
555
556 return cvalue
557
558 def parameter(self, name, withunits=False, tempounits=False):
559 """
560 Return the parameter given by name.
561
562 Args:
563 name (str): the name of the parameter to return
564 withunits (bool): if True return the parameter in a form with its appropriate units
565 tempounits (bool): of True return the parameter converted into the units required in a
566 TEMPO-style parameter file
567 """
568
569 value = self[name]
570
571 if value is None:
572 return None
573
574 if name.upper()[-4:] == "_ERR":
575 uname = name.upper()[:-4]
576 iserr = True
577 else:
578 uname = name
579 iserr = False
580
581 if tempounits:
582 ovalue = self.convert_to_tempo_units(uname, value, iserr=iserr)
583 elif withunits:
584 ovalue = self.convert_to_units(uname, value)
585 else:
586 ovalue = value
587
588 return ovalue
589
590 def keys(self):
591 """
592 Return a list of the parameter names stored in the PulsarParameters structure
593 """
594
595 thisitem = self._pulsarparameters.head
596 self.keynameskeynames = [] # clear any previous key names
598 while thisitem:
599 tname = thisitem.name
600 self.keynameskeynames.append(tname)
601 self.lengthlength += 1
602
603 thisitem = thisitem.next # move on to next value
604
605 # reverse the order of the names, so they're in the same order as read from a par file
606 self.keynameskeynames = self.keynameskeynames[::-1]
607
608 return self.keynameskeynames
609
610 def values(self):
611 """
612 Return the values of each parameter in the structure
613 """
614
615 tvalues = []
616
617 keys = self.keys()
618 for key in keys:
619 if lalpulsar.PulsarCheckParam(self._pulsarparameters, key):
620 # get type of parameter
621 ptype = lalpulsar.PulsarGetParamType(self._pulsarparameters, key)
622
623 if ptype == lalpulsar.PULSARTYPE_REAL8_t:
624 value = lalpulsar.PulsarGetREAL8Param(self._pulsarparameters, key)
625 elif ptype == lalpulsar.PULSARTYPE_REAL8Vector_t:
626 tmpvalue = lalpulsar.PulsarGetREAL8VectorParam(
627 self._pulsarparameters, key
628 )
629 value = (
630 tmpvalue.data
631 ) # 'data' in a REAL8Vector gets returned as a numpy array
632 elif ptype == lalpulsar.PULSARTYPE_string_t:
633 value = lalpulsar.PulsarGetStringParam(self._pulsarparameters, key)
634 elif ptype == lalpulsar.PULSARTYPE_UINT4_t:
635 value = lalpulsar.PulsarGetUINT4Param(self._pulsarparameters, key)
636 else:
637 raise ValueError("UINT4-type cannot have an error")
638 else:
639 raise ValueError("Could not find {} in strcuture".format(key))
640
641 tvalues.append(value)
642
643 return tvalues
644
645 def as_dict(self):
646 """
647 Return the contents (not error at the moment) of the structure as a dictionary
648 """
649
650 tdict = {}
651
652 for tpair in self.items():
653 tdict[tpair[0]] = tpair[1]
654
655 return tdict
656
657 def items(self):
658 """
659 Return list of item tuples for each parameter in the structure
660 """
661
662 tkeys = self.keys()
663 tvalues = self.values()
664
665 titems = []
666
667 for tk, tv in zip(tkeys, tvalues):
668 titems.append((tk, tv))
669
670 return titems
671
672 def read(self, filename):
673 """
674 Read a TEMPO-style parameter file into a PulsarParameters structure
675
676 Args:
677 filename (str): the path to the pulsar `.par` file.
678 """
679
680 # remove existing pulsarparameters
681 if self._pulsarparameters is not None:
682 del self._pulsarparameters
683
684 pp = lalpulsar.ReadTEMPOParFile(filename)
685
686 if pp is None:
687 raise IOError(
688 "Problem reading in pulsar parameter file '{}'".format(filename)
689 )
690
691 self._pulsarparameters = pp
692
693 def get_error(self, name):
694 """
695 Return the error value for a particular parameter
696
697 Args:
698 name (str): the name of the parameter
699 """
700
701 try:
702 if name.upper()[-4:] == "_ERR":
703 return self[name.upper()]
704 else:
705 uname = name.upper() + "_ERR"
706 return self[uname]
707 except ValueError:
708 return None
709
710 def get_fitflag(self, name):
711 """
712 Return the "fit flag" (a 1 or 0 depending whether the parameter with fit for by TEMPO(2)
713 in the `.par` file).
714
715 Args:
716 name (str): the name of the parameter
717 """
718
719 if name.upper() not in self.keys():
720 return 0.0
721
722 fitflag = lalpulsar.PulsarGetParamFitFlagAsVector(self._pulsarparameters, name)
723
724 if len(fitflag.data) > 1 or isinstance(self[name.upper()], (list, np.ndarray)):
725 return fitflag.data
726 else:
727 return fitflag.data[0]
728
729 def pp_to_str(self, precision=19):
730 """
731 Convert the PulsarParameter structure to a string in the format of a
732 TEMPO-style `.par` file.
733
734 Args:
735 precision (int): the number of decimal places for an output value
736
737 Returns:
738 str: the contents in TEMPO-format
739 """
740
741 # output string format (set so that values should line up)
742 mkl = (
743 max([len(kn) for kn in self.keys()]) + 5
744 ) # max key length for output alignment
745 vlb = precision + 10 # allow extra space for minus sign/exponents
746 outputstr = "{{name: <{0}}}{{value: <{1}}}{{fitflag}}\t{{error}}".format(
747 mkl, vlb
748 )
749
750 parstr = ""
751
752 parsedwaves = False
753
754 # get names of parameters in file
755 for item in self.items():
756 key = item[0]
757 value = item[1]
758
759 if key in ["WAVESIN", "WAVECOS"] and parsedwaves:
760 # already added the FITWAVES values
761 continue
762
763 # get error
764 evalue = self.get_error(key)
765
766 # check for required conversion back to TEMPO-par file units
767 if key in TEMPOUNITS:
768 uvalue = self.convert_to_tempo_units(key, value)
769 if evalue is not None:
770 uevalue = self.convert_to_tempo_units(key, evalue, iserr=True)
771
772 # just get values without units
773 if isinstance(uvalue, list):
774 tvalue = []
775 tevalue = []
776 for tv, te in zip(uvalue, uevalue):
777 tvalue.append(tv.value)
778 tevalue.append(te.value)
779 elif isinstance(uvalue, str):
780 tvalue = uvalue
781 tevalue = uevalue.value
782 else:
783 tvalue = uvalue.value
784 tevalue = uevalue.value
785 else:
786 tvalue = value
787 tevalue = evalue
788
789 fitflag = None
790 if evalue is not None:
791 fitflag = self.get_fitflag(key)
792
793 oevalue = "" # default output error value
794 ofitflag = " " # default output fit flag value (a single space for alignment purposes)
795 if isinstance(tvalue, list) or isinstance(tvalue, np.ndarray):
796 idxoffset = 0
797 idxsep = ""
798 if key in [
799 "WAVESIN",
800 "WAVECOS",
801 "GLEP",
802 "GLPH",
803 "GLF0",
804 "GLF1",
805 "GLF2",
806 "GLF0D",
807 "GLTD",
808 ]:
809 # the TEMPO variable name for these parameter start with an index a 1
810 idxoffset = 1
811
812 if key[:2] == "GL":
813 # glitch parameters have an "_" seperating the name from the index
814 idxsep = "_"
815
816 if key in ["WAVESIN", "WAVECOS"]:
817 # do things differently for FITWAVES parameters
818 parsedwaves = True
819
820 if key == "WAVESIN":
821 wavesin = tvalue
822 wavecos = self["WAVECOS"]
823 else:
824 wavesin = self["WAVESIN"]
825 wavecos = tvalue
826
827 for ws, wc in zip(wavesin, wavecos):
828 precstrs = "{{0:.{}f}}".format(precision) # print out float
829 if ws < 1e-6 or ws > 1e6:
830 # print out float in scientific notation
831 precstrs = "{{0:.{}e}}".format(precision)
832
833 precstrc = "{{0:.{}f}}".format(precision) # print out float
834 if wc < 1e-6 or wc > 1e6:
835 # print out float in scientific notation
836 precstrc = "{{0:.{}e}}".format(precision)
837
838 outputdic = {}
839 outputdic["name"] = "WAVE{}{}".format(idxsep, idxoffset)
840 idxoffset += 1
841 outputdic["value"] = "{}\t{}".format(
842 precstrs.format(ws), precstrc.format(wc)
843 )
844 outputdic["fitflag"] = ""
845 outputdic["error"] = ""
846
847 parstr += outputstr.format(**outputdic).strip() + "\n"
848 else:
849 for tv, te, tf in zip(tvalue, tevalue, fitflag):
850 precstr = "{{0:.{}f}}".format(precision) # print out float
851 if tv < 1e-6 or tv > 1e6:
852 # print out float in scientific notation
853 precstr = "{{0:.{}e}}".format(precision)
854
855 precstre = "{{0:.{}f}}".format(precision) # print out float
856 if te < 1e-6 or te > 1e6:
857 # print out float in scientific notation
858 precstre = "{{0:.{}e}}".format(precision)
859
860 outputdic = {}
861 outputdic["name"] = "{}{}{}".format(key, idxsep, idxoffset)
862 idxoffset += 1
863 outputdic["value"] = precstr.format(tv)
864 outputdic["fitflag"] = "1" if tf == 1 else ""
865 outputdic["error"] = precstre.format(te) if te != 0.0 else ""
866
867 parstr += outputstr.format(**outputdic).strip() + "\n"
868 else:
869 if isinstance(tvalue, float) or key in ["RA", "RAJ", "DEC", "DECJ"]:
870 if isinstance(tvalue, float):
871 precstr = "{{0:.{}f}}".format(precision) # print out float
872 if tvalue < 1e-6 or tvalue > 1e6:
873 # print out float in scientific notation
874 precstr = "{{0:.{}e}}".format(precision)
875
876 ovalue = precstr.format(tvalue)
877 else:
878 ovalue = tvalue
879
880 precstre = "{{0:.{}f}}".format(precision) # print out float
881 oevalue = ""
882 if tevalue is not None:
883 if tevalue != 0.0:
884 if tevalue < 1e-6 or tevalue > 1e6:
885 # print out float in scientific notation
886 precstre = "{{0:.{}e}}".format(precision)
887
888 oevalue = precstre.format(tevalue)
889
890 if fitflag is not None:
891 if fitflag == 1:
892 ofitflag = "1"
893 else:
894 ovalue = tvalue
895
896 outputdic = {}
897 outputdic["name"] = key
898 outputdic["value"] = ovalue
899 outputdic["fitflag"] = ofitflag
900 outputdic["error"] = oevalue
901
902 parstr += outputstr.format(**outputdic).strip() + "\n"
903
904 return parstr
905
906 def pp_to_par(self, filename, precision=19):
907 """
908 Output the PulsarParameter structure to a `.par` file.
909
910 Args:
911 filename (str): the path to the output file
912 precision (int): the number of decimal places for an output value
913 """
914
915 try:
916 fp = open(filename, "w")
917 except IOError:
918 raise IOError("Could not open file '{}' for writing".format(filename))
919
920 fp.write(self.pp_to_str(precision=precision))
921 fp.close()
922
923 def PulsarParameters(self):
924 """
925 Return the PulsarParameters structure
926 """
927
928 return self._pulsarparameters
929
930 def __deepcopy__(self, memo):
931 """
932 Create a copy of the parameters.
933 """
934
935 newpar = PulsarParametersPy()
936 memo[id(self)] = newpar
937 lalpulsar.PulsarCopyParams(self.PulsarParameters(), newpar.PulsarParameters())
938 return newpar
A class to wrap the SWIG-wrapped lalpulsar.PulsarParameters structure.
def read(self, filename)
Read a TEMPO-style parameter file into a PulsarParameters structure.
def pp_to_par(self, filename, precision=19)
Output the PulsarParameter structure to a .par file.
def get_fitflag(self, name)
Return the "fit flag" (a 1 or 0 depending whether the parameter with fit for by TEMPO(2) in the ....
def convert_to_tempo_units(self, name, value, iserr=False)
Convert from PulsarParameter units to TEMPO-par-file units.
def __deepcopy__(self, memo)
Create a copy of the parameters.
def __setitem__(self, key, value)
Set the value of a key.
def __getitem__(self, key)
Get value from pulsar parameters.
def parameter(self, name, withunits=False, tempounits=False)
Return the parameter given by name.
def pp_to_str(self, precision=19)
Convert the PulsarParameter structure to a string in the format of a TEMPO-style ....
def as_dict(self)
Return the contents (not error at the moment) of the structure as a dictionary.
def keys(self)
Return a list of the parameter names stored in the PulsarParameters structure.
def PulsarParameters(self)
Return the PulsarParameters structure.
def get_error(self, name)
Return the error value for a particular parameter.
def convert_to_units(self, name, value)
Convert parameter values to equivalent dimensional versions.
def items(self)
Return list of item tuples for each parameter in the structure.
def values(self)
Return the values of each parameter in the structure.
static const INT4 a
#define max(a, b)