mirror of
https://github.com/fjosw/pyerrors.git
synced 2025-05-14 19:43:41 +02:00
refactor: computation of inverse cholesky composition for correlated
fits simplified. Threshold for ill conditioned correlation matrix warning raised to 1e8.
This commit is contained in:
parent
b14314b424
commit
b6eed60ec1
1 changed files with 5 additions and 9 deletions
|
@ -458,15 +458,11 @@ def _standard_fit(x, y, func, silent=False, **kwargs):
|
||||||
x0 = [0.1] * n_parms
|
x0 = [0.1] * n_parms
|
||||||
|
|
||||||
if kwargs.get('correlated_fit') is True:
|
if kwargs.get('correlated_fit') is True:
|
||||||
cov = covariance(y)
|
corr = covariance(y, correlation=True)
|
||||||
covdiag = np.diag(1. / np.sqrt(np.diag(cov)))
|
covdiag = np.diag(1 / np.asarray(dy_f))
|
||||||
corr = np.copy(cov)
|
|
||||||
for i in range(len(y)):
|
|
||||||
for j in range(len(y)):
|
|
||||||
corr[i][j] = cov[i][j] / np.sqrt(cov[i][i] * cov[j][j])
|
|
||||||
condn = np.linalg.cond(corr)
|
condn = np.linalg.cond(corr)
|
||||||
if condn > 1e4:
|
if condn > 1e8:
|
||||||
warnings.warn("Correlation matrix may be ill-conditioned! condition number: %1.2e" % (condn), RuntimeWarning)
|
warnings.warn("Correlation matrix may be ill-conditioned, condition number: %1.2e" % (condn), RuntimeWarning)
|
||||||
chol = np.linalg.cholesky(corr)
|
chol = np.linalg.cholesky(corr)
|
||||||
chol_inv = np.linalg.inv(chol)
|
chol_inv = np.linalg.inv(chol)
|
||||||
chol_inv = np.dot(chol_inv, covdiag)
|
chol_inv = np.dot(chol_inv, covdiag)
|
||||||
|
@ -487,7 +483,7 @@ def _standard_fit(x, y, func, silent=False, **kwargs):
|
||||||
|
|
||||||
if output.method != 'Levenberg-Marquardt':
|
if output.method != 'Levenberg-Marquardt':
|
||||||
if output.method == 'migrad':
|
if output.method == 'migrad':
|
||||||
fit_result = iminuit.minimize(chisqfunc, x0, tol=1e-4) # Stopping crieterion 0.002 * tol * errordef
|
fit_result = iminuit.minimize(chisqfunc, x0, tol=1e-4) # Stopping criterion 0.002 * tol * errordef
|
||||||
output.iterations = fit_result.nfev
|
output.iterations = fit_result.nfev
|
||||||
else:
|
else:
|
||||||
fit_result = scipy.optimize.minimize(chisqfunc, x0, method=kwargs.get('method'), tol=1e-12)
|
fit_result = scipy.optimize.minimize(chisqfunc, x0, method=kwargs.get('method'), tol=1e-12)
|
||||||
|
|
Loading…
Add table
Add a link
Reference in a new issue