pint package

PINT Is Not TEMPO3!

Subpackages

Submodules

pint.config module

pint.config.datapath(fname)[source]

Returns the full path to the requested data file. Will first search the appdirs user_data_dir (typically $HOME/.local/share/pint on linux) then the installed data files dir (__file__/datafiles). If the file is not found, returns None.

pint.erfautils module

pint.erfautils.astropy_topo_posvels(loc, toas, obsname='obs')[source]
pint.erfautils.topo_posvels(loc, toas, obsname='obs')[source]

Return a list of PosVel instances for the observatory at the TOA times.

Observatory location should be given in the loc argument as an astropy EarthLocation object.

The optional obsname argument will be used as label in the returned PosVel instance.

This routine returns a list of PosVel instances, containing the positions (m) and velocities (m / s) at the times of the toas and referenced to the ITRF geocentric coordinates. This routine is basically SOFA’s pvtob() with an extra rotation from c2ixys.

pint.eventstats module

A collection of various routines for calculating pulsation test test statistics on event data and helper functions.

author: M. Kerr <matthew.kerr@gmail.com>

pint.eventstats.best_m(phases, weights=None, m=100)[source]
pint.eventstats.em_four(phases, m=2, weights=None)[source]

Return the empirical Fourier coefficients up to the mth harmonic. These are derived from the empirical trignometric moments.

pint.eventstats.em_lc(coeffs, dom)[source]

Evaluate the light curve at the provided phases (0 to 1) for the provided coeffs, e.g., as estimated by em_four.

pint.eventstats.from_array(x)[source]
pint.eventstats.h2sig(h)[source]

Shortcut function for calculating sigma from the H-test.

pint.eventstats.hm(phases, m=20, c=4)[source]

Calculate the H statistic (de Jager et al. 1989) for given phases. H_m = max(Z^2_k - c*(k-1)), 1 <= k <= m m == maximum search harmonic c == offset for each successive harmonic

pint.eventstats.hmw(phases, weights, m=20, c=4)[source]

Calculate the H statistic (de Jager et al. 1989) and weight each sine/cosine with the weights in the argument. The distribution is corrected such that the CLT still applies, i.e., it maintains the same calibration as the unweighted version.

pint.eventstats.sf_h20_dj2010(h)[source]

Use the approximate asymptotic calibration from de Jager et al. 2010.

pint.eventstats.sf_hm(h, m=20, c=4, logprob=False)[source]

Return (analytic, asymptotic) survival function (1-F(h)) for the generalized H-test. For more details see:

docstrings for hm, hmw M. Kerr dissertation (arXiv:1101.6072) Kerr, ApJ 732, 38, 2011 (arXiv:1103.2128)

logprob [False] return natural logarithm of probability

pint.eventstats.sf_stackedh(k, h, l=0.398405)[source]

Return the chance probability for a stacked H test assuming the null df for H is exponentially distributed with scale l and that there are k sub-integrations yielding a total TS of h. See, e.g. de Jager & Busching 2010.

pint.eventstats.sf_z2m(ts, m=2)[source]

Return the survival function (chance probability) according to the asymptotic calibration for the Z^2_m test.

ts result of the Z^2_m test

pint.eventstats.sig2h20(sig)[source]

Use approximate (de Jager 2010) relation to invert.

pint.eventstats.sig2sigma(sig, two_tailed=True, logprob=False)[source]

Convert tail probability to “sigma” units, i.e., find the value of the argument for the normal distribution beyond which the integrated tail probability is sig. Note that the default is to interpret this number as the two-tailed value, as this is the quantity that goes to 0 when sig goes to 1.

sig the chance probability

two_tailed [True] interpret sig as two-tailed or one-tailed
probability in converting
logprob [False] if True, the argument is the natural logarithm
of the probability
pint.eventstats.sigma2sig(sigma, two_tailed=True)[source]

Convert “sigma” units to chance probability, i.e., return the integral of the normal distribution from sigma to infinity, or twice that quantity if two_tailed. args —- sigma argument for the normal distribution

