pint.models.timing_model.TimingModel

class pint.models.timing_model.TimingModel(name='', components=[])[source]

Bases: object

Timing model object built from Components.

This object is the primary object to represent a timing model in PINT. It is normally constructed with get_model(), and it contains a variety of Component objects, each representing a physical process that either introduces delays in the pulse arrival time or introduces shifts in the pulse arrival phase. These components have parameters, described by Parameter objects, and methods. Both the parameters and the methods are accessible through this object using attribute access, for example as model.F0 or model.coords_as_GAL().

Components in a TimingModel objects are accessible through the model.components property, using their class name to index the TimingModel, as model.components["Spindown"]. They can be added and removed with methods on this object, and for many of them additional parameters in families (DMXEP_1234) can be added.

Parameters in a TimingModel object are listed in the model.params object. Each Parameter can be set as free or frozen using its .frozen attribute, and a list of the free parameters is available through the model.free_params property; this can also be used to set which parameters are free. Several methods are available to get and set some or all parameters in the forms of dictionaries.

TimingModel objects also support a number of functions for computing various things like orbital phase, and barycentric versions of TOAs, as well as the various derivatives and matrices needed to support fitting.

TimingModel objects forward attribute lookups to their components, so that you can access any method or attribute (in particular Parameters) of any Component directly on the TimingModel object, for example as model.F0.

TimingModel objects can be written out to .par files using pint.models.timing_model.TimingModel.write_parfile() or . pint.models.timing_model.TimingModel.as_parfile():

>>> model.write_parfile("output.par")

PINT Parameters supported (here, rather than in any Component):

Name / Aliases

Description

Kind

PSR / PSRJ, PSRB

Source name

string

TRACK

Tracking Information

string

EPHEM

Ephemeris to use

string

CLOCK / CLK

Timescale to use

string

UNITS

Units (TDB assumed)

string

START

Start MJD for fitting

d

FINISH

End MJD for fitting

d

RM

Rotation measure

rad / m2

INFO

Tells TEMPO to write some extra information about frontend/backend combinations; -f is recommended

string

TIMEEPH

Time ephemeris to use for TDB conversion; for PINT, always FB90

string

T2CMETHOD

Method for transforming from terrestrial to celestial frame (IAU2000B/TEMPO; PINT only supports ????)

string

BINARY

Pulsar System/Binary model

string

DILATEFREQ

Whether or not TEMPO2 should apply gravitational redshift and time dilation to observing frequency (Y/N; PINT only supports N)

boolean

DMDATA

Was the fit done using per-TOA DM information?

boolean

NTOA

Number of TOAs used in the fitting

integer

CHI2

Chi-squared value obtained during fitting

number

CHI2R

Reduced chi-squared value obtained during fitting

number

TRES

TOA residual after fitting

us

DMRES

DM residual after fitting (wideband only)

pc / cm3

Parameters:
  • name (str, optional) – The name of the timing model.

  • components (list of Component, optional) – The model components for timing model.

Notes

PINT models pulsar pulse time of arrival at observer from its emission process and propagation to observer. Emission generally modeled as pulse ‘Phase’ and propagation. ‘time delay’. In pulsar timing different astrophysics phenomenons are separated to time model components for handling a specific emission or propagation effect.

Each timing model component generally requires the following parts:

  • Timing Parameters

  • Delay/Phase functions which implements the time delay and phase.

  • Derivatives of delay and phase respect to parameter for fitting toas.

Each timing parameters are stored as TimingModel attribute in the type of Parameter delay or phase and its derivatives are implemented as TimingModel Methods.

Variables:
  • name (str) – The name of the timing model

  • component_types (list) – A list of the distinct categories of component. For example, delay components will be register as ‘DelayComponent’.

  • top_level_params (list) – Names of parameters belonging to the TimingModel as a whole rather than to any particular component.

Methods

add_component(component[, order, force, ...])

Add a component into TimingModel.

add_param_from_top(param, target_component)

Add a parameter to a timing model component.

add_tzr_toa(toas)

Create a TZR TOA for the given TOAs object and add it to the timing model.

as_ECL([epoch, ecl])

Return TimingModel in PulsarEcliptic frame.

