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