diff --git a/pyerrors/fits.py b/pyerrors/fits.py index 6db0b1dc..17da5e71 100644 --- a/pyerrors/fits.py +++ b/pyerrors/fits.py @@ -464,8 +464,10 @@ def _standard_fit(x, y, func, silent=False, **kwargs): corr = covariance(y, correlation=True) covdiag = np.diag(1 / np.asarray(dy_f)) condn = np.linalg.cond(corr) - if condn > 1e8: - warnings.warn("Correlation matrix may be ill-conditioned, condition number: %1.2e" % (condn), RuntimeWarning) + if condn > 0.1 / np.finfo(float).eps: + raise Exception(f"Cannot invert correlation matrix as its condition number exceeds machine precision ({condn:1.2e})") + if condn > 1 / np.sqrt(np.finfo(float).eps): + warnings.warn("Correlation matrix may be ill-conditioned, condition number: {%1.2e}" % (condn), RuntimeWarning) chol = np.linalg.cholesky(corr) chol_inv = np.linalg.inv(chol) chol_inv = np.dot(chol_inv, covdiag) diff --git a/tests/fits_test.py b/tests/fits_test.py index cc958dbc..b6236fb9 100644 --- a/tests/fits_test.py +++ b/tests/fits_test.py @@ -375,6 +375,12 @@ def test_fit_no_autograd(): pe.total_least_squares(oy, oy, func) +def test_singular_correlated_fit(): + obs1 = pe.pseudo_Obs(1.0, 0.1, 'test') + with pytest.raises(Exception): + pe.fits.fit_lin([0, 1], [obs1, obs1], correlated_fit=True) + + def test_ks_test(): def f(a, x): y = a[0] + a[1] * x