as_ICRS([epoch])

Return TimingModel in ICRS frame.

as_parfile([start_order, last_order, ...])

Represent the entire model as a parfile string.

companion_radial_velocity(barytimes, massratio)

Return line-of-sight velocity of the companion relative to the system barycenter at barycentric MJD times.

compare(othermodel[, nodmx, ...])

Print comparison with another model

conjunction(baryMJD)

Return the time(s) of the first superior conjunction(s) after baryMJD.

d_delay_d_param(toas, param[, acc_delay])

Return the derivative of delay with respect to the parameter.

d_delay_d_param_num(toas, param[, step])

Return the derivative of delay with respect to the parameter.

d_dm_d_param(data, param)

Return the derivative of DM with respect to the parameter.

d_phase_d_param(toas, delay, param)

Return the derivative of phase with respect to the parameter.

d_phase_d_param_num(toas, param[, step])

Return the derivative of phase with respect to the parameter.

d_phase_d_toa(toas[, sample_step])

Return the finite-difference derivative of phase wrt TOA.

d_phase_d_tpulsar(toas)

Return the derivative of phase wrt time at the pulsar.

d_toasigma_d_param(data, param)

Return the derivative of the scaled TOA uncertainty with respect to the parameter.

delay(toas[, cutoff_component, include_last])

Total delay for the TOAs.

delete_jump_and_flags(toa_table, jump_num)

Delete jump object from PhaseJump and remove its flags from TOA table.

designmatrix(toas[, acc_delay, incfrozen, ...])

Return the design matrix.

dm_covariance_matrix(toas)

Get the DM covariance matrix for noise models.

find_empty_masks(toas[, freeze])

Find unfrozen mask parameters with no TOAs before trying to fit

get_barycentric_toas(toas[, cutoff_component])

Conveniently calculate the barycentric TOAs.

get_component_type(component)

Identify the component object's type.

get_components_by_category()

Return a dict of this model's component objects keyed by the category name.

get_deriv_funcs(component_type[, ...])

Return a dictionary of derivative functions.

get_derived_params([rms, ntoas, returndict])

Return a string with various derived parameters from the fitted model

get_params_dict([which, kind])

Return a dict mapping parameter names to values.

get_params_mapping()

Report which component each parameter name comes from.

get_params_of_component_type(component_type)

Get a list of parameters belonging to a component type.

get_params_of_type_top(param_type)

get_prefix_list(prefix[, start_index])

Return the Quantities associated with a sequence of prefix parameters.

get_prefix_mapping(prefix)

Get the index mapping for the prefix parameters.

items()

jump_flags_to_params(toas)

Add JUMP parameters corresponding to tim_jump flags.

keys()

map_component(component)

Get the location of component.

match_param_aliases(alias)

Return PINT parameter name corresponding to this alias.

noise_model_basis_weight(toas)

noise_model_designmatrix(toas)

noise_model_dimensions(toas)

Number of basis functions for each noise model component.

orbital_phase(barytimes[, anom, radians])

Return orbital phase (in radians) at barycentric MJD times.

param_help()

Print help lines for all available parameters in model.

phase(toas[, abs_phase])

Return the model-predicted pulse phase for the given TOAs.

pulsar_radial_velocity(barytimes)

Return line-of-sight velocity of the pulsar relative to the system barycenter at barycentric MJD times.

remove_component(component)

Remove one component from the timing model.

remove_param(param)

Remove a parameter from timing model.

scaled_dm_uncertainty(toas)

Get the scaled DM data uncertainties noise models.

scaled_toa_uncertainty(toas)

Get the scaled TOA data uncertainties noise models.

search_cmp_attr(name)

Search for an attribute in all components.

set_param_uncertainties(fitp)

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

set_param_values(fitp)

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

setup()

Run setup methods on all components.

toa_covariance_matrix(toas)

Get the TOA covariance matrix for noise models.

total_dispersion_slope(toas)

Calculate the dispersion slope from all the dispersion-type components.

total_dm(toas)

Calculate dispersion measure from all the dispersion type of components.

use_aliases([reset_to_default, ...])

Set the parameters to use aliases as specified upon writing.