two_tailed [True] if True, return twice the integral, else once

pint.eventstats.sigma_trials(sigma, trials)[source]
pint.eventstats.to_array(x)[source]
pint.eventstats.vec(func)[source]
pint.eventstats.z2m(phases, m=2)[source]

Return the Z^2_m test for each harmonic up to the specified m. See de Jager et al. 1989 for definition.

pint.eventstats.z2mw(phases, weights, m=2)[source]

Return the Z^2_m test for each harmonic up to the specified m.

The user provides a list of weights. In the case that they are well-distributed or assumed to be fixed, the CLT applies and the statistic remains calibrated. Nice!

pint.fermi_toas module

pint.fermi_toas.calc_lat_weights(energies, angseps, logeref=4.1, logesig=0.5)[source]

This function computes photon weights based on the energy-dependent PSF, as defined in Philippe Bruel’s SearchPulsation code. It was built by David Smith, based on some code from Lucas Guillemot. This computation uses only the PSF as a function of energy, not a full spectral model of the region, so is less exact than gtsrcprob.

The input values are: energies : Array of photon energies in MeV angseps : Angular separations between photon direction and target

This should be astropy Angle array, such as returned from SkyCoord_photons.separation(SkyCoord_target)

logeref : Parameter from SearchPulsation optimization logesig : Parameter from SearchPulsation optimization

Returns a numpy array of weights (probabilities that the photons came from the target, based on the PSF).

pint.fermi_toas.load_Fermi_TOAs(ft1name, ft2name=None, weightcolumn=None, targetcoord=None, logeref=4.1, logesig=0.5, minweight=0.0)[source]
TOAlist = load_Fermi_TOAs(ft1name,ft2name=None)

Read photon event times out of a Fermi FT1 file and return a list of PINT TOA objects. Correctly handles raw FT1 files, or ones processed with gtbary to have barycentered or geocentered TOAs.

weightcolumn specifies the FITS column name to read the photon weights from. The special value ‘CALC’ causes the weights to be computed empirically as in Philippe Bruel’s SearchPulsation code. logeref and logesig are parameters for the weight computation and are only used when weightcolumn=’CALC’.

When weights are loaded, or computed, events are filtered by weight >= minweight

pint.fermi_toas.phaseogram(mjds, phases, weights=None, title=None, bins=100, rotate=0.0, size=5, alpha=0.25, width=6, maxphs=2.0, file=False)[source]

Make a nice 2-panel phaseogram

pint.fitter module

class pint.fitter.Fitter(toas, model)[source]

Bases: object

Base class for fitter.

The fitting function should be defined as the fit_toas() method.

toas : a pint TOAs instance
The input toas.
model : a pint timing model instance
The initial timing model for fitting.
fit_toas(maxiter=None)[source]
get_allparams()[source]

Return a dict of all param names and values.

get_designmatrix()[source]
get_fitparams()[source]

Return a dict of fittable param names and quantity.

get_fitparams_num()[source]

Return a dict of fittable param names and numeric values.

get_fitparams_uncertainty()[source]
minimize_func(x, *args)[source]

Wrapper function for the residual class, meant to be passed to scipy.optimize.minimize. The function must take a single list of input values, x, and a second optional tuple of input arguments. It returns a quantity to be minimized (in this case chi^2).

reset_model()[source]

Reset the current model to the initial model.

set_fitparams(*params)[source]

Update the “frozen” attribute of model parameters.

Ex. fitter.set_fitparams(‘F0’,’F1’)

set_param_uncertainties(fitp)[source]
set_params(fitp)[source]

Set the model parameters to the value contained in the input dict.

Ex. fitter.set_params({‘F0’:60.1,’F1’:-1.3e-15})

update_resids()[source]

Update the residuals. Run after updating a model parameter.

class pint.fitter.PowellFitter(toas, model)[source]

Bases: pint.fitter.Fitter

