pint.fitter.GLSFitter

class pint.fitter.GLSFitter(toas=None, model=None, track_mode=None, residuals=None)[source]

Bases: Fitter

Generalized least-squares fitting.

This fitter extends the pint.fitter.WLSFitter to permit the errors on the data points to be correlated. The fitting proceeds by decomposing the data covariance matrix.

Methods

auto(toas, model[, downhill, track_mode, ...])

Automatically return the proper pint.fitter.Fitter object depending on the TOAs and model.

fit_toas([maxiter, threshold, full_cov, debug])

Run a generalized least-squares fitting method.

ftest(parameter, component[, remove, ...])

Compare the significance of adding/removing parameters to a timing model.

get_allparams()

Return a dict of all param names and values.

get_derived_params([returndict])

Return a string with various derived parameters from the fitted model

get_designmatrix()

Return the model's design matrix for these TOAs.

get_fitparams()

Return a dict of fittable param names and quantity.

get_fitparams_num()

Return a dict of fittable param names and numeric values.

get_fitparams_uncertainty()

Return a dict of fittable param names and numeric values.

get_parameter_correlation_matrix([...])

Show the parameter correlation matrix post-fit.

get_parameter_covariance_matrix([...])

Show the parameter covariance matrix post-fit.

get_params_dict([which, kind])

Return a dict mapping parameter names to values.

get_summary([nodmx])

Return a human-readable summary of the Fitter results.

make_resids(model)

minimize_func(x, *args)

Wrapper function for the residual class.

plot()

Make residuals plot.

print_summary()

Write a summary of the TOAs to stdout.

reset_model()

Reset the current model to the initial model.

set_fitparams(*params)

Update the "frozen" attribute of model parameters.

set_param_uncertainties(fitp)

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

set_params(fitp)

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

update_model([chi2])

Update the model to reflect fit results and TOA properties.

update_resids()

Update the residuals.

Attributes

covariance_matrix

fit_toas(maxiter=1, threshold=0, full_cov=False, debug=False)[source]

Run a generalized least-squares fitting method.

A first attempt is made to solve the fitting problem by Cholesky decomposition, but if this fails singular value decomposition is used instead. In this case singular values below threshold are removed.

Parameters:
  • maxiter (int) – How many times to run the linear least-squares fit, re-evaluating the derivatives at each step. If maxiter is less than one, no fitting is done, just the chi-squared computation. In this case, you must provide the residuals argument when constructing the class. If maxiter is one or more, so fitting is actually done, the chi-squared value returned is only approximately the chi-squared of the improved(?) model. In fact it is the chi-squared of the solution to the linear fitting problem, and the full non-linear model should be evaluated and new residuals produced if an accurate chi-squared is desired.

  • threshold (float) – When to start discarding singular values. Typical values are about 1e-14 - a singular value smaller than this indicates parameters that are so degenerate the numerical precision cannot distinguish them. Such highly degenerate parameter sets are reported to the user with a DegeneracyWarning. Negative values force a faster but less stable Cholesky decomposition method.

  • full_cov (bool) – full_cov determines which calculation is used. If true, the full covariance matrix is constructed and the calculation is relatively straightforward but the full covariance matrix may be enormous. If false, an algorithm is used that takes advantage of the structure of the covariance matrix, based on information provided by the noise model. The two algorithms should give the same result to numerical accuracy where they both can be applied.

classmethod auto(toas, model, downhill=True, track_mode=None, residuals=None, **kwargs)

Automatically return the proper pint.fitter.Fitter object depending on the TOAs and model.

In general the downhill fitters are to be preferred. See https://github.com/nanograv/PINT/wiki/How-To#choose-a-fitter for the logic used.

Parameters:
  • toas (a pint TOAs instance) – The input toas.

  • model (a pint timing model instance) – The initial timing model for fitting.

  • downhill (bool, optional) – Whether or not to use the downhill fitter variant

  • track_mode (str, optional) – How to handle phase wrapping. This is used when creating pint.residuals.Residuals objects, and its meaning is defined there.

  • residuals (pint.residuals.Residuals) – Initial residuals. This argument exists to support an optimization, where GLSFitter is used to compute chi2 for appropriate Residuals objects.

Returns:

Returns appropriate subclass

Return type:

pint.fitter.Fitter