validate([allow_tcb])

Validate component setup.

validate_component_types()

Physically motivated validation of a timing model.

validate_toas(toas)

Sanity check to verify that this model is compatible with these toas.

write_parfile(filename[, start_order, ...])

Write the entire model as a parfile.

Attributes

basis_funcs

List of scaled uncertainty functions.

components

All the components in a dictionary indexed by name.

covariance_matrix_funcs

List of covariance matrix functions.

d_phase_d_delay_funcs

List of d_phase_d_delay functions.

delay_deriv_funcs

List of derivative functions for delay components.

delay_funcs

List of all delay functions.

dm_covariance_matrix_funcs

List of covariance matrix functions.

dm_derivs

List of DM derivative functions.

dm_funcs

List of all dm value functions.

fittable_params

List of parameters that are fittable, i.e., the parameters which have a derivative implemented.

free_params

List of all the free parameters in the timing model.

has_correlated_errors

Whether or not this model has correlated errors.

has_time_correlated_errors

Whether or not this model has time-correlated errors.

is_binary

Does the model describe a binary pulsar?

params

List of all parameter names in this model and all its components, in a sensible order.

params_ordered

List of all parameter names in this model and all its components.

phase_deriv_funcs

List of derivative functions for phase components.

phase_funcs

List of all phase functions.

scaled_dm_uncertainty_funcs

List of scaled dm uncertainty functions.

scaled_toa_uncertainty_funcs

List of scaled toa uncertainty functions.

toasigma_derivs

List of scaled TOA uncertainty derivative functions

validate(allow_tcb=False)[source]

Validate component setup.

The checks include required parameters and parameter values, and component types. See also: pint.models.timing_model.TimingModel.validate_component_types().

validate_component_types()[source]

Physically motivated validation of a timing model. This method checks the compatibility of different model components when used together.

This function throws an error if multiple deterministic components that model the same effect are used together (e.g. pint.models.astrometry.AstrometryEquatorial and pint.models.astrometry.AstrometryEcliptic). It emits a warning if a deterministic component and a stochastic component that model the same effect are used together (e.g. pint.models.noise_model.PLDMNoise and pint.models.dispersion_model.DispersionDMX). It also requires that one and only one pint.models.spindown.SpindownBase component is present in a timing model.

property params_ordered

List of all parameter names in this model and all its components. This is the same as params.

property params

List of all parameter names in this model and all its components, in a sensible order.

property free_params

List of all the free parameters in the timing model. Can be set to change which are free.

These are ordered as self.params does.

Upon setting, order does not matter, and aliases are accepted. ValueError is raised if a parameter is not recognized.

On setting, parameter aliases are converted with pint.models.timing_model.TimingModel.match_param_aliases().

property fittable_params

List of parameters that are fittable, i.e., the parameters which have a derivative implemented. These derivatives are usually accessed via the d_delay_d_param and d_phase_d_param methods.

match_param_aliases(alias)[source]

Return PINT parameter name corresponding to this alias.

Parameters:

alias (str) – Parameter’s alias.

Returns:

PINT parameter name corresponding to the input alias.

Return type:

str

get_params_dict(which='free', kind='quantity')[source]

Return a dict mapping parameter names to values.

This can return only the free parameters or all; and it can return the parameter objects, the floating-point values, or the uncertainties.

Parameters:
  • which ("free", "all") –

  • kind ("quantity", "value", "uncertainty") –

Return type:

OrderedDict

get_params_of_component_type(component_type)[source]

Get a list of parameters belonging to a component type.

Parameters:

component_type ("PhaseComponent", "DelayComponent", "NoiseComponent") –

Return type:

list

set_param_values(fitp)[source]

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

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

set_param_uncertainties(fitp)[source]

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

property components

All the components in a dictionary indexed by name.

property delay_funcs

List of all delay functions.

property phase_funcs

List of all phase functions.

property is_binary

Does the model describe a binary pulsar?

orbital_phase(barytimes, anom='mean', radians=True)[source]

Return orbital phase (in radians) at barycentric MJD times.

