pint.residuals.Residuals
- class pint.residuals.Residuals(toas: TOAs | None = None, model: TimingModel | None = None, residual_type: Literal['toa', 'dm'] = 'toa', unit: Unit | str = Unit('s'), subtract_mean: bool = True, use_weighted_mean: bool = True, track_mode: Literal['use_pulse_numbers', 'nearest'] | None = None, use_abs_phase: bool = True)[source]
Bases:
objectClass 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_residsand the.time_residsattributes.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:
phase_resids (
astropy.units.Quantity) – Residuals in phase unitstime_resids (
astropy.units.Quantity) – Residuals in time units
- Parameters:
toas (
pint.toa.TOAs, optional) – The input TOAs object. Default: Nonemodel (
pint.models.timing_model.TimingModel, optional) – Input model object. Default: Noneresidual_type (str, optional) – The type of the residuals. Default: ‘toa’
unit (
astropy.units.Unit, optional) – The default unit of the residuals. Default: u.ssubtract_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 thepulse_numbercolumn 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).
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.
Compute the log-likelihood for the model and TOAs.
Compute weighted RMS of the residuals in time.
update()Recalculate everything in residuals class after changing model or TOAs
Checks if the whitened residuals are unit-normal distributed using an Anderson-Darling test.
Checks if the whitened residuals are unit-normal distributed using a KS test.
Attributes
Compute chi-squared as needed and cache the result.
Reduced chi-squared.
Return number of degrees of freedom for the model.
noise_residsReturn the weighted reduced chi-squared for the model and toas.
Residuals in time units.
Residuals in seconds, with the units stripped.
- get_data_error(scaled: bool = True) Quantity[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.
- get_PSR_freq(calctype='modelF0') Quantity[source]
Return pulsar rotational frequency in Hz.
- Parameters:
calctype ({'modelF0', 'numerical', 'taylor'}) – Type of calculation. If calctype == “modelF0”, then simply the
F0parameter 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
F0in the model or the spin frequency at the moment of each TOA.- Return type:
- calc_phase_resids(subtract_mean: bool | None = None, use_weighted_mean: bool | None = None, use_abs_phase: bool | None = None) Quantity[source]
Compute timing model residuals in pulse phase.
if
subtract_meanoruse_weighted_meanis 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:
- calc_phase_mean(weighted: bool = True) Quantity[source]
Calculate mean phase of residuals, optionally weighted
- Parameters:
weighted (bool, optional) –
- Return type:
- calc_time_mean(calctype: Literal['taylor', 'modelF0', 'numerical'] = 'taylor', weighted: bool = True) Quantity[source]
Calculate mean time of residuals, optionally weighted
- Parameters:
calctype (str, optional) – Calculation time for phase to time conversion. See
pint.residuals.Residuals.calc_time_resids()for details.weighted (bool, optional) –
- Return type:
- calc_time_resids(calctype: Literal['taylor', 'modelF0', 'numerical'] = 'taylor', subtract_mean: bool | None = None, use_weighted_mean: bool | None = None, use_abs_phase: bool | None = None) Quantity[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_meanoruse_weighted_meanis None, will use the values set for the object itself- Parameters:
calctype ({'taylor', 'modelF0', 'numerical'}) – Type of calculation. If calctype == “modelF0”, then simply the
F0parameter from the model. If calctype == “numerical”, then try a numerical derivative If calctype == “taylor”, evaluate the frequency with a Taylor seriessubtract_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:
- calc_whitened_resids() Quantity[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()
- 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
See also
- 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
See also
Notes
Calculation done in
pint.utils.anderson_darling().
- calc_chi2(lognorm: bool = 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 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
- d_Ndiag_d_param(param: str) Quantity[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: bool = True) Quantity[source]
Uses the ECORR noise model time-binning to compute “epoch-averaged” residuals.
Requires ECORR be used in the timing model. If
use_noise_modelis 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