ftest(parameter, component, remove=False, full_output=False, maxiter=1)

Compare the significance of adding/removing parameters to a timing model.

Parameters:
  • parameter (PINT parameter object) – (may be a list of parameter objects)

  • component (String) – Name of component of timing model that the parameter should be added to (may be a list) The number of components must equal number of parameters.

  • remove (Bool) – If False, will add the listed parameters to the model. If True will remove the input parameters from the timing model.

  • full_output (Bool) – If False, just returns the result of the F-Test. If True, will also return the new model’s residual RMS (us), chi-squared, and number of degrees of freedom of new model.

  • maxiter (int) – How many times to run the linear least-squares fit, re-evaluating the derivatives at each step for the F-tested model. Default is one.

Returns:

ftFloat

F-test significance value for the model with the larger number of components over the other. Computed with pint.utils.FTest().

resid_rms_testFloat (Quantity)

If full_output is True, returns the RMS of the residuals of the tested model fit. Will be in units of microseconds as an astropy quantity. If wideband fitter this will be the time residuals.

resid_wrms_testFloat (Quantity)

If full_output is True, returns the Weighted RMS of the residuals of the tested model fit. Will be in units of microseconds as an astropy quantity. If wideband fitter this will be the time residuals.

chi2_testFloat

If full_output is True, returns the chi-squared of the tested model. If wideband fitter this will be the total chi-squared of the combined residual.

dof_testInt

If full_output is True, returns the degrees of freedom of the tested model. If wideband fitter this will be the total chi-squared of the combined residual.

dm_resid_rms_testFloat (Quantity)

If full_output is True and a wideband timing fitter is used, returns the RMS of the DM residuals of the tested model fit. Will be in units of pc/cm^3 as an astropy quantity.

dm_resid_wrms_testFloat (Quantity)

If full_output is True and a wideband timing fitter is used, returns the Weighted RMS of the DM residuals of the tested model fit. Will be in units of pc/cm^3 as an astropy quantity.

Return type:

dictionary

get_allparams()

Return a dict of all param names and values. Deprecated.

get_derived_params(returndict=False)

Return a string with various derived parameters from the fitted model

Parameters:

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

Returns:

  • results (str)

  • parameters (dict, optional)

get_designmatrix()

Return the model’s design matrix for these TOAs.

get_fitparams()

Return a dict of fittable param names and quantity. Deprecated.

get_fitparams_num()

Return a dict of fittable param names and numeric values. Deprecated.

get_fitparams_uncertainty()

Return a dict of fittable param names and numeric values. Deprecated.

get_parameter_correlation_matrix(with_phase=False, pretty_print=False, prec=3, usecolor=True)

Show the parameter correlation matrix post-fit.

If with_phase, then show and return the phase column as well. If pretty_print, then also pretty-print on stdout the matrix. prec is the precision of the floating point results. If usecolor is True, then pretty printing will have color.

get_parameter_covariance_matrix(with_phase=False, pretty_print=False, prec=3)

Show the parameter covariance matrix post-fit.

If with_phase, then show and return the phase column as well. If pretty_print, then also pretty-print on stdout the matrix. prec is the precision of the floating point results.

get_params_dict(which='free', kind='quantity')

Return a dict mapping parameter names to values.

See pint.models.timing_model.TimingModel.get_params_dict().

get_summary(nodmx=False)

Return a human-readable summary of the Fitter results.

Parameters:

nodmx (bool) – Set to True to suppress printing DMX parameters in summary

minimize_func(x, *args)

Wrapper function for the residual class.

This is 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).

plot()

Make residuals plot.

This produces a time residual plot.

print_summary()

Write a summary of the TOAs to stdout.

reset_model()

Reset the current model to the initial model.

set_fitparams(*params)

Update the “frozen” attribute of model parameters. Deprecated.

set_param_uncertainties(fitp)

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

See pint.models.timing_model.TimingModel.set_param_uncertainties().

set_params(fitp)

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

See pint.models.timing_model.TimingModel.set_param_values().

update_model(chi2=None)

Update the model to reflect fit results and TOA properties.

This is called by fit_toas to ensure that parameters like START, FINISH, EPHEM, and DMDATA are set in the model to reflect the TOAs in actual use.

update_resids()

Update the residuals.

Run after updating a model parameter.