pint.residuals.Residuals

class pint.residuals.Residuals(toas=None, model=None, residual_type='toa', unit=Unit('s'), subtract_mean=True, use_weighted_mean=True, track_mode=None, use_abs_phase=True)[source]

Bases: object

Class to compute residuals between TOAs and a TimingModel.

This class serves to store the results of a comparison between TOAs and a model. It also implements certain basic statistical calculations. This class also serves as a base class providing some infrastructure to support residuals from other kinds of data/model comparison.

This class provides access to the residuals in both phase (turns) and time (seconds) form through the .phase_resids and the .time_resids attributes.

Uncertainties on these residuals are available in time units using .get_data_error(); this can include or not include any rescaling of the uncertainties implied by the model’s EFAC or EQUAD.

Variables:
Parameters:
  • toas (pint.toa.TOAs, optional) – The input TOAs object. Default: None

  • model (pint.models.timing_model.TimingModel, optional) – Input model object. Default: None

  • residual_type (str, optional) – The type of the residuals. Default: ‘toa’

  • unit (astropy.units.Unit, optional) – The default unit of the residuals. Default: u.s

  • subtract_mean (bool) – Controls whether mean will be subtracted from the residuals. This option will be ignored if a PhaseOffset is present in the timing model.

  • use_weighted_mean (bool) – Controls whether mean computation is weighted (by errors) or not.

  • track_mode (None, "nearest", "use_pulse_numbers") – Controls how pulse numbers are assigned. "nearest" assigns each TOA to the nearest integer pulse. "use_pulse_numbers" uses the pulse_number column of the TOAs table to assign pulse numbers. If the default, None, is passed, use the pulse numbers if the model has the parameter TRACK == “-2” and not if it has TRACK == “0”. If neither of the above is set, use pulse numbers if there are pulse numbers present and not if there aren’t.

Methods

calc_chi2([lognorm])

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

calc_phase_mean([weighted])

Calculate mean phase of residuals, optionally weighted

calc_phase_resids([subtract_mean, ...])

Compute timing model residuals in pulse phase.

calc_time_mean([calctype, weighted])

Calculate mean time of residuals, optionally weighted

calc_time_resids([calctype, subtract_mean, ...])

Compute timing model residuals in time (seconds).

calc_whitened_resids()

Compute whitened timing residuals (dimensionless).

d_Ndiag_d_param(param)

Derivative of the white noise covariance matrix diagonal elements w.r.t.

d_lnlikelihood_d_ECORR(param)

d_lnlikelihood_d_Ndiag()

d_lnlikelihood_d_param(param)

ecorr_average([use_noise_model])

Uses the ECORR noise model time-binning to compute "epoch-averaged" residuals.

get_PSR_freq([calctype])

Return pulsar rotational frequency in Hz.

get_data_error([scaled])

Get errors on time residuals.

lnlikelihood()

Compute the log-likelihood for the model and TOAs.

rms_weighted()

Compute weighted RMS of the residuals in time.

update()

Recalculate everything in residuals class after changing model or TOAs

Attributes

chi2

Compute chi-squared as needed and cache the result.

chi2_reduced

Reduced chi-squared.

dof

Return number of degrees of freedom for the model.

reduced_chi2

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

resids

Residuals in time units.

resids_value

Residuals in seconds, with the units stripped.

property resids

Residuals in time units.

property resids_value

Residuals in seconds, with the units stripped.

update()[source]

Recalculate everything in residuals class after changing model or TOAs

property chi2

Compute chi-squared as needed and cache the result.

property dof

Return number of degrees of freedom for the model.

property reduced_chi2

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

property chi2_reduced

Reduced chi-squared.

get_data_error(scaled=True)[source]

Get errors on time residuals.

This returns the uncertainties on the time residuals, optionally scaled by the noise model.

Parameters:

scaled (bool, optional) – If errors get scaled by the noise model.

rms_weighted()[source]

Compute weighted RMS of the residuals in time.

get_PSR_freq(calctype='modelF0')[source]

Return pulsar rotational frequency in Hz.

Parameters:

calctype ({'modelF0', 'numerical', 'taylor'}) – Type of calculation. If calctype == “modelF0”, then simply the F0 parameter from the model. If calctype == “numerical”, then try a numerical derivative If calctype == “taylor”, evaluate the frequency with a Taylor series

Returns:

freq – Either the single F0 in the model or the spin frequency at the moment of each TOA.

Return type:

