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.

comb

n choose r combination function

corr

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

demo_hessian

demo_jacobian

demo_stderr_hilbert

demo_stderr_perturbed

gradient

hessian

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

hessian_cov

Given Hessian H, return the covariance matrix inv(H).

hilbert

Generate ill-conditioned Hilbert matrix of size n x n

hilbertinv

Analytical inverse for ill-conditioned Hilbert matrix of size n x n

jacobian

Returns the derivative wrt the fit parameters at point p.

jacobian_cov

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

max_correlation

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

perturbed_hessian

DEPRECATED Numerical testing has shown that the perturbed Hessian is too aggressive with its perturbation, and it is distorting the error too much, so use hessian_cov(H) instead.

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.

jacobian_cov() takes the Jacobian and computes the covariance matrix \(C\). hessian_cov() takes the Hessian 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}\)

Warning: assumes H = L@L.T (numpy default) not H = U.T@U (scipy default).

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.

Note that this calls chol_cov to compute the inverse from the Cholesky decomposition, so use stderr(C) if you are already computing C = chol_cov().

Warning: assumes H = L@L.T (numpy default) not H = U.T@U (scipy default).

bumps.lsqerror.comb(n, r)[source]

n choose r combination function

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.demo_hessian()[source]
bumps.lsqerror.demo_jacobian()[source]
bumps.lsqerror.demo_stderr_hilbert(n=5)[source]
bumps.lsqerror.demo_stderr_perturbed()[source]
bumps.lsqerror.gradient(problem, p=None, step=None)[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.hessian_cov(H, tol=1e-15)[source]

Given Hessian H, return the covariance matrix inv(H).

We provide some protection against singular matrices by setting singular values smaller than tolerance tol (relative to the largest singular value) to zero (see np.linalg.pinv for details).

bumps.lsqerror.hilbert(n)[source]

Generate ill-conditioned Hilbert matrix of size n x n

bumps.lsqerror.hilbertinv(n)[source]

Analytical inverse for ill-conditioned Hilbert matrix of size n x n

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.

Note that the problem.residuals() method should not reuse memory for the returned value otherwise the derivative calculation (f(x+dx) - f(x))/dx will always be zero. The returned value need not be 1D, but it should be contiguous so that it can be reshaped to 1D without an extra copy. This will only be an issue for very large datasets.

bumps.lsqerror.jacobian_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.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]

DEPRECATED Numerical testing has shown that the perturbed Hessian is too aggressive with its perturbation, and it is distorting the error too much, so use hessian_cov(H) instead.

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.