A class for Scipy Powell fitting method. This method searches over parameter space. It is a relative basic method.

fit_toas(maxiter=20)[source]
class pint.fitter.WlsFitter(toas=None, model=None)[source]

Bases: pint.fitter.Fitter

A class for weighted least square fitting method. The design matrix is required.

fit_toas(maxiter=1, thershold=False)[source]

Run a linear weighted least-squared fitting method

pint.ltinterface module

An interface for pint compatible to the interface of libstempo

class pint.ltinterface.pintpar(par, parname, *args, **kwargs)[source]

Bases: object

Similar to the parameter class defined in libstempo, this class gives a nice interface to the timing model parameters

Initialize parameter object

Parameters:
  • par – The PINT parameter object
  • parname – The name (key) of the parameter
err
fit
frozen
set
val
class pint.ltinterface.pintpulsar(parfile, timfile=None, warnings=False, fixprefiterrors=True, dofit=False, maxobs=None, units=False)[source]

Bases: object

Pulsar object class with an interface similar to tempopulsar of libstempo

The same init function as used in libstempo

Parameters:
  • parfile – Name of the parfile
  • timfile – Name of the timfile, if we want to load it
  • warnings – Whether we are shoing warnings
  • fixprefiterrors – TODO: check what this should do
  • maxobs – PINT has no need for a maxobs parameter. Included here for compatibility
  • units – Whether or not we are using the ‘units’ interface of libstempo
add_phasejump(mjd, phasejump)[source]

Add a phase jump of value phasejump at time mjd.

Note: due to the comparison observation_SAT > phasejump_SAT in tempo2, the exact MJD itself where the jump was added is not affected.

binarydelay()[source]

Return a long-double numpy array of the delay introduced by the binary model. Does not reform residuals.

chisq(removemean='weighted')[source]

Computes the chisq of current residuals vs errors, removing the noise-weighted residual, unless specified otherwise.

deletedmask()[source]

Returns a numpy.bool array of the delection station of observations. You get a copy of the current values.

designmatrix(updatebats=True, fixunits=True, incoffset=True)[source]

Returns the design matrix [nobs x (ndim+1)] as a numpy.longdouble array for current fit-parameter values. If fixunits=True, adjust the units of the design-matrix columns so that they match the tempo2 parameter units. If fixsigns=True, adjust the sign of the columns corresponding to FX (F0, F1, ...) and JUMP parameters, so that they match finite-difference derivatives. If incoffset=False, the constant phaseoffset column is not included in the designmatrix.

elevation()[source]

Return a numpy double array of the elevation of the pulsar at the time of the observations.

errs(values=None, which='fit')[source]

Same as vals(), but for parameter errors.

fit(iters=1)[source]

Runs iters iterations of the tempo2 fit, recomputing barycentric TOAs and residuals each time.

flags()[source]

Returns the list of flags defined in this dataset (for at least some observations).

flagvals(flagname, values=None)[source]

Returns (or sets, if values are given) a numpy unicode-string array containing the values of flag flagname for every observation.

formresiduals()[source]

Form the residuals

loadparfile(parfile)[source]

Load a parfile with pint

Parameters:parfile – Name of the parfile
loadtimfile(timfile)[source]

Load a pulsar with pint

Parameters:
  • parfile – Name of the parfile
  • timfile – Name of the timfile, if we want to load it
pars(which='fit')[source]

Return parameter keys

phasejumps()[source]

Return an array of phase-jump tuples: (MJD, phase). These are copies.

NOTE: As in tempo2, we refer the phase-jumps to site arrival times. The tempo2 definition of a phasejump is such that it is applied when observation_SAT > phasejump_SAT

pulsenumbers(updatebats=True, formresiduals=True, removemean=True)[source]

Return the pulse number relative to PEPOCH, as detected by tempo2

Returns the pulse numbers as a numpy array. Will update the TOAs/recompute residuals if updatebats/formresiduals is True (default for both). If that is requested, the residual mean is removed removemean is True. All this just like in residuals.

