pint.residuals.WidebandTOAResiduals

class pint.residuals.WidebandTOAResiduals(toas: TOAs, model: TimingModel, toa_resid_args: dict = {}, dm_resid_args: dict = {})[source]

Bases: CombinedResiduals

A class for handling the wideband toa residuals.

Wideband TOAs have independent measurement of DM values. The residuals for wideband TOAs have two parts, the TOA residuals and DM residuals. Both residuals will be used for fitting one timing model. Currently, the DM values are stored at the TOA object.

The TOA and DM residuals are probably best accessed using the .toa and .dm properties.

This class inherits the .chi2 property from pint.residuals.CombinedResiduals.

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

  • model (pint.models.timing_model.TimingModel, optional) – The input timing model. Default: None

  • toa_resid_args (dict, optional) – The additional arguments(not including toas and model) for TOA residuals. Default: {}

  • dm_resid_args (dict, optional) – The additional arguments(not including toas and model) for DM residuals. Default: {}

Methods

calc_chi2([full_cov])

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

calc_whitened_dm_resids()

Compute the whitened wideband DM residuals.

calc_whitened_resids()

Compute the whitened wideband TOA residuals.

calc_wideband_resids()

Returns the combined TOA and DM residuals as a numpy array.

calc_wideband_whitened_resids()

Compute the whitened wideband residuals.

rms_weighted()

Compute weighted RMS of the residuals in time.

whitened_resids_adtest()

Checks if the whitened residuals are unit-normal distributed using an Anderson-Darling test.

whitened_resids_kstest()

Checks if the whitened residuals are unit-normal distributed using a KS test.

Attributes

chi2

Compute chi-squared as needed and cache the result.

data_error

dm

WidebandDMResiduals object containing the DM residuals.

dof

The number of degrees of freedom for the wideband residuals.

model

The model used to construct the residuals.

noise_resids

reduced_chi2

Return the weighted reduced chi-squared.

toa

Residuals object containing the TOA residuals.

unit

property toa: Residuals

Residuals object containing the TOA residuals.

property dm: WidebandDMResiduals

WidebandDMResiduals object containing the DM residuals.

property chi2: float

Compute chi-squared as needed and cache the result.

calc_chi2(full_cov=False) float[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, here, by constructing a GLSFitter and asking it to do the chi-squared computation but not a fit.

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.

property model: TimingModel

The model used to construct the residuals.

Modifying this model, even changing its parameters, may have confusing effects. It is probably best to use copy.deepcopy() to duplicate it before making any changes.

property dof: int

The number of degrees of freedom for the wideband residuals.

property reduced_chi2: float

Return the weighted reduced chi-squared.

calc_wideband_resids() ndarray[source]

Returns the combined TOA and DM residuals as a numpy array. The TOA residuals are in s and the DM residuals are in dmu.

Use pint.residuals.Residuals.calc_time_resids() and pint.residuals.WidebandDMResiduals.calc_dm_resids() to get time and DM residuals respectively as quantities.

rms_weighted() ndarray

Compute weighted RMS of the residuals in time.

calc_whitened_resids() ndarray[source]

Compute the whitened wideband TOA residuals. The whitened TOA residuals are dimensionless and should be unit-normal distributed if the timing model is correctly fit.

calc_whitened_dm_resids() ndarray[source]

Compute the whitened wideband DM residuals. The whitened DM residuals are dimensionless and should be unit-normal distributed if the timing model is correctly fit.

calc_wideband_whitened_resids() ndarray[source]

Compute the whitened wideband residuals. The first Ntoas elements are the whitened TOA residuals and the rest are DM residuals. The whitened residuals are dimensionless and should be unit-normal distributed if the timing model is correctly fit.

whitened_resids_kstest() Tuple[float, float][source]

Checks if the whitened residuals are unit-normal distributed using a KS test. A small p-value indicates a significant departure from unit-Gaussianity.

Returns:

  • float – A-D statistic

  • float – p-value

whitened_resids_adtest() Tuple[float, float][source]

Checks if the whitened residuals are unit-normal distributed using an Anderson-Darling test. A small p-value indicates a significant departure from unit-Gaussianity.

The Anderson-Darling test has more discriminatory power than the K-S test.

Returns:

  • float – A-D statistic

  • float – p-value

Notes

Calculation done in pint.utils.anderson_darling().