astropy.units.Quantity

calc_phase_resids(subtract_mean=None, use_weighted_mean=None, use_abs_phase=None)[source]

Compute timing model residuals in pulse phase.

if subtract_mean or use_weighted_mean is None, will use the values set for the object itself

Parameters:
  • subtract_mean (bool or None, optional) – Subtract the mean of the residuals. This is ignored if the PhaseOffset component is present in the model. Default is to use the class attribute.

  • use_weighted_mean (bool or None, optional) – Whether to use weighted mean for mean subtraction. Default is to use the class attribute.

  • use_abs_phase (bool or None, optional) – Whether to use absolute phase (w.r.t. the TZR TOA). Default is to use the class attribute.

Return type:

Phase

calc_phase_mean(weighted=True)[source]

Calculate mean phase of residuals, optionally weighted

Parameters:

weighted (bool, optional) –

Return type:

astropy.units.Quantity

calc_time_mean(calctype='taylor', weighted=True)[source]

Calculate mean time of residuals, optionally weighted

Parameters:
Return type:

astropy.units.Quantity

calc_time_resids(calctype='taylor', subtract_mean=None, use_weighted_mean=None, use_abs_phase=None)[source]

Compute timing model residuals in time (seconds).

Converts from phase residuals to time residuals using several possible ways to calculate the frequency.

If subtract_mean or use_weighted_mean is None, will use the values set for the object itself

Parameters:
  • calctype ({'taylor', 'modelF0', 'numerical'}) – Type of calculation. If calctype == “modelF0”, then simply the F0 parameter from the model. If calctype == “numerical”, then try a numerical derivative If calctype == “taylor”, evaluate the frequency with a Taylor series

  • subtract_mean (bool or None, optional) – Subtract the mean of the residuals. This is ignored if the PhaseOffset component is present in the model. Default is to use the class attribute.

  • use_weighted_mean (bool or None, optional) – Whether to use weighted mean for mean subtraction. Default is to use the class attribute.

  • use_abs_phase (bool or None, optional) – Whether to use absolute phase (w.r.t. the TZR TOA). Default is to use the class attribute.

Returns:

residuals

Return type:

astropy.units.Quantity

calc_whitened_resids()[source]

Compute whitened timing residuals (dimensionless).

Whitened residuals are computed by subtracting the correlated noise realization from the time residuals and normalizes the result using scaled TOA uncertainties.

This requires the noise_resids attribute to be set. This is usually available in a post-fit residuals object.

Example usage:

>>> ftr = Fitter.auto(toas, model)
>>> ftr.fit_toas()
>>> res = ftr.resids
>>> white_res = res.calc_whitened_resids()
calc_chi2(lognorm=False)[source]

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

If the errors on the TOAs are independent this is a straightforward calculation, but if the noise model introduces correlated errors then obtaining a meaningful chi-squared value requires a Cholesky decomposition. This is carried out using the :method:`~pint.residuals.Residuals._calc_gls_chi2` helper function.

The return value here is available as self.chi2, which will not redo the computation unless necessary.

The chi-squared value calculated here is suitable for use in downhill minimization algorithms and Bayesian approaches.

Handling of problematic results - degenerate conditions explored by a minimizer for example - may need to be checked to confirm that they correctly return infinity.

If lognorm=True is given, the log-normalization-factor of the likelihood function will also be returned.

Parameters:

lognorm (bool) – If True, return the the log-normalization-factor of the likelihood function along with the chi2 value.

Returns:

  • chi2 if lognorm is False

  • (chi2, log_norm) if lognorm is True

lnlikelihood()[source]

Compute the log-likelihood for the model and TOAs.

d_Ndiag_d_param(param)[source]

Derivative of the white noise covariance matrix diagonal elements w.r.t. a white noise parameter (EFAC or EQUAD).

ecorr_average(use_noise_model=True)[source]

Uses the ECORR noise model time-binning to compute “epoch-averaged” residuals.

Requires ECORR be used in the timing model. If use_noise_model is true, the noise model terms (EFAC, EQUAD, ECORR) will be applied to the TOA uncertainties, otherwise only the raw uncertainties will be used.

Returns a dictionary with the following entries:

mjds Average MJD for each segment

freqs Average topocentric frequency for each segment

time_resids Average residual for each segment, time units

noise_resids Dictionary of per-noise-component average residual

errors Uncertainty on averaged residuals

indices List of lists giving the indices of TOAs in the original TOA table for each segment