remove_phasejumps()[source]

Remove all phase jumps.

residuals(updatebats=True, formresiduals=True, removemean=True)[source]

Returns residuals

rms(removemean='weighted')[source]

Computes the current residual rms, removing the noise-weighted residual, unless specified otherwise.

savepar(parfile)[source]

Save current par file (calls tempo2’s textOutput(...)).

savetim(timfile)[source]

Save current par file (calls tempo2’s writeTim(...)).

ssbfreqs()[source]

Return SSB frequencies

telescope()[source]

Returns a numpy character array of the telescope for each observation, mapping tempo2 telID values to names by way of the tempo2 runtime file observatory/aliases.

toas(updatebats=False)[source]

Return TDB arrival times in MJDs

vals(values=None, which='fit')[source]

Get (if no values provided) or set the parameter values, depending on which:

  • if which is ‘fit’ (default), fitted parameters;
  • if which is ‘set’, all parameters with a defined value;
  • if which is ‘all’, all parameters;
  • if which is a sequence, all parameters listed there.

Parameter values are returned as a numpy longdouble array.

Values to be set can be passed as a numpy array, sequence (in which case they are taken to correspond to parameters in the order given by pars(which=which)), or dict (in which case which will be ignored).

Notes:

  • Passing values as anything else than numpy longdoubles may result in loss of precision.
  • Not all parameters in the selection need to be set.
  • Setting an unset parameter sets its set flag (obviously).
  • Unlike in earlier libstempo versions, setting a parameter does not set its error to zero.
binarymodel

Return the binary model

freqs

Returns observation frequencies in units of MHz as a numpy.double array.

name

Get or set pulsar name.

ndim

Number of dimensions.

nobs

Number of observations

nphasejumps

Return the number of phase jumps.

pulse_number

Return the pulse number relative to PEPOCH, as detected by tempo2

WARNING: Will be deprecated in the future. Use pulsenumbers.

stoas

Return site arrival times

toaerrs

Return the TOA uncertainty in microseconds

pint.nicer_toas module

pint.nicer_toas.load_NICER_TOAs(eventname)[source]
TOAlist = load_NICER_TOAs(eventname)
Read photon event times out of a NICER event FITS file and return a list of PINT TOA objects. Correctly handles raw NICER files, or ones processed with axBary to have barycentered TOAs.
pint.nicer_toas.nicer_phaseogram(mjds, phases, weights=None, title=None, bins=64, rotate=0.0, size=5, alpha=0.25, width=6, maxphs=2.0, file=False)[source]

Make a nice 2-panel phaseogram

pint.phase module

class pint.phase.Phase[source]

Bases: pint.phase.Phase

Phase class array version

pint.pulsar_ecliptic module

class pint.pulsar_ecliptic.PulsarEcliptic(*args, **kwargs)[source]

Bases: astropy.coordinates.baseframe.BaseCoordinateFrame

A Pulsar Ecliptic coordinate system is defined by rotating ICRS coordinate about x-axis by obliquity angle. Historical, This coordinate is used by tempo/tempo2 for a better fitting error treatment. The obliquity angle values respect to time are given in the file named “ecliptic.dat” in the pint/datafile directory. Parameters ———- representation : BaseRepresentation or None

A representation object or None to have no data (or use the other keywords)
Lambda : Angle, optional, must be keyword
The longitude-like angle corresponding to Sagittarius’ orbit.
Beta : Angle, optional, must be keyword
The latitude-like angle corresponding to Sagittarius’ orbit.
distance : Quantity, optional, must be keyword
The Distance for this object along the line-of-sight.
default_representation
frame_specific_representation_info
name = 'pulsarecliptic'
obliquity = <Quantity 84381.406 arcsec>
pint.pulsar_ecliptic.icrs_to_pulsarecliptic(from_coo, to_frame)[source]
pint.pulsar_ecliptic.load_obliquity_file(filename)[source]
pint.pulsar_ecliptic.pulsarecliptic_to_icrs(from_coo, to_frame)[source]

