lsqerror - Least squares eorror analysis

chol_cov Given the cholesky decomposition of the Hessian matrix H, compute the covariance matrix \(C = H^{-1}\)
chol_stderr Return parameter uncertainty from the Cholesky decomposition of the Hessian matrix, as returned, e.g., from the quasi-Newton optimizer BFGS or as calculated from perturbed_hessian() on the output of hessian() applied to the cost function problem.nllf.
corr Convert covariance matrix \(C\) to correlation matrix \(R^2\).
cov Given Jacobian J, return the covariance matrix inv(J’J).
demo_hessian
demo_jacobian
hessian Returns the derivative wrt to the fit parameters at point p.
jacobian Returns the derivative wrt the fit parameters at point p.
max_correlation Return the maximum correlation coefficient for any pair of variables in correlation matrix Rsq.
perturbed_hessian Adjust Hessian matrix to be positive definite.
stderr Return parameter uncertainty from the covariance matrix C.

Least squares error analysis.

Given a data set with gaussian uncertainty on the points, and a model which is differentiable at the minimum, the parameter uncertainty can be estimated from the covariance matrix at the minimum. The model and data are wrapped in a problem object, which must define the following methods:

getp() get the current value of the model
setp(p) set a new value in the model
nllf(p) negative log likelihood function
residuals(p) residuals around its current value
bounds() get the bounds on the parameter p [optional]

jacobian() computes the Jacobian matrix \(J\) using numerical differentiation on residuals. Derivatives are computed using the center point formula, with two evaluations per dimension. If the problem has analytic derivatives with respect to the fitting parameters available, then these should be used to compute the Jacobian instead.

hessian() computes the Hessian matrix \(H\) using numerical differentiation on nllf. This uses the center point formula, with two evaluations for each (i,j) combination.

cov() takes the Jacobian and computes the covariance matrix \(C\).

corr() uses the off-diagonal elements of \(C\) to compute correlation coefficients \(R^2_{ij}\) between the parameters.

stderr() computes the uncertain \(\sigma_i\) from covariance matrix \(C\), assuming that the \(C_\text{diag}\) contains \(\sigma_i^2\), which should be the case for functions which are approximately linear near the minimum.

max_correlation() takes \(R^2\) and returns the maximum correlation.

The user should be shown the uncertainty \(\sigma_i\) for each parameter, and if there are strong parameter correlations (e.g., \(R^2_\text{max} > 0.2\)), the correlation matrix as well.

The bounds method for the problem is optional, and is used only to determine the step size needed for the numerical derivative. If bounds are not present and finite, the current value for the parameter is used as a basis to estimate step size.

bumps.lsqerror.chol_cov(L)[source]

Given the cholesky decomposition of the Hessian matrix H, compute the covariance matrix \(C = H^{-1}\)

bumps.lsqerror.chol_stderr(L)[source]

Return parameter uncertainty from the Cholesky decomposition of the Hessian matrix, as returned, e.g., from the quasi-Newton optimizer BFGS or as calculated from perturbed_hessian() on the output of hessian() applied to the cost function problem.nllf.

bumps.lsqerror.corr(C)[source]

Convert covariance matrix \(C\) to correlation matrix \(R^2\).

Uses \(R = D^{-1} C D^{-1}\) where \(D\) is the square root of the diagonal of the covariance matrix, or the standard error of each variable.

bumps.lsqerror.cov(J, tol=1e-08)[source]

Given Jacobian J, return the covariance matrix inv(J’J).

We provide some protection against singular matrices by setting singular values smaller than tolerance tol to the tolerance value.

bumps.lsqerror.demo_hessian()[source]
bumps.lsqerror.demo_jacobian()[source]
bumps.lsqerror.hessian(problem, p=None, step=None)[source]

Returns the derivative wrt to the fit parameters at point p.

The current point is preserved.

bumps.lsqerror.jacobian(problem, p=None, step=None)[source]

Returns the derivative wrt the fit parameters at point p.

Numeric derivatives are calculated based on step, where step is the portion of the total range for parameter j, or the portion of point value p_j if the range on parameter j is infinite.

The current point is preserved.

bumps.lsqerror.max_correlation(Rsq)[source]

Return the maximum correlation coefficient for any pair of variables in correlation matrix Rsq.

bumps.lsqerror.perturbed_hessian(H, scale=None)[source]

Adjust Hessian matrix to be positive definite.

Returns the adjusted Hessian and its Cholesky decomposition.

bumps.lsqerror.stderr(C)[source]

Return parameter uncertainty from the covariance matrix C.

This is just the square root of the diagonal, without any correction for covariance.

If measurement uncertainty is unknown, scale the returned uncertainties by \(\sqrt{\chi^2_N}\), where \(\chi^2_N\) is the sum squared residuals divided by the degrees of freedom. This will match the uncertainty on the parameters to the observed scatter assuming the model is correct and the fit is optimal. This will also be appropriate for weighted fits when the true measurement uncertainty dy_i is known up to a scaling constant for all y_i.

Standard error on scipy.optimize.curve_fit always includes the chisq correction, whereas scipy.optimize.leastsq never does.