Parameters:
  • barytimes (Time, TOAs, array-like, or float) – MJD barycentric time(s). The times to compute the orbital phases. Needs to be a barycentric time in TDB. If a TOAs instance is passed, the barycentering will happen automatically. If an astropy Time object is passed, it must be in scale=’tdb’. If an array-like object is passed or a simple float, the time must be in MJD format.

  • anom (str, optional) – Type of phase/anomaly. Defaults to “mean”. Other options are “eccentric” or “true”

  • radians (bool, optional) – Units to return. Defaults to True. If False, will return unitless phases in cycles (i.e. 0-1).

Raises:

ValueError – If anom.lower() is not “mean”, “ecc*”, or “true”, or if an astropy Time object is passed with scale!=”tdb”.

Returns:

The specified anomaly in radians (with unit), unless radians=False, which return unitless cycles (0-1).

Return type:

array

pulsar_radial_velocity(barytimes)[source]

Return line-of-sight velocity of the pulsar relative to the system barycenter at barycentric MJD times.

Parameters:

barytimes (Time, TOAs, array-like, or float) – MJD barycentric time(s). The times to compute the orbital phases. Needs to be a barycentric time in TDB. If a TOAs instance is passed, the barycentering will happen automatically. If an astropy Time object is passed, it must be in scale=’tdb’. If an array-like object is passed or a simple float, the time must be in MJD format.

Raises:

ValueError – If an astropy Time object is passed with scale!=”tdb”.

Returns:

The line-of-sight velocity

Return type:

array

Notes

This is the radial velocity of the pulsar.

See [1]_

companion_radial_velocity(barytimes, massratio)[source]

Return line-of-sight velocity of the companion relative to the system barycenter at barycentric MJD times.

Parameters:
  • barytimes (Time, TOAs, array-like, or float) – MJD barycentric time(s). The times to compute the orbital phases. Needs to be a barycentric time in TDB. If a TOAs instance is passed, the barycentering will happen automatically. If an astropy Time object is passed, it must be in scale=’tdb’. If an array-like object is passed or a simple float, the time must be in MJD format.

  • massratio (float) – Ratio of pulsar mass to companion mass

Raises:

ValueError – If an astropy Time object is passed with scale!=”tdb”.

Returns:

The line-of-sight velocity

Return type:

array

Notes

This is the radial velocity of the companion.

See [1]_

conjunction(baryMJD)[source]

Return the time(s) of the first superior conjunction(s) after baryMJD.

Parameters:

baryMJD (floats or Time) – barycentric (tdb) MJD(s) prior to the conjunction we are looking for. Can be an array.

Raises:

ValueError – If baryMJD is an astropy Time object with scale!=”tdb”.

Returns:

The barycentric MJD(tdb) time(s) of the next superior conjunction(s) after baryMJD

Return type:

float or array

property dm_funcs

List of all dm value functions.

property has_correlated_errors

Whether or not this model has correlated errors.

property has_time_correlated_errors

Whether or not this model has time-correlated errors.

property covariance_matrix_funcs

List of covariance matrix functions.

property dm_covariance_matrix_funcs

List of covariance matrix functions.

property scaled_toa_uncertainty_funcs

List of scaled toa uncertainty functions.

property scaled_dm_uncertainty_funcs

List of scaled dm uncertainty functions.

property basis_funcs

List of scaled uncertainty functions.

property phase_deriv_funcs

List of derivative functions for phase components.

property delay_deriv_funcs

List of derivative functions for delay components.

property dm_derivs

List of DM derivative functions.

property toasigma_derivs

List of scaled TOA uncertainty derivative functions

property d_phase_d_delay_funcs

List of d_phase_d_delay functions.

get_deriv_funcs(component_type, derivative_type='')[source]

Return a dictionary of derivative functions.

Parameters:
  • component_type (str) – Type of component to look for derivatives (“PhaseComponent”, “DelayComponent”, or “NoiseComponent”)

  • derivative_type (str) – Derivative type (“”, “dm”, or “toasigma”. Empty string denotes delay and phase derivatives.)

search_cmp_attr(name)[source]

Search for an attribute in all components.

Return the component, or None.

If multiple components have same attribute, it will return the first component.

get_component_type(component)[source]

Identify the component object’s type.

Parameters:

component (component instance) – The component object need to be inspected.

Note

Since a component can be an inheritance from other component We inspect all the component object bases. “inspect getmro” method returns the base classes (including ‘object’) in method resolution order. The third level of inheritance class name is what we want. Object –> component –> TypeComponent. (i.e. DelayComponent) This class type is in the third to the last of the getmro returned result.

map_component(component)[source]

Get the location of component.

Parameters:

component (str or Component object) – Component name or component object.

Returns:

  • comp (Component object) – Component object.

  • order (int) – The index/order of the component in the component list

  • host_list (List) – The host list of the component.

  • comp_type (str) – The component type (e.g., Delay or Phase)

add_component(component, order=['astrometry', 'jump_delay', 'troposphere', 'solar_system_shapiro', 'solar_wind', 'dispersion_constant', 'dispersion_dmx', 'dispersion_jump', 'pulsar_system', 'frequency_dependent', 'absolute_phase', 'spindown', 'phase_jump', 'wave', 'wavex'], force=False, setup=True, validate=True)[source]

Add a component into TimingModel.

Parameters:
  • component (Component) – The component to be added to the timing model.

  • order (list, optional) – The component category order list. Default is the DEFAULT_ORDER.

  • force (bool, optional) – If true, add a duplicate component. Default is False.

remove_component(component)[source]

Remove one component from the timing model.

Parameters:

component (str or Component object) – Component name or component object.

get_components_by_category()[source]

Return a dict of this model’s component objects keyed by the category name.

add_param_from_top(param, target_component, setup=False)[source]

Add a parameter to a timing model component.

Parameters:
  • param (pint.models.parameter.Parameter) – Parameter instance

  • target_component (str) – Parameter host component name. If given as “” it would add parameter to the top level TimingModel class

  • setup (bool, optional) – Flag to run setup() function.

remove_param(param)[source]

Remove a parameter from timing model.

Parameters:

param (str) – The name of parameter to be removed.

get_params_mapping()[source]

Report which component each parameter name comes from.

get_prefix_mapping(prefix)[source]

Get the index mapping for the prefix parameters.

Parameters:

prefix (str) – Name of prefix.

Returns:

A dictionary with prefix parameter real index as key and parameter name as value.

Return type:

dict

get_prefix_list(prefix, start_index=0)[source]

Return the Quantities associated with a sequence of prefix parameters.

Parameters:
  • prefix (str) – Name of prefix.

  • start_index (int) – The index to start the sequence at (DM1, DM2, … vs F0, F1, …)

Returns:

The .quantity associated with parameter prefix + start_index, prefix + (start_index+1), … up to the last that exists and is set.

Return type:

list of astropy.units.Quantity

Raises:

ValueError – If any prefix parameters exist outside the sequence that would be returned (for example if there are DM1 and DM3 but not DM2, or F0 exists but start_index was given as 1).

param_help()[source]

Print help lines for all available parameters in model.

delay(toas, cutoff_component='', include_last=True)[source]

Total delay for the TOAs.

Return the total delay which will be subtracted from the given TOA to get time of emission at the pulsar.

Parameters:
  • toas (pint.toa.TOAs) – The toas for analysis delays.

  • cutoff_component (str) – The delay component name that a user wants the calculation to stop at.

  • include_last (bool) – If the cutoff delay component is included.

phase(toas, abs_phase=None)[source]

Return the model-predicted pulse phase for the given TOAs.

This is the phase as observed at the observatory at the exact moment specified in each TOA. The result is a pint.phase.Phase object.

add_tzr_toa(toas)[source]

Create a TZR TOA for the given TOAs object and add it to the timing model. This corresponds to TOA closest to the PEPOCH.

total_dm(toas)[source]

Calculate dispersion measure from all the dispersion type of components.

total_dispersion_slope(toas)[source]

Calculate the dispersion slope from all the dispersion-type components.

toa_covariance_matrix(toas)[source]

Get the TOA covariance matrix for noise models.

If there is no noise model component provided, a diagonal matrix with TOAs error as diagonal element will be returned.

dm_covariance_matrix(toas)[source]

Get the DM covariance matrix for noise models.