pint.pulsar_mjd module

class pint.pulsar_mjd.TimePulsarMJD(val1, val2, scale, precision, in_subfmt, out_subfmt, from_jd=False)[source]

Bases: astropy.time.formats.TimeFormat

MJD using tempo/tempo2 convention for time within leap second days. This is only relevant if scale=’utc’, otherwise will act like the standard astropy MJD time format.

set_jds(val1, val2)[source]
name = 'pulsar_mjd'
value

pint.residuals module

class pint.residuals.resids(toa=None, model=None)[source]

Bases: object

calc_chi2()[source]

Return the weighted chi-squared for the model and toas.

calc_phase_resids(weighted_mean=True)[source]

Return timing model residuals in pulse phase.

calc_time_resids(weighted_mean=True)[source]

Return timing model residuals in time (seconds).

get_PSR_freq()[source]

Return pulsar rotational frequency in Hz. model.F0 must be defined.

get_dof()[source]

Return number of degrees of freedom for the model.

get_reduced_chi2()[source]

Return the weighted reduced chi-squared for the model and toas.

pint.solar_system_ephemerides module

pint.solar_system_ephemerides.objPosVel(obj1, obj2, t, ephem)[source]

Compute the position and velocity for solar system obj2 referenced at obj1. This function uses astropy solar system Ephemerides module. Parameters ———- obj1: str

The name of reference solar system object
obj2: str
The name of target solar system object
tdb: Astropy.time.Time object
TDB time in Astropy.time.Time object format
ephem: str
The ephem to for computing solar system object position and velocity
PosVel object.
solar system obj1’s position and velocity with respect to obj2 in the J2000 cartesian coordinate.
pint.solar_system_ephemerides.objPosVel_wrt_SSB(objname, t, ephem, path=None, link=None)[source]

This function computes a solar system object position and velocity respect to solar system barycenter using astropy coordinates get_body_barycentric() method.

The coordinate frame is that of the underlying solar system ephemeris, which has been the ICRF (J2000) since the DE4XX series.

objname: str
Solar system object name. Current support solar system bodies are listed in astropy.coordinates.solar_system_ephemeris.bodies attribution.
t: Astropy.time.Time object
Observation time in Astropy.time.Time object format.
ephem: str
The ephem to for computing solar system object position and velocity
path: str optional
The data directory point to a local ephemeris.
link: str optional
The link where to download the ephemeris.

PosVel object with 3-vectors for the position and velocity of the object

If both path and link are provided. Path will be first to try.

pint.str2ld module

pint.str2ld.str2ldarr1(*args, **kwargs)[source]

pint.tmtestfuncs module

pint.tmtestfuncs.F0(toa, model)[source]

pint.toa module

class pint.toa.TOA(MJD, error=0.0, obs='Barycenter', freq=inf, scale=None, **kwargs)[source]

Bases: object

A time of arrival (TOA) class.

MJD will be stored in astropy.time.Time format, and can be
passed as a double (not recommended), a string, a tuple of component parts (day and fraction of day).

error is the TOA uncertainty in microseconds obs is the observatory name as defined in XXX freq is the observatory-centric frequency in MHz other keyword/value pairs can be specified as needed

# SUGGESTION(paulr): Here, or in a higher level document, the time system # philosophy should be specified. It looks like TOAs are assumed to be # input as UTC(observatory) and then are converted to UTC(???) using the # observatory clock correction file.

Example:
>>> a = TOA((54567, 0.876876876876876), 4.5, freq=1400.0,
...         obs="GBT", backend="GUPPI")
>>> print a
54567.876876876876876:  4.500 us error from 'GBT' at 1400.0000 MHz {'backend': 'GUPPI'}
What happens if IERS data is not available for the date:
>>> a = TOA((154567, 0.876876876876876), 4.5, freq=1400.0,
...         obs="GBT", backend="GUPPI")
Traceback (most recent call last):
omitted

IndexError: (some) times are outside of range covered by IERS table.

Construct a TOA object

MJD : astropy Time, float, or tuple of floats
The time of the TOA, which can be expressed as an astropy Time, a floating point MJD (64 or 128 bit precision), or a tuple of (MJD1,MJD2) whose sum is the full precision MJD (usually the integer and fractional part of the MJD)
obs : string
The observatory code for the TOA
freq : float or astropy Quantity
Frequency corresponding to the TOA. Either a Quantity with frequency units, or a number for which MHz is assumed.
scale : string
Time scale for the TOA time. Defaults to the timescale appropriate to the site, but can be overridden

It is VERY important that all astropy.Time() objects are created with precision=9. This is ensured in the code and is checked for any Time object passed to the TOA constructor.

class pint.toa.TOAs(toafile=None, toalist=None)[source]

Bases: object

A class of multiple TOAs, loaded from zero or more files.

adjust_TOAs(delta)[source]

Apply a time delta to TOAs

Adjusts the time (MJD) of the TOAs by applying delta, which should be a numpy.time.TimeDelta instance with the same shape as self.table[‘mjd’]

delta : astropy.time.TimeDelta
The time difference to add to the MJD of each TOA
apply_clock_corrections(include_bipm=True, include_gps=True)[source]

Apply observatory clock corrections and TIME statments.

Apply clock corrections to all the TOAs where corrections are available. This routine actually changes the value of the TOA, although the correction is also listed as a new flag for the TOA called ‘clkcorr’ so that it can be reversed if necessary. This routine also applies all ‘TIME’ commands and treats them exactly as if they were a part of the observatory clock corrections.

Options to include GPS or BIPM clock corrections are set to True by default in order to give the most accurate clock corrections.

# SUGGESTION(paulr): Somewhere in this docstring, or in a higher level # documentation, the assumptions about the timescales should be specified. # The docstring says apply “correction” but does not say what it is correcting. # Be more specific.

compute_TDBs()[source]

Compute and add TDB and TDB long double columns to the TOA table.

This routine creates new columns ‘tdb’ and ‘tdbld’ in a TOA table for TDB times, using the Observatory locations and IERS A Earth rotation corrections for UT1.

compute_posvels(ephem='DE421', planets=False)[source]

Compute positions and velocities of the observatories and Earth.

Compute the positions and velocities of the observatory (wrt the Geocenter) and the center of the Earth (referenced to the SSB) for each TOA. The JPL solar system ephemeris can be set using the ‘ephem’ parameter. The positions and velocities are set with PosVel class instances which have astropy units.

get_errors()[source]

Return a numpy array of the TOA errors in us

get_flags()[source]

Return a numpy array of the TOA flags

get_freqs()[source]

Return a numpy array of the observing frequencies in MHz for the TOAs

get_mjds(high_precision=False)[source]

With high_precision is True Return an array of the astropy.times (UTC) of the TOAs

With high_precision is False Return an array of toas in mjd as double precision floats

WARNING: Depending on the situation, you may get MJDs in a different scales (e.g. UTC, TT, or TDB) or even a mixture of scales if some TOAs are barycentred and some are not (a perfectly valid situation when fitting both Fermi and radio TOAs)

get_obss()[source]

Return a numpy array of the observatories for each TOA

get_summary()[source]

Return a short ASCII summary of the TOAs.

pickle(filename=None)[source]

Write the TOAs to a .pickle file with optional filename.

print_summary()[source]

Write a summary of the TOAs to stdout.

read_pickle_file(filename)[source]

Read the TOAs from the pickle file specified in filename. Note the filename should include any pickle-specific extensions (ie ”.pickle.gz” or similar), these will not be added automatically.

read_toa_file(filename, process_includes=True, top=True)[source]

Read the given filename and return a list of TOA objects.

Will process INCLUDEd files unless process_includes is False.

write_TOA_file(filename, name='pint', format='Princeton')[source]