If there is no noise model component provided, a diagonal matrix with TOAs error as diagonal element will be returned.

scaled_toa_uncertainty(toas)[source]

Get the scaled TOA data uncertainties noise models.

If there is no noise model component provided, a vector with TOAs error as values will be returned.

Parameters:

toas (pint.toa.TOAs) – The input data object for TOAs uncertainty.

scaled_dm_uncertainty(toas)[source]

Get the scaled DM data uncertainties noise models.

If there is no noise model component provided, a vector with DM error as values will be returned.

Parameters:

toas (pint.toa.TOAs) – The input data object for DM uncertainty.

noise_model_dimensions(toas)[source]

Number of basis functions for each noise model component.

Returns a dictionary of correlated-noise components in the noise model. Each entry contains a tuple (offset, size) where size is the number of basis functions for the component, and offset is their starting location in the design matrix and weights vector.

jump_flags_to_params(toas)[source]

Add JUMP parameters corresponding to tim_jump flags.

When a .tim file contains pairs of JUMP lines, the user’s expectation is that the TOAs between each pair of flags will be affected by a JUMP, even if that JUMP does not appear in the .par file. (This is how TEMPO works.) In PINT, those TOAs have a flag attached, -tim_jump N, where N is a number that is different for each JUMPed set of TOAs. The goal of this function is to add JUMP parameters to the model corresponding to these.

Some complexities arise: TOAs may also have -tim_jump flags associated with them, just as flags, for example if such a .tim file were exported in PINT-native format and then reloaded. And models may already have JUMPs associated with some or all tim_jump values.

This function looks at all the tim_jump values and adds JUMP parameters for any that do not have any. It does not change the TOAs object it is passed.

delete_jump_and_flags(toa_table, jump_num)[source]

Delete jump object from PhaseJump and remove its flags from TOA table.

This is a helper function for pintk.

Parameters:
  • toa_table (list or None) – The TOA table which must be modified. In pintk (pulsar.py), for the prefit model, this will be all_toas.table[“flags”]. For the postfit model, it will be None (one set of TOA tables for both models).

  • jump_num (int) – Specifies the index of the jump to be deleted.

get_barycentric_toas(toas, cutoff_component='')[source]

Conveniently calculate the barycentric TOAs.

Parameters:
  • toas (TOAs object) – The TOAs the barycentric corrections are applied on

  • cutoff_delay (str, optional) – The cutoff delay component name. If it is not provided, it will search for binary delay and apply all the delay before binary.

Returns:

Barycentered TOAs.

Return type:

astropy.units.Quantity

d_phase_d_toa(toas, sample_step=None)[source]

Return the finite-difference derivative of phase wrt TOA.

Parameters:
  • toas (pint.toa.TOAs) – The toas when the derivative of phase will be evaluated at.

  • sample_step (float, optional) – Finite difference steps. If not specified, it will take 1000 times the spin period.

d_phase_d_tpulsar(toas)[source]

Return the derivative of phase wrt time at the pulsar.

NOT implemented yet.

d_phase_d_param(toas, delay, param)[source]

Return the derivative of phase with respect to the parameter.

This is the derivative of the phase observed at each TOA with respect to each parameter. This is closely related to the derivative of residuals with respect to each parameter, differing only by a factor of the spin frequency and possibly a minus sign. See pint.models.timing_model.TimingModel.designmatrix() for a way of evaluating many derivatives at once.

The calculation is done by combining the analytical derivatives reported by all the components in the model.

Parameters:
  • toas (pint.toa.TOAs) – The TOAs at which the derivative should be evaluated.

  • delay (astropy.units.Quantity or None) – The delay at the TOAs where the derivatives should be evaluated. This permits certain optimizations in the derivative calculations; the value should be self.delay(toas).

  • param (str) – The name of the parameter to differentiate with respect to.

Returns:

The derivative of observed phase with respect to the model parameter.

Return type:

astropy.units.Quantity

d_delay_d_param(toas, param, acc_delay=None)[source]

Return the derivative of delay with respect to the parameter.

d_phase_d_param_num(toas, param, step=0.01)[source]

Return the derivative of phase with respect to the parameter.