Dump current TOA table out as a TOA file

filename : str
File name to write to
format : str
Format specifier for file (‘TEMPO’ or ‘Princeton’) or (‘Tempo2’ or ‘1’)

Currently does not undo any clock corrections that were applied, so TOA file won’t match the input TOA file if any were applied.

pint.toa.format_toa_line(toatime, toaerr, freq, obs, dm=<Quantity 0.0 pc / cm3>, name='unk', flags={}, format='Princeton')[source]

Format TOA line for writing

toatime Time object containing TOA arrival time toaerr TOA error as a Quantity with units freq Frequency as a Quantity with units (NB: value of np.inf is allowed) obs Observatory object

dm DM for the TOA as a Quantity with units (not printed if 0.0 pc/cm^3) name Name to embed in TOA line (conventionally the data file name) format (Princeton | Tempo2) flags Any Tempo2 flags to append to the TOA line

This implementation is currently incomplete in that it will not undo things like TIME statements and probably other things.

columns item 1-1 Observatory (one-character code) ‘@’ is barycenter 2-2 must be blank 16-24 Observing frequency (MHz) 25-44 TOA (decimal point must be in column 30 or column 31) 45-53 TOA uncertainty (microseconds) 69-78 DM correction (pc cm^-3)

First line of file should be “FORMAT 1” TOA format is “file freq sat satErr siteID <flags>”

out : string
Formatted TOA line
pint.toa.get_TOAs(timfile, ephem='DE421', include_bipm=True, include_gps=True, planets=False, usepickle=False)[source]

Convenience function to load and prepare TOAs for PINT use.

Loads TOAs from a ‘.tim’ file, applies clock corrections, computes key values (like TDB), computes the observatory position and velocity vectors, and pickles the file for later use (if requested).

Includes options to specify solar system ephemeris [default DE421], gps clock corrections [default=True], and BIPM clock corrections [default=True].

pint.toa.get_TOAs_list(toa_list, ephem='DE421', include_bipm=True, include_gps=True, planets=False)[source]

Load TOAs from a list of TOA objects.

Compute the TDB time and observatory positions and velocity vectors.

Includes options to specify solar system ephemeris [default DE421], gps clock corrections [default=True], and BIPM clock corrections [default=True].

pint.toa.get_obs(obscode)[source]

Return the standard name for the given code.

pint.toa.parse_TOA_line(line, fmt='Unknown')[source]

Parse a one-line ASCII time-of-arrival.

Return an MJD tuple and a dictionary of other TOA information. The format can be one of: Comment, Command, Blank, Tempo2, Princeton, ITOA, Parkes, or Unknown.

pint.toa.toa_format(line, fmt='Unknown')[source]

Determine the type of a TOA line.

Identifies a TOA line as one of the following types: Comment, Command, Blank, Tempo2, Princeton, ITOA, Parkes, Unknown.

pint.toa_select module

class pint.toa_select.TOASelect(key, key_value)[source]

Bases: object

This class is designed for select toas based on a key word and key value. It will check toa table and do table operation.

check_table_keys(toas)[source]
get_key_section(toas)[source]
get_toa_key_mask(toas)[source]

This is a function to return toas mask depend on key word and key value.

pint.utils module

Miscellaneous potentially-helpful functions.

class pint.utils.PosVel(pos, vel, obj=None, origin=None)[source]

Bases: object

Position/Velocity class.

The class is used to represent the 6 values describing position and velocity vectors. Instances have ‘pos’ and ‘vel’ attributes that are numpy arrays of floats (and can have attached astropy units). The ‘pos’ and ‘vel’ params are 3-vectors of the positions and velocities respectively.

The coordinates are generally assumed to be aligned with ITRF (J2000)

The ‘obj’ and ‘origin’ components are strings that can optionally be used to specify names for endpoints of the vectors. If present, addition/subtraction will check that vectors are being combined in a consistent way.

pint.utils.GEO_WGS84_to_ITRF(lon, lat, hgt)[source]

Convert lat/long/height to rectangular.

Convert WGS-84 references lon, lat, height (using astropy units) to ITRF x,y,z rectangular coords (m) .

pint.utils.MJD_string2longdouble(s)[source]
MJD_string2longdouble(s):
Convert a MJD string to a numpy longdouble
pint.utils.check_all_partials(f, args, delta=1e-06, atol=0.0001, rtol=0.0001)[source]

Check the partial derivatives of a function that returns derivatives.

The function is assumed to return a pair (values, partials), where partials is supposed to be a matrix of the partial derivatives of f with respect to all its arguments. These values are checked against numerical partial derivatives.

pint.utils.data2longdouble(data)[source]

Return a numpy long double scalar form different type of data Parameters ——— data : str, ndarray, or a number Return ——— numpy long double type of data.

pint.utils.ddouble2ldouble(t1, t2, format='jd')[source]
ddouble2ldouble(t1, t2, format=’jd’):
inputs two double-precision numbers representing JD times, converts them to a single long double MJD value
pint.utils.fortran_float(x)[source]

Convert Fortran-format floating-point strings.

returns a copy of the input string with all ‘D’ or ‘d’ turned into ‘e’ characters. Intended for dealing with exponential notation in tempo1-generated parfiles.

pint.utils.has_astropy_unit(x)[source]

has_astropy_unit(x):

Return True/False if x has an astropy unit type associated with it. This is useful, because different data types can still have units associated with them.

pint.utils.is_number(s)[source]

Check if it is a number string.

pint.utils.longdouble2string(x)[source]

Convert numpy longdouble to string

pint.utils.numeric_partial(f, args, ix=0, delta=1e-06)[source]

Compute the partial derivative of f numerically.

This uses symmetric differences to estimate the partial derivative of a function (that takes some number of numeric arguments and may return an array) with respect to one of its arguments.

pint.utils.numeric_partials(f, args, delta=1e-06)[source]

Compute all the partial derivatives of f numerically.

Returns a matrix of the partial derivative of every return value with respect to every input argument. f is assumed to take a flat list of numeric arguments and return a list or array of values.

pint.utils.split_prefixed_name(name)[source]

A utility function that splits a prefixed name. Parameter ———- name : str

Prefixed name
prefixPart : str
The prefix part of the name
indexPart : str
The index part from the name
indexValue : int
The absolute index valeu
pint.utils.str2longdouble(str_data)[source]

Return a numpy long double scalar from the input string, using strtold()

pint.utils.taylor_horner(x, coeffs)[source]

Evaluate a Taylor series of coefficients at x via the Horner scheme. For example, if we want: 10 + 3*x/1! + 4*x^2/2! + 12*x^3/3! with x evaluated at 2.0, we would do: In [1]: taylor_horner(2.0, [10, 3, 4, 12]) Out[1]: 40.0

pint.utils.taylor_horner_deriv(x, coeffs)[source]

Evaluate a Taylor series of coefficients at x via the Horner scheme. For example, if we want: 3/1! + 4*x/2! + 12*x^2/3! with x evaluated at 2.0, we would do: In [1]: taylor_horner_deriv(2.0, [10, 3, 4, 12]) Out[1]: 15.0

pint.utils.time_from_longdouble(t, scale='utc')[source]
pint.utils.time_from_mjd_string(s, scale='utc')[source]

Returns an astropy Time object generated from a MJD string input.

pint.utils.time_to_longdouble(t)[source]

Return an astropy Time value as MJD in longdouble

## SUGGESTION(paulr): This function is at least partly redundant with ## ddouble2ldouble() below...

## Also, is it certain that this calculation retains the full precision?

pint.utils.time_to_mjd_string(t, prec=15)[source]

Print an MJD time with lots of digits (number is ‘prec’).

astropy does not seem to provide this capability (yet?).

pint.utils.time_to_mjd_string_array(t, prec=15)[source]

Print and MJD time array from an astropy time object as array in time.