Compute the value numerically, using a symmetric finite difference.

d_delay_d_param_num(toas, param, step=0.01)[source]

Return the derivative of delay with respect to the parameter.

Compute the value numerically, using a symmetric finite difference.

d_dm_d_param(data, param)[source]

Return the derivative of DM with respect to the parameter.

d_toasigma_d_param(data, param)[source]

Return the derivative of the scaled TOA uncertainty with respect to the parameter.

designmatrix(toas, acc_delay=None, incfrozen=False, incoffset=True)[source]

Return the design matrix.

The design matrix is the matrix with columns of d_phase_d_param/F0 or d_toa_d_param; it is used in fitting and calculating parameter covariances.

The value of F0 used here is the parameter value in the model.

The order of parameters that are included is that returned by self.params.

Parameters:
  • toas (pint.toa.TOAs) – The TOAs at which to compute the design matrix.

  • acc_delay

    ???

  • incfrozen (bool) – Whether to include frozen parameters in the design matrix

  • incoffset (bool) – Whether to include the constant offset in the design matrix This option is ignored if a PhaseOffset component is present.

Returns:

  • M (array) – The design matrix, with shape (len(toas), len(self.free_params)+1)

  • names (list of str) – The names of parameters in the corresponding parts of the design matrix

  • units (astropy.units.Unit) – The units of the corresponding parts of the design matrix

Notes

1. We have negative sign here. Since the residuals are calculated as (Phase - int(Phase)) in pulsar timing, which is different from the conventional definition of least square definition (Data - model), we have decided to add a minus sign here in the design matrix so that the fitter keeps the conventional sign.

2. Design matrix entries can be computed only for parameters for which the derivatives are implemented. If a parameter without a derivative is unfrozen while calling this method, it will raise an informative error, except in the case of unfrozen noise parameters, which are simply ignored.

compare(othermodel, nodmx=True, convertcoordinates=True, threshold_sigma=3.0, unc_rat_threshold=1.05, verbosity='max', usecolor=True, format='text')[source]

Print comparison with another model

Parameters:
  • othermodel – TimingModel object to compare to

  • nodmx (bool, optional) – If True, don’t print the DMX parameters in the comparison

  • convertcoordinates (bool, optional) – Convert coordinates from ICRS<->ECL to make models consistent

  • threshold_sigma (float, optional) – Pulsar parameters for which diff_sigma > threshold will be printed with an exclamation point at the end of the line

  • unc_rat_threshold (float, optional) – Pulsar parameters for which the uncertainty has increased by a factor of unc_rat_threshold will be printed with an asterisk at the end of the line

  • verbosity (string, optional) –

    Dictates amount of information returned. Options include “max”, “med”, and “min”, which have the following results:

    ”max” - print all lines from both models whether they are fit or not (note that nodmx will override this); DEFAULT “med” - only print lines for parameters that are fit “min” - only print lines for fit parameters for which diff_sigma > threshold “check” - only print significant changes with logging.warning, not as string (note that all other modes will still print this)

  • usecolor (bool, optional) – Use colors on the output to complement use of “!” and “*”

  • format (string, optional) – One of “text” or “markdown”

Returns:

Human readable comparison, for printing. Formatted as a five column table with titles of PARAMETER NAME | Model1 | Model2 | Diff_Sigma1 | Diff_Sigma2 where ModelX refer to self and othermodel Timing Model objects, and Diff_SigmaX is the difference in a given parameter as reported by the two models, normalized by the uncertainty in model X. If model X has no reported uncertainty, nothing will be printed.

If format="text", when either Diff_SigmaX value is greater than threshold_sigma, an exclamation point (!) will be appended to the line and color will be added if usecolor=True. If the uncertainty in the first model if smaller than the second, an asterisk (*) will be appended to the line and color will be added if usecolor=True.

If format="markdown" then will be formatted as a markdown table with bold, colored, and highlighted text as appropriate.

For both output formats, warnings and info statements will be printed.

Return type:

str

Note

Prints logging warnings for parameters that have changed significantly and/or have increased in uncertainty.

Examples

To use this in a Jupyter notebook with and without markdown:

>>> from pint.models import get_model
>>> import pint.logging
>>> from IPython.display import display_markdown
>>> pint.logging.setup(level="WARNING")
>>> m1 = get_model(<file1>)
>>> m2 = get_model(<file2>)
>>> print(m1.compare(m2))
>>> display_markdown(m1.compare(m2, format="markdown"), raw=True)

Make sure to use raw=True to get the markdown output in a notebook.

use_aliases(reset_to_default=True, alias_translation=None)[source]

Set the parameters to use aliases as specified upon writing.

Parameters:
  • reset_to_default (bool) – If True, forget what name was used for each parameter in the input par file.

  • alias_translation (dict or None) – If not None, use this to map PINT parameter names to output names. This overrides input names even if they are not otherwise being reset to default. This is to allow compatibility with TEMPO/TEMPO2. The dictionary pint.toa.tempo_aliases should provide a reasonable selection.

as_parfile(start_order=['astrometry', 'spindown', 'dispersion'], last_order=['jump_delay'], *, include_info=True, comment=None, format='pint')[source]

Represent the entire model as a parfile string.

See also pint.models.TimingModel.write_parfile().

Parameters:
  • start_order (list) – Categories to include at the beginning

  • last_order (list) – Categories to include at the end

  • include_info (bool, optional) – Include information string if True

  • comment (str, optional) – Additional comment string to include in parfile

  • format (str, optional) – Parfile output format. PINT outputs in ‘tempo’, ‘tempo2’ and ‘pint’ formats. The defaul format is pint.

write_parfile(filename, start_order=['astrometry', 'spindown', 'dispersion'], last_order=['jump_delay'], *, include_info=True, comment=None, format='pint')[source]

Write the entire model as a parfile.

See also pint.models.TimingModel.as_parfile().

Parameters:
  • filename (str or Path or file-like) – The destination to write the parfile to

  • start_order (list) – Categories to include at the beginning

  • last_order (list) – Categories to include at the end

  • include_info (bool, optional) – Include information string if True

  • comment (str, optional) – Additional comment string to include in parfile

  • format (str, optional) – Parfile output format. PINT outputs in ‘tempo’, ‘tempo2’ and ‘pint’ formats. The defaul format is pint.

validate_toas(toas)[source]

Sanity check to verify that this model is compatible with these toas.

This checks that where this model needs TOAs to constrain parameters, that there is at least one TOA. This includes checking that every DMX range for with the DMX is free has at least one TOA, and it verifies that each “mask parameter” (for example JUMP) corresponds to at least one TOA.

Individual components can implement a validate_toas method; this method will automatically call such a method on each component that has one.

If some TOAs are missing, this method will raise a MissingTOAError that lists some (at least one) of the problem parameters.

find_empty_masks(toas, freeze=False)[source]

Find unfrozen mask parameters with no TOAs before trying to fit

Parameters:
  • toas (pint.toa.TOAs) –

  • freeze (bool, optional) – Should the parameters with on TOAs be frozen

Returns:

Parameters with no TOAs

Return type:

list

setup()[source]

Run setup methods on all components.

as_ECL(epoch=None, ecl='IERS2010')[source]

Return TimingModel in PulsarEcliptic frame.

Parameters:
Returns:

In PulsarEcliptic frame

Return type:

pint.models.timing_model.TimingModel

Notes

For the DDK model, the KOM vector is also transformed

as_ICRS(epoch=None)[source]

Return TimingModel in ICRS frame.

Parameters:

epoch (float or astropy.time.Time or astropy.units.Quantity, optional) – If float or Quantity, MJD(TDB) is assumed. New epoch for position.

Returns:

In ICRS frame

Return type:

pint.models.timing_model.TimingModel

Notes

For the DDK model, the KOM vector is also transformed

get_derived_params(rms=None, ntoas=None, returndict=False)[source]

Return a string with various derived parameters from the fitted model

Parameters:
  • rms (astropy.units.Quantity, optional) – RMS of fit for checking ELL1 validity

  • ntoas (int, optional) – Number of TOAs for checking ELL1 validity

  • returndict (bool, optional) – Whether to only return the string of results or also a dictionary

Returns:

  • results (str)

  • parameters